38 use comm,
only: neko_comm
39 use utils,
only: neko_error
40 use num_types,
only: rp
41 use,
intrinsic :: iso_fortran_env, only: error_unit
42 use rhs_maker,
only : rhs_maker_bdf_t, rhs_maker_ext_t, rhs_maker_oifs_t, &
43 rhs_maker_ext_fctry, rhs_maker_bdf_fctry, rhs_maker_oifs_fctry
45 use checkpoint,
only : chkp_t
46 use field,
only : field_t
47 use bc_list,
only : bc_list_t
48 use mesh,
only : mesh_t
49 use coefs,
only : coef_t
50 use device,
only : host_to_device, device_memcpy
51 use gather_scatter,
only : gs_t, gs_op_add
52 use scalar_residual,
only : scalar_residual_t, scalar_residual_factory
53 use ax_product,
only : ax_t, ax_helm_factory
54 use field_series,
only: field_series_t
55 use facet_normal,
only : facet_normal_t
56 use krylov,
only : ksp_monitor_t
57 use device_math,
only : device_add2s2, device_col2
58 use time_scheme_controller,
only : time_scheme_controller_t
59 use projection,
only : projection_t
60 use math,
only : glsc2, col2, add2s2
61 use logger,
only : neko_log, log_size, neko_log_debug
63 use profiler,
only : profiler_start_region, profiler_end_region
64 use json_utils,
only : json_get, json_get_or_default, json_extract_item
65 use json_module,
only : json_file, json_core, json_value
66 use user_intf,
only : user_t
67 use neko_config,
only : neko_bcknd_device
68 use zero_dirichlet,
only : zero_dirichlet_t
69 use time_step_controller,
only : time_step_controller_t
70 use scratch_registry,
only : neko_scratch_registry
71 use time_state,
only : time_state_t
73 use mpi_f08,
only: mpi_integer, mpi_sum, mpi_max
81 type(field_t) :: s_adj_res
84 type(field_t) :: ds_adj
87 class(ax_t),
allocatable :: ax
90 type(projection_t) :: proj_s
95 type(zero_dirichlet_t) :: bc_res
101 type(bc_list_t) :: bclst_ds
110 type(field_t) :: advs
113 class(scalar_residual_t),
allocatable :: res
116 class(rhs_maker_ext_t),
allocatable :: makeext
119 class(rhs_maker_bdf_t),
allocatable :: makebdf
122 class(rhs_maker_oifs_t),
allocatable :: makeoifs
126 procedure, pass(this) :: init => adjoint_scalar_pnpn_init
128 procedure, pass(this) :: restart => adjoint_scalar_pnpn_restart
130 procedure, pass(this) :: free => adjoint_scalar_pnpn_free
132 procedure, pass(this) :: step => adjoint_scalar_pnpn_step
134 procedure, pass(this) :: setup_bcs_ => adjoint_scalar_pnpn_setup_bcs_
144 module subroutine adjoint_bc_factory(object, scheme, json, coef, user)
145 class(bc_t),
pointer,
intent(inout) :: object
146 type(adjoint_scalar_pnpn_t),
intent(in) :: scheme
147 type(json_file),
intent(inout) :: json
148 type(coef_t),
intent(in) :: coef
149 type(user_t),
intent(in) :: user
150 end subroutine adjoint_bc_factory
151 end interface adjoint_bc_factory
171 subroutine adjoint_scalar_pnpn_init(this, msh, coef, gs, params_adjoint, &
172 params_primal, numerics_params, user, chkp, ulag, vlag, wlag, &
174 class(adjoint_scalar_pnpn_t),
target,
intent(inout) :: this
175 type(mesh_t),
target,
intent(in) :: msh
176 type(coef_t),
target,
intent(in) :: coef
177 type(gs_t),
target,
intent(inout) :: gs
178 type(json_file),
target,
intent(inout) :: params_adjoint
179 type(json_file),
target,
intent(inout) :: params_primal
180 type(json_file),
target,
intent(inout) :: numerics_params
181 type(user_t),
target,
intent(in) :: user
182 type(chkp_t),
target,
intent(inout) :: chkp
183 type(field_series_t),
target,
intent(in) :: ulag, vlag, wlag
184 type(time_scheme_controller_t),
target,
intent(in) :: time_scheme
185 type(field_t),
target,
intent(in) :: rho
187 class(bc_t),
pointer :: bc_i
188 character(len=15),
parameter :: scheme =
'Modular (Pn/Pn)'
194 call this%scheme_init(msh, coef, gs, params_adjoint, params_primal, &
198 call ax_helm_factory(this%ax, full_formulation = .false.)
201 call scalar_residual_factory(this%res)
204 call rhs_maker_ext_fctry(this%makeext)
207 call rhs_maker_bdf_fctry(this%makebdf)
210 call rhs_maker_oifs_fctry(this%makeoifs)
213 associate(xh_lx => this%Xh%lx, xh_ly => this%Xh%ly, xh_lz => this%Xh%lz, &
214 dm_xh => this%dm_Xh, nelv => this%msh%nelv)
216 call this%s_adj_res%init(dm_xh,
"s_adj_res")
218 call this%abx1%init(dm_xh,
"abx1")
220 call this%abx2%init(dm_xh,
"abx2")
222 call this%advs%init(dm_xh,
"advs")
224 call this%ds_adj%init(dm_xh,
'ds_adj')
229 call this%setup_bcs_(user)
232 call this%bc_res%init(this%c_Xh, params_adjoint)
233 do i = 1, this%bcs%size()
234 if (this%bcs%strong(i))
then
235 bc_i => this%bcs%get(i)
236 call this%bc_res%mark_facets(bc_i%marked_facet)
241 call this%bc_res%finalize()
243 call this%bclst_ds%init()
244 call this%bclst_ds%append(this%bc_res)
248 call this%proj_s%init(this%dm_Xh%size(), this%projection_dim, &
249 this%projection_activ_step)
252 call json_get_or_default(numerics_params,
'oifs', this%oifs, .false.)
258 call json_get_or_default(params_adjoint,
'advection', advection, .true.)
260 call advection_adjoint_factory(this%adv, numerics_params, this%c_Xh, &
261 ulag, vlag, wlag, this%chkp%dtlag, &
262 this%chkp%tlag, time_scheme, .not. advection, &
273 end subroutine adjoint_scalar_pnpn_init
276 subroutine adjoint_scalar_pnpn_restart(this, chkp)
277 class(adjoint_scalar_pnpn_t),
target,
intent(inout) :: this
278 type(chkp_t),
intent(inout) :: chkp
279 real(kind=rp) :: dtlag(10), tlag(10)
284 n = this%s_adj%dof%size()
286 call col2(this%s_adj%x, this%c_Xh%mult, n)
287 call col2(this%s_adj_lag%lf(1)%x, this%c_Xh%mult, n)
288 call col2(this%s_adj_lag%lf(2)%x, this%c_Xh%mult, n)
289 if (neko_bcknd_device .eq. 1)
then
290 call device_memcpy(this%s_adj%x, this%s_adj%x_d, &
291 n, host_to_device, sync = .false.)
292 call device_memcpy(this%s_adj_lag%lf(1)%x, this%s_adj_lag%lf(1)%x_d, &
293 n, host_to_device, sync = .false.)
294 call device_memcpy(this%s_adj_lag%lf(2)%x, this%s_adj_lag%lf(2)%x_d, &
295 n, host_to_device, sync = .false.)
296 call device_memcpy(this%abx1%x, this%abx1%x_d, &
297 n, host_to_device, sync = .false.)
298 call device_memcpy(this%abx2%x, this%abx2%x_d, &
299 n, host_to_device, sync = .false.)
300 call device_memcpy(this%advs%x, this%advs%x_d, &
301 n, host_to_device, sync = .false.)
304 call this%gs_Xh%op(this%s_adj, gs_op_add)
305 call this%gs_Xh%op(this%s_adj_lag%lf(1), gs_op_add)
306 call this%gs_Xh%op(this%s_adj_lag%lf(2), gs_op_add)
308 end subroutine adjoint_scalar_pnpn_restart
310 subroutine adjoint_scalar_pnpn_free(this)
311 class(adjoint_scalar_pnpn_t),
intent(inout) :: this
314 call this%scheme_free()
316 call this%bclst_ds%free()
317 call this%bc_res%free()
318 call this%proj_s%free()
320 call this%s_adj_res%free()
322 call this%ds_adj%free()
324 call this%abx1%free()
325 call this%abx2%free()
327 call this%advs%free()
329 if (
allocated(this%Ax))
then
333 if (
allocated(this%res))
then
337 if (
allocated(this%makeext))
then
338 deallocate(this%makeext)
341 if (
allocated(this%makebdf))
then
342 deallocate(this%makebdf)
345 if (
allocated(this%makeoifs))
then
346 deallocate(this%makeoifs)
349 end subroutine adjoint_scalar_pnpn_free
351 subroutine adjoint_scalar_pnpn_step(this, time, ext_bdf, dt_controller, &
353 class(adjoint_scalar_pnpn_t),
intent(inout) :: this
354 type(time_state_t),
intent(in) :: time
355 type(time_scheme_controller_t),
intent(in) :: ext_bdf
356 type(time_step_controller_t),
intent(in) :: dt_controller
357 type(ksp_monitor_t),
intent(inout) :: ksp_results
361 if (this%freeze)
return
363 n = this%dm_Xh%size()
365 call profiler_start_region(
'Adjoint Scalar')
366 associate(u => this%u, v => this%v, w => this%w, s_adj => this%s_adj, &
367 cp => this%cp, rho => this%rho, lambda => this%lambda, &
368 ds_adj => this%ds_adj, &
369 s_adj_res => this%s_adj_res, &
370 ax => this%Ax, f_xh => this%f_Xh, xh => this%Xh, &
371 c_xh => this%c_Xh, dm_xh => this%dm_Xh, gs_xh => this%gs_Xh, &
372 s_adj_lag => this%s_adj_lag, oifs => this%oifs, &
373 projection_dim => this%projection_dim, &
374 msh => this%msh, res => this%res, makeoifs => this%makeoifs, &
375 makeext => this%makeext, makebdf => this%makebdf, &
376 t => time%t, tstep => time%tstep, dt => time%dt)
379 call print_debug(this)
381 call this%source_term%compute(time)
384 call this%bcs%apply_scalar(this%f_Xh%x, dm_xh%size(), time, .false.)
398 call this%adv%compute_adjoint_scalar(u, v, w, s_adj, f_xh, &
399 xh, this%c_Xh, dm_xh%size())
406 call makeext%compute_scalar(this%abx1, this%abx2, f_xh%x, &
407 rho%x(1,1,1,1), ext_bdf%advection_coeffs%x, n)
410 call makebdf%compute_scalar(s_adj_lag, f_xh%x, s_adj, c_xh%B, &
411 rho%x(1,1,1,1), dt, ext_bdf%diffusion_coeffs%x, ext_bdf%ndiff, n)
414 call s_adj_lag%update()
417 call this%bcs%apply_scalar(this%s_adj%x, this%dm_Xh%size(), time, &
421 call this%update_material_properties(time)
424 call profiler_start_region(
'Adjoint_scalar_residual')
425 call res%compute(ax, s_adj, s_adj_res, f_xh, c_xh, msh, xh, &
426 lambda, rho%x(1,1,1,1)*cp%x(1,1,1,1), ext_bdf%diffusion_coeffs%x(1), &
429 call gs_xh%op(s_adj_res, gs_op_add)
433 call this%bclst_ds%apply_scalar(s_adj_res%x, dm_xh%size())
435 call profiler_end_region(
'Adjoint_scalar_residual')
437 call this%proj_s%pre_solving(s_adj_res%x, tstep, c_xh, n, dt_controller)
439 call this%pc%update()
440 call profiler_start_region(
'Adjoint_scalar_solve')
441 ksp_results = this%ksp%solve(ax, ds_adj, s_adj_res%x, n, &
442 c_xh, this%bclst_ds, gs_xh)
443 call profiler_end_region(
'Adjoint_scalar_solve')
445 ksp_results%name =
'Adjoint Scalar'
447 call this%proj_s%post_solving(ds_adj%x, ax, c_xh, this%bclst_ds, gs_xh, &
448 n, tstep, dt_controller)
451 if (neko_bcknd_device .eq. 1)
then
452 call device_add2s2(s_adj%x_d, ds_adj%x_d, 1.0_rp, n)
454 call add2s2(s_adj%x, ds_adj%x, 1.0_rp, n)
458 call profiler_end_region(
'Adjoint Scalar')
459 end subroutine adjoint_scalar_pnpn_step
461 subroutine print_debug(this)
462 class(adjoint_scalar_pnpn_t),
intent(inout) :: this
466 n = this%dm_Xh%size()
477 end subroutine print_debug
482 subroutine adjoint_scalar_pnpn_setup_bcs_(this, user)
483 class(adjoint_scalar_pnpn_t),
intent(inout) :: this
484 type(user_t),
target,
intent(in) :: user
485 integer :: i, j, n_bcs, zone_size, global_zone_size, ierr
486 type(json_core) :: core
487 type(json_value),
pointer :: bc_object
488 type(json_file) :: bc_subdict
489 class(bc_t),
pointer :: bc_i
492 logical,
allocatable :: marked_zones(:)
493 integer,
allocatable :: zone_indices(:)
495 if (this%params%valid_path(
'boundary_conditions'))
then
496 call this%params%info(
'boundary_conditions', &
498 call this%params%get_core(core)
499 call this%params%get(
'boundary_conditions', bc_object, found)
501 call this%bcs%init(n_bcs)
503 allocate(marked_zones(
size(this%msh%labeled_zones)))
504 marked_zones = .false.
508 call json_extract_item(core, bc_object, i, bc_subdict)
513 call json_get(bc_subdict,
"zone_indices", zone_indices)
515 do j = 1,
size(zone_indices)
516 zone_size = this%msh%labeled_zones(zone_indices(j))%size
517 call mpi_allreduce(zone_size, global_zone_size, 1, &
518 mpi_integer, mpi_max, neko_comm, ierr)
520 if (global_zone_size .eq. 0)
then
521 write(error_unit,
'(A, A, I0, A, A, I0, A)')
"*** ERROR ***: ",&
522 "Zone index ", zone_indices(j), &
523 " is invalid as this zone has 0 size, meaning it ", &
524 "does not exist in the mesh. Check adjoint scalar BC ", &
529 if (marked_zones(zone_indices(j)) .eqv. .true.)
then
530 write(error_unit,
'(A, A, I0, A, A, A, A)')
"*** ERROR ***: ", &
531 "Zone with index ", zone_indices(j), &
532 " has already been assigned a boundary condition. ", &
533 "Please check your boundary_conditions entry for the ", &
534 "adjoint scalar and make sure that each zone index ", &
535 "appears only in a single boundary condition."
538 marked_zones(zone_indices(j)) = .true.
544 call adjoint_bc_factory(bc_i, this, bc_subdict, this%c_Xh, user)
545 call this%bcs%append(bc_i)
549 do i = 1,
size(this%msh%labeled_zones)
550 if ((this%msh%labeled_zones(i)%size .gt. 0) .and. &
551 (marked_zones(i) .eqv. .false.))
then
552 write(error_unit,
'(A, A, I0)')
"*** ERROR ***: ", &
553 "No adjoint scalar boundary condition assigned to zone ", i
559 do i = 1,
size(this%msh%labeled_zones)
560 if (this%msh%labeled_zones(i)%size .gt. 0)
then
561 write(error_unit,
'(A, A, A)')
"*** ERROR ***: ", &
562 "No boundary_conditions entry in the case file for " // &
563 " adjoint scalar ", &
569 end subroutine adjoint_scalar_pnpn_setup_bcs_
Boundary condition factory. Both constructs and initializes the object.
Contains the adjoint_scalar_pnpn_t type.
Contains the adjoint_scalar_scheme_t type.
Subroutines to add advection terms to the RHS of a transport equation.
Base type for a scalar advection-diffusion solver.
Base abstract type for computing the advection operator.