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 json_utils, only : json_get, json_get_or_default, json_extract_item
61 use json_module, only : json_file
62 use user_intf, only : user_t, dummy_user_material_properties, &
63 user_material_properties_intf
64 use utils, only : neko_error, neko_warning
65 use comm, only: neko_comm
66 use scalar_source_term, only : scalar_source_term_t
67 use field_series, only : field_series_t
68 use math, only : cfill, add2s2
69 use field_math, only : field_cmult, field_col3, field_cfill, field_add2, &
70 field_col2
71 use device_math, only : device_cfill, device_add2s2
72 use neko_config, only : neko_bcknd_device
73 use field_series, only : field_series_t
74 use time_step_controller, only : time_step_controller_t
75 use json_utils_ext, only: json_key_fallback
76 use scalar_scheme, only: scalar_scheme_precon_factory, &
77 scalar_scheme_solver_factory
78 use scratch_registry, only : neko_scratch_registry
79 use time_state, only : time_state_t
80 use device, only : device_memcpy, device_to_host
81 use field_math, only : field_col3, field_cmult2, field_add2, field_cfill
82 use mpi_f08, only: mpi_integer, mpi_sum
83 implicit none
84
86 type, abstract :: adjoint_scalar_scheme_t
88 character(len=:), allocatable :: name
90 character(len=:), allocatable :: primal_name
92 type(field_t), pointer :: u
94 type(field_t), pointer :: v
96 type(field_t), pointer :: w
98 type(field_t), pointer :: s
100 type(field_t), pointer :: s_adj
102 type(field_series_t) :: s_adj_lag
104 type(space_t), pointer :: xh
106 type(dofmap_t), pointer :: dm_xh
108 type(gs_t), pointer :: gs_xh
110 type(coef_t), pointer :: c_xh
112 type(field_t), pointer :: f_xh => null()
114 type(scalar_source_term_t) :: source_term
116 class(ksp_t), allocatable :: ksp
118 integer :: ksp_maxiter
120 integer :: projection_dim
121
122 integer :: projection_activ_step
124 class(pc_t), allocatable :: pc
126 type(bc_list_t) :: bcs
128 type(json_file), pointer :: params
130 type(mesh_t), pointer :: msh => null()
132 type(chkp_t), pointer :: chkp => null()
134 character(len=:), allocatable :: nut_field_name
136 type(field_t), pointer :: rho => null()
138 type(field_t) :: lambda
140 type(field_t) :: cp
142 real(kind=rp) :: pr_turb
144 type(field_list_t) :: material_properties
146 logical :: variable_material_properties = .false.
147 ! Lag arrays for the RHS.
148 type(field_t) :: abx1, abx2
149 procedure(user_material_properties_intf), nopass, pointer :: &
150 user_material_properties => null()
151 contains
153 procedure, pass(this) :: scheme_init => adjoint_scalar_scheme_init
155 procedure, pass(this) :: scheme_free => adjoint_scalar_scheme_free
157 procedure, pass(this) :: validate => adjoint_scalar_scheme_validate
159 procedure, pass(this) :: set_material_properties => &
160 adjoint_scalar_scheme_set_material_properties
162 procedure, pass(this) :: update_material_properties => &
163 adjoint_scalar_scheme_update_material_properties
165 procedure(adjoint_scalar_scheme_init_intrf), pass(this), deferred :: init
167 procedure(adjoint_scalar_scheme_free_intrf), pass(this), deferred :: free
169 procedure(adjoint_scalar_scheme_step_intrf), pass(this), deferred :: step
171 procedure(adjoint_scalar_scheme_restart_intrf), pass(this), deferred :: &
172 restart
174
176 abstract interface
177 subroutine adjoint_scalar_scheme_init_intrf(this, msh, coef, gs, &
178 params_adjoint, params_primal, numerics_params, user, chkp, ulag, &
179 vlag, wlag, time_scheme, rho)
181 import json_file
182 import coef_t
183 import gs_t
184 import mesh_t
185 import user_t
186 import field_series_t, field_t
187 import time_scheme_controller_t
188 import rp
189 import chkp_t
190 class(adjoint_scalar_scheme_t), target, intent(inout) :: this
191 type(mesh_t), target, intent(in) :: msh
192 type(coef_t), target, intent(in) :: coef
193 type(gs_t), target, intent(inout) :: gs
194 type(json_file), target, intent(inout) :: params_adjoint
195 type(json_file), target, intent(inout) :: params_primal
196 type(json_file), target, intent(inout) :: numerics_params
197 type(user_t), target, intent(in) :: user
198 type(chkp_t), target, intent(inout) :: chkp
199 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
200 type(time_scheme_controller_t), target, intent(in) :: time_scheme
201 type(field_t), target, intent(in) :: rho
203 end interface
204
206 abstract interface
209 import chkp_t
210 import rp
211 class(adjoint_scalar_scheme_t), target, intent(inout) :: this
212 type(chkp_t), intent(inout) :: chkp
214 end interface
215
217 abstract interface
220 class(adjoint_scalar_scheme_t), intent(inout) :: this
222 end interface
223
225 abstract interface
226 subroutine adjoint_scalar_scheme_step_intrf(this, time, ext_bdf, &
227 dt_controller)
229 import time_state_t
230 import time_scheme_controller_t
231 import time_step_controller_t
232 class(adjoint_scalar_scheme_t), intent(inout) :: this
233 type(time_state_t), intent(in) :: time
234 type(time_scheme_controller_t), intent(in) :: ext_bdf
235 type(time_step_controller_t), intent(in) :: dt_controller
237 end interface
238
239contains
240
251 subroutine adjoint_scalar_scheme_init(this, msh, c_Xh, gs_Xh, &
252 params_adjoint, params_primal, scheme, user, rho)
253 class(adjoint_scalar_scheme_t), target, intent(inout) :: this
254 type(mesh_t), target, intent(in) :: msh
255 type(coef_t), target, intent(in) :: c_Xh
256 type(gs_t), target, intent(inout) :: gs_Xh
257 type(json_file), target, intent(inout) :: params_primal, params_adjoint
258 character(len=*), intent(in) :: scheme
259 type(user_t), target, intent(in) :: user
260 type(field_t), target, intent(in) :: rho
261 type(json_file), pointer :: params_selected
262 ! IO buffer for log output
263 character(len=LOG_SIZE) :: log_buf
264 ! Variables for retrieving json parameters
265 logical :: logical_val
266 real(kind=rp) :: solver_abstol
267 integer :: integer_val
268 character(len=:), allocatable :: solver_type, solver_precon
269 type(json_file) :: precon_params
270
271 this%u => neko_field_registry%get_field('u')
272 this%v => neko_field_registry%get_field('v')
273 this%w => neko_field_registry%get_field('w')
274 this%rho => rho
275
276 ! get the primal adjoint's name
277 call json_get(params_adjoint, 'primal_name', this%primal_name)
278 ! Assign a name
279 call json_get_or_default(params_adjoint, 'name', this%name, &
280 'scalar adjoint')
281
282 call neko_log%section('Adjoint scalar')
283 params_selected => json_key_fallback(params_adjoint, params_primal, &
284 'solver.type')
285 call json_get(params_selected, 'solver.type', solver_type)
286
287 params_selected => json_key_fallback(params_adjoint, params_primal, &
288 'solver.preconditioner.type')
289 call json_get(params_selected, 'solver.preconditioner.type', solver_precon)
290
291 params_selected => json_key_fallback(params_adjoint, params_primal, &
292 'solver.preconditioner')
293 call json_get(params_selected, 'solver.preconditioner', &
294 precon_params)
295
296 params_selected => json_key_fallback(params_adjoint, params_primal, &
297 'solver.absolute_tolerance')
298 call json_get(params_selected, 'solver.absolute_tolerance', &
299 solver_abstol)
300
301 params_selected => json_key_fallback(params_adjoint, params_primal, &
302 'solver.projection_space_size')
303 call json_get_or_default(params_selected, &
304 'solver.projection_space_size', &
305 this%projection_dim, 0)
306
307 params_selected => json_key_fallback(params_adjoint, params_primal, &
308 'solver.projection_hold_steps')
309 call json_get_or_default(params_selected, &
310 'solver.projection_hold_steps', &
311 this%projection_activ_step, 5)
312
313
314 write(log_buf, '(A, A)') 'Type : ', trim(scheme)
315 call neko_log%message(log_buf)
316 write(log_buf, '(A, A)') 'Name : ', trim(this%name)
317 call neko_log%message(log_buf)
318 call neko_log%message('Ksp adjoint scalar : ('// trim(solver_type) // &
319 ', ' // trim(solver_precon) // ')')
320 write(log_buf, '(A,ES13.6)') ' `-abs tol :', solver_abstol
321 call neko_log%message(log_buf)
322
323 this%Xh => this%u%Xh
324 this%dm_Xh => this%u%dof
325 this%params => params_adjoint
326 this%msh => msh
327
328 if (.not. neko_field_registry%field_exists(this%name)) then
329 call neko_field_registry%add_field(this%dm_Xh, this%name)
330 end if
331
332 this%s_adj => neko_field_registry%get_field(this%name)
333
334 call this%s_adj_lag%init(this%s_adj, 2)
335
336 this%s => neko_field_registry%get_field(this%primal_name)
337
338 this%gs_Xh => gs_xh
339 this%c_Xh => c_xh
340
341 !
342 ! Material properties
343 !
344 call this%set_material_properties(params_primal, user)
345
346
347 !
348 ! Turbulence modelling and variable material properties
349 !
350 params_selected => json_key_fallback(params_adjoint, params_primal, &
351 'variable_material_properties')
352 if (params_selected%valid_path('variable_material_properties')) then
353 call neko_error('variable material properties no supported for adjoint')
354 end if
355
356 write(log_buf, '(A,L1)') 'LES : ', this%variable_material_properties
357 call neko_log%message(log_buf)
358
359 !
360 ! Setup right-hand side field.
361 !
362 allocate(this%f_Xh)
363 call this%f_Xh%init(this%dm_Xh, fld_name = "adjoint_scalar_rhs")
364
365 ! Initialize the source term
366 call this%source_term%init(this%f_Xh, this%c_Xh, user, this%name)
367 ! We should ONLY read the adjoint
368 call this%source_term%add(params_primal, 'source_terms')
369
370 ! todo parameter file ksp tol should be added
371 params_selected => json_key_fallback(params_adjoint, params_primal, &
372 'solver.max_iterations')
373 call json_get_or_default(params_selected, &
374 'solver.max_iterations', &
375 integer_val, ksp_max_iter)
376 params_selected => json_key_fallback(params_adjoint, params_primal, &
377 'solver.monitor')
378 call json_get_or_default(params_selected, &
379 'solver.monitor', &
380 logical_val, .false.)
381 call adjoint_scalar_scheme_solver_factory(this%ksp, this%dm_Xh%size(), &
382 solver_type, integer_val, solver_abstol, logical_val)
383 call scalar_scheme_precon_factory(this%pc, this%ksp, &
384 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs, &
385 solver_precon, precon_params)
386
387 call neko_log%end_section()
388
389 end subroutine adjoint_scalar_scheme_init
390
391
393 subroutine adjoint_scalar_scheme_free(this)
394 class(adjoint_scalar_scheme_t), intent(inout) :: this
395
396 nullify(this%Xh)
397 nullify(this%dm_Xh)
398 nullify(this%gs_Xh)
399 nullify(this%c_Xh)
400 nullify(this%params)
401
402 if (allocated(this%ksp)) then
403 call this%ksp%free()
404 deallocate(this%ksp)
405 end if
406
407 if (allocated(this%pc)) then
408 call precon_destroy(this%pc)
409 deallocate(this%pc)
410 end if
411
412 call this%source_term%free()
413
414 call this%bcs%free()
415
416 call this%cp%free()
417 call this%lambda%free()
418 call this%s_adj_lag%free()
419
420 end subroutine adjoint_scalar_scheme_free
421
424 subroutine adjoint_scalar_scheme_validate(this)
425 class(adjoint_scalar_scheme_t), target, intent(inout) :: this
426
427 if ( (.not. allocated(this%u%x)) .or. &
428 (.not. allocated(this%v%x)) .or. &
429 (.not. allocated(this%w%x)) .or. &
430 (.not. allocated(this%s%x)) .or. &
431 (.not. allocated(this%s_adj%x))) then
432 call neko_error('Fields are not allocated')
433 end if
434
435 if (.not. allocated(this%ksp)) then
436 call neko_error('No Krylov solver for velocity defined')
437 end if
438
439 if (.not. associated(this%Xh)) then
440 call neko_error('No function space defined')
441 end if
442
443 if (.not. associated(this%dm_Xh)) then
444 call neko_error('No dofmap defined')
445 end if
446
447 if (.not. associated(this%c_Xh)) then
448 call neko_error('No coefficients defined')
449 end if
450
451 if (.not. associated(this%f_Xh)) then
452 call neko_error('No rhs allocated')
453 end if
454
455 if (.not. associated(this%params)) then
456 call neko_error('No parameters defined')
457 end if
458
459 if (.not. associated(this%rho)) then
460 call neko_error('No density field defined')
461 end if
462
463 end subroutine adjoint_scalar_scheme_validate
464
467 subroutine adjoint_scalar_scheme_solver_factory(ksp, n, solver, max_iter, &
468 abstol, monitor)
469 class(ksp_t), allocatable, target, intent(inout) :: ksp
470 integer, intent(in), value :: n
471 integer, intent(in) :: max_iter
472 character(len=*), intent(in) :: solver
473 real(kind=rp) :: abstol
474 logical, intent(in) :: monitor
475
476 call krylov_solver_factory(ksp, n, solver, max_iter, &
477 abstol, monitor = monitor)
478
479 end subroutine adjoint_scalar_scheme_solver_factory
480
482 subroutine adjoint_scalar_scheme_precon_factory(pc, ksp, coef, dof, gs, &
483 bclst, pctype, pcparams)
484 class(pc_t), allocatable, target, intent(inout) :: pc
485 class(ksp_t), target, intent(inout) :: ksp
486 type(coef_t), target, intent(in) :: coef
487 type(dofmap_t), target, intent(in) :: dof
488 type(gs_t), target, intent(inout) :: gs
489 type(bc_list_t), target, intent(inout) :: bclst
490 character(len=*) :: pctype
491 type(json_file), intent(inout) :: pcparams
492
493 call precon_factory(pc, pctype)
494
495 select type (pcp => pc)
496 type is (jacobi_t)
497 call pcp%init(coef, dof, gs)
498 type is (sx_jacobi_t)
499 call pcp%init(coef, dof, gs)
500 type is (device_jacobi_t)
501 call pcp%init(coef, dof, gs)
502 type is (hsmg_t)
503 call pcp%init(coef, bclst, pcparams)
504 end select
505
506 call ksp%set_pc(pc)
507
508 end subroutine adjoint_scalar_scheme_precon_factory
509
514 subroutine adjoint_scalar_scheme_update_material_properties(this, time)
515 class(adjoint_scalar_scheme_t), intent(inout) :: this
516 type(time_state_t), intent(in) :: time
517 type(field_t), pointer :: nut
518 integer :: index
519 ! Factor to transform nu_t to lambda_t
520 type(field_t), pointer :: lambda_factor
521
522 call this%user_material_properties(this%name, this%material_properties, &
523 time)
524
525 ! factor = rho * cp / pr_turb
526 if (this%variable_material_properties .and. &
527 len(trim(this%nut_field_name)) > 0) then
528 nut => neko_field_registry%get_field(this%nut_field_name)
529
530 ! lambda = lambda + rho * cp * nut / pr_turb
531 call neko_scratch_registry%request_field(lambda_factor, index)
532
533 call field_col3(lambda_factor, this%cp, this%rho)
534 call field_col2(lambda_factor, nut)
535 call field_cmult(lambda_factor, 1.0_rp / this%pr_turb)
536 call field_add2(this%lambda, lambda_factor)
537 call neko_scratch_registry%relinquish_field(index)
538 end if
539
540 ! Since cp is a fields and we use the %x(1,1,1,1) of the
541 ! host array data to pass constant material properties
542 ! to some routines, we need to make sure that the host
543 ! values are also filled
544 if (neko_bcknd_device .eq. 1) then
545 call device_memcpy(this%cp%x, this%cp%x_d, this%cp%size(), &
546 device_to_host, sync=.false.)
547 end if
548
549 end subroutine adjoint_scalar_scheme_update_material_properties
550
555 subroutine adjoint_scalar_scheme_set_material_properties(this, &
556 params_primal, user)
557 class(adjoint_scalar_scheme_t), intent(inout) :: this
558 type(json_file), intent(inout) :: params_primal
559 type(user_t), target, intent(in) :: user
560 character(len=LOG_SIZE) :: log_buf
561 ! A local pointer that is needed to make Intel happy
562 procedure(user_material_properties_intf), pointer :: dummy_mp_ptr
563 real(kind=rp) :: const_cp, const_lambda
564 ! Dummy time state set to 0
565 type(time_state_t) :: time
566
567 dummy_mp_ptr => dummy_user_material_properties
568
569 ! Fill lambda field with the physical value
570 call this%lambda%init(this%dm_Xh, "lambda")
571 call this%cp%init(this%dm_Xh, "cp")
572
573 call this%material_properties%init(2)
574 call this%material_properties%assign_to_field(1, this%cp)
575 call this%material_properties%assign_to_field(2, this%lambda)
576
577 if (.not. associated(user%material_properties, dummy_mp_ptr)) then
578
579 write(log_buf, '(A)') "Material properties must be set in the user " // &
580 "file!"
581 call neko_log%message(log_buf)
582 this%user_material_properties => user%material_properties
583 call user%material_properties(this%name, this%material_properties, time)
584 else
585 this%user_material_properties => dummy_user_material_properties
586 if (params_primal%valid_path('Pe') .and. &
587 (params_primal%valid_path('lambda') .or. &
588 params_primal%valid_path('cp'))) then
589 call neko_error("To set the material properties for the scalar, " // &
590 "either provide Pe OR lambda and cp in the case file.")
591 ! Non-dimensional case
592 else if (params_primal%valid_path('Pe')) then
593 write(log_buf, '(A)') 'Non-dimensional scalar material properties' //&
594 ' input.'
595 call neko_log%message(log_buf, lvl = neko_log_verbose)
596 write(log_buf, '(A)') 'Specific heat capacity will be set to 1,'
597 call neko_log%message(log_buf, lvl = neko_log_verbose)
598 write(log_buf, '(A)') 'conductivity to 1/Pe. Assumes density is 1.'
599 call neko_log%message(log_buf, lvl = neko_log_verbose)
600
601 ! Read Pe into lambda for further manipulation.
602 call json_get(params_primal, 'Pe', const_lambda)
603 write(log_buf, '(A,ES13.6)') 'Pe :', const_lambda
604 call neko_log%message(log_buf)
605
606 ! Set cp and rho to 1 since the setup is non-dimensional.
607 const_cp = 1.0_rp
608 ! Invert the Pe to get conductivity
609 const_lambda = 1.0_rp/const_lambda
610 ! Dimensional case
611 else
612 call json_get(params_primal, 'lambda', const_lambda)
613 call json_get(params_primal, 'cp', const_cp)
614 end if
615 end if
616 ! We need to fill the fields based on the parsed const values
617 ! if the user routine is not used.
618 if (associated(user%material_properties, dummy_mp_ptr)) then
619 ! Fill mu and rho field with the physical value
620 call field_cfill(this%lambda, const_lambda)
621 call field_cfill(this%cp, const_cp)
622
623 write(log_buf, '(A,ES13.6)') 'lambda :', const_lambda
624 call neko_log%message(log_buf)
625 write(log_buf, '(A,ES13.6)') 'cp :', const_cp
626 call neko_log%message(log_buf)
627 end if
628
629 ! Since cp is a field and we use the %x(1,1,1,1) of the
630 ! host array data to pass constant material properties
631 ! to some routines, we need to make sure that the host
632 ! values are also filled
633 if (neko_bcknd_device .eq. 1) then
634 call device_memcpy(this%cp%x, this%cp%x_d, this%cp%size(), &
635 device_to_host, sync=.false.)
636 end if
637 end subroutine adjoint_scalar_scheme_set_material_properties
638
639end 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.