Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_scalar_scheme.f90
1! Copyright (c) 2022-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!
34
36 use gather_scatter, only : gs_t
37 use checkpoint, only : chkp_t
38 use num_types, only: rp
39 use field, only : field_t
40 use field_list, only: field_list_t
41 use space, only : space_t
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 dirichlet, only : dirichlet_t
46 use neumann, only : neumann_t
47 use jacobi, only : jacobi_t
48 use device_jacobi, only : device_jacobi_t
49 use sx_jacobi, only : sx_jacobi_t
50 use hsmg, only : hsmg_t
51 use bc, only : bc_t
52 use bc_list, only : bc_list_t
53 use precon, only : pc_t, precon_factory, precon_destroy
54 use field_dirichlet, only: field_dirichlet_t, field_dirichlet_update
55 use mesh, only : mesh_t, neko_msh_max_zlbls, neko_msh_max_zlbl_len
56 use facet_zone, only : facet_zone_t
57 use time_scheme_controller, only : time_scheme_controller_t
58 use logger, only : neko_log, log_size, neko_log_verbose
59 use field_registry, only : neko_field_registry
60 use usr_scalar, only : usr_scalar_t, usr_scalar_bc_eval
61 use json_utils, only : json_get, json_get_or_default, json_extract_item, &
62 json_extract_object
63 use json_module, only : json_file
64 use user_intf, only : user_t, dummy_user_material_properties, &
65 user_material_properties
66 use utils, only : neko_error, neko_warning
67 use comm, only: neko_comm, mpi_integer, mpi_sum
68 use scalar_source_term, only : scalar_source_term_t
69 use field_series, only : field_series_t
70 use math, only : cfill, add2s2
71 use field_math, only : field_cmult, field_col3, field_cfill, field_add2, &
72 field_col2
73 use device_math, only : device_cfill, device_add2s2
74 use neko_config, only : neko_bcknd_device
75 use field_series, only : field_series_t
76 use time_step_controller, only : time_step_controller_t
77 use json_utils_ext, only: json_key_fallback
78 use scalar_scheme, only: scalar_scheme_precon_factory, &
79 scalar_scheme_solver_factory
80 use scratch_registry, only : neko_scratch_registry
81 use time_state, only : time_state_t
82 use device, only : device_memcpy, device_to_host
83 use field_math, only : field_col3, field_cmult2, field_add2, field_cfill
84 implicit none
85
87 type, abstract :: adjoint_scalar_scheme_t
89 character(len=:), allocatable :: name
91 character(len=:), allocatable :: primal_name
93 type(field_t), pointer :: u
95 type(field_t), pointer :: v
97 type(field_t), pointer :: w
99 type(field_t), pointer :: s
101 type(field_t), pointer :: s_adj
103 type(field_series_t) :: s_adj_lag
105 type(space_t), pointer :: xh
107 type(dofmap_t), pointer :: dm_xh
109 type(gs_t), pointer :: gs_xh
111 type(coef_t), pointer :: c_xh
113 type(field_t), pointer :: f_xh => null()
115 type(scalar_source_term_t) :: source_term
117 class(ksp_t), allocatable :: ksp
119 integer :: ksp_maxiter
121 integer :: projection_dim
122
123 integer :: projection_activ_step
125 class(pc_t), allocatable :: pc
127 type(bc_list_t) :: bcs
129 type(json_file), pointer :: params
131 type(mesh_t), pointer :: msh => null()
133 type(chkp_t), pointer :: chkp => null()
135 character(len=:), allocatable :: nut_field_name
137 type(field_t), pointer :: rho => null()
139 type(field_t) :: lambda
141 type(field_t) :: cp
143 real(kind=rp) :: pr_turb
145 type(field_list_t) :: material_properties
147 logical :: variable_material_properties = .false.
148 ! Lag arrays for the RHS.
149 type(field_t) :: abx1, abx2
150 procedure(user_material_properties), nopass, pointer :: &
151 user_material_properties => null()
152 contains
154 procedure, pass(this) :: scheme_init => adjoint_scalar_scheme_init
156 procedure, pass(this) :: scheme_free => adjoint_scalar_scheme_free
158 procedure, pass(this) :: validate => adjoint_scalar_scheme_validate
160 procedure, pass(this) :: set_material_properties => &
161 adjoint_scalar_scheme_set_material_properties
163 procedure, pass(this) :: update_material_properties => &
164 adjoint_scalar_scheme_update_material_properties
166 procedure(adjoint_scalar_scheme_init_intrf), pass(this), deferred :: init
168 procedure(adjoint_scalar_scheme_free_intrf), pass(this), deferred :: free
170 procedure(adjoint_scalar_scheme_step_intrf), pass(this), deferred :: step
172 procedure(adjoint_scalar_scheme_restart_intrf), pass(this), deferred :: &
173 restart
175
177 abstract interface
178 subroutine adjoint_scalar_scheme_init_intrf(this, msh, coef, gs, &
179 params_adjoint, params_primal, numerics_params, user, chkp, ulag, &
180 vlag, wlag, time_scheme, rho)
182 import json_file
183 import coef_t
184 import gs_t
185 import mesh_t
186 import user_t
187 import field_series_t, field_t
188 import time_scheme_controller_t
189 import rp
190 import chkp_t
191 class(adjoint_scalar_scheme_t), target, intent(inout) :: this
192 type(mesh_t), target, intent(in) :: msh
193 type(coef_t), target, intent(in) :: coef
194 type(gs_t), target, intent(inout) :: gs
195 type(json_file), target, intent(inout) :: params_adjoint
196 type(json_file), target, intent(inout) :: params_primal
197 type(json_file), target, intent(inout) :: numerics_params
198 type(user_t), target, intent(in) :: user
199 type(chkp_t), target, intent(inout) :: chkp
200 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
201 type(time_scheme_controller_t), target, intent(in) :: time_scheme
202 type(field_t), target, intent(in) :: rho
204 end interface
205
207 abstract interface
210 import chkp_t
211 import rp
212 class(adjoint_scalar_scheme_t), target, intent(inout) :: this
213 type(chkp_t), intent(inout) :: chkp
215 end interface
216
218 abstract interface
221 class(adjoint_scalar_scheme_t), intent(inout) :: this
223 end interface
224
226 abstract interface
227 subroutine adjoint_scalar_scheme_step_intrf(this, time, ext_bdf, &
228 dt_controller)
230 import time_state_t
231 import time_scheme_controller_t
232 import time_step_controller_t
233 class(adjoint_scalar_scheme_t), intent(inout) :: this
234 type(time_state_t), intent(in) :: time
235 type(time_scheme_controller_t), intent(in) :: ext_bdf
236 type(time_step_controller_t), intent(in) :: dt_controller
238 end interface
239
240contains
241
252 subroutine adjoint_scalar_scheme_init(this, msh, c_Xh, gs_Xh, &
253 params_adjoint, params_primal, scheme, user, rho)
254 class(adjoint_scalar_scheme_t), target, intent(inout) :: this
255 type(mesh_t), target, intent(in) :: msh
256 type(coef_t), target, intent(in) :: c_Xh
257 type(gs_t), target, intent(inout) :: gs_Xh
258 type(json_file), target, intent(inout) :: params_primal, params_adjoint
259 character(len=*), intent(in) :: scheme
260 type(user_t), target, intent(in) :: user
261 type(field_t), target, intent(in) :: rho
262 type(json_file), pointer :: params_selected
263 ! IO buffer for log output
264 character(len=LOG_SIZE) :: log_buf
265 ! Variables for retrieving json parameters
266 logical :: logical_val
267 real(kind=rp) :: solver_abstol
268 integer :: integer_val
269 character(len=:), allocatable :: solver_type, solver_precon
270 type(json_file) :: precon_params
271
272 this%u => neko_field_registry%get_field('u')
273 this%v => neko_field_registry%get_field('v')
274 this%w => neko_field_registry%get_field('w')
275 this%rho => rho
276
277 ! get the primal adjoint's name
278 call json_get(params_adjoint, 'primal_name', this%primal_name)
279 ! Assign a name
280 call json_get_or_default(params_adjoint, 'name', this%name, &
281 'scalar adjoint')
282
283 call neko_log%section('Adjoint scalar')
284 params_selected => json_key_fallback(params_adjoint, params_primal, &
285 'solver.type')
286 call json_get(params_selected, 'solver.type', solver_type)
287
288 params_selected => json_key_fallback(params_adjoint, params_primal, &
289 'solver.preconditioner.type')
290 call json_get(params_selected, 'solver.preconditioner.type', solver_precon)
291
292 params_selected => json_key_fallback(params_adjoint, params_primal, &
293 'solver.preconditioner')
294 call json_extract_object(params_selected, 'solver.preconditioner', &
295 precon_params)
296
297 params_selected => json_key_fallback(params_adjoint, params_primal, &
298 'solver.absolute_tolerance')
299 call json_get(params_selected, 'solver.absolute_tolerance', &
300 solver_abstol)
301
302 params_selected => json_key_fallback(params_adjoint, params_primal, &
303 'solver.projection_space_size')
304 call json_get_or_default(params_selected, &
305 'solver.projection_space_size', &
306 this%projection_dim, 0)
307
308 params_selected => json_key_fallback(params_adjoint, params_primal, &
309 'solver.projection_hold_steps')
310 call json_get_or_default(params_selected, &
311 'solver.projection_hold_steps', &
312 this%projection_activ_step, 5)
313
314
315 write(log_buf, '(A, A)') 'Type : ', trim(scheme)
316 call neko_log%message(log_buf)
317 write(log_buf, '(A, A)') 'Name : ', trim(this%name)
318 call neko_log%message(log_buf)
319 call neko_log%message('Ksp adjoint scalar : ('// trim(solver_type) // &
320 ', ' // trim(solver_precon) // ')')
321 write(log_buf, '(A,ES13.6)') ' `-abs tol :', solver_abstol
322 call neko_log%message(log_buf)
323
324 this%Xh => this%u%Xh
325 this%dm_Xh => this%u%dof
326 this%params => params_adjoint
327 this%msh => msh
328
329 if (.not. neko_field_registry%field_exists(this%name)) then
330 call neko_field_registry%add_field(this%dm_Xh, this%name)
331 end if
332
333 this%s_adj => neko_field_registry%get_field(this%name)
334
335 call this%s_adj_lag%init(this%s_adj, 2)
336
337 this%s => neko_field_registry%get_field(this%primal_name)
338
339 this%gs_Xh => gs_xh
340 this%c_Xh => c_xh
341
342 !
343 ! Material properties
344 !
345 call this%set_material_properties(params_primal, user)
346
347
348 !
349 ! Turbulence modelling and variable material properties
350 !
351 params_selected => json_key_fallback(params_adjoint, params_primal, &
352 'variable_material_properties')
353 if (params_selected%valid_path('variable_material_properties')) then
354 call neko_error('variable material properties no supported for adjoint')
355 end if
356
357 write(log_buf, '(A,L1)') 'LES : ', this%variable_material_properties
358 call neko_log%message(log_buf)
359
360 !
361 ! Setup right-hand side field.
362 !
363 allocate(this%f_Xh)
364 call this%f_Xh%init(this%dm_Xh, fld_name = "adjoint_scalar_rhs")
365
366 ! Initialize the source term
367 call this%source_term%init(this%f_Xh, this%c_Xh, user, this%name)
368 ! We should ONLY read the adjoint
369 call this%source_term%add(params_primal, 'source_terms')
370
371 ! todo parameter file ksp tol should be added
372 params_selected => json_key_fallback(params_adjoint, params_primal, &
373 'solver.max_iterations')
374 call json_get_or_default(params_selected, &
375 'solver.max_iterations', &
376 integer_val, ksp_max_iter)
377 params_selected => json_key_fallback(params_adjoint, params_primal, &
378 'solver.monitor')
379 call json_get_or_default(params_selected, &
380 'solver.monitor', &
381 logical_val, .false.)
382 call adjoint_scalar_scheme_solver_factory(this%ksp, this%dm_Xh%size(), &
383 solver_type, integer_val, solver_abstol, logical_val)
384 call scalar_scheme_precon_factory(this%pc, this%ksp, &
385 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs, &
386 solver_precon, precon_params)
387
388 call neko_log%end_section()
389
390 end subroutine adjoint_scalar_scheme_init
391
392
394 subroutine adjoint_scalar_scheme_free(this)
395 class(adjoint_scalar_scheme_t), intent(inout) :: this
396
397 nullify(this%Xh)
398 nullify(this%dm_Xh)
399 nullify(this%gs_Xh)
400 nullify(this%c_Xh)
401 nullify(this%params)
402
403 if (allocated(this%ksp)) then
404 call this%ksp%free()
405 deallocate(this%ksp)
406 end if
407
408 if (allocated(this%pc)) then
409 call precon_destroy(this%pc)
410 deallocate(this%pc)
411 end if
412
413 call this%source_term%free()
414
415 call this%bcs%free()
416
417 call this%cp%free()
418 call this%lambda%free()
419 call this%s_adj_lag%free()
420
421 end subroutine adjoint_scalar_scheme_free
422
425 subroutine adjoint_scalar_scheme_validate(this)
426 class(adjoint_scalar_scheme_t), target, intent(inout) :: this
427
428 if ( (.not. allocated(this%u%x)) .or. &
429 (.not. allocated(this%v%x)) .or. &
430 (.not. allocated(this%w%x)) .or. &
431 (.not. allocated(this%s%x)) .or. &
432 (.not. allocated(this%s_adj%x))) then
433 call neko_error('Fields are not allocated')
434 end if
435
436 if (.not. allocated(this%ksp)) then
437 call neko_error('No Krylov solver for velocity defined')
438 end if
439
440 if (.not. associated(this%Xh)) then
441 call neko_error('No function space defined')
442 end if
443
444 if (.not. associated(this%dm_Xh)) then
445 call neko_error('No dofmap defined')
446 end if
447
448 if (.not. associated(this%c_Xh)) then
449 call neko_error('No coefficients defined')
450 end if
451
452 if (.not. associated(this%f_Xh)) then
453 call neko_error('No rhs allocated')
454 end if
455
456 if (.not. associated(this%params)) then
457 call neko_error('No parameters defined')
458 end if
459
460 if (.not. associated(this%rho)) then
461 call neko_error('No density field defined')
462 end if
463
464 end subroutine adjoint_scalar_scheme_validate
465
468 subroutine adjoint_scalar_scheme_solver_factory(ksp, n, solver, max_iter, &
469 abstol, monitor)
470 class(ksp_t), allocatable, target, intent(inout) :: ksp
471 integer, intent(in), value :: n
472 integer, intent(in) :: max_iter
473 character(len=*), intent(in) :: solver
474 real(kind=rp) :: abstol
475 logical, intent(in) :: monitor
476
477 call krylov_solver_factory(ksp, n, solver, max_iter, &
478 abstol, monitor = monitor)
479
480 end subroutine adjoint_scalar_scheme_solver_factory
481
483 subroutine adjoint_scalar_scheme_precon_factory(pc, ksp, coef, dof, gs, &
484 bclst, pctype, pcparams)
485 class(pc_t), allocatable, target, intent(inout) :: pc
486 class(ksp_t), target, intent(inout) :: ksp
487 type(coef_t), target, intent(in) :: coef
488 type(dofmap_t), target, intent(in) :: dof
489 type(gs_t), target, intent(inout) :: gs
490 type(bc_list_t), target, intent(inout) :: bclst
491 character(len=*) :: pctype
492 type(json_file), intent(inout) :: pcparams
493
494 call precon_factory(pc, pctype)
495
496 select type (pcp => pc)
497 type is (jacobi_t)
498 call pcp%init(coef, dof, gs)
499 type is (sx_jacobi_t)
500 call pcp%init(coef, dof, gs)
501 type is (device_jacobi_t)
502 call pcp%init(coef, dof, gs)
503 type is (hsmg_t)
504 call pcp%init(coef, bclst, pcparams)
505 end select
506
507 call ksp%set_pc(pc)
508
509 end subroutine adjoint_scalar_scheme_precon_factory
510
516 subroutine adjoint_scalar_scheme_update_material_properties(t, tstep, this)
517 class(adjoint_scalar_scheme_t), intent(inout) :: this
518 real(kind=rp),intent(in) :: t
519 integer, intent(in) :: tstep
520 type(field_t), pointer :: nut
521 integer :: index
522 ! Factor to transform nu_t to lambda_t
523 type(field_t), pointer :: lambda_factor
524
525 call this%user_material_properties(t, tstep, this%name, &
526 this%material_properties)
527
528 ! factor = rho * cp / pr_turb
529 if (this%variable_material_properties .and. &
530 len(trim(this%nut_field_name)) > 0) then
531 nut => neko_field_registry%get_field(this%nut_field_name)
532
533 ! lambda = lambda + rho * cp * nut / pr_turb
534 call neko_scratch_registry%request_field(lambda_factor, index)
535
536 call field_col3(lambda_factor, this%cp, this%rho)
537 call field_col2(lambda_factor, nut)
538 call field_cmult(lambda_factor, 1.0_rp / this%pr_turb)
539 call field_add2(this%lambda, lambda_factor)
540 call neko_scratch_registry%relinquish_field(index)
541 end if
542
543 ! Since cp is a fields and we use the %x(1,1,1,1) of the
544 ! host array data to pass constant material properties
545 ! to some routines, we need to make sure that the host
546 ! values are also filled
547 if (neko_bcknd_device .eq. 1) then
548 call device_memcpy(this%cp%x, this%cp%x_d, this%cp%size(), &
549 device_to_host, sync=.false.)
550 end if
551
552 end subroutine adjoint_scalar_scheme_update_material_properties
553
558 subroutine adjoint_scalar_scheme_set_material_properties(this, &
559 params_primal, user)
560 class(adjoint_scalar_scheme_t), intent(inout) :: this
561 type(json_file), intent(inout) :: params_primal
562 type(user_t), target, intent(in) :: user
563 character(len=LOG_SIZE) :: log_buf
564 ! A local pointer that is needed to make Intel happy
565 procedure(user_material_properties), pointer :: dummy_mp_ptr
566 real(kind=rp) :: const_cp, const_lambda
567
568 dummy_mp_ptr => dummy_user_material_properties
569
570 ! Fill lambda field with the physical value
571 call this%lambda%init(this%dm_Xh, "lambda")
572 call this%cp%init(this%dm_Xh, "cp")
573
574 call this%material_properties%init(2)
575 call this%material_properties%assign_to_field(1, this%cp)
576 call this%material_properties%assign_to_field(2, this%lambda)
577
578 if (.not. associated(user%material_properties, dummy_mp_ptr)) then
579
580 write(log_buf, '(A)') "Material properties must be set in the user " // &
581 "file!"
582 call neko_log%message(log_buf)
583 this%user_material_properties => user%material_properties
584 call user%material_properties(0.0_rp, 0, this%name, &
585 this%material_properties)
586 else
587 this%user_material_properties => dummy_user_material_properties
588 if (params_primal%valid_path('Pe') .and. &
589 (params_primal%valid_path('lambda') .or. &
590 params_primal%valid_path('cp'))) then
591 call neko_error("To set the material properties for the scalar, " // &
592 "either provide Pe OR lambda and cp in the case file.")
593 ! Non-dimensional case
594 else if (params_primal%valid_path('Pe')) then
595 write(log_buf, '(A)') 'Non-dimensional scalar material properties' //&
596 ' input.'
597 call neko_log%message(log_buf, lvl = neko_log_verbose)
598 write(log_buf, '(A)') 'Specific heat capacity will be set to 1,'
599 call neko_log%message(log_buf, lvl = neko_log_verbose)
600 write(log_buf, '(A)') 'conductivity to 1/Pe. Assumes density is 1.'
601 call neko_log%message(log_buf, lvl = neko_log_verbose)
602
603 ! Read Pe into lambda for further manipulation.
604 call json_get(params_primal, 'Pe', const_lambda)
605 write(log_buf, '(A,ES13.6)') 'Pe :', const_lambda
606 call neko_log%message(log_buf)
607
608 ! Set cp and rho to 1 since the setup is non-dimensional.
609 const_cp = 1.0_rp
610 ! Invert the Pe to get conductivity
611 const_lambda = 1.0_rp/const_lambda
612 ! Dimensional case
613 else
614 call json_get(params_primal, 'lambda', const_lambda)
615 call json_get(params_primal, 'cp', const_cp)
616 end if
617 end if
618 ! We need to fill the fields based on the parsed const values
619 ! if the user routine is not used.
620 if (associated(user%material_properties, dummy_mp_ptr)) then
621 ! Fill mu and rho field with the physical value
622 call field_cfill(this%lambda, const_lambda)
623 call field_cfill(this%cp, const_cp)
624
625 write(log_buf, '(A,ES13.6)') 'lambda :', const_lambda
626 call neko_log%message(log_buf)
627 write(log_buf, '(A,ES13.6)') 'cp :', const_cp
628 call neko_log%message(log_buf)
629 end if
630
631 ! Since cp is a field and we use the %x(1,1,1,1) of the
632 ! host array data to pass constant material properties
633 ! to some routines, we need to make sure that the host
634 ! values are also filled
635 if (neko_bcknd_device .eq. 1) then
636 call device_memcpy(this%cp%x, this%cp%x_d, this%cp%size(), &
637 device_to_host, sync=.false.)
638 end if
639 end subroutine adjoint_scalar_scheme_set_material_properties
640
641end module adjoint_scalar_scheme
Abstract interface to dealocate a scalar formulation.
Abstract interface to initialize a scalar formulation.
Abstract interface to restart a scalar formulation.
Contains the adjoint_scalar_scheme_t type.
Base type for a scalar advection-diffusion solver.