Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_fluid_scheme.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 gather_scatter, only: gs_t
36 use mean_sqr_flow, only: mean_sqr_flow_t
37 use neko_config, only: neko_bcknd_device
38 use checkpoint, only: chkp_t
39 use mean_flow, only: mean_flow_t
40 use num_types, only: rp, i8
41 use comm, only: neko_comm
43 use field, only: field_t
44 use space, only: space_t, gll
45 use dofmap, only: dofmap_t
46 use zero_dirichlet, only: zero_dirichlet_t
47 use krylov, only: ksp_t, krylov_solver_factory, ksp_max_iter
48 use coefs, only: coef_t
49 use usr_inflow, only: usr_inflow_t, usr_inflow_eval
50 use dirichlet, only: dirichlet_t
51 use field_dirichlet, only: field_dirichlet_t
52 use field_dirichlet_vector, only: field_dirichlet_vector_t
53 use jacobi, only: jacobi_t
54 use sx_jacobi, only: sx_jacobi_t
55 use device_jacobi, only: device_jacobi_t
56 use hsmg, only: hsmg_t
57 use phmg, only: phmg_t
58 use precon, only: pc_t, precon_factory, precon_destroy
59 use fluid_stats, only: fluid_stats_t
60 use bc, only: bc_t
61 use bc_list, only: bc_list_t
62 use mesh, only: mesh_t, neko_msh_max_zlbls, neko_msh_max_zlbl_len
63 use math, only: cfill, add2s2, glsum
64 use device_math, only: device_cfill, device_add2s2
65 use time_scheme_controller, only: time_scheme_controller_t
66 use operators, only: cfl
67 use logger, only: neko_log, log_size, neko_log_verbose
68 use field_registry, only: neko_field_registry
69 use json_utils, only: json_get, json_get_or_default, json_extract_object, &
70 json_extract_item
71 use json_module, only: json_file, json_core, json_value
72 use scratch_registry, only: scratch_registry_t
73 use user_intf, only: user_t, dummy_user_material_properties, &
74 user_material_properties
75 use utils, only: neko_error, neko_warning
76 use field_series, only: field_series_t
77 use time_step_controller, only: time_step_controller_t
78 use field_math, only: field_cfill, field_add2s2
79 use wall_model_bc, only: wall_model_bc_t
80 use shear_stress, only: shear_stress_t
81 use field_list, only : field_list_t
82 use gradient_jump_penalty, only: gradient_jump_penalty_t
83 use field_math, only: field_addcol3
84
85 use mpi_f08, only: mpi_integer, mpi_sum, mpi_allreduce
86 use json_utils_ext, only: json_key_fallback
87 use device, only : device_event_sync, glb_cmd_event, device_to_host, &
88 device_memcpy
89 implicit none
90 private
91
93 type, abstract :: adjoint_fluid_scheme_t
95 character(len=:), allocatable :: name
97 type(field_t), pointer :: u_adj => null()
99 type(field_t), pointer :: v_adj => null()
101 type(field_t), pointer :: w_adj => null()
103 type(field_t), pointer :: p_adj => null()
105 type(field_series_t) :: ulag, vlag, wlag
107 type(space_t) :: xh
109 type(dofmap_t) :: dm_xh
111 type(gs_t) :: gs_xh
113 type(coef_t) :: c_xh
115 type(adjoint_source_term_t) :: source_term
117 type(field_t), pointer :: f_adj_x => null()
119 type(field_t), pointer :: f_adj_y => null()
121 type(field_t), pointer :: f_adj_z => null()
122
123 ! Krylov solvers and settings
125 class(ksp_t), allocatable :: ksp_vel
127 class(ksp_t), allocatable :: ksp_prs
129 class(pc_t), allocatable :: pc_vel
131 class(pc_t), allocatable :: pc_prs
133 integer :: vel_projection_dim
135 integer :: pr_projection_dim
137 integer :: vel_projection_activ_step
139 integer :: pr_projection_activ_step
141 logical :: strict_convergence
142
143 ! List of boundary conditions for pressure
144 type(bc_list_t) :: bcs_prs
145 ! List of boundary conditions for velocity
146 type(bc_list_t) :: bcs_vel
148 type(mesh_t), pointer :: msh => null()
150 type(chkp_t) :: chkp
152 type(mean_flow_t) :: mean
154 type(fluid_stats_t) :: stats
156 type(mean_sqr_flow_t) :: mean_sqr
158 logical :: forced_flow_rate = .false.
160 logical :: freeze = .false.
161
163 character(len=:), allocatable :: nut_field_name
165 logical :: variable_material_properties = .false.
166
168 integer(kind=i8) :: glb_n_points
170 integer(kind=i8) :: glb_unique_points
172 type(scratch_registry_t) :: scratch
173
175 type(field_t) :: rho
176
178 type(field_t) :: mu
179
181 type(field_list_t) :: material_properties
182
184 procedure(user_material_properties), nopass, pointer :: &
185 user_material_properties => null()
186 contains
188 procedure, pass(this) :: init_base => adjoint_fluid_scheme_init_base
190 procedure, pass(this) :: scheme_free => adjoint_fluid_scheme_free
192 procedure, pass(this) :: validate => adjoint_fluid_scheme_validate
194 procedure, pass(this) :: bc_apply_vel => adjoint_fluid_scheme_bc_apply_vel
196 procedure, pass(this) :: bc_apply_prs => adjoint_fluid_scheme_bc_apply_prs
198 procedure, pass(this) :: compute_cfl => adjoint_compute_cfl
200 procedure, pass(this) :: set_material_properties => &
201 adjoint_fluid_scheme_set_material_properties
203 procedure(adjoint_fluid_scheme_init_intrf), pass(this), deferred :: init
205 procedure(adjoint_fluid_scheme_free_intrf), pass(this), deferred :: free
207 procedure(adjoint_fluid_scheme_step_intrf), pass(this), deferred :: step
209 procedure(adjoint_fluid_scheme_restart_intrf), &
210 pass(this), deferred :: restart
212 procedure(adjoint_fluid_scheme_setup_bcs_intrf), pass(this), deferred :: &
213 setup_bcs
215 procedure, pass(this) :: update_material_properties => &
216 adjoint_fluid_scheme_update_material_properties
218 procedure, nopass :: solver_factory => adjoint_fluid_scheme_solver_factory
220 procedure, pass(this) :: precon_factory_ => &
221 adjoint_fluid_scheme_precon_factory
223
224
226 abstract interface
227 subroutine adjoint_fluid_scheme_init_intrf(this, msh, lx, params, user, &
228 time_scheme)
230 import json_file
231 import mesh_t
232 import user_t
233 import time_scheme_controller_t
234 class(adjoint_fluid_scheme_t), target, intent(inout) :: this
235 type(mesh_t), target, intent(inout) :: msh
236 integer, intent(in) :: lx
237 type(json_file), target, intent(inout) :: params
238 type(user_t), target, intent(in) :: user
239 type(time_scheme_controller_t), target, intent(in) :: time_scheme
241 end interface
242
244 abstract interface
245 subroutine adjoint_fluid_scheme_free_intrf(this)
247 class(adjoint_fluid_scheme_t), intent(inout) :: this
248 end subroutine adjoint_fluid_scheme_free_intrf
249 end interface
250
252 abstract interface
253 subroutine adjoint_fluid_scheme_step_intrf(this, t, tstep, dt, ext_bdf, &
254 dt_controller)
256 import time_scheme_controller_t
257 import time_step_controller_t
258 import rp
259 class(adjoint_fluid_scheme_t), target, intent(inout) :: this
260 real(kind=rp), intent(in) :: t
261 integer, intent(in) :: tstep
262 real(kind=rp), intent(in) :: dt
263 type(time_scheme_controller_t), intent(in) :: ext_bdf
264 type(time_step_controller_t), intent(in) :: dt_controller
265 end subroutine adjoint_fluid_scheme_step_intrf
266 end interface
267
269 abstract interface
270 subroutine adjoint_fluid_scheme_restart_intrf(this, dtlag, tlag)
272 import rp
273 class(adjoint_fluid_scheme_t), target, intent(inout) :: this
274 real(kind=rp) :: dtlag(10), tlag(10)
275
276 end subroutine adjoint_fluid_scheme_restart_intrf
277 end interface
278
280 abstract interface
281 subroutine adjoint_fluid_scheme_setup_bcs_intrf(this, user, params)
282 import adjoint_fluid_scheme_t, user_t, json_file
283 class(adjoint_fluid_scheme_t), intent(inout) :: this
284 type(user_t), target, intent(in) :: user
285 type(json_file), intent(inout) :: params
286 end subroutine adjoint_fluid_scheme_setup_bcs_intrf
287 end interface
288
289 interface
290
291 module subroutine adjoint_fluid_scheme_factory(object, type_name)
292 class(adjoint_fluid_scheme_t), intent(inout), allocatable :: object
293 character(len=*) :: type_name
294 end subroutine adjoint_fluid_scheme_factory
295 end interface
296
297 public :: adjoint_fluid_scheme_t, adjoint_fluid_scheme_factory
298
299contains
300
302 subroutine adjoint_fluid_scheme_init_base(this, msh, lx, params, scheme, &
303 user, kspv_init)
304 class(adjoint_fluid_scheme_t), target, intent(inout) :: this
305 type(mesh_t), target, intent(inout) :: msh
306 integer, intent(in) :: lx
307 character(len=*), intent(in) :: scheme
308 type(json_file), target, intent(inout) :: params
309 type(user_t), target, intent(in) :: user
310 logical, intent(in) :: kspv_init
311 character(len=LOG_SIZE) :: log_buf
312 real(kind=rp) :: real_val
313 logical :: logical_val
314 integer :: integer_val
315 character(len=:), allocatable :: string_val1, string_val2
316 character(len=:), allocatable :: json_key
317 type(json_file) :: json_subdict
318
319 !
320 ! SEM simulation fundamentals
321 !
322
323 this%msh => msh
324
325 if (msh%gdim .eq. 2) then
326 call this%Xh%init(gll, lx, lx)
327 else
328 call this%Xh%init(gll, lx, lx, lx)
329 end if
330
331 call this%dm_Xh%init(msh, this%Xh)
332
333 call this%gs_Xh%init(this%dm_Xh)
334
335 call this%c_Xh%init(this%gs_Xh)
336
337 ! Local scratch registry
338 this%scratch = scratch_registry_t(this%dm_Xh, 10, 2)
339
340 !
341 ! First section of fluid log
342 !
343
344 call neko_log%section('Adjoint fluid')
345 write(log_buf, '(A, A)') 'Type : ', trim(scheme)
346 call neko_log%message(log_buf)
347
348 !
349 ! Material properties
350 !
351 call this%set_material_properties(params, user)
352
353 ! Projection spaces
354 json_key = json_key_fallback(params, &
355 'case.adjoint_fluid.velocity_solver.projection_space_size', &
356 'case.fluid.velocity_solver.projection_space_size')
357 call json_get_or_default(params, json_key, this%vel_projection_dim, 20)
358 json_key = json_key_fallback(params, &
359 'case.adjoint_fluid.pressure_solver.projection_space_size', &
360 'case.fluid.pressure_solver.projection_space_size')
361 call json_get_or_default(params, json_key, this%pr_projection_dim, 20)
362
363 json_key = json_key_fallback(params, &
364 'case.adjoint_fluid.velocity_solver.projection_hold_steps', &
365 'case.fluid.velocity_solver.projection_hold_steps')
366 call json_get_or_default(params, json_key, &
367 this%vel_projection_activ_step, 5)
368 json_key = json_key_fallback(params, &
369 'case.adjoint_fluid.pressure_solver.projection_hold_steps', &
370 'case.fluid.pressure_solver.projection_hold_steps')
371 call json_get_or_default(params, json_key, &
372 this%pr_projection_activ_step, 5)
373
374 json_key = json_key_fallback(params, 'case.adjoint_fluid.freeze', &
375 'case.fluid.freeze')
376 call json_get_or_default(params, json_key, this%freeze, .false.)
377
378 ! TODO
379 ! This needs to be discussed... In pricipal, I think if the forward is
380 ! forced to a
381 ! certain flow rate, then the adjoint should be forced to zero flow rate,
382 ! we had a derivation in
383 ! https://www.sciencedirect.com/science/article/pii/S0045782522006764
384 ! for now I'm commenting this out
385 if (params%valid_path("case.fluid.flow_rate_force")) then
386 call neko_error('Flow rate forcing not yet implemented')
387 this%forced_flow_rate = .true.
388 end if
389
390
391 if (lx .lt. 10) then
392 write(log_buf, '(A, I1)') 'Poly order : ', lx-1
393 else if (lx .ge. 10) then
394 write(log_buf, '(A, I2)') 'Poly order : ', lx-1
395 else
396 write(log_buf, '(A, I3)') 'Poly order : ', lx-1
397 end if
398 call neko_log%message(log_buf)
399 this%glb_n_points = int(this%msh%glb_nelv, i8)*int(this%Xh%lxyz, i8)
400 this%glb_unique_points = int(glsum(this%c_Xh%mult, this%dm_Xh%size()), i8)
401
402 write(log_buf, '(A, I0)') 'GLL points : ', this%glb_n_points
403 call neko_log%message(log_buf)
404 write(log_buf, '(A, I0)') 'Unique pts.: ', this%glb_unique_points
405 call neko_log%message(log_buf)
406
407 call json_get(params, 'case.numerics.dealias', logical_val)
408 write(log_buf, '(A, L1)') 'Dealias : ', logical_val
409 call neko_log%message(log_buf)
410
411 write(log_buf, '(A, L1)') 'LES : ', this%variable_material_properties
412 call neko_log%message(log_buf)
413
414 call json_get_or_default(params, 'case.output_boundary', logical_val, &
415 .false.)
416 write(log_buf, '(A, L1)') 'Save bdry : ', logical_val
417 call neko_log%message(log_buf)
418
419 !
420 ! Setup right-hand side fields.
421 !
422 allocate(this%f_adj_x)
423 allocate(this%f_adj_y)
424 allocate(this%f_adj_z)
425 call this%f_adj_x%init(this%dm_Xh, fld_name = "adjoint_rhs_x")
426 call this%f_adj_y%init(this%dm_Xh, fld_name = "adjoint_rhs_y")
427 call this%f_adj_z%init(this%dm_Xh, fld_name = "adjoint_rhs_z")
428
429 ! Initialize the source term
430 ! TODO
431 ! HARRY
432 ! Ok, this one I dislike the most... I think the CORRECT solution here is to
433 ! modify
434 ! the fluid_source_term.f90 so that case.fluid is not hardcoded
435 ! whoever coded it originally is probably much smarter than me and has their
436 ! reasons
437 !
438 ! I would envision somehow merging fluid_source_term and scalar_source_term
439 ! (maybe they're different because because it's a vector forcing vs a scalar
440 ! forcing,
441 ! but still...
442 !
443 ! this is my attempt...
444 ! params replaced by adjoint_json
445 !call this%source_term%init(params, this%f_adj_x, this%f_adj_y, &
446 ! this%f_adj_z, this%c_Xh,&
447 ! user)
448
449 call this%source_term%init(this%f_adj_x, this%f_adj_y, this%f_adj_z, &
450 this%c_Xh, user)
451 call this%source_term%add(params, 'case.adjoint_fluid.source_term')
452
453 ! Todo: HARRY
454 ! --------------------------------------------------------------------------
455 ! ok here is a chance that we can maybe steal the krylov solvers from the
456 ! forward??
457 !
458 ! something along the lines of
459 !
460 ! if (.not.params%valid_path('case.adjoint_fluid.velocity_solver')) then
461 ! this%ksp_vel => case%fluid%ksp_vel
462 ! this%pc_vel => case%fluid%ksp_vel
463 ! else
464 ! initialize everything
465 ! end if
466 !
467 ! but I don't know how we can steal the krylov solvers out of the forward,
468 ! so for now we'll just initialize two of them...
469 ! Initialize velocity solver
470 if (kspv_init) then
471 call neko_log%section("Adjoint Velocity solver")
472
473 json_key = json_key_fallback(params, &
474 'case.adjoint_fluid.velocity_solver.max_iterations', &
475 'case.fluid.velocity_solver.max_iterations')
476 call json_get_or_default(params, json_key, integer_val, ksp_max_iter)
477
478 json_key = json_key_fallback(params, &
479 'case.adjoint_fluid.velocity_solver.type', &
480 'case.fluid.velocity_solver.type')
481 call json_get(params, json_key, string_val1)
482
483 json_key = json_key_fallback(params, &
484 'case.adjoint_fluid.velocity_solver.preconditioner.type', &
485 'case.fluid.velocity_solver.preconditioner.type')
486 call json_get(params, json_key, string_val2)
487
488 json_key = json_key_fallback(params, &
489 'case.adjoint_fluid.velocity_solver.preconditioner', &
490 'case.fluid.velocity_solver.preconditioner')
491 call json_extract_object(params, json_key, json_subdict)
492
493 json_key = json_key_fallback(params, &
494 'case.adjoint_fluid.velocity_solver.absolute_tolerance', &
495 'case.fluid.velocity_solver.absolute_tolerance')
496 call json_get(params, json_key, real_val)
497
498 json_key = json_key_fallback(params, &
499 'case.adjoint_fluid.velocity_solver.monitor', &
500 'case.fluid.velocity_solver.monitor')
501 call json_get_or_default(params, json_key, logical_val, .false.)
502
503 call neko_log%message('Type : ('// trim(string_val1) // &
504 ', ' // trim(string_val2) // ')')
505
506 write(log_buf, '(A,ES13.6)') 'Abs tol :', real_val
507 call neko_log%message(log_buf)
508 call this%solver_factory(this%ksp_vel, this%dm_Xh%size(), &
509 string_val1, integer_val, real_val, logical_val)
510 call this%precon_factory_(this%pc_vel, this%ksp_vel, &
511 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs_vel, &
512 string_val2, json_subdict)
513 call neko_log%end_section()
514 end if
515
516 ! Strict convergence for the velocity solver
517 call json_get_or_default(params, 'case.fluid.strict_convergence', &
518 this%strict_convergence, .false.)
519
520 ! Assign velocity fields
521 call neko_field_registry%add_field(this%dm_Xh, 'u_adj')
522 call neko_field_registry%add_field(this%dm_Xh, 'v_adj')
523 call neko_field_registry%add_field(this%dm_Xh, 'w_adj')
524 this%u_adj => neko_field_registry%get_field('u_adj')
525 this%v_adj => neko_field_registry%get_field('v_adj')
526 this%w_adj => neko_field_registry%get_field('w_adj')
527
528 !! Initialize time-lag fields
529 call this%ulag%init(this%u_adj, 2)
530 call this%vlag%init(this%v_adj, 2)
531 call this%wlag%init(this%w_adj, 2)
532
533 call neko_log%end_section()
534
535 end subroutine adjoint_fluid_scheme_init_base
536
537 ! ========================================================================== !
538 ! Todo: This section need to be moved
539
540 ! !> Initialize all components of the current scheme
541 ! subroutine adjoint_fluid_scheme_init_all(this, msh, lx, params, kspv_init, &
542 ! kspp_init, scheme, user)
543 ! implicit none
544 ! class(adjoint_fluid_scheme_t), target, intent(inout) :: this
545 ! type(mesh_t), target, intent(inout) :: msh
546 ! integer, intent(inout) :: lx
547 ! type(json_file), target, intent(inout) :: params
548 ! type(user_t), target, intent(in) :: user
549 ! logical :: kspv_init
550 ! logical :: kspp_init
551 ! character(len=*), intent(in) :: scheme
552 ! real(kind=rp) :: abs_tol
553 ! integer :: integer_val, ierr
554 ! character(len=:), allocatable :: solver_type, precon_type, json_key
555 ! character(len=LOG_SIZE) :: log_buf
556
557 ! call adjoint_fluid_scheme_init_base(this, msh, lx, params, scheme, user, &
558 ! kspv_init)
559
560 ! call neko_field_registry%add_field(this%dm_Xh, 'p_adj')
561 ! this%p_adj => neko_field_registry%get_field('p_adj')
562
563 ! ! HARRY
564 ! ! Holy shit pressure BC's confuse me on a good day...
565 ! ! but they're especially confusing for the adjoint.
566 ! ! for now... I'm keeping them the same as the primal
567 ! ! but I wouldn't be surprised if this is incorrect for the adjoint
568 ! !
569 ! ! Setup pressure boundary conditions
570 ! !
571 ! call this%bcs_prs%init()
572 ! call this%bc_prs%init_base(this%c_Xh)
573 ! call this%bc_prs%mark_zones_from_list(msh%labeled_zones,&
574 ! 'o', this%bc_labels)
575 ! call this%bc_prs%mark_zones_from_list(msh%labeled_zones,&
576 ! 'on', this%bc_labels)
577
578 ! ! Field dirichlet pressure bc
579 ! call this%user_field_bc_prs%init_base(this%c_Xh)
580 ! call this%user_field_bc_prs%mark_zones_from_list(msh%labeled_zones,&
581 ! 'd_pres', this%bc_labels)
582 ! call this%user_field_bc_prs%finalize()
583 ! call MPI_Allreduce(this%user_field_bc_prs%msk(0), integer_val, 1, &
584 ! MPI_INTEGER, MPI_SUM, NEKO_COMM, ierr)
585
586 ! if (integer_val .gt. 0) call this%user_field_bc_prs%init_field('d_pres')
587 ! call this%bcs_prs%append(this%user_field_bc_prs)
588 ! call this%user_field_bc_vel%bc_list%append(this%user_field_bc_prs)
589
590 ! if (msh%outlet%size .gt. 0) then
591 ! call this%bc_prs%mark_zone(msh%outlet)
592 ! end if
593 ! if (msh%outlet_normal%size .gt. 0) then
594 ! call this%bc_prs%mark_zone(msh%outlet_normal)
595 ! end if
596
597 ! call this%bc_prs%finalize()
598 ! call this%bc_prs%set_g(0.0_rp)
599 ! call this%bcs_prs%append(this%bc_prs)
600 ! call this%bc_dong%init_base(this%c_Xh)
601 ! call this%bc_dong%mark_zones_from_list(msh%labeled_zones,&
602 ! 'o+dong', this%bc_labels)
603 ! call this%bc_dong%mark_zones_from_list(msh%labeled_zones,&
604 ! 'on+dong', this%bc_labels)
605 ! call this%bc_dong%finalize()
606
607
608 ! call this%bc_dong%init(this%c_Xh, params)
609
610 ! call this%bcs_prs%append(this%bc_dong)
611
612 ! ! Pressure solver
613 ! ! hopefully we can do the same trick here as we will eventually do for the
614 ! ! velocity solver
615 ! if (kspp_init) then
616 ! call neko_log%section("Pressure solver")
617
618 ! json_key = json_key_fallback(params, &
619 ! 'case.adjoint_fluid.pressure_solver.max_iterations', &
620 ! 'case.fluid.pressure_solver.max_iterations')
621 ! call json_get_or_default(params, json_key, integer_val, 800)
622
623 ! json_key = json_key_fallback(params, &
624 ! 'case.adjoint_fluid.pressure_solver.type', &
625 ! 'case.fluid.pressure_solver.type')
626 ! call json_get(params, json_key, solver_type)
627
628 ! json_key = json_key_fallback(params, &
629 ! 'case.adjoint_fluid.pressure_solver.preconditioner', &
630 ! 'case.fluid.pressure_solver.preconditioner')
631 ! call json_get(params, json_key, precon_type)
632
633 ! json_key = json_key_fallback(params, &
634 ! 'case.adjoint_fluid.pressure_solver.absolute_tolerance', &
635 ! 'case.fluid.pressure_solver.absolute_tolerance')
636 ! call json_get(params, json_key, abs_tol)
637
638 ! call neko_log%message('Type : ('// trim(solver_type) // &
639 ! ', ' // trim(precon_type) // ')')
640 ! write(log_buf, '(A,ES13.6)') 'Abs tol :', abs_tol
641 ! call neko_log%message(log_buf)
642
643 ! call adjoint_fluid_scheme_solver_factory(this%ksp_prs, &
644 ! this%dm_Xh%size(), &
645 ! solver_type, integer_val, abs_tol)
646 ! call adjoint_fluid_scheme_precon_factory(this%pc_prs, this%ksp_prs, &
647 ! this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs_prs, precon_type)
648
649 ! call neko_log%end_section()
650
651 ! end if
652
653
654 ! call neko_log%end_section()
655
656 ! end subroutine adjoint_fluid_scheme_init_all
657
658 ! End of section to be moved
659 ! ========================================================================== !
660
662 subroutine adjoint_fluid_scheme_free(this)
663 class(adjoint_fluid_scheme_t), intent(inout) :: this
664
665 !
666 ! Free everything related to field_dirichlet BCs
667 !
668
669 call this%Xh%free()
670
671 if (allocated(this%ksp_vel)) then
672 call this%ksp_vel%free()
673 deallocate(this%ksp_vel)
674 end if
675
676 if (allocated(this%ksp_prs)) then
677 call this%ksp_prs%free()
678 deallocate(this%ksp_prs)
679 end if
680
681 if (allocated(this%pc_vel)) then
682 call precon_destroy(this%pc_vel)
683 deallocate(this%pc_vel)
684 end if
685
686 if (allocated(this%pc_prs)) then
687 call precon_destroy(this%pc_prs)
688 deallocate(this%pc_prs)
689 end if
690
691 call this%source_term%free()
692
693 call this%gs_Xh%free()
694
695 call this%c_Xh%free()
696
697 call this%scratch%free()
698
699 nullify(this%u_adj)
700 nullify(this%v_adj)
701 nullify(this%w_adj)
702 nullify(this%p_adj)
703
704 call this%ulag%free()
705 call this%vlag%free()
706 call this%wlag%free()
707
708
709 if (associated(this%f_adj_x)) then
710 call this%f_adj_x%free()
711 end if
712
713 if (associated(this%f_adj_y)) then
714 call this%f_adj_y%free()
715 end if
716
717 if (associated(this%f_adj_z)) then
718 call this%f_adj_z%free()
719 end if
720
721 nullify(this%f_adj_x)
722 nullify(this%f_adj_y)
723 nullify(this%f_adj_z)
724
725 call this%mu%free()
726
727 end subroutine adjoint_fluid_scheme_free
728
731 subroutine adjoint_fluid_scheme_validate(this)
732 class(adjoint_fluid_scheme_t), target, intent(inout) :: this
733 ! Variables for retrieving json parameters
734
735 if ((.not. associated(this%u_adj)) .or. &
736 (.not. associated(this%v_adj)) .or. &
737 (.not. associated(this%w_adj)) .or. &
738 (.not. associated(this%p_adj))) then
739 call neko_error('Fields are not registered')
740 end if
741
742 if ((.not. allocated(this%u_adj%x)) .or. &
743 (.not. allocated(this%v_adj%x)) .or. &
744 (.not. allocated(this%w_adj%x)) .or. &
745 (.not. allocated(this%p_adj%x))) then
746 call neko_error('Fields are not allocated')
747 end if
748
749 if (.not. allocated(this%ksp_vel)) then
750 call neko_error('No Krylov solver for velocity defined')
751 end if
752
753 if (.not. allocated(this%ksp_prs)) then
754 call neko_error('No Krylov solver for pressure defined')
755 end if
756
757 !
758 ! Setup checkpoint structure (if everything is fine)
759 !
760 call this%chkp%init(this%u_adj, this%v_adj, this%w_adj, this%p_adj)
761
762 end subroutine adjoint_fluid_scheme_validate
763
768 subroutine adjoint_fluid_scheme_bc_apply_vel(this, t, tstep, strong)
769 class(adjoint_fluid_scheme_t), intent(inout) :: this
770 real(kind=rp), intent(in) :: t
771 integer, intent(in) :: tstep
772 logical, intent(in) :: strong
773
774 call this%bcs_vel%apply_vector(&
775 this%u_adj%x, this%v_adj%x, this%w_adj%x, this%dm_Xh%size(), &
776 t, tstep, strong)
777
778 end subroutine adjoint_fluid_scheme_bc_apply_vel
779
782 subroutine adjoint_fluid_scheme_bc_apply_prs(this, t, tstep)
783 class(adjoint_fluid_scheme_t), intent(inout) :: this
784 real(kind=rp), intent(in) :: t
785 integer, intent(in) :: tstep
786
787 call this%bcs_prs%apply(this%p_adj, t, tstep)
788
789 end subroutine adjoint_fluid_scheme_bc_apply_prs
790
793 subroutine adjoint_fluid_scheme_solver_factory(ksp, n, solver, &
794 max_iter, abstol, monitor)
795 class(ksp_t), allocatable, target, intent(inout) :: ksp
796 integer, intent(in), value :: n
797 character(len=*), intent(in) :: solver
798 integer, intent(in) :: max_iter
799 real(kind=rp), intent(in) :: abstol
800 logical, intent(in) :: monitor
801
802 call krylov_solver_factory(ksp, n, solver, max_iter, abstol, &
803 monitor = monitor)
804
805 end subroutine adjoint_fluid_scheme_solver_factory
806
808 subroutine adjoint_fluid_scheme_precon_factory(this, pc, ksp, coef, dof, gs, &
809 bclst, pctype, pcparams)
810 class(adjoint_fluid_scheme_t), intent(inout) :: this
811 class(pc_t), allocatable, target, intent(inout) :: pc
812 class(ksp_t), target, intent(inout) :: ksp
813 type(coef_t), target, intent(in) :: coef
814 type(dofmap_t), target, intent(in) :: dof
815 type(gs_t), target, intent(inout) :: gs
816 type(bc_list_t), target, intent(inout) :: bclst
817 character(len=*) :: pctype
818 type(json_file), intent(inout) :: pcparams
819
820 call precon_factory(pc, pctype)
821
822 select type (pcp => pc)
823 type is (jacobi_t)
824 call pcp%init(coef, dof, gs)
825 type is (sx_jacobi_t)
826 call pcp%init(coef, dof, gs)
827 type is (device_jacobi_t)
828 call pcp%init(coef, dof, gs)
829 type is (hsmg_t)
830 call pcp%init(coef, bclst, pcparams)
831 type is (phmg_t)
832 call pcp%init(dof%msh, dof%Xh, coef, dof, gs, bclst)
833 end select
834
835 call ksp%set_pc(pc)
836
837 end subroutine adjoint_fluid_scheme_precon_factory
838
840 ! TODO
841 ! HARRY
842 ! ok this needs to be discussed.
843 ! To me, the CFL is based on the convecting velocity, so I would argue
844 ! that the CFL of the adjoint is based on
845 ! u, v, w
846 ! not
847 ! u_adj, v_adj, w_adj
848 !
849 ! Then again, this function is really only used for varaible timestep,
850 ! which we can't use if we're checkpointing...
851 !
852 ! hmmmmm.... requires more thought. Because people optimal ICs or Arnouldi
853 ! may use variable timesteps... or us in a steady case..
854 !
855 ! for now.... let's ignore it
856 function adjoint_compute_cfl(this, dt) result(c)
857 class(adjoint_fluid_scheme_t), intent(in) :: this
858 real(kind=rp), intent(in) :: dt
859 real(kind=rp) :: c
860
861 c = cfl(dt, this%u_adj%x, this%v_adj%x, this%w_adj%x, &
862 this%Xh, this%c_Xh, this%msh%nelv, this%msh%gdim)
863
864 end function adjoint_compute_cfl
865
866 ! ========================================================================== !
867 ! Todo: This section need to be moved
868
869 ! !> Set boundary types for the diagnostic output.
870 ! !! @param params The JSON case file.
871 ! subroutine adjoint_fluid_scheme_set_bc_type_output(this, params)
872 ! class(adjoint_fluid_scheme_t), intent(inout) :: this
873 ! type(json_file), intent(inout) :: params
874 ! type(dirichlet_t) :: bdry_mask
875 ! logical :: found
876
877 ! !
878 ! ! Check if we need to output boundaries
879 ! !
880 ! call json_get_or_default(params, 'case.output_boundary', found, .false.)
881
882 ! if (found) then
883 ! call this%bdry%init(this%dm_Xh, 'bdry')
884 ! this%bdry = 0.0_rp
885
886 ! call bdry_mask%init_base(this%c_Xh)
887 ! call bdry_mask%mark_zone(this%msh%wall)
888 ! call bdry_mask%mark_zones_from_list(this%msh%labeled_zones,&
889 ! 'w', this%bc_labels)
890 ! call bdry_mask%finalize()
891 ! call bdry_mask%set_g(1.0_rp)
892 ! call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
893 ! call bdry_mask%free()
894
895 ! call bdry_mask%init_base(this%c_Xh)
896 ! call bdry_mask%mark_zone(this%msh%inlet)
897 ! call bdry_mask%mark_zones_from_list(this%msh%labeled_zones,&
898 ! 'v', this%bc_labels)
899 ! call bdry_mask%finalize()
900 ! call bdry_mask%set_g(2.0_rp)
901 ! call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
902 ! call bdry_mask%free()
903
904 ! call bdry_mask%init_base(this%c_Xh)
905 ! call bdry_mask%mark_zone(this%msh%outlet)
906 ! call bdry_mask%mark_zones_from_list(this%msh%labeled_zones,&
907 ! 'o', this%bc_labels)
908 ! call bdry_mask%finalize()
909 ! call bdry_mask%set_g(3.0_rp)
910 ! call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
911 ! call bdry_mask%free()
912
913 ! call bdry_mask%init_base(this%c_Xh)
914 ! call bdry_mask%mark_zone(this%msh%sympln)
915 ! call bdry_mask%mark_zones_from_list(this%msh%labeled_zones,&
916 ! 'sym', this%bc_labels)
917 ! call bdry_mask%finalize()
918 ! call bdry_mask%set_g(4.0_rp)
919 ! call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
920 ! call bdry_mask%free()
921
922 ! call bdry_mask%init_base(this%c_Xh)
923 ! call bdry_mask%mark_zone(this%msh%outlet_normal)
924 ! call bdry_mask%mark_zones_from_list(this%msh%labeled_zones,&
925 ! 'on', this%bc_labels)
926 ! call bdry_mask%finalize()
927 ! call bdry_mask%set_g(5.0_rp)
928 ! call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
929 ! call bdry_mask%free()
930
931 ! call bdry_mask%init_base(this%c_Xh)
932 ! call bdry_mask%mark_zone(this%msh%periodic)
933 ! call bdry_mask%finalize()
934 ! call bdry_mask%set_g(6.0_rp)
935 ! call bdry_mask%apply_scalar(this%bdry%x, this%dm_Xh%size())
936 ! call bdry_mask%free()
937 ! end if
938
939 ! end subroutine adjoint_fluid_scheme_set_bc_type_output
940
941 ! End of section to be moved
942 ! ========================================================================== !
943
944
945
951 subroutine adjoint_fluid_scheme_update_material_properties(this, t, tstep)
952 class(adjoint_fluid_scheme_t), intent(inout) :: this
953 real(kind=rp),intent(in) :: t
954 integer, intent(in) :: tstep
955 type(field_t), pointer :: nut
956
957 call this%user_material_properties(t, tstep, this%name, &
958 this%material_properties)
959
960 if (len(trim(this%nut_field_name)) > 0) then
961 nut => neko_field_registry%get_field(this%nut_field_name)
962 call field_addcol3(this%mu, nut, this%rho)
963 end if
964
965 ! Since mu, rho is a field_t, and we use the %x(1,1,1,1)
966 ! host array data to pass constant density and viscosity
967 ! to some routines, we need to make sure that the host
968 ! values are also filled
969 if (neko_bcknd_device .eq. 1) then
970 call device_memcpy(this%rho%x, this%rho%x_d, this%rho%size(), &
971 device_to_host, sync=.false.)
972 call device_memcpy(this%mu%x, this%mu%x_d, this%mu%size(), &
973 device_to_host, sync=.false.)
974 end if
975 end subroutine adjoint_fluid_scheme_update_material_properties
976
981 subroutine adjoint_fluid_scheme_set_material_properties(this, params, user)
982 class(adjoint_fluid_scheme_t), intent(inout) :: this
983 type(json_file), intent(inout) :: params
984 type(user_t), target, intent(in) :: user
985 character(len=LOG_SIZE) :: log_buf
986 ! A local pointer that is needed to make Intel happy
987 procedure(user_material_properties), pointer :: dummy_mp_ptr
988 real(kind=rp) :: const_mu, const_rho
989
990
991 dummy_mp_ptr => dummy_user_material_properties
992
993 call this%mu%init(this%dm_Xh, "mu")
994 call this%rho%init(this%dm_Xh, "rho")
995 call this%material_properties%init(2)
996 call this%material_properties%assign_to_field(1, this%rho)
997 call this%material_properties%assign_to_field(2, this%mu)
998
999 if (.not. associated(user%material_properties, dummy_mp_ptr)) then
1000
1001 write(log_buf, '(A)') "Material properties must be set in the user&
1002 & file!"
1003 call neko_log%message(log_buf)
1004 this%user_material_properties => user%material_properties
1005
1006 call user%material_properties(0.0_rp, 0, this%name, &
1007 this%material_properties)
1008
1009 else
1010 this%user_material_properties => dummy_user_material_properties
1011 ! Incorrect user input
1012 if (params%valid_path('case.fluid.Re') .and. &
1013 (params%valid_path('case.fluid.mu') .or. &
1014 params%valid_path('case.fluid.rho'))) then
1015 call neko_error("To set the material properties for the fluid, " // &
1016 "either provide Re OR mu and rho in the case file.")
1017
1018 else if (params%valid_path('case.fluid.Re')) then
1019 ! Non-dimensional case
1020 write(log_buf, '(A)') 'Non-dimensional fluid material properties &
1021 & input.'
1022 call neko_log%message(log_buf, lvl = neko_log_verbose)
1023 write(log_buf, '(A)') 'Density will be set to 1, dynamic viscosity to&
1024 & 1/Re.'
1025 call neko_log%message(log_buf, lvl = neko_log_verbose)
1026
1027 ! Read Re into mu for further manipulation.
1028 call json_get(params, 'case.fluid.Re', const_mu)
1029 write(log_buf, '(A)') 'Read non-dimensional material properties'
1030 call neko_log%message(log_buf)
1031 write(log_buf, '(A,ES13.6)') 'Re :', const_mu
1032 call neko_log%message(log_buf)
1033
1034 ! Set rho to 1 since the setup is non-dimensional.
1035 const_rho = 1.0_rp
1036 ! Invert the Re to get viscosity.
1037 const_mu = 1.0_rp/const_mu
1038 else
1039 ! Dimensional case
1040 call json_get(params, 'case.fluid.mu', const_mu)
1041 call json_get(params, 'case.fluid.rho', const_rho)
1042 end if
1043 end if
1044
1045 ! We need to fill the fields based on the parsed const values
1046 ! if the user routine is not used.
1047 if (associated(user%material_properties, dummy_mp_ptr)) then
1048 ! Fill mu and rho field with the physical value
1049 call field_cfill(this%mu, const_mu)
1050 call field_cfill(this%rho, const_rho)
1051
1052
1053 write(log_buf, '(A,ES13.6)') 'rho :', const_rho
1054 call neko_log%message(log_buf)
1055 write(log_buf, '(A,ES13.6)') 'mu :', const_mu
1056 call neko_log%message(log_buf)
1057 end if
1058
1059 ! Since mu, rho is a field_t, and we use the %x(1,1,1,1)
1060 ! host array data to pass constant density and viscosity
1061 ! to some routines, we need to make sure that the host
1062 ! values are also filled
1063 if (neko_bcknd_device .eq. 1) then
1064 call device_memcpy(this%rho%x, this%rho%x_d, this%rho%size(), &
1065 device_to_host, sync=.false.)
1066 call device_memcpy(this%mu%x, this%mu%x_d, this%mu%size(), &
1067 device_to_host, sync=.false.)
1068 end if
1069 end subroutine adjoint_fluid_scheme_set_material_properties
1070
1071end module adjoint_fluid_scheme
Abstract interface to initialize an adjoint formulation.
Implements the adjoint_source_term_t type.
Wrapper contaning and executing the adjoint source terms.