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