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