Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_fluid_scheme_incompressible.f90
1! Copyright (c) 2020-2024, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
35 use adjoint_fluid_scheme, only: adjoint_fluid_scheme_t
36 use gather_scatter, only: gs_t, gs_op_min, gs_op_max
37 use neko_config, only: neko_bcknd_device
38 use num_types, only: rp, i8
40 use field, only: field_t
41 use space, only: space_t, gll
42 use dofmap, only: dofmap_t
43 use krylov, only: ksp_t, krylov_solver_factory, ksp_max_iter
44 use coefs, only: coef_t
45 use jacobi, only: jacobi_t
46 use sx_jacobi, only: sx_jacobi_t
47 use device_jacobi, only: device_jacobi_t
48 use hsmg, only: hsmg_t
49 use phmg, only: phmg_t
50 use precon, only: pc_t, precon_factory, precon_destroy
51 use fluid_stats, only: fluid_stats_t
52 use bc, only: bc_t
53 use bc_list, only: bc_list_t
54 use mesh, only: mesh_t, neko_msh_max_zlbls, neko_msh_max_zlbl_len
55 use operators, only: cfl
56 use logger, only: neko_log, log_size, neko_log_verbose
57 use field_registry, only: neko_field_registry
58 use json_utils, only: json_get, json_get_or_default, json_extract_object
59 use json_module, only: json_file
60 use scratch_registry, only: scratch_registry_t
61 use user_intf, only: user_t, dummy_user_material_properties, &
62 user_material_properties_intf
63 use utils, only: neko_error
64 use time_state, only: time_state_t
65
66 use math, only: glsum
67 use device_math, only: device_cfill, device_add2s2
68 use field_math, only: field_cfill, field_add2s2, field_addcol3
69
70 use json_utils_ext, only: json_key_fallback
71 use device, only : device_event_sync, glb_cmd_event, device_to_host, &
72 device_memcpy
73 implicit none
74 private
75
77 type, abstract, extends(adjoint_fluid_scheme_t) :: &
80 type(adjoint_source_term_t) :: source_term
81 class(ksp_t), allocatable :: ksp_vel
82 class(ksp_t), allocatable :: ksp_prs
83 class(pc_t), allocatable :: pc_vel
84 class(pc_t), allocatable :: pc_prs
85 integer :: vel_projection_dim
86 integer :: pr_projection_dim
87 integer :: vel_projection_activ_step
88 integer :: pr_projection_activ_step
89 logical :: strict_convergence
90
92 type(field_t), pointer :: u_adj_e => null()
93 type(field_t), pointer :: v_adj_e => null()
94 type(field_t), pointer :: w_adj_e => null()
95
96 type(fluid_stats_t) :: stats
97 logical :: forced_flow_rate = .false.
98
100 character(len=:), allocatable :: nut_field_name
101
103 integer(kind=i8) :: glb_n_points
105 integer(kind=i8) :: glb_unique_points
106 type(scratch_registry_t) :: scratch
107 contains
109 procedure, pass(this) :: init_base => adjoint_fluid_scheme_init_base
110 procedure, pass(this) :: scheme_free => adjoint_fluid_scheme_free
112 procedure, pass(this) :: validate => adjoint_fluid_scheme_validate
114 procedure, pass(this) :: bc_apply_vel => adjoint_fluid_scheme_bc_apply_vel
116 procedure, pass(this) :: bc_apply_prs => adjoint_fluid_scheme_bc_apply_prs
118 procedure, pass(this) :: compute_cfl => adjoint_compute_cfl
120 procedure, pass(this) :: set_material_properties => &
121 adjoint_fluid_scheme_set_material_properties
122
124 procedure, pass(this) :: update_material_properties => &
125 adjoint_fluid_scheme_update_material_properties
127 procedure, nopass :: solver_factory => adjoint_fluid_scheme_solver_factory
129 procedure, pass(this) :: precon_factory_ => &
130 adjoint_fluid_scheme_precon_factory
132
133 interface
134
135 module subroutine adjoint_fluid_scheme_factory(object, type_name)
136 class(adjoint_fluid_scheme_t), intent(inout), allocatable :: object
137 character(len=*) :: type_name
138 end subroutine adjoint_fluid_scheme_factory
139 end interface
140
141 public :: adjoint_fluid_scheme_incompressible_t, adjoint_fluid_scheme_factory
142
143contains
144
146 subroutine adjoint_fluid_scheme_init_base(this, msh, lx, params, scheme, &
147 user, kspv_init)
148 class(adjoint_fluid_scheme_incompressible_t), target, intent(inout) :: this
149 type(mesh_t), target, intent(inout) :: msh
150 integer, intent(in) :: lx
151 character(len=*), intent(in) :: scheme
152 type(json_file), target, intent(inout) :: params
153 type(user_t), target, intent(in) :: user
154 logical, intent(in) :: kspv_init
155 character(len=LOG_SIZE) :: log_buf
156 real(kind=rp) :: real_val
157 logical :: logical_val
158 integer :: integer_val
159 character(len=:), allocatable :: string_val1, string_val2
160 character(len=:), allocatable :: json_key
161 type(json_file) :: json_subdict
162
163 !
164 ! SEM simulation fundamentals
165 !
166
167 this%msh => msh
168
169 if (msh%gdim .eq. 2) then
170 call this%Xh%init(gll, lx, lx)
171 else
172 call this%Xh%init(gll, lx, lx, lx)
173 end if
174
175 call this%dm_Xh%init(msh, this%Xh)
176
177 call this%gs_Xh%init(this%dm_Xh)
178
179 call this%c_Xh%init(this%gs_Xh)
180
181 ! Local scratch registry
182 call this%scratch%init(this%dm_Xh, 10, 2)
183
184 ! Assign a name
185 call json_get_or_default(params, 'case.fluid.name', this%name, "fluid")
186
187 !
188 ! First section of fluid log
189 !
190
191 call neko_log%section('Adjoint fluid')
192 write(log_buf, '(A, A)') 'Type : ', trim(scheme)
193 call neko_log%message(log_buf)
194 write(log_buf, '(A, A)') 'Name : ', trim(this%name)
195 call neko_log%message(log_buf)
196
197 ! Assign velocity fields
198 call neko_field_registry%add_field(this%dm_Xh, 'u_adj')
199 call neko_field_registry%add_field(this%dm_Xh, 'v_adj')
200 call neko_field_registry%add_field(this%dm_Xh, 'w_adj')
201 this%u_adj => neko_field_registry%get_field('u_adj')
202 this%v_adj => neko_field_registry%get_field('v_adj')
203 this%w_adj => neko_field_registry%get_field('w_adj')
204
205 !
206 ! Material properties
207 !
208 call this%set_material_properties(params, user)
209
210 ! Projection spaces
211 json_key = json_key_fallback(params, &
212 'case.adjoint_fluid.velocity_solver.projection_space_size', &
213 'case.fluid.velocity_solver.projection_space_size')
214 call json_get_or_default(params, json_key, this%vel_projection_dim, 0)
215 json_key = json_key_fallback(params, &
216 'case.adjoint_fluid.pressure_solver.projection_space_size', &
217 'case.fluid.pressure_solver.projection_space_size')
218 call json_get_or_default(params, json_key, this%pr_projection_dim, 0)
219
220 json_key = json_key_fallback(params, &
221 'case.adjoint_fluid.velocity_solver.projection_hold_steps', &
222 'case.fluid.velocity_solver.projection_hold_steps')
223 call json_get_or_default(params, json_key, &
224 this%vel_projection_activ_step, 5)
225 json_key = json_key_fallback(params, &
226 'case.adjoint_fluid.pressure_solver.projection_hold_steps', &
227 'case.fluid.pressure_solver.projection_hold_steps')
228 call json_get_or_default(params, json_key, &
229 this%pr_projection_activ_step, 5)
230
231 json_key = json_key_fallback(params, 'case.adjoint_fluid.freeze', &
232 'case.fluid.freeze')
233 call json_get_or_default(params, json_key, this%freeze, .false.)
234
235 ! TODO
236 ! This needs to be discussed... In pricipal, I think if the forward is
237 ! forced to a
238 ! certain flow rate, then the adjoint should be forced to zero flow rate,
239 ! we had a derivation in
240 ! https://www.sciencedirect.com/science/article/pii/S0045782522006764
241 ! for now I'm commenting this out
242 if (params%valid_path("case.fluid.flow_rate_force")) then
243 call neko_error('Flow rate forcing not yet implemented')
244 this%forced_flow_rate = .true.
245 end if
246
247
248 if (lx .lt. 10) then
249 write(log_buf, '(A, I1)') 'Poly order : ', lx-1
250 else if (lx .ge. 10) then
251 write(log_buf, '(A, I2)') 'Poly order : ', lx-1
252 else
253 write(log_buf, '(A, I3)') 'Poly order : ', lx-1
254 end if
255 call neko_log%message(log_buf)
256 this%glb_n_points = int(this%msh%glb_nelv, i8)*int(this%Xh%lxyz, i8)
257 this%glb_unique_points = int(glsum(this%c_Xh%mult, this%dm_Xh%size()), i8)
258
259 write(log_buf, '(A, I0)') 'GLL points : ', this%glb_n_points
260 call neko_log%message(log_buf)
261 write(log_buf, '(A, I0)') 'Unique pts.: ', this%glb_unique_points
262 call neko_log%message(log_buf)
263
264 call json_get(params, 'case.numerics.dealias', logical_val)
265 write(log_buf, '(A, L1)') 'Dealias : ', logical_val
266 call neko_log%message(log_buf)
267
268 call json_get_or_default(params, 'case.output_boundary', logical_val, &
269 .false.)
270 write(log_buf, '(A, L1)') 'Save bdry : ', logical_val
271 call neko_log%message(log_buf)
272
273 call json_get_or_default(params, "case.fluid.full_stress_formulation", &
274 logical_val, .false.)
275 write(log_buf, '(A, L1)') 'Full stress: ', logical_val
276 call neko_log%message(log_buf)
277
278 !
279 ! Setup right-hand side fields.
280 !
281 allocate(this%f_adj_x)
282 allocate(this%f_adj_y)
283 allocate(this%f_adj_z)
284 call this%f_adj_x%init(this%dm_Xh, fld_name = "adjoint_rhs_x")
285 call this%f_adj_y%init(this%dm_Xh, fld_name = "adjoint_rhs_y")
286 call this%f_adj_z%init(this%dm_Xh, fld_name = "adjoint_rhs_z")
287
288 ! Todo: HARRY
289 ! --------------------------------------------------------------------------
290 ! ok here is a chance that we can maybe steal the krylov solvers from the
291 ! forward??
292 !
293 ! something along the lines of
294 !
295 ! if (.not.params%valid_path('case.adjoint_fluid.velocity_solver')) then
296 ! this%ksp_vel => case%fluid%ksp_vel
297 ! this%pc_vel => case%fluid%ksp_vel
298 ! else
299 ! initialize everything
300 ! end if
301 !
302 ! but I don't know how we can steal the krylov solvers out of the forward,
303 ! so for now we'll just initialize two of them...
304 ! Initialize velocity solver
305 if (kspv_init) then
306 call neko_log%section("Adjoint Velocity solver")
307
308 json_key = json_key_fallback(params, &
309 'case.adjoint_fluid.velocity_solver.max_iterations', &
310 'case.fluid.velocity_solver.max_iterations')
311 call json_get_or_default(params, json_key, integer_val, ksp_max_iter)
312
313 json_key = json_key_fallback(params, &
314 'case.adjoint_fluid.velocity_solver.type', &
315 'case.fluid.velocity_solver.type')
316 call json_get(params, json_key, string_val1)
317
318 json_key = json_key_fallback(params, &
319 'case.adjoint_fluid.velocity_solver.preconditioner.type', &
320 'case.fluid.velocity_solver.preconditioner.type')
321 call json_get(params, json_key, string_val2)
322
323 json_key = json_key_fallback(params, &
324 'case.adjoint_fluid.velocity_solver.preconditioner', &
325 'case.fluid.velocity_solver.preconditioner')
326 call json_extract_object(params, json_key, json_subdict)
327
328 json_key = json_key_fallback(params, &
329 'case.adjoint_fluid.velocity_solver.absolute_tolerance', &
330 'case.fluid.velocity_solver.absolute_tolerance')
331 call json_get(params, json_key, real_val)
332
333 json_key = json_key_fallback(params, &
334 'case.adjoint_fluid.velocity_solver.monitor', &
335 'case.fluid.velocity_solver.monitor')
336 call json_get_or_default(params, json_key, logical_val, .false.)
337
338 call neko_log%message('Type : (' // trim(string_val1) // &
339 ', ' // trim(string_val2) // ')')
340
341 write(log_buf, '(A,ES13.6)') 'Abs tol :', real_val
342 call neko_log%message(log_buf)
343 call this%solver_factory(this%ksp_vel, this%dm_Xh%size(), &
344 string_val1, integer_val, real_val, logical_val)
345 call this%precon_factory_(this%pc_vel, this%ksp_vel, &
346 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs_vel, &
347 string_val2, json_subdict)
348 call neko_log%end_section()
349 end if
350
351 ! Strict convergence for the velocity solver
352 call json_get_or_default(params, 'case.fluid.strict_convergence', &
353 this%strict_convergence, .false.)
354
355 !! Initialize time-lag fields
356 call this%ulag%init(this%u_adj, 2)
357 call this%vlag%init(this%v_adj, 2)
358 call this%wlag%init(this%w_adj, 2)
359
360 call neko_field_registry%add_field(this%dm_Xh, 'u_adj_e')
361 call neko_field_registry%add_field(this%dm_Xh, 'v_adj_e')
362 call neko_field_registry%add_field(this%dm_Xh, 'w_adj_e')
363 this%u_adj_e => neko_field_registry%get_field('u_adj_e')
364 this%v_adj_e => neko_field_registry%get_field('v_adj_e')
365 this%w_adj_e => neko_field_registry%get_field('w_adj_e')
366
367 ! Initialize the source term
368 call this%source_term%init(this%f_adj_x, this%f_adj_y, this%f_adj_z, &
369 this%c_Xh, user, this%name)
370 call this%source_term%add(params, 'case.adjoint_fluid.source_term')
371
372 call neko_log%end_section()
373
374 end subroutine adjoint_fluid_scheme_init_base
375
377 subroutine adjoint_fluid_scheme_free(this)
378 class(adjoint_fluid_scheme_incompressible_t), intent(inout) :: this
379
380 call this%Xh%free()
381
382 if (allocated(this%ksp_vel)) then
383 call this%ksp_vel%free()
384 deallocate(this%ksp_vel)
385 end if
386
387 if (allocated(this%ksp_prs)) then
388 call this%ksp_prs%free()
389 deallocate(this%ksp_prs)
390 end if
391
392 if (allocated(this%pc_vel)) then
393 call precon_destroy(this%pc_vel)
394 deallocate(this%pc_vel)
395 end if
396
397 if (allocated(this%pc_prs)) then
398 call precon_destroy(this%pc_prs)
399 deallocate(this%pc_prs)
400 end if
401
402 call this%source_term%free()
403
404 call this%gs_Xh%free()
405
406 call this%c_Xh%free()
407
408 call this%scratch%free()
409
410 nullify(this%u_adj)
411 nullify(this%v_adj)
412 nullify(this%w_adj)
413 nullify(this%p_adj)
414
415 nullify(this%u_adj_e)
416 nullify(this%v_adj_e)
417 nullify(this%w_adj_e)
418
419 call this%ulag%free()
420 call this%vlag%free()
421 call this%wlag%free()
422
423
424 if (associated(this%f_adj_x)) then
425 call this%f_adj_x%free()
426 end if
427
428 if (associated(this%f_adj_y)) then
429 call this%f_adj_y%free()
430 end if
431
432 if (associated(this%f_adj_z)) then
433 call this%f_adj_z%free()
434 end if
435
436 nullify(this%f_adj_x)
437 nullify(this%f_adj_y)
438 nullify(this%f_adj_z)
439
440 call this%rho%free()
441 call this%mu%free()
442
443 end subroutine adjoint_fluid_scheme_free
444
447 subroutine adjoint_fluid_scheme_validate(this)
448 class(adjoint_fluid_scheme_incompressible_t), target, intent(inout) :: this
449 ! Variables for retrieving json parameters
450
451 if ((.not. associated(this%u_adj)) .or. &
452 (.not. associated(this%v_adj)) .or. &
453 (.not. associated(this%w_adj)) .or. &
454 (.not. associated(this%p_adj))) then
455 call neko_error('Fields are not registered')
456 end if
457
458 if ((.not. allocated(this%u_adj%x)) .or. &
459 (.not. allocated(this%v_adj%x)) .or. &
460 (.not. allocated(this%w_adj%x)) .or. &
461 (.not. allocated(this%p_adj%x))) then
462 call neko_error('Fields are not allocated')
463 end if
464
465 if (.not. allocated(this%ksp_vel)) then
466 call neko_error('No Krylov solver for velocity defined')
467 end if
468
469 if (.not. allocated(this%ksp_prs)) then
470 call neko_error('No Krylov solver for pressure defined')
471 end if
472
473 !
474 ! Setup checkpoint structure (if everything is fine)
475 !
476 call this%chkp%init(this%u_adj, this%v_adj, this%w_adj, this%p_adj)
477
478 end subroutine adjoint_fluid_scheme_validate
479
484 subroutine adjoint_fluid_scheme_bc_apply_vel(this, time, strong)
485 class(adjoint_fluid_scheme_incompressible_t), intent(inout) :: this
486 type(time_state_t), intent(in) :: time
487 logical, intent(in) :: strong
488
489 integer :: i
490 class(bc_t), pointer :: b => null()
491
492 call this%bcs_vel%apply_vector(&
493 this%u_adj%x, this%v_adj%x, this%w_adj%x, this%dm_Xh%size(), &
494 time, strong)
495 call this%gs_Xh%op(this%u_adj, gs_op_min, glb_cmd_event)
496 call this%gs_Xh%op(this%v_adj, gs_op_min, glb_cmd_event)
497 call this%gs_Xh%op(this%w_adj, gs_op_min, glb_cmd_event)
498 call device_event_sync(glb_cmd_event)
499
500 call this%bcs_vel%apply_vector(&
501 this%u_adj%x, this%v_adj%x, this%w_adj%x, this%dm_Xh%size(), &
502 time, strong)
503 call this%gs_Xh%op(this%u_adj, gs_op_max, glb_cmd_event)
504 call this%gs_Xh%op(this%v_adj, gs_op_max, glb_cmd_event)
505 call this%gs_Xh%op(this%w_adj, gs_op_max, glb_cmd_event)
506 call device_event_sync(glb_cmd_event)
507
508 do i = 1, this%bcs_vel%size()
509 b => this%bcs_vel%get(i)
510 b%updated = .false.
511 end do
512 nullify(b)
513
514
515 end subroutine adjoint_fluid_scheme_bc_apply_vel
516
519 subroutine adjoint_fluid_scheme_bc_apply_prs(this, time)
520 class(adjoint_fluid_scheme_incompressible_t), intent(inout) :: this
521 type(time_state_t), intent(in) :: time
522
523 integer :: i
524 class(bc_t), pointer :: b => null()
525
526 call this%bcs_prs%apply(this%p_adj, time)
527 call this%gs_Xh%op(this%p_adj, gs_op_min, glb_cmd_event)
528 call device_event_sync(glb_cmd_event)
529
530 call this%bcs_prs%apply(this%p_adj, time)
531 call this%gs_Xh%op(this%p_adj, gs_op_max, glb_cmd_event)
532 call device_event_sync(glb_cmd_event)
533
534 do i = 1, this%bcs_prs%size()
535 b => this%bcs_prs%get(i)
536 b%updated = .false.
537 end do
538 nullify(b)
539
540 end subroutine adjoint_fluid_scheme_bc_apply_prs
541
544 subroutine adjoint_fluid_scheme_solver_factory(ksp, n, solver, &
545 max_iter, abstol, monitor)
546 class(ksp_t), allocatable, target, intent(inout) :: ksp
547 integer, intent(in), value :: n
548 character(len=*), intent(in) :: solver
549 integer, intent(in) :: max_iter
550 real(kind=rp), intent(in) :: abstol
551 logical, intent(in) :: monitor
552
553 call krylov_solver_factory(ksp, n, solver, max_iter, abstol, &
554 monitor = monitor)
555
556 end subroutine adjoint_fluid_scheme_solver_factory
557
559 subroutine adjoint_fluid_scheme_precon_factory(this, pc, ksp, coef, dof, gs, &
560 bclst, pctype, pcparams)
561 class(adjoint_fluid_scheme_incompressible_t), intent(inout) :: this
562 class(pc_t), allocatable, target, intent(inout) :: pc
563 class(ksp_t), target, intent(inout) :: ksp
564 type(coef_t), target, intent(in) :: coef
565 type(dofmap_t), target, intent(in) :: dof
566 type(gs_t), target, intent(inout) :: gs
567 type(bc_list_t), target, intent(inout) :: bclst
568 character(len=*) :: pctype
569 type(json_file), intent(inout) :: pcparams
570
571 call precon_factory(pc, pctype)
572
573 select type (pcp => pc)
574 type is (jacobi_t)
575 call pcp%init(coef, dof, gs)
576 type is (sx_jacobi_t)
577 call pcp%init(coef, dof, gs)
578 type is (device_jacobi_t)
579 call pcp%init(coef, dof, gs)
580 type is (hsmg_t)
581 call pcp%init(coef, bclst, pcparams)
582 type is (phmg_t)
583 call pcp%init(coef, bclst, pcparams)
584 end select
585
586 call ksp%set_pc(pc)
587
588 end subroutine adjoint_fluid_scheme_precon_factory
589
591 ! TODO
592 ! HARRY
593 ! ok this needs to be discussed.
594 ! To me, the CFL is based on the convecting velocity, so I would argue
595 ! that the CFL of the adjoint is based on
596 ! u, v, w
597 ! not
598 ! u_adj, v_adj, w_adj
599 !
600 ! Then again, this function is really only used for varaible timestep,
601 ! which we can't use if we're checkpointing...
602 !
603 ! hmmmmm.... requires more thought. Because people optimal ICs or Arnouldi
604 ! may use variable timesteps... or us in a steady case..
605 !
606 ! for now.... let's ignore it
607 function adjoint_compute_cfl(this, dt) result(c)
608 class(adjoint_fluid_scheme_incompressible_t), intent(in) :: this
609 real(kind=rp), intent(in) :: dt
610 real(kind=rp) :: c
611
612 c = cfl(dt, this%u_adj%x, this%v_adj%x, this%w_adj%x, &
613 this%Xh, this%c_Xh, this%msh%nelv, this%msh%gdim)
614
615 end function adjoint_compute_cfl
616
621 subroutine adjoint_fluid_scheme_update_material_properties(this, time)
622 class(adjoint_fluid_scheme_incompressible_t), intent(inout) :: this
623 type(time_state_t), intent(in) :: time
624 type(field_t), pointer :: nut
625
626 call this%user_material_properties(this%name, &
627 this%material_properties, time)
628
629 if (len(trim(this%nut_field_name)) > 0) then
630 nut => neko_field_registry%get_field(this%nut_field_name)
631 call field_addcol3(this%mu, nut, this%rho)
632 end if
633
634 ! Since mu, rho is a field_t, and we use the %x(1,1,1,1)
635 ! host array data to pass constant density and viscosity
636 ! to some routines, we need to make sure that the host
637 ! values are also filled
638 if (neko_bcknd_device .eq. 1) then
639 call device_memcpy(this%rho%x, this%rho%x_d, this%rho%size(), &
640 device_to_host, sync=.false.)
641 call device_memcpy(this%mu%x, this%mu%x_d, this%mu%size(), &
642 device_to_host, sync=.false.)
643 end if
644 end subroutine adjoint_fluid_scheme_update_material_properties
645
650 subroutine adjoint_fluid_scheme_set_material_properties(this, params, user)
651 class(adjoint_fluid_scheme_incompressible_t), intent(inout) :: this
652 type(json_file), intent(inout) :: params
653 type(user_t), target, intent(in) :: user
654 character(len=LOG_SIZE) :: log_buf
655 ! A local pointer that is needed to make Intel happy
656 procedure(user_material_properties_intf), pointer :: dummy_mp_ptr
657 real(kind=rp) :: const_mu, const_rho
658 type(time_state_t) :: time
659
660
661 dummy_mp_ptr => dummy_user_material_properties
662
663 call this%mu%init(this%dm_Xh, "mu")
664 call this%rho%init(this%dm_Xh, "rho")
665 call this%material_properties%init(2)
666 call this%material_properties%assign_to_field(1, this%rho)
667 call this%material_properties%assign_to_field(2, this%mu)
668
669 if (.not. associated(user%material_properties, dummy_mp_ptr)) then
670
671 write(log_buf, '(A)') "Material properties must be set in the user&
672 & file!"
673 call neko_log%message(log_buf)
674 this%user_material_properties => user%material_properties
675
676 call user%material_properties(this%name, &
677 this%material_properties, time)
678
679 else
680 this%user_material_properties => dummy_user_material_properties
681 ! Incorrect user input
682 if (params%valid_path('case.fluid.Re') .and. &
683 (params%valid_path('case.fluid.mu') .or. &
684 params%valid_path('case.fluid.rho'))) then
685 call neko_error("To set the material properties for the fluid, " // &
686 "either provide Re OR mu and rho in the case file.")
687
688 else if (params%valid_path('case.fluid.Re')) then
689 ! Non-dimensional case
690 write(log_buf, '(A)') 'Non-dimensional fluid material properties &
691 & input.'
692 call neko_log%message(log_buf, lvl = neko_log_verbose)
693 write(log_buf, '(A)') 'Density will be set to 1, dynamic viscosity to&
694 & 1/Re.'
695 call neko_log%message(log_buf, lvl = neko_log_verbose)
696
697 ! Read Re into mu for further manipulation.
698 call json_get(params, 'case.fluid.Re', const_mu)
699 write(log_buf, '(A)') 'Read non-dimensional material properties'
700 call neko_log%message(log_buf)
701 write(log_buf, '(A,ES13.6)') 'Re :', const_mu
702 call neko_log%message(log_buf)
703
704 ! Set rho to 1 since the setup is non-dimensional.
705 const_rho = 1.0_rp
706 ! Invert the Re to get viscosity.
707 const_mu = 1.0_rp/const_mu
708 else
709 ! Dimensional case
710 call json_get(params, 'case.fluid.mu', const_mu)
711 call json_get(params, 'case.fluid.rho', const_rho)
712 end if
713 end if
714
715 ! We need to fill the fields based on the parsed const values
716 ! if the user routine is not used.
717 if (associated(user%material_properties, dummy_mp_ptr)) then
718 ! Fill mu and rho field with the physical value
719 call field_cfill(this%mu, const_mu)
720 call field_cfill(this%rho, const_rho)
721
722
723 write(log_buf, '(A,ES13.6)') 'rho :', const_rho
724 call neko_log%message(log_buf)
725 write(log_buf, '(A,ES13.6)') 'mu :', const_mu
726 call neko_log%message(log_buf)
727 end if
728
729 ! Since mu, rho is a field_t, and we use the %x(1,1,1,1)
730 ! host array data to pass constant density and viscosity
731 ! to some routines, we need to make sure that the host
732 ! values are also filled
733 if (neko_bcknd_device .eq. 1) then
734 call device_memcpy(this%rho%x, this%rho%x_d, this%rho%size(), &
735 device_to_host, sync=.false.)
736 call device_memcpy(this%mu%x, this%mu%x_d, this%mu%size(), &
737 device_to_host, sync=.false.)
738 end if
739 end subroutine adjoint_fluid_scheme_set_material_properties
740
Implements the adjoint_source_term_t type.
Wrapper contaning and executing the adjoint source terms.