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 field,
only : field_t
43 use bc_list,
only : bc_list_t
44 use mesh,
only : mesh_t
45 use checkpoint,
only : chkp_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
62 use profiler,
only : profiler_start_region, profiler_end_region
63 use json_utils,
only : json_get, json_get_or_default, json_extract_item
64 use json_module,
only : json_file, json_core, json_value
65 use user_intf,
only : user_t
66 use neko_config,
only : neko_bcknd_device
67 use zero_dirichlet,
only : zero_dirichlet_t
68 use time_step_controller,
only : time_step_controller_t
69 use scratch_registry,
only : neko_scratch_registry
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
96 type(bc_list_t) :: bclst_ds
105 type(field_t) :: abx1, abx2
108 type(field_t) :: advs
111 class(scalar_residual_t),
allocatable :: res
114 class(rhs_maker_ext_t),
allocatable :: makeext
117 class(rhs_maker_bdf_t),
allocatable :: makebdf
120 class(rhs_maker_oifs_t),
allocatable :: makeoifs
124 procedure, pass(this) :: init => adjoint_scalar_pnpn_init
126 procedure, pass(this) :: restart => adjoint_scalar_pnpn_restart
128 procedure, pass(this) :: free => adjoint_scalar_pnpn_free
130 procedure, pass(this) :: step => adjoint_scalar_pnpn_step
132 procedure, pass(this) :: setup_bcs_ => adjoint_scalar_pnpn_setup_bcs_
142 module subroutine adjoint_bc_factory(object, scheme, json, coef, user)
143 class(bc_t),
pointer,
intent(inout) :: object
144 type(adjoint_scalar_pnpn_t),
intent(in) :: scheme
145 type(json_file),
intent(inout) :: json
146 type(coef_t),
intent(in) :: coef
147 type(user_t),
intent(in) :: user
148 end subroutine adjoint_bc_factory
149 end interface adjoint_bc_factory
166 subroutine adjoint_scalar_pnpn_init(this, msh, coef, gs, params, user, &
167 ulag, vlag, wlag, time_scheme, rho)
168 class(adjoint_scalar_pnpn_t),
target,
intent(inout) :: this
169 type(mesh_t),
target,
intent(in) :: msh
170 type(coef_t),
target,
intent(in) :: coef
171 type(gs_t),
target,
intent(inout) :: gs
172 type(json_file),
target,
intent(inout) :: params
173 type(user_t),
target,
intent(in) :: user
174 type(field_series_t),
target,
intent(in) :: ulag, vlag, wlag
175 type(time_scheme_controller_t),
target,
intent(in) :: time_scheme
176 type(field_t),
target,
intent(in) :: rho
178 class(bc_t),
pointer :: bc_i
179 character(len=15),
parameter :: scheme =
'Modular (Pn/Pn)'
184 call this%scheme_init(msh, coef, gs, params, scheme, user, rho)
187 call ax_helm_factory(this%ax, full_formulation = .false.)
190 call scalar_residual_factory(this%res)
193 call rhs_maker_ext_fctry(this%makeext)
196 call rhs_maker_bdf_fctry(this%makebdf)
199 call rhs_maker_oifs_fctry(this%makeoifs)
202 associate(xh_lx => this%Xh%lx, xh_ly => this%Xh%ly, xh_lz => this%Xh%lz, &
203 dm_xh => this%dm_Xh, nelv => this%msh%nelv)
205 call this%s_adj_res%init(dm_xh,
"s_adj_res")
207 call this%abx1%init(dm_xh,
"abx1")
209 call this%abx2%init(dm_xh,
"abx2")
211 call this%advs%init(dm_xh,
"advs")
213 call this%ds_adj%init(dm_xh,
'ds_adj')
218 call this%setup_bcs_(user)
221 call this%bc_res%init(this%c_Xh, params)
222 do i = 1, this%bcs%size()
223 if (this%bcs%strong(i))
then
224 bc_i => this%bcs%get(i)
225 call this%bc_res%mark_facets(bc_i%marked_facet)
230 call this%bc_res%finalize()
232 call this%bclst_ds%init()
233 call this%bclst_ds%append(this%bc_res)
237 call this%proj_s%init(this%dm_Xh%size(), this%projection_dim, &
238 this%projection_activ_step)
245 call json_get_or_default(params,
'case.numerics.oifs', this%oifs, .false.)
256 call advection_adjoint_factory(this%adv, params, this%c_Xh)
257 end subroutine adjoint_scalar_pnpn_init
260 subroutine adjoint_scalar_pnpn_restart(this, dtlag, tlag)
261 class(adjoint_scalar_pnpn_t),
target,
intent(inout) :: this
262 real(kind=rp) :: dtlag(10), tlag(10)
266 n = this%s_adj%dof%size()
268 call col2(this%s_adj%x, this%c_Xh%mult, n)
269 call col2(this%s_adj_lag%lf(1)%x, this%c_Xh%mult, n)
270 call col2(this%s_adj_lag%lf(2)%x, this%c_Xh%mult, n)
271 if (neko_bcknd_device .eq. 1)
then
272 call device_memcpy(this%s_adj%x, this%s_adj%x_d, &
273 n, host_to_device, sync = .false.)
274 call device_memcpy(this%s_adj_lag%lf(1)%x, this%s_adj_lag%lf(1)%x_d, &
275 n, host_to_device, sync = .false.)
276 call device_memcpy(this%s_adj_lag%lf(2)%x, this%s_adj_lag%lf(2)%x_d, &
277 n, host_to_device, sync = .false.)
278 call device_memcpy(this%abx1%x, this%abx1%x_d, &
279 n, host_to_device, sync = .false.)
280 call device_memcpy(this%abx2%x, this%abx2%x_d, &
281 n, host_to_device, sync = .false.)
282 call device_memcpy(this%advs%x, this%advs%x_d, &
283 n, host_to_device, sync = .false.)
286 call this%gs_Xh%op(this%s_adj, gs_op_add)
287 call this%gs_Xh%op(this%s_adj_lag%lf(1), gs_op_add)
288 call this%gs_Xh%op(this%s_adj_lag%lf(2), gs_op_add)
290 end subroutine adjoint_scalar_pnpn_restart
292 subroutine adjoint_scalar_pnpn_free(this)
293 class(adjoint_scalar_pnpn_t),
intent(inout) :: this
296 call this%scheme_free()
298 call this%bclst_ds%free()
299 call this%proj_s%free()
301 call this%s_adj_res%free()
303 call this%ds_adj%free()
305 call this%abx1%free()
306 call this%abx2%free()
308 call this%advs%free()
310 if (
allocated(this%Ax))
then
314 if (
allocated(this%res))
then
318 if (
allocated(this%makeext))
then
319 deallocate(this%makeext)
322 if (
allocated(this%makebdf))
then
323 deallocate(this%makebdf)
326 if (
allocated(this%makeoifs))
then
327 deallocate(this%makeoifs)
330 end subroutine adjoint_scalar_pnpn_free
332 subroutine adjoint_scalar_pnpn_step(this, t, tstep, dt, ext_bdf, &
334 class(adjoint_scalar_pnpn_t),
intent(inout) :: this
335 real(kind=rp),
intent(in) :: t
336 integer,
intent(in) :: tstep
337 real(kind=rp),
intent(in) :: dt
338 type(time_scheme_controller_t),
intent(in) :: ext_bdf
339 type(time_step_controller_t),
intent(in) :: dt_controller
343 type(ksp_monitor_t) :: ksp_results(1)
345 n = this%dm_Xh%size()
347 call profiler_start_region(
'Scalar', 2)
348 associate(u => this%u, v => this%v, w => this%w, s_adj => this%s_adj, &
349 cp => this%cp, rho => this%rho, &
350 ds_adj => this%ds_adj, &
351 s_adj_res => this%s_adj_res, &
352 ax => this%Ax, f_xh => this%f_Xh, xh => this%Xh, &
353 c_xh => this%c_Xh, dm_xh => this%dm_Xh, gs_xh => this%gs_Xh, &
354 s_adj_lag => this%s_adj_lag, oifs => this%oifs, &
355 lambda_field => this%lambda, &
356 projection_dim => this%projection_dim, &
357 msh => this%msh, res => this%res, makeoifs => this%makeoifs, &
358 makeext => this%makeext, makebdf => this%makebdf, &
359 if_variable_dt => dt_controller%if_variable_dt, &
360 dt_last_change => dt_controller%dt_last_change)
363 call print_debug(this)
365 call this%source_term%compute(t, tstep)
368 if (this%if_gradient_jump_penalty .eqv. .true.)
then
369 call neko_error(
"gradient jump not implemented for adjoint scalar")
375 call this%bcs%apply_scalar(this%f_Xh%x, dm_xh%size(), t, tstep, .false.)
389 call this%adv%compute_adjoint_scalar(u, v, w, s_adj, f_xh, &
390 xh, this%c_Xh, dm_xh%size())
397 call makeext%compute_scalar(this%abx1, this%abx2, f_xh%x, rho%x(1,1,1,1), &
398 ext_bdf%advection_coeffs, n)
401 call makebdf%compute_scalar(s_adj_lag, f_xh%x, s_adj, c_xh%B, rho%x(1,1,1,1), &
402 dt, ext_bdf%diffusion_coeffs, ext_bdf%ndiff, n)
405 call s_adj_lag%update()
408 call this%bcs%apply_scalar(this%s_adj%x, this%dm_Xh%size(), t, tstep, &
412 call this%update_material_properties(t, tstep)
415 call profiler_start_region(
'Adjoint_scalar_residual', 20)
416 call res%compute(ax, s_adj, s_adj_res, f_xh, c_xh, msh, xh, &
417 lambda_field, rho%x(1,1,1,1)*cp%x(1,1,1,1), ext_bdf%diffusion_coeffs(1), dt, dm_xh%size())
419 call gs_xh%op(s_adj_res, gs_op_add)
423 call this%bclst_ds%apply_scalar(s_adj_res%x, dm_xh%size())
425 call profiler_end_region(
'Adjoint_scalar_residual', 20)
427 call this%proj_s%pre_solving(s_adj_res%x, tstep, c_xh, n, dt_controller)
429 call this%pc%update()
430 call profiler_start_region(
'Adjoint_scalar_solve', 21)
431 ksp_results(1) = this%ksp%solve(ax, ds_adj, s_adj_res%x, n, &
432 c_xh, this%bclst_ds, gs_xh)
433 call profiler_end_region(
'Adjoint_scalar_solve', 21)
435 call this%proj_s%post_solving(ds_adj%x, ax, c_xh, this%bclst_ds, gs_xh, &
436 n, tstep, dt_controller)
439 if (neko_bcknd_device .eq. 1)
then
440 call device_add2s2(s_adj%x_d, ds_adj%x_d, 1.0_rp, n)
442 call add2s2(s_adj%x, ds_adj%x, 1.0_rp, n)
445 call scalar_step_info(tstep, t, dt, ksp_results)
448 call profiler_end_region(
'Scalar', 2)
449 end subroutine adjoint_scalar_pnpn_step
451 subroutine print_debug(this)
452 class(adjoint_scalar_pnpn_t),
intent(inout) :: this
455 n = this%dm_Xh%size()
466 end subroutine print_debug
471 subroutine adjoint_scalar_pnpn_setup_bcs_(this, user)
472 class(adjoint_scalar_pnpn_t),
intent(inout) :: this
473 type(user_t),
target,
intent(in) :: user
474 integer :: i, j, n_bcs, zone_size, global_zone_size, ierr
475 type(json_core) :: core
476 type(json_value),
pointer :: bc_object
477 type(json_file) :: bc_subdict
478 class(bc_t),
pointer :: bc_i
481 logical,
allocatable :: marked_zones(:)
482 integer,
allocatable :: zone_indices(:)
485 if (this%params%valid_path(
'case.adjoint_scalar.boundary_conditions'))
then
486 call this%params%info(
'case.adjoint_scalar.boundary_conditions', &
488 call this%params%get_core(core)
489 call this%params%get(
'case.adjoint_scalar.boundary_conditions', &
492 call this%bcs%init(n_bcs)
494 allocate(marked_zones(
size(this%msh%labeled_zones)))
495 marked_zones = .false.
499 call json_extract_item(core, bc_object, i, bc_subdict)
504 call json_get(bc_subdict,
"zone_indices", zone_indices)
506 do j = 1,
size(zone_indices)
507 zone_size = this%msh%labeled_zones(zone_indices(j))%size
508 call mpi_allreduce(zone_size, global_zone_size, 1, &
509 mpi_integer, mpi_max, neko_comm, ierr)
511 if (global_zone_size .eq. 0)
then
512 write(error_unit,
'(A, A, I0, A, A, I0, A)')
"*** ERROR ***: ",&
513 "Zone index ", zone_indices(j), &
514 " is invalid as this zone has 0 size, meaning it ", &
515 "does not in the mesh. Check adjoint scalar BC ", &
520 if (marked_zones(zone_indices(j)) .eqv. .true.)
then
521 write(error_unit,
'(A, A, I0, A, A, A, A)')
"*** ERROR ***: ", &
522 "Zone with index ", zone_indices(j), &
523 " has already been assigned a boundary condition. ", &
524 "Please check your boundary_conditions entry for the ", &
525 "adjoint scalar and make sure that each zone index ", &
526 "appears only in a single boundary condition."
529 marked_zones(zone_indices(j)) = .true.
535 call adjoint_bc_factory(bc_i, this, bc_subdict, this%c_Xh, user)
536 call this%bcs%append(bc_i)
540 do i = 1,
size(this%msh%labeled_zones)
541 if ((this%msh%labeled_zones(i)%size .gt. 0) .and. &
542 (marked_zones(i) .eqv. .false.))
then
543 write(error_unit,
'(A, A, I0)')
"*** ERROR ***: ", &
544 "No adjoint scalar boundary condition assigned to zone ", i
549 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.
Contains the factory routine for advection_t children.
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.