Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_fluid_scheme_incompressible.f90
Go to the documentation of this file.
1
34!
37 use adjoint_fluid_scheme, only: adjoint_fluid_scheme_t
38 use gather_scatter, only: gs_t, gs_op_min, gs_op_max
39 use neko_config, only: neko_bcknd_device
40 use num_types, only: rp, i8
42 use field, only: field_t
43 use space, only: space_t, gll, gl
44 use dofmap, only: dofmap_t
45 use krylov, only: ksp_t, krylov_solver_factory, ksp_max_iter
46 use coefs, only: coef_t
47 use jacobi, only: jacobi_t
48 use sx_jacobi, only: sx_jacobi_t
49 use device_jacobi, only: device_jacobi_t
50 use hsmg, only: hsmg_t
51 use phmg, only: phmg_t
52 use precon, only: pc_t, precon_factory, precon_destroy
53 use fluid_stats, only: fluid_stats_t
54 use bc, only: bc_t
55 use bc_list, only: bc_list_t
56 use mesh, only: mesh_t, neko_msh_max_zlbls, neko_msh_max_zlbl_len
57 use operators, only: cfl
58 use logger, only: neko_log, log_size, neko_log_verbose
59 use registry, only: neko_registry
60 use json_utils, only: json_get, json_get_or_default
61 use json_module, only: json_file
62 use user_intf, only: user_t, dummy_user_material_properties, &
63 user_material_properties_intf
64 use utils, only: neko_error
65 use time_state, only: time_state_t
66
67 use math, only: glsum
68 use device_math, only: device_cfill, device_add2s2
69 use field_math, only: field_cfill, field_add2s2, field_addcol3
70
71 use json_utils_ext, only: json_key_fallback
72 use device, only : device_event_sync, glb_cmd_event, device_to_host, &
73 device_memcpy
74 implicit none
75 private
76
78 type, abstract, extends(adjoint_fluid_scheme_t) :: &
81 type(adjoint_source_term_t) :: source_term
82 class(ksp_t), allocatable :: ksp_vel
83 class(ksp_t), allocatable :: ksp_prs
84 class(pc_t), allocatable :: pc_vel
85 class(pc_t), allocatable :: pc_prs
86 integer :: vel_projection_dim
87 integer :: pr_projection_dim
88 integer :: vel_projection_activ_step
89 integer :: pr_projection_activ_step
90 logical :: strict_convergence
91
93 type(field_t), pointer :: u_adj_e => null()
94 type(field_t), pointer :: v_adj_e => null()
95 type(field_t), pointer :: w_adj_e => null()
96
97 type(fluid_stats_t) :: stats
98 logical :: forced_flow_rate = .false.
99
101 character(len=:), allocatable :: nut_field_name
102
104 integer(kind=i8) :: glb_n_points
106 integer(kind=i8) :: glb_unique_points
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 integer :: lxd
163
164 !
165 ! SEM simulation fundamentals
166 !
167
168 this%msh => msh
169
170 ! over intergration order (hard coded now, should be optional)
171 lxd = (3 * (lx + 1)) / 2
172 if (msh%gdim .eq. 2) then
173 call this%Xh%init(gll, lx, lx)
174 else
175 call this%Xh%init(gll, lx, lx, lx)
176 end if
177 call this%Xh_GL%init(gl, lxd, lxd, lxd)
178
179 ! NOTE. This shouldn't require remaking all this stuff. What should be
180 ! changed on the neko side is a way of initializing a coef FULLY (ie, Bs)
181 ! with just a space.
182 call this%dm_Xh%init(msh, this%Xh)
183 call this%dm_Xh_GL%init(msh, this%Xh_GL)
184
185 call this%gs_Xh%init(this%dm_Xh)
186 call this%gs_Xh_GL%init(this%dm_Xh_GL)
187
188 call this%c_Xh%init(this%gs_Xh)
189 call this%c_Xh_GL%init(this%gs_Xh_GL)
190
191 call this%GLL_to_GL%init(this%Xh_GL, this%Xh)
192
193 ! Overintegration scratch registry (5 should be sufficient)
194 call this%scratch_GL%init(5, 2, this%dm_Xh_GL)
195
196 ! Assign a name
197 call json_get_or_default(params, 'case.fluid.name', this%name, "fluid")
198
199 !
200 ! First section of fluid log
201 !
202
203 call neko_log%section('Adjoint fluid')
204 write(log_buf, '(A, A)') 'Type : ', trim(scheme)
205 call neko_log%message(log_buf)
206 write(log_buf, '(A, A)') 'Name : ', trim(this%name)
207 call neko_log%message(log_buf)
208
209 ! Assign velocity fields
210 call neko_registry%add_field(this%dm_Xh, 'u_adj')
211 call neko_registry%add_field(this%dm_Xh, 'v_adj')
212 call neko_registry%add_field(this%dm_Xh, 'w_adj')
213 this%u_adj => neko_registry%get_field('u_adj')
214 this%v_adj => neko_registry%get_field('v_adj')
215 this%w_adj => neko_registry%get_field('w_adj')
216
217 !
218 ! Material properties
219 !
220 call this%set_material_properties(params, user)
221
222 ! Projection spaces
223 json_key = json_key_fallback(params, &
224 'case.adjoint_fluid.velocity_solver.projection_space_size', &
225 'case.fluid.velocity_solver.projection_space_size')
226 call json_get_or_default(params, json_key, this%vel_projection_dim, 0)
227 json_key = json_key_fallback(params, &
228 'case.adjoint_fluid.pressure_solver.projection_space_size', &
229 'case.fluid.pressure_solver.projection_space_size')
230 call json_get_or_default(params, json_key, this%pr_projection_dim, 0)
231
232 json_key = json_key_fallback(params, &
233 'case.adjoint_fluid.velocity_solver.projection_hold_steps', &
234 'case.fluid.velocity_solver.projection_hold_steps')
235 call json_get_or_default(params, json_key, &
236 this%vel_projection_activ_step, 5)
237 json_key = json_key_fallback(params, &
238 'case.adjoint_fluid.pressure_solver.projection_hold_steps', &
239 'case.fluid.pressure_solver.projection_hold_steps')
240 call json_get_or_default(params, json_key, &
241 this%pr_projection_activ_step, 5)
242
243 json_key = json_key_fallback(params, 'case.adjoint_fluid.freeze', &
244 'case.fluid.freeze')
245 call json_get_or_default(params, json_key, this%freeze, .false.)
246
247 ! TODO
248 ! This needs to be discussed... In pricipal, I think if the forward is
249 ! forced to a
250 ! certain flow rate, then the adjoint should be forced to zero flow rate,
251 ! we had a derivation in
252 ! https://www.sciencedirect.com/science/article/pii/S0045782522006764
253 ! for now I'm commenting this out
254 if (params%valid_path("case.fluid.flow_rate_force")) then
255 call neko_error('Flow rate forcing not yet implemented')
256 this%forced_flow_rate = .true.
257 end if
258
259
260 if (lx .lt. 10) then
261 write(log_buf, '(A, I1)') 'Poly order : ', lx-1
262 else if (lx .ge. 10) then
263 write(log_buf, '(A, I2)') 'Poly order : ', lx-1
264 else
265 write(log_buf, '(A, I3)') 'Poly order : ', lx-1
266 end if
267 call neko_log%message(log_buf)
268 this%glb_n_points = int(this%msh%glb_nelv, i8)*int(this%Xh%lxyz, i8)
269 this%glb_unique_points = int(glsum(this%c_Xh%mult, this%dm_Xh%size()), i8)
270
271 write(log_buf, '(A, I0)') 'GLL points : ', this%glb_n_points
272 call neko_log%message(log_buf)
273 write(log_buf, '(A, I0)') 'Unique pts.: ', this%glb_unique_points
274 call neko_log%message(log_buf)
275
276 call json_get(params, 'case.numerics.dealias', logical_val)
277 write(log_buf, '(A, L1)') 'Dealias : ', logical_val
278 call neko_log%message(log_buf)
279
280 call json_get_or_default(params, 'case.output_boundary', logical_val, &
281 .false.)
282 write(log_buf, '(A, L1)') 'Save bdry : ', logical_val
283 call neko_log%message(log_buf)
284
285 call json_get_or_default(params, "case.fluid.full_stress_formulation", &
286 logical_val, .false.)
287 write(log_buf, '(A, L1)') 'Full stress: ', logical_val
288 call neko_log%message(log_buf)
289
290 !
291 ! Setup right-hand side fields.
292 !
293 allocate(this%f_adj_x)
294 allocate(this%f_adj_y)
295 allocate(this%f_adj_z)
296 call this%f_adj_x%init(this%dm_Xh, fld_name = "adjoint_rhs_x")
297 call this%f_adj_y%init(this%dm_Xh, fld_name = "adjoint_rhs_y")
298 call this%f_adj_z%init(this%dm_Xh, fld_name = "adjoint_rhs_z")
299
300 ! Todo: HARRY
301 ! --------------------------------------------------------------------------
302 ! ok here is a chance that we can maybe steal the krylov solvers from the
303 ! forward??
304 !
305 ! something along the lines of
306 !
307 ! if (.not.params%valid_path('case.adjoint_fluid.velocity_solver')) then
308 ! this%ksp_vel => case%fluid%ksp_vel
309 ! this%pc_vel => case%fluid%ksp_vel
310 ! else
311 ! initialize everything
312 ! end if
313 !
314 ! but I don't know how we can steal the krylov solvers out of the forward,
315 ! so for now we'll just initialize two of them...
316 ! Initialize velocity solver
317 if (kspv_init) then
318 call neko_log%section("Adjoint Velocity solver")
319
320 json_key = json_key_fallback(params, &
321 'case.adjoint_fluid.velocity_solver.max_iterations', &
322 'case.fluid.velocity_solver.max_iterations')
323 call json_get_or_default(params, json_key, integer_val, ksp_max_iter)
324
325 json_key = json_key_fallback(params, &
326 'case.adjoint_fluid.velocity_solver.type', &
327 'case.fluid.velocity_solver.type')
328 call json_get(params, json_key, string_val1)
329
330 json_key = json_key_fallback(params, &
331 'case.adjoint_fluid.velocity_solver.preconditioner.type', &
332 'case.fluid.velocity_solver.preconditioner.type')
333 call json_get(params, json_key, string_val2)
334
335 json_key = json_key_fallback(params, &
336 'case.adjoint_fluid.velocity_solver.preconditioner', &
337 'case.fluid.velocity_solver.preconditioner')
338 call json_get(params, json_key, json_subdict)
339
340 json_key = json_key_fallback(params, &
341 'case.adjoint_fluid.velocity_solver.absolute_tolerance', &
342 'case.fluid.velocity_solver.absolute_tolerance')
343 call json_get(params, json_key, real_val)
344
345 json_key = json_key_fallback(params, &
346 'case.adjoint_fluid.velocity_solver.monitor', &
347 'case.fluid.velocity_solver.monitor')
348 call json_get_or_default(params, json_key, logical_val, .false.)
349
350 call neko_log%message('Type : (' // trim(string_val1) // &
351 ', ' // trim(string_val2) // ')')
352
353 write(log_buf, '(A,ES13.6)') 'Abs tol :', real_val
354 call neko_log%message(log_buf)
355 call this%solver_factory(this%ksp_vel, this%dm_Xh%size(), &
356 string_val1, integer_val, real_val, logical_val)
357 call this%precon_factory_(this%pc_vel, this%ksp_vel, &
358 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs_vel, &
359 string_val2, json_subdict)
360 call neko_log%end_section()
361 end if
362
363 ! Strict convergence for the velocity solver
364 call json_get_or_default(params, 'case.fluid.strict_convergence', &
365 this%strict_convergence, .false.)
366
367 !! Initialize time-lag fields
368 call this%ulag%init(this%u_adj, 2)
369 call this%vlag%init(this%v_adj, 2)
370 call this%wlag%init(this%w_adj, 2)
371
372 call neko_registry%add_field(this%dm_Xh, 'u_adj_e')
373 call neko_registry%add_field(this%dm_Xh, 'v_adj_e')
374 call neko_registry%add_field(this%dm_Xh, 'w_adj_e')
375 this%u_adj_e => neko_registry%get_field('u_adj_e')
376 this%v_adj_e => neko_registry%get_field('v_adj_e')
377 this%w_adj_e => neko_registry%get_field('w_adj_e')
378
379 ! Initialize the source term
380 call this%source_term%init(this%f_adj_x, this%f_adj_y, this%f_adj_z, &
381 this%c_Xh, user, this%name)
382 call this%source_term%add(params, 'case.adjoint_fluid.source_term')
383
384 call neko_log%end_section()
385
386 end subroutine adjoint_fluid_scheme_init_base
387
389 subroutine adjoint_fluid_scheme_free(this)
390 class(adjoint_fluid_scheme_incompressible_t), intent(inout) :: this
391
392 call this%Xh%free()
393
394 if (allocated(this%ksp_vel)) then
395 call this%ksp_vel%free()
396 deallocate(this%ksp_vel)
397 end if
398
399 if (allocated(this%ksp_prs)) then
400 call this%ksp_prs%free()
401 deallocate(this%ksp_prs)
402 end if
403
404 if (allocated(this%pc_vel)) then
405 call precon_destroy(this%pc_vel)
406 deallocate(this%pc_vel)
407 end if
408
409 if (allocated(this%pc_prs)) then
410 call precon_destroy(this%pc_prs)
411 deallocate(this%pc_prs)
412 end if
413
414 call this%source_term%free()
415
416 call this%gs_Xh%free()
417
418 call this%c_Xh%free()
419
420 call this%scratch_GL%free()
421
422 nullify(this%u_adj)
423 nullify(this%v_adj)
424 nullify(this%w_adj)
425 nullify(this%p_adj)
426
427 nullify(this%u_adj_e)
428 nullify(this%v_adj_e)
429 nullify(this%w_adj_e)
430
431 call this%ulag%free()
432 call this%vlag%free()
433 call this%wlag%free()
434
435
436 if (associated(this%f_adj_x)) then
437 call this%f_adj_x%free()
438 end if
439
440 if (associated(this%f_adj_y)) then
441 call this%f_adj_y%free()
442 end if
443
444 if (associated(this%f_adj_z)) then
445 call this%f_adj_z%free()
446 end if
447
448 nullify(this%f_adj_x)
449 nullify(this%f_adj_y)
450 nullify(this%f_adj_z)
451
452 call this%rho%free()
453 call this%mu%free()
454
455 end subroutine adjoint_fluid_scheme_free
456
459 subroutine adjoint_fluid_scheme_validate(this)
460 class(adjoint_fluid_scheme_incompressible_t), target, intent(inout) :: this
461 ! Variables for retrieving json parameters
462
463 if ((.not. associated(this%u_adj)) .or. &
464 (.not. associated(this%v_adj)) .or. &
465 (.not. associated(this%w_adj)) .or. &
466 (.not. associated(this%p_adj))) then
467 call neko_error('Fields are not registered')
468 end if
469
470 if ((.not. allocated(this%u_adj%x)) .or. &
471 (.not. allocated(this%v_adj%x)) .or. &
472 (.not. allocated(this%w_adj%x)) .or. &
473 (.not. allocated(this%p_adj%x))) then
474 call neko_error('Fields are not allocated')
475 end if
476
477 if (.not. allocated(this%ksp_vel)) then
478 call neko_error('No Krylov solver for velocity defined')
479 end if
480
481 if (.not. allocated(this%ksp_prs)) then
482 call neko_error('No Krylov solver for pressure defined')
483 end if
484
485 !
486 ! Setup checkpoint structure (if everything is fine)
487 !
488 call this%chkp%init(this%u_adj, this%v_adj, this%w_adj, this%p_adj)
489
490 end subroutine adjoint_fluid_scheme_validate
491
496 subroutine adjoint_fluid_scheme_bc_apply_vel(this, time, strong)
497 class(adjoint_fluid_scheme_incompressible_t), intent(inout) :: this
498 type(time_state_t), intent(in) :: time
499 logical, intent(in) :: strong
500
501 integer :: i
502 class(bc_t), pointer :: b => null()
503
504 call this%bcs_vel%apply_vector(&
505 this%u_adj%x, this%v_adj%x, this%w_adj%x, this%dm_Xh%size(), &
506 time, strong)
507 call this%gs_Xh%op(this%u_adj, gs_op_min, glb_cmd_event)
508 call this%gs_Xh%op(this%v_adj, gs_op_min, glb_cmd_event)
509 call this%gs_Xh%op(this%w_adj, gs_op_min, glb_cmd_event)
510 call device_event_sync(glb_cmd_event)
511
512 call this%bcs_vel%apply_vector(&
513 this%u_adj%x, this%v_adj%x, this%w_adj%x, this%dm_Xh%size(), &
514 time, strong)
515 call this%gs_Xh%op(this%u_adj, gs_op_max, glb_cmd_event)
516 call this%gs_Xh%op(this%v_adj, gs_op_max, glb_cmd_event)
517 call this%gs_Xh%op(this%w_adj, gs_op_max, glb_cmd_event)
518 call device_event_sync(glb_cmd_event)
519
520 do i = 1, this%bcs_vel%size()
521 b => this%bcs_vel%get(i)
522 b%updated = .false.
523 end do
524 nullify(b)
525
526
527 end subroutine adjoint_fluid_scheme_bc_apply_vel
528
531 subroutine adjoint_fluid_scheme_bc_apply_prs(this, time)
532 class(adjoint_fluid_scheme_incompressible_t), intent(inout) :: this
533 type(time_state_t), intent(in) :: time
534
535 integer :: i
536 class(bc_t), pointer :: b => null()
537
538 call this%bcs_prs%apply(this%p_adj, time)
539 call this%gs_Xh%op(this%p_adj, gs_op_min, glb_cmd_event)
540 call device_event_sync(glb_cmd_event)
541
542 call this%bcs_prs%apply(this%p_adj, time)
543 call this%gs_Xh%op(this%p_adj, gs_op_max, glb_cmd_event)
544 call device_event_sync(glb_cmd_event)
545
546 do i = 1, this%bcs_prs%size()
547 b => this%bcs_prs%get(i)
548 b%updated = .false.
549 end do
550 nullify(b)
551
552 end subroutine adjoint_fluid_scheme_bc_apply_prs
553
556 subroutine adjoint_fluid_scheme_solver_factory(ksp, n, solver, &
557 max_iter, abstol, monitor)
558 class(ksp_t), allocatable, target, intent(inout) :: ksp
559 integer, intent(in), value :: n
560 character(len=*), intent(in) :: solver
561 integer, intent(in) :: max_iter
562 real(kind=rp), intent(in) :: abstol
563 logical, intent(in) :: monitor
564
565 call krylov_solver_factory(ksp, n, solver, max_iter, abstol, &
566 monitor = monitor)
567
568 end subroutine adjoint_fluid_scheme_solver_factory
569
571 subroutine adjoint_fluid_scheme_precon_factory(this, pc, ksp, coef, dof, gs, &
572 bclst, pctype, pcparams)
573 class(adjoint_fluid_scheme_incompressible_t), intent(inout) :: this
574 class(pc_t), allocatable, target, intent(inout) :: pc
575 class(ksp_t), target, intent(inout) :: ksp
576 type(coef_t), target, intent(in) :: coef
577 type(dofmap_t), target, intent(in) :: dof
578 type(gs_t), target, intent(inout) :: gs
579 type(bc_list_t), target, intent(inout) :: bclst
580 character(len=*) :: pctype
581 type(json_file), intent(inout) :: pcparams
582
583 call precon_factory(pc, pctype)
584
585 select type (pcp => pc)
586 type is (jacobi_t)
587 call pcp%init(coef, dof, gs)
588 type is (sx_jacobi_t)
589 call pcp%init(coef, dof, gs)
590 type is (device_jacobi_t)
591 call pcp%init(coef, dof, gs)
592 type is (hsmg_t)
593 call pcp%init(coef, bclst, pcparams)
594 type is (phmg_t)
595 call pcp%init(coef, bclst, pcparams)
596 end select
597
598 call ksp%set_pc(pc)
599
600 end subroutine adjoint_fluid_scheme_precon_factory
601
603 ! TODO
604 ! HARRY
605 ! ok this needs to be discussed.
606 ! To me, the CFL is based on the convecting velocity, so I would argue
607 ! that the CFL of the adjoint is based on
608 ! u, v, w
609 ! not
610 ! u_adj, v_adj, w_adj
611 !
612 ! Then again, this function is really only used for varaible timestep,
613 ! which we can't use if we're checkpointing...
614 !
615 ! hmmmmm.... requires more thought. Because people optimal ICs or Arnouldi
616 ! may use variable timesteps... or us in a steady case..
617 !
618 ! for now.... let's ignore it
619 function adjoint_compute_cfl(this, dt) result(c)
620 class(adjoint_fluid_scheme_incompressible_t), intent(in) :: this
621 real(kind=rp), intent(in) :: dt
622 real(kind=rp) :: c
623
624 c = cfl(dt, this%u_adj%x, this%v_adj%x, this%w_adj%x, &
625 this%Xh, this%c_Xh, this%msh%nelv, this%msh%gdim)
626
627 end function adjoint_compute_cfl
628
633 subroutine adjoint_fluid_scheme_update_material_properties(this, time)
634 class(adjoint_fluid_scheme_incompressible_t), intent(inout) :: this
635 type(time_state_t), intent(in) :: time
636 type(field_t), pointer :: nut
637
638 call this%user_material_properties(this%name, &
639 this%material_properties, time)
640
641 if (len(trim(this%nut_field_name)) > 0) then
642 nut => neko_registry%get_field(this%nut_field_name)
643 call field_addcol3(this%mu, nut, this%rho)
644 end if
645
646 ! Since mu, rho is a field_t, and we use the %x(1,1,1,1)
647 ! host array data to pass constant density and viscosity
648 ! to some routines, we need to make sure that the host
649 ! values are also filled
650 if (neko_bcknd_device .eq. 1) then
651 call device_memcpy(this%rho%x, this%rho%x_d, this%rho%size(), &
652 device_to_host, sync=.false.)
653 call device_memcpy(this%mu%x, this%mu%x_d, this%mu%size(), &
654 device_to_host, sync=.false.)
655 end if
656 end subroutine adjoint_fluid_scheme_update_material_properties
657
662 subroutine adjoint_fluid_scheme_set_material_properties(this, params, user)
663 class(adjoint_fluid_scheme_incompressible_t), intent(inout) :: this
664 type(json_file), intent(inout) :: params
665 type(user_t), target, intent(in) :: user
666 character(len=LOG_SIZE) :: log_buf
667 ! A local pointer that is needed to make Intel happy
668 procedure(user_material_properties_intf), pointer :: dummy_mp_ptr
669 real(kind=rp) :: const_mu, const_rho
670 type(time_state_t) :: time
671
672
673 dummy_mp_ptr => dummy_user_material_properties
674
675 call this%mu%init(this%dm_Xh, "mu")
676 call this%rho%init(this%dm_Xh, "rho")
677 call this%material_properties%init(2)
678 call this%material_properties%assign_to_field(1, this%rho)
679 call this%material_properties%assign_to_field(2, this%mu)
680
681 if (.not. associated(user%material_properties, dummy_mp_ptr)) then
682
683 write(log_buf, '(A)') "Material properties must be set in the user&
684 & file!"
685 call neko_log%message(log_buf)
686 this%user_material_properties => user%material_properties
687
688 call user%material_properties(this%name, &
689 this%material_properties, time)
690
691 else
692 this%user_material_properties => dummy_user_material_properties
693 ! Incorrect user input
694 if (params%valid_path('case.fluid.Re') .and. &
695 (params%valid_path('case.fluid.mu') .or. &
696 params%valid_path('case.fluid.rho'))) then
697 call neko_error("To set the material properties for the fluid, " // &
698 "either provide Re OR mu and rho in the case file.")
699
700 else if (params%valid_path('case.fluid.Re')) then
701 ! Non-dimensional case
702 write(log_buf, '(A)') 'Non-dimensional fluid material properties &
703 & input.'
704 call neko_log%message(log_buf, lvl = neko_log_verbose)
705 write(log_buf, '(A)') 'Density will be set to 1, dynamic viscosity to&
706 & 1/Re.'
707 call neko_log%message(log_buf, lvl = neko_log_verbose)
708
709 ! Read Re into mu for further manipulation.
710 call json_get(params, 'case.fluid.Re', const_mu)
711 write(log_buf, '(A)') 'Read non-dimensional material properties'
712 call neko_log%message(log_buf)
713 write(log_buf, '(A,ES13.6)') 'Re :', const_mu
714 call neko_log%message(log_buf)
715
716 ! Set rho to 1 since the setup is non-dimensional.
717 const_rho = 1.0_rp
718 ! Invert the Re to get viscosity.
719 const_mu = 1.0_rp/const_mu
720 else
721 ! Dimensional case
722 call json_get(params, 'case.fluid.mu', const_mu)
723 call json_get(params, 'case.fluid.rho', const_rho)
724 end if
725 end if
726
727 ! We need to fill the fields based on the parsed const values
728 ! if the user routine is not used.
729 if (associated(user%material_properties, dummy_mp_ptr)) then
730 ! Fill mu and rho field with the physical value
731 call field_cfill(this%mu, const_mu)
732 call field_cfill(this%rho, const_rho)
733
734
735 write(log_buf, '(A,ES13.6)') 'rho :', const_rho
736 call neko_log%message(log_buf)
737 write(log_buf, '(A,ES13.6)') 'mu :', const_mu
738 call neko_log%message(log_buf)
739 end if
740
741 ! Since mu, rho is a field_t, and we use the %x(1,1,1,1)
742 ! host array data to pass constant density and viscosity
743 ! to some routines, we need to make sure that the host
744 ! values are also filled
745 if (neko_bcknd_device .eq. 1) then
746 call device_memcpy(this%rho%x, this%rho%x_d, this%rho%size(), &
747 device_to_host, sync=.false.)
748 call device_memcpy(this%mu%x, this%mu%x_d, this%mu%size(), &
749 device_to_host, sync=.false.)
750 end if
751 end subroutine adjoint_fluid_scheme_set_material_properties
752
subroutine adjoint_fluid_scheme_init_base(this, msh, lx, params, scheme, user, kspv_init)
Initialise an adjoint scheme.
Implements the adjoint_source_term_t type.
Wrapper contaning and executing the adjoint source terms.