36 use comm,
only: mpi_integer, mpi_max, neko_comm, neko_error
37 use num_types,
only: rp
38 use,
intrinsic :: iso_fortran_env, only: error_unit
39 use rhs_maker,
only : rhs_maker_bdf_t, rhs_maker_ext_t, rhs_maker_oifs_t, &
40 rhs_maker_ext_fctry, rhs_maker_bdf_fctry, rhs_maker_oifs_fctry
42 use checkpoint,
only : chkp_t
43 use field,
only : field_t
44 use bc_list,
only : bc_list_t
45 use mesh,
only : mesh_t
46 use coefs,
only : coef_t
47 use device,
only : host_to_device, device_memcpy
48 use gather_scatter,
only : gs_t, gs_op_add
49 use scalar_residual,
only : scalar_residual_t, scalar_residual_factory
50 use ax_product,
only : ax_t, ax_helm_factory
51 use field_series,
only: field_series_t
52 use facet_normal,
only : facet_normal_t
53 use krylov,
only : ksp_monitor_t
54 use device_math,
only : device_add2s2, device_col2
55 use scalar_aux,
only : scalar_step_info
56 use time_scheme_controller,
only : time_scheme_controller_t
57 use projection,
only : projection_t
58 use math,
only : glsc2, col2, add2s2
59 use logger,
only : neko_log, log_size, neko_log_debug
61 use profiler,
only : profiler_start_region, profiler_end_region
62 use json_utils,
only : json_get, json_get_or_default, json_extract_item
63 use json_module,
only : json_file, json_core, json_value
64 use user_intf,
only : user_t
65 use neko_config,
only : neko_bcknd_device
66 use zero_dirichlet,
only : zero_dirichlet_t
67 use time_step_controller,
only : time_step_controller_t
68 use scratch_registry,
only : neko_scratch_registry
69 use time_state,
only : time_state_t
78 type(field_t) :: s_adj_res
81 type(field_t) :: ds_adj
84 class(ax_t),
allocatable :: ax
87 type(projection_t) :: proj_s
92 type(zero_dirichlet_t) :: bc_res
98 type(bc_list_t) :: bclst_ds
107 type(field_t) :: advs
110 class(scalar_residual_t),
allocatable :: res
113 class(rhs_maker_ext_t),
allocatable :: makeext
116 class(rhs_maker_bdf_t),
allocatable :: makebdf
119 class(rhs_maker_oifs_t),
allocatable :: makeoifs
123 procedure, pass(this) :: init => adjoint_scalar_pnpn_init
125 procedure, pass(this) :: restart => adjoint_scalar_pnpn_restart
127 procedure, pass(this) :: free => adjoint_scalar_pnpn_free
129 procedure, pass(this) :: step => adjoint_scalar_pnpn_step
131 procedure, pass(this) :: setup_bcs_ => adjoint_scalar_pnpn_setup_bcs_
141 module subroutine adjoint_bc_factory(object, scheme, json, coef, user)
142 class(bc_t),
pointer,
intent(inout) :: object
143 type(adjoint_scalar_pnpn_t),
intent(in) :: scheme
144 type(json_file),
intent(inout) :: json
145 type(coef_t),
intent(in) :: coef
146 type(user_t),
intent(in) :: user
147 end subroutine adjoint_bc_factory
148 end interface adjoint_bc_factory
168 subroutine adjoint_scalar_pnpn_init(this, msh, coef, gs, params_adjoint, &
169 params_primal, numerics_params, user, chkp, ulag, vlag, wlag, &
171 class(adjoint_scalar_pnpn_t),
target,
intent(inout) :: this
172 type(mesh_t),
target,
intent(in) :: msh
173 type(coef_t),
target,
intent(in) :: coef
174 type(gs_t),
target,
intent(inout) :: gs
175 type(json_file),
target,
intent(inout) :: params_adjoint
176 type(json_file),
target,
intent(inout) :: params_primal
177 type(json_file),
target,
intent(inout) :: numerics_params
178 type(user_t),
target,
intent(in) :: user
179 type(chkp_t),
target,
intent(inout) :: chkp
180 type(field_series_t),
target,
intent(in) :: ulag, vlag, wlag
181 type(time_scheme_controller_t),
target,
intent(in) :: time_scheme
182 type(field_t),
target,
intent(in) :: rho
184 class(bc_t),
pointer :: bc_i
185 character(len=15),
parameter :: scheme =
'Modular (Pn/Pn)'
191 call this%scheme_init(msh, coef, gs, params_adjoint, params_primal, &
195 call ax_helm_factory(this%ax, full_formulation = .false.)
198 call scalar_residual_factory(this%res)
201 call rhs_maker_ext_fctry(this%makeext)
204 call rhs_maker_bdf_fctry(this%makebdf)
207 call rhs_maker_oifs_fctry(this%makeoifs)
210 associate(xh_lx => this%Xh%lx, xh_ly => this%Xh%ly, xh_lz => this%Xh%lz, &
211 dm_xh => this%dm_Xh, nelv => this%msh%nelv)
213 call this%s_adj_res%init(dm_xh,
"s_adj_res")
215 call this%abx1%init(dm_xh,
"abx1")
217 call this%abx2%init(dm_xh,
"abx2")
219 call this%advs%init(dm_xh,
"advs")
221 call this%ds_adj%init(dm_xh,
'ds_adj')
226 call this%setup_bcs_(user)
229 call this%bc_res%init(this%c_Xh, params_adjoint)
230 do i = 1, this%bcs%size()
231 if (this%bcs%strong(i))
then
232 bc_i => this%bcs%get(i)
233 call this%bc_res%mark_facets(bc_i%marked_facet)
238 call this%bc_res%finalize()
240 call this%bclst_ds%init()
241 call this%bclst_ds%append(this%bc_res)
245 call this%proj_s%init(this%dm_Xh%size(), this%projection_dim, &
246 this%projection_activ_step)
249 call json_get_or_default(numerics_params,
'oifs', this%oifs, .false.)
255 call json_get_or_default(params_adjoint,
'advection', advection, .true.)
257 call advection_adjoint_factory(this%adv, numerics_params, this%c_Xh, &
258 ulag, vlag, wlag, this%chkp%dtlag, &
259 this%chkp%tlag, time_scheme, .not. advection, &
270 end subroutine adjoint_scalar_pnpn_init
273 subroutine adjoint_scalar_pnpn_restart(this, chkp)
274 class(adjoint_scalar_pnpn_t),
target,
intent(inout) :: this
275 type(chkp_t),
intent(inout) :: chkp
276 real(kind=rp) :: dtlag(10), tlag(10)
281 n = this%s_adj%dof%size()
283 call col2(this%s_adj%x, this%c_Xh%mult, n)
284 call col2(this%s_adj_lag%lf(1)%x, this%c_Xh%mult, n)
285 call col2(this%s_adj_lag%lf(2)%x, this%c_Xh%mult, n)
286 if (neko_bcknd_device .eq. 1)
then
287 call device_memcpy(this%s_adj%x, this%s_adj%x_d, &
288 n, host_to_device, sync = .false.)
289 call device_memcpy(this%s_adj_lag%lf(1)%x, this%s_adj_lag%lf(1)%x_d, &
290 n, host_to_device, sync = .false.)
291 call device_memcpy(this%s_adj_lag%lf(2)%x, this%s_adj_lag%lf(2)%x_d, &
292 n, host_to_device, sync = .false.)
293 call device_memcpy(this%abx1%x, this%abx1%x_d, &
294 n, host_to_device, sync = .false.)
295 call device_memcpy(this%abx2%x, this%abx2%x_d, &
296 n, host_to_device, sync = .false.)
297 call device_memcpy(this%advs%x, this%advs%x_d, &
298 n, host_to_device, sync = .false.)
301 call this%gs_Xh%op(this%s_adj, gs_op_add)
302 call this%gs_Xh%op(this%s_adj_lag%lf(1), gs_op_add)
303 call this%gs_Xh%op(this%s_adj_lag%lf(2), gs_op_add)
305 end subroutine adjoint_scalar_pnpn_restart
307 subroutine adjoint_scalar_pnpn_free(this)
308 class(adjoint_scalar_pnpn_t),
intent(inout) :: this
311 call this%scheme_free()
313 call this%bclst_ds%free()
314 call this%proj_s%free()
316 call this%s_adj_res%free()
318 call this%ds_adj%free()
320 call this%abx1%free()
321 call this%abx2%free()
323 call this%advs%free()
325 if (
allocated(this%Ax))
then
329 if (
allocated(this%res))
then
333 if (
allocated(this%makeext))
then
334 deallocate(this%makeext)
337 if (
allocated(this%makebdf))
then
338 deallocate(this%makebdf)
341 if (
allocated(this%makeoifs))
then
342 deallocate(this%makeoifs)
345 end subroutine adjoint_scalar_pnpn_free
347 subroutine adjoint_scalar_pnpn_step(this, time, ext_bdf, dt_controller)
348 class(adjoint_scalar_pnpn_t),
intent(inout) :: this
349 type(time_state_t),
intent(in) :: time
350 type(time_scheme_controller_t),
intent(in) :: ext_bdf
351 type(time_step_controller_t),
intent(in) :: dt_controller
355 type(ksp_monitor_t) :: ksp_results(1)
357 n = this%dm_Xh%size()
359 call profiler_start_region(
'Scalar', 2)
360 associate(u => this%u, v => this%v, w => this%w, s_adj => this%s_adj, &
361 cp => this%cp, rho => this%rho, lambda => this%lambda, &
362 ds_adj => this%ds_adj, &
363 s_adj_res => this%s_adj_res, &
364 ax => this%Ax, f_xh => this%f_Xh, xh => this%Xh, &
365 c_xh => this%c_Xh, dm_xh => this%dm_Xh, gs_xh => this%gs_Xh, &
366 s_adj_lag => this%s_adj_lag, oifs => this%oifs, &
367 projection_dim => this%projection_dim, &
368 msh => this%msh, res => this%res, makeoifs => this%makeoifs, &
369 makeext => this%makeext, makebdf => this%makebdf, &
370 t => time%t, tstep => time%tstep, dt => time%dt)
373 call print_debug(this)
375 call this%source_term%compute(t, tstep)
378 call this%bcs%apply_scalar(this%f_Xh%x, dm_xh%size(), t, tstep, .false.)
392 call this%adv%compute_adjoint_scalar(u, v, w, s_adj, f_xh, &
393 xh, this%c_Xh, dm_xh%size())
400 call makeext%compute_scalar(this%abx1, this%abx2, f_xh%x, &
401 rho%x(1,1,1,1), ext_bdf%advection_coeffs, n)
404 call makebdf%compute_scalar(s_adj_lag, f_xh%x, s_adj, c_xh%B, &
405 rho%x(1,1,1,1), dt, ext_bdf%diffusion_coeffs, ext_bdf%ndiff, n)
408 call s_adj_lag%update()
411 call this%bcs%apply_scalar(this%s_adj%x, this%dm_Xh%size(), t, tstep, &
415 call this%update_material_properties(t, tstep)
418 call profiler_start_region(
'Adjoint_scalar_residual', 20)
419 call res%compute(ax, s_adj, s_adj_res, f_xh, c_xh, msh, xh, &
420 lambda, rho%x(1,1,1,1)*cp%x(1,1,1,1), ext_bdf%diffusion_coeffs(1), &
423 call gs_xh%op(s_adj_res, gs_op_add)
427 call this%bclst_ds%apply_scalar(s_adj_res%x, dm_xh%size())
429 call profiler_end_region(
'Adjoint_scalar_residual', 20)
431 call this%proj_s%pre_solving(s_adj_res%x, tstep, c_xh, n, dt_controller)
433 call this%pc%update()
434 call profiler_start_region(
'Adjoint_scalar_solve', 21)
435 ksp_results(1) = this%ksp%solve(ax, ds_adj, s_adj_res%x, n, &
436 c_xh, this%bclst_ds, gs_xh)
437 call profiler_end_region(
'Adjoint_scalar_solve', 21)
439 call this%proj_s%post_solving(ds_adj%x, ax, c_xh, this%bclst_ds, gs_xh, &
440 n, tstep, dt_controller)
443 if (neko_bcknd_device .eq. 1)
then
444 call device_add2s2(s_adj%x_d, ds_adj%x_d, 1.0_rp, n)
446 call add2s2(s_adj%x, ds_adj%x, 1.0_rp, n)
449 call scalar_step_info(time, ksp_results)
452 call profiler_end_region(
'Scalar', 2)
453 end subroutine adjoint_scalar_pnpn_step
455 subroutine print_debug(this)
456 class(adjoint_scalar_pnpn_t),
intent(inout) :: this
460 n = this%dm_Xh%size()
471 end subroutine print_debug
476 subroutine adjoint_scalar_pnpn_setup_bcs_(this, user)
477 class(adjoint_scalar_pnpn_t),
intent(inout) :: this
478 type(user_t),
target,
intent(in) :: user
479 integer :: i, j, n_bcs, zone_size, global_zone_size, ierr
480 type(json_core) :: core
481 type(json_value),
pointer :: bc_object
482 type(json_file) :: bc_subdict
483 class(bc_t),
pointer :: bc_i
486 logical,
allocatable :: marked_zones(:)
487 integer,
allocatable :: zone_indices(:)
489 if (this%params%valid_path(
'boundary_conditions'))
then
490 call this%params%info(
'boundary_conditions', &
492 call this%params%get_core(core)
493 call this%params%get(
'boundary_conditions', bc_object, found)
495 call this%bcs%init(n_bcs)
497 allocate(marked_zones(
size(this%msh%labeled_zones)))
498 marked_zones = .false.
502 call json_extract_item(core, bc_object, i, bc_subdict)
507 call json_get(bc_subdict,
"zone_indices", zone_indices)
509 do j = 1,
size(zone_indices)
510 zone_size = this%msh%labeled_zones(zone_indices(j))%size
511 call mpi_allreduce(zone_size, global_zone_size, 1, &
512 mpi_integer, mpi_max, neko_comm, ierr)
514 if (global_zone_size .eq. 0)
then
515 write(error_unit,
'(A, A, I0, A, A, I0, A)')
"*** ERROR ***: ",&
516 "Zone index ", zone_indices(j), &
517 " is invalid as this zone has 0 size, meaning it ", &
518 "does not exist in the mesh. Check adjoint scalar BC ", &
523 if (marked_zones(zone_indices(j)) .eqv. .true.)
then
524 write(error_unit,
'(A, A, I0, A, A, A, A)')
"*** ERROR ***: ", &
525 "Zone with index ", zone_indices(j), &
526 " has already been assigned a boundary condition. ", &
527 "Please check your boundary_conditions entry for the ", &
528 "adjoint scalar and make sure that each zone index ", &
529 "appears only in a single boundary condition."
532 marked_zones(zone_indices(j)) = .true.
538 call adjoint_bc_factory(bc_i, this, bc_subdict, this%c_Xh, user)
539 call this%bcs%append(bc_i)
543 do i = 1,
size(this%msh%labeled_zones)
544 if ((this%msh%labeled_zones(i)%size .gt. 0) .and. &
545 (marked_zones(i) .eqv. .false.))
then
546 write(error_unit,
'(A, A, I0)')
"*** ERROR ***: ", &
547 "No adjoint scalar boundary condition assigned to zone ", i
553 do i = 1,
size(this%msh%labeled_zones)
554 if (this%msh%labeled_zones(i)%size .gt. 0)
then
555 write(error_unit,
'(A, A, A)')
"*** ERROR ***: ", &
556 "No boundary_conditions entry in the case file for " // &
557 " adjoint scalar ", &
563 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.