Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_scalar_scheme.f90
Go to the documentation of this file.
1
34!
36
38 use gather_scatter, only : gs_t
39 use checkpoint, only : chkp_t
40 use num_types, only: rp
41 use field, only : field_t
42 use field_list, only: field_list_t
43 use space, only : space_t
44 use dofmap, only : dofmap_t
45 use krylov, only : ksp_t, krylov_solver_factory, ksp_max_iter, ksp_monitor_t
46 use coefs, only : coef_t
47 use dirichlet, only : dirichlet_t
48 use neumann, only : neumann_t
49 use jacobi, only : jacobi_t
50 use device_jacobi, only : device_jacobi_t
51 use sx_jacobi, only : sx_jacobi_t
52 use hsmg, only : hsmg_t
53 use bc, only : bc_t
54 use bc_list, only : bc_list_t
55 use precon, only : pc_t, precon_factory, precon_destroy
56 use field_dirichlet, only: field_dirichlet_t, field_dirichlet_update
57 use mesh, only : mesh_t, neko_msh_max_zlbls, neko_msh_max_zlbl_len
58 use facet_zone, only : facet_zone_t
59 use time_scheme_controller, only : time_scheme_controller_t
60 use logger, only : neko_log, log_size, neko_log_verbose
61 use registry, only : neko_registry
62 use json_utils, only : json_get, json_get_or_default, json_extract_item
63 use json_module, only : json_file
64 use user_intf, only : user_t, dummy_user_material_properties, &
65 user_material_properties_intf
66 use utils, only : neko_error, neko_warning
67 use comm, only: neko_comm
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 use mpi_f08, only: mpi_integer, mpi_sum
85 implicit none
86
88 type, abstract :: adjoint_scalar_scheme_t
90 character(len=:), allocatable :: name
92 character(len=:), allocatable :: primal_name
94 type(field_t), pointer :: u
96 type(field_t), pointer :: v
98 type(field_t), pointer :: w
100 type(field_t), pointer :: s
102 type(field_t), pointer :: s_adj
104 type(field_series_t) :: s_adj_lag
106 type(space_t), pointer :: xh
108 type(dofmap_t), pointer :: dm_xh
110 type(gs_t), pointer :: gs_xh
112 type(coef_t), pointer :: c_xh
114 type(field_t), pointer :: f_xh => null()
116 type(scalar_source_term_t) :: source_term
118 class(ksp_t), allocatable :: ksp
120 integer :: ksp_maxiter
122 integer :: projection_dim
123
124 integer :: projection_activ_step
126 class(pc_t), allocatable :: pc
128 type(bc_list_t) :: bcs
130 type(json_file), pointer :: params
132 type(mesh_t), pointer :: msh => null()
134 type(chkp_t), pointer :: chkp => null()
136 character(len=:), allocatable :: nut_field_name
138 type(field_t), pointer :: rho => null()
140 type(field_t) :: lambda
142 type(field_t) :: cp
144 real(kind=rp) :: pr_turb
146 type(field_list_t) :: material_properties
148 logical :: variable_material_properties = .false.
149 ! Lag arrays for the RHS.
150 type(field_t) :: abx1, abx2
151 procedure(user_material_properties_intf), nopass, pointer :: &
152 user_material_properties => null()
154 logical :: freeze = .false.
155 contains
157 procedure, pass(this) :: scheme_init => adjoint_scalar_scheme_init
159 procedure, pass(this) :: scheme_free => adjoint_scalar_scheme_free
161 procedure, pass(this) :: validate => adjoint_scalar_scheme_validate
163 procedure, pass(this) :: set_material_properties => &
166 procedure, pass(this) :: update_material_properties => &
169 procedure(adjoint_scalar_scheme_init_intrf), pass(this), deferred :: init
171 procedure(adjoint_scalar_scheme_free_intrf), pass(this), deferred :: free
173 procedure(adjoint_scalar_scheme_step_intrf), pass(this), deferred :: step
175 procedure(adjoint_scalar_scheme_restart_intrf), pass(this), deferred :: &
176 restart
178
180 abstract interface
181 subroutine adjoint_scalar_scheme_init_intrf(this, msh, coef, gs, &
182 params_adjoint, params_primal, numerics_params, user, chkp, ulag, &
183 vlag, wlag, time_scheme, rho)
185 import json_file
186 import coef_t
187 import gs_t
188 import mesh_t
189 import user_t
190 import field_series_t, field_t
191 import time_scheme_controller_t
192 import rp
193 import chkp_t
194 class(adjoint_scalar_scheme_t), target, intent(inout) :: this
195 type(mesh_t), target, intent(in) :: msh
196 type(coef_t), target, intent(in) :: coef
197 type(gs_t), target, intent(inout) :: gs
198 type(json_file), target, intent(inout) :: params_adjoint
199 type(json_file), target, intent(inout) :: params_primal
200 type(json_file), target, intent(inout) :: numerics_params
201 type(user_t), target, intent(in) :: user
202 type(chkp_t), target, intent(inout) :: chkp
203 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
204 type(time_scheme_controller_t), target, intent(in) :: time_scheme
205 type(field_t), target, intent(in) :: rho
207 end interface
208
210 abstract interface
213 import chkp_t
214 import rp
215 class(adjoint_scalar_scheme_t), target, intent(inout) :: this
216 type(chkp_t), intent(inout) :: chkp
218 end interface
219
221 abstract interface
224 class(adjoint_scalar_scheme_t), intent(inout) :: this
226 end interface
227
229 abstract interface
230 subroutine adjoint_scalar_scheme_step_intrf(this, time, ext_bdf, &
231 dt_controller, ksp_results)
233 import time_state_t
234 import time_scheme_controller_t
235 import time_step_controller_t
236 import ksp_monitor_t
237 class(adjoint_scalar_scheme_t), intent(inout) :: this
238 type(time_state_t), intent(in) :: time
239 type(time_scheme_controller_t), intent(in) :: ext_bdf
240 type(time_step_controller_t), intent(in) :: dt_controller
241 type(ksp_monitor_t), intent(inout) :: ksp_results
243 end interface
244
245contains
246
257 subroutine adjoint_scalar_scheme_init(this, msh, c_Xh, gs_Xh, &
258 params_adjoint, params_primal, scheme, user, rho)
259 class(adjoint_scalar_scheme_t), target, intent(inout) :: this
260 type(mesh_t), target, intent(in) :: msh
261 type(coef_t), target, intent(in) :: c_Xh
262 type(gs_t), target, intent(inout) :: gs_Xh
263 type(json_file), target, intent(inout) :: params_primal, params_adjoint
264 character(len=*), intent(in) :: scheme
265 type(user_t), target, intent(in) :: user
266 type(field_t), target, intent(in) :: rho
267 type(json_file), pointer :: params_selected
268 ! IO buffer for log output
269 character(len=LOG_SIZE) :: log_buf
270 ! Variables for retrieving json parameters
271 logical :: logical_val
272 real(kind=rp) :: solver_abstol
273 integer :: integer_val
274 character(len=:), allocatable :: solver_type, solver_precon
275 type(json_file) :: precon_params
276
277 this%u => neko_registry%get_field('u')
278 this%v => neko_registry%get_field('v')
279 this%w => neko_registry%get_field('w')
280 this%rho => rho
281
282 ! get the primal adjoint's name
283 call json_get_or_default(params_adjoint, 'primal_name', this%primal_name, &
284 's')
285 ! Assign a name
286 call json_get_or_default(params_adjoint, 'name', this%name, &
287 this%primal_name // '_adj')
288
289 call neko_log%section('Adjoint scalar')
290 params_selected => json_key_fallback(params_adjoint, params_primal, &
291 'solver.type')
292 call json_get(params_selected, 'solver.type', solver_type)
293
294 params_selected => json_key_fallback(params_adjoint, params_primal, &
295 'solver.preconditioner.type')
296 call json_get(params_selected, 'solver.preconditioner.type', solver_precon)
297
298 params_selected => json_key_fallback(params_adjoint, params_primal, &
299 'solver.preconditioner')
300 call json_get(params_selected, 'solver.preconditioner', &
301 precon_params)
302
303 params_selected => json_key_fallback(params_adjoint, params_primal, &
304 'solver.absolute_tolerance')
305 call json_get(params_selected, 'solver.absolute_tolerance', &
306 solver_abstol)
307
308 params_selected => json_key_fallback(params_adjoint, params_primal, &
309 'solver.projection_space_size')
310 call json_get_or_default(params_selected, &
311 'solver.projection_space_size', &
312 this%projection_dim, 0)
313
314 params_selected => json_key_fallback(params_adjoint, params_primal, &
315 'solver.projection_hold_steps')
316 call json_get_or_default(params_selected, &
317 'solver.projection_hold_steps', &
318 this%projection_activ_step, 5)
319
320
321 write(log_buf, '(A, A)') 'Type : ', trim(scheme)
322 call neko_log%message(log_buf)
323 write(log_buf, '(A, A)') 'Name : ', trim(this%name)
324 call neko_log%message(log_buf)
325 call neko_log%message('Ksp adjoint scalar : ('// trim(solver_type) // &
326 ', ' // trim(solver_precon) // ')')
327 write(log_buf, '(A,ES13.6)') ' `-abs tol :', solver_abstol
328 call neko_log%message(log_buf)
329
330 this%Xh => this%u%Xh
331 this%dm_Xh => this%u%dof
332 this%params => params_adjoint
333 this%msh => msh
334
335 if (.not. neko_registry%field_exists(this%name)) then
336 call neko_registry%add_field(this%dm_Xh, this%name)
337 end if
338
339 this%s_adj => neko_registry%get_field(this%name)
340
341 call this%s_adj_lag%init(this%s_adj, 2)
342
343 this%s => neko_registry%get_field(this%primal_name)
344
345 this%gs_Xh => gs_xh
346 this%c_Xh => c_xh
347
348 !
349 ! Material properties
350 !
351 call this%set_material_properties(params_primal, user)
352
353
354 !
355 ! Turbulence modelling and variable material properties
356 !
357 params_selected => json_key_fallback(params_adjoint, params_primal, &
358 'variable_material_properties')
359 if (params_selected%valid_path('variable_material_properties')) then
360 call neko_error('variable material properties no supported for adjoint')
361 end if
362
363 write(log_buf, '(A,L1)') 'LES : ', this%variable_material_properties
364 call neko_log%message(log_buf)
365
366 !
367 ! Setup right-hand side field.
368 !
369 allocate(this%f_Xh)
370 call this%f_Xh%init(this%dm_Xh, fld_name = "adjoint_scalar_rhs")
371
372 ! Initialize the source term
373 call this%source_term%init(this%f_Xh, this%c_Xh, user, this%name)
374 ! We should ONLY read the adjoint
375 call this%source_term%add(params_primal, 'source_terms')
376
377 ! todo parameter file ksp tol should be added
378 params_selected => json_key_fallback(params_adjoint, params_primal, &
379 'solver.max_iterations')
380 call json_get_or_default(params_selected, &
381 'solver.max_iterations', &
382 integer_val, ksp_max_iter)
383 params_selected => json_key_fallback(params_adjoint, params_primal, &
384 'solver.monitor')
385 call json_get_or_default(params_selected, &
386 'solver.monitor', &
387 logical_val, .false.)
388 call adjoint_scalar_scheme_solver_factory(this%ksp, this%dm_Xh%size(), &
389 solver_type, integer_val, solver_abstol, logical_val)
390 call scalar_scheme_precon_factory(this%pc, this%ksp, &
391 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs, &
392 solver_precon, precon_params)
393
394 call neko_log%end_section()
395
396 end subroutine adjoint_scalar_scheme_init
397
398
401 class(adjoint_scalar_scheme_t), intent(inout) :: this
402
403 nullify(this%Xh)
404 nullify(this%dm_Xh)
405 nullify(this%gs_Xh)
406 nullify(this%c_Xh)
407 nullify(this%params)
408
409 if (allocated(this%ksp)) then
410 call this%ksp%free()
411 deallocate(this%ksp)
412 end if
413
414 if (allocated(this%pc)) then
415 call precon_destroy(this%pc)
416 deallocate(this%pc)
417 end if
418
419 call this%source_term%free()
420
421 call this%bcs%free()
422
423 call this%cp%free()
424 call this%lambda%free()
425 call this%s_adj_lag%free()
426
427 end subroutine adjoint_scalar_scheme_free
428
432 class(adjoint_scalar_scheme_t), target, intent(inout) :: this
433
434 if ( (.not. allocated(this%u%x)) .or. &
435 (.not. allocated(this%v%x)) .or. &
436 (.not. allocated(this%w%x)) .or. &
437 (.not. allocated(this%s%x)) .or. &
438 (.not. allocated(this%s_adj%x))) then
439 call neko_error('Fields are not allocated')
440 end if
441
442 if (.not. allocated(this%ksp)) then
443 call neko_error('No Krylov solver for velocity defined')
444 end if
445
446 if (.not. associated(this%Xh)) then
447 call neko_error('No function space defined')
448 end if
449
450 if (.not. associated(this%dm_Xh)) then
451 call neko_error('No dofmap defined')
452 end if
453
454 if (.not. associated(this%c_Xh)) then
455 call neko_error('No coefficients defined')
456 end if
457
458 if (.not. associated(this%f_Xh)) then
459 call neko_error('No rhs allocated')
460 end if
461
462 if (.not. associated(this%params)) then
463 call neko_error('No parameters defined')
464 end if
465
466 if (.not. associated(this%rho)) then
467 call neko_error('No density field defined')
468 end if
469
470 end subroutine adjoint_scalar_scheme_validate
471
474 subroutine adjoint_scalar_scheme_solver_factory(ksp, n, solver, max_iter, &
475 abstol, monitor)
476 class(ksp_t), allocatable, target, intent(inout) :: ksp
477 integer, intent(in), value :: n
478 integer, intent(in) :: max_iter
479 character(len=*), intent(in) :: solver
480 real(kind=rp) :: abstol
481 logical, intent(in) :: monitor
482
483 call krylov_solver_factory(ksp, n, solver, max_iter, &
484 abstol, monitor = monitor)
485
487
489 subroutine adjoint_scalar_scheme_precon_factory(pc, ksp, coef, dof, gs, &
490 bclst, pctype, pcparams)
491 class(pc_t), allocatable, target, intent(inout) :: pc
492 class(ksp_t), target, intent(inout) :: ksp
493 type(coef_t), target, intent(in) :: coef
494 type(dofmap_t), target, intent(in) :: dof
495 type(gs_t), target, intent(inout) :: gs
496 type(bc_list_t), target, intent(inout) :: bclst
497 character(len=*) :: pctype
498 type(json_file), intent(inout) :: pcparams
499
500 call precon_factory(pc, pctype)
501
502 select type (pcp => pc)
503 type is (jacobi_t)
504 call pcp%init(coef, dof, gs)
505 type is (sx_jacobi_t)
506 call pcp%init(coef, dof, gs)
507 type is (device_jacobi_t)
508 call pcp%init(coef, dof, gs)
509 type is (hsmg_t)
510 call pcp%init(coef, bclst, pcparams)
511 end select
512
513 call ksp%set_pc(pc)
514
516
522 class(adjoint_scalar_scheme_t), intent(inout) :: this
523 type(time_state_t), intent(in) :: time
524 type(field_t), pointer :: nut
525 integer :: index
526 ! Factor to transform nu_t to lambda_t
527 type(field_t), pointer :: lambda_factor
528
529 call this%user_material_properties(this%name, this%material_properties, &
530 time)
531
532 ! factor = rho * cp / pr_turb
533 if (this%variable_material_properties .and. &
534 len(trim(this%nut_field_name)) > 0) then
535 nut => neko_registry%get_field(this%nut_field_name)
536
537 ! lambda = lambda + rho * cp * nut / pr_turb
538 call neko_scratch_registry%request_field(lambda_factor, index, .false.)
539
540 call field_col3(lambda_factor, this%cp, this%rho)
541 call field_col2(lambda_factor, nut)
542 call field_cmult(lambda_factor, 1.0_rp / this%pr_turb)
543 call field_add2(this%lambda, lambda_factor)
544 call neko_scratch_registry%relinquish_field(index)
545 end if
546
547 ! Since cp is a fields and we use the %x(1,1,1,1) of the
548 ! host array data to pass constant material properties
549 ! to some routines, we need to make sure that the host
550 ! values are also filled
551 if (neko_bcknd_device .eq. 1) then
552 call device_memcpy(this%cp%x, this%cp%x_d, this%cp%size(), &
553 device_to_host, sync=.false.)
554 end if
555
557
563 params_primal, user)
564 class(adjoint_scalar_scheme_t), intent(inout) :: this
565 type(json_file), intent(inout) :: params_primal
566 type(user_t), target, intent(in) :: user
567 character(len=LOG_SIZE) :: log_buf
568 ! A local pointer that is needed to make Intel happy
569 procedure(user_material_properties_intf), pointer :: dummy_mp_ptr
570 real(kind=rp) :: const_cp, const_lambda
571 ! Dummy time state set to 0
572 type(time_state_t) :: time
573
574 dummy_mp_ptr => dummy_user_material_properties
575
576 ! Fill lambda field with the physical value
577 call this%lambda%init(this%dm_Xh, "lambda")
578 call this%cp%init(this%dm_Xh, "cp")
579
580 call this%material_properties%init(2)
581 call this%material_properties%assign_to_field(1, this%cp)
582 call this%material_properties%assign_to_field(2, this%lambda)
583
584 if (.not. associated(user%material_properties, dummy_mp_ptr)) then
585
586 write(log_buf, '(A)') "Material properties must be set in the user " // &
587 "file!"
588 call neko_log%message(log_buf)
589 this%user_material_properties => user%material_properties
590 call user%material_properties(this%name, this%material_properties, time)
591 else
592 this%user_material_properties => dummy_user_material_properties
593 if (params_primal%valid_path('Pe') .and. &
594 (params_primal%valid_path('lambda') .or. &
595 params_primal%valid_path('cp'))) then
596 call neko_error("To set the material properties for the scalar, " // &
597 "either provide Pe OR lambda and cp in the case file.")
598 ! Non-dimensional case
599 else if (params_primal%valid_path('Pe')) then
600 write(log_buf, '(A)') 'Non-dimensional scalar material properties' //&
601 ' input.'
602 call neko_log%message(log_buf, lvl = neko_log_verbose)
603 write(log_buf, '(A)') 'Specific heat capacity will be set to 1,'
604 call neko_log%message(log_buf, lvl = neko_log_verbose)
605 write(log_buf, '(A)') 'conductivity to 1/Pe. Assumes density is 1.'
606 call neko_log%message(log_buf, lvl = neko_log_verbose)
607
608 ! Read Pe into lambda for further manipulation.
609 call json_get(params_primal, 'Pe', const_lambda)
610 write(log_buf, '(A,ES13.6)') 'Pe :', const_lambda
611 call neko_log%message(log_buf)
612
613 ! Set cp and rho to 1 since the setup is non-dimensional.
614 const_cp = 1.0_rp
615 ! Invert the Pe to get conductivity
616 const_lambda = 1.0_rp/const_lambda
617 ! Dimensional case
618 else
619 call json_get(params_primal, 'lambda', const_lambda)
620 call json_get(params_primal, 'cp', const_cp)
621 end if
622 end if
623 ! We need to fill the fields based on the parsed const values
624 ! if the user routine is not used.
625 if (associated(user%material_properties, dummy_mp_ptr)) then
626 ! Fill mu and rho field with the physical value
627 call field_cfill(this%lambda, const_lambda)
628 call field_cfill(this%cp, const_cp)
629
630 write(log_buf, '(A,ES13.6)') 'lambda :', const_lambda
631 call neko_log%message(log_buf)
632 write(log_buf, '(A,ES13.6)') 'cp :', const_cp
633 call neko_log%message(log_buf)
634 end if
635
636 ! Since cp is a field and we use the %x(1,1,1,1) of the
637 ! host array data to pass constant material properties
638 ! to some routines, we need to make sure that the host
639 ! values are also filled
640 if (neko_bcknd_device .eq. 1) then
641 call device_memcpy(this%cp%x, this%cp%x_d, this%cp%size(), &
642 device_to_host, sync=.false.)
643 end if
645
646end 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.
subroutine adjoint_scalar_scheme_set_material_properties(this, params_primal, user)
Set lamdba and cp.
subroutine adjoint_scalar_scheme_init(this, msh, c_xh, gs_xh, params_adjoint, params_primal, scheme, user, rho)
Initialize all related components of the current scheme.
subroutine adjoint_scalar_scheme_solver_factory(ksp, n, solver, max_iter, abstol, monitor)
Initialize a linear solver.
subroutine adjoint_scalar_scheme_precon_factory(pc, ksp, coef, dof, gs, bclst, pctype, pcparams)
Initialize a Krylov preconditioner.
subroutine adjoint_scalar_scheme_update_material_properties(this, time)
Call user material properties routine and update the values of lambda if necessary.
subroutine adjoint_scalar_scheme_free(this)
Deallocate a scalar formulation.
subroutine adjoint_scalar_scheme_validate(this)
Validate that all fields, solvers etc necessary for performing time-stepping are defined.
Base type for a scalar advection-diffusion solver.