37 use,
intrinsic :: iso_fortran_env, only: error_unit
39 rhs_maker_ext_fctry, rhs_maker_bdf_fctry, rhs_maker_oifs_fctry
63 use json_module,
only : json_file, json_core, json_value
71 use mpi_f08,
only : mpi_allreduce, mpi_integer, mpi_max
85 class(
ax_t),
allocatable :: ax
148 module subroutine bc_factory(object, scheme, json, coef,
user)
149 class(
bc_t),
pointer,
intent(inout) :: object
150 type(scalar_pnpn_t),
intent(in) :: scheme
151 type(json_file),
intent(inout) :: json
152 type(
coef_t),
intent(in) :: coef
154 end subroutine bc_factory
171 subroutine scalar_pnpn_init(this, msh, coef, gs, params, numerics_params, &
172 user, chkp, ulag, vlag, wlag, time_scheme, rho)
173 class(scalar_pnpn_t),
target,
intent(inout) :: this
174 type(
mesh_t),
target,
intent(in) :: msh
175 type(
coef_t),
target,
intent(in) :: coef
176 type(
gs_t),
target,
intent(inout) :: gs
177 type(json_file),
target,
intent(inout) :: params
178 type(json_file),
target,
intent(inout) :: numerics_params
179 type(
user_t),
target,
intent(in) :: user
180 type(
chkp_t),
target,
intent(inout) :: chkp
183 type(
field_t),
target,
intent(in) :: rho
185 class(
bc_t),
pointer :: bc_i
186 character(len=15),
parameter :: scheme =
'Modular (Pn/Pn)'
192 call this%scheme_init(msh, coef, gs, params, scheme,
user, rho)
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_res%init(dm_xh,
"s_res")
215 call this%abx1%init(dm_xh, trim(this%name)//
"_abx1")
216 call neko_registry%add_field(dm_xh, trim(this%name)//
"_abx1", ignore_existing = .true.)
218 call this%abx2%init(dm_xh, trim(this%name)//
"_abx2")
219 call neko_registry%add_field(dm_xh, trim(this%name)//
"_abx2", ignore_existing = .true.)
221 call this%advs%init(dm_xh,
"advs")
223 call this%ds%init(dm_xh,
'ds')
228 call this%setup_bcs_(
user)
231 call this%bc_res%init(this%c_Xh, params)
232 do i = 1, this%bcs%size()
233 if (this%bcs%strong(i))
then
234 bc_i => this%bcs%get(i)
235 call this%bc_res%mark_facets(bc_i%marked_facet)
240 call this%bc_res%finalize()
242 call this%bclst_ds%init()
243 call this%bclst_ds%append(this%bc_res)
247 call this%proj_s%init(this%dm_Xh%size(), this%projection_dim, &
248 this%projection_activ_step)
257 call advection_factory(this%adv, numerics_params, this%c_Xh, &
258 ulag, vlag, wlag, this%chkp%dtlag, &
261 end subroutine scalar_pnpn_init
264 subroutine scalar_pnpn_restart(this, chkp)
265 class(scalar_pnpn_t),
target,
intent(inout) :: this
266 type(chkp_t),
intent(inout) :: chkp
267 real(kind=rp) :: dtlag(10), tlag(10)
269 type(field_t),
pointer :: temp_field
273 n = this%s%dof%size()
277 call col2(this%s%x, this%c_Xh%mult, n)
278 call col2(this%slag%lf(1)%x, this%c_Xh%mult, n)
279 call col2(this%slag%lf(2)%x, this%c_Xh%mult, n)
280 if (neko_bcknd_device .eq. 1)
then
281 call device_memcpy(this%s%x, this%s%x_d, &
282 n, host_to_device, sync = .false.)
283 call device_memcpy(this%slag%lf(1)%x, this%slag%lf(1)%x_d, &
284 n, host_to_device, sync = .false.)
285 call device_memcpy(this%slag%lf(2)%x, this%slag%lf(2)%x_d, &
286 n, host_to_device, sync = .false.)
287 call device_memcpy(this%abx1%x, this%abx1%x_d, &
288 n, host_to_device, sync = .false.)
289 call device_memcpy(this%abx2%x, this%abx2%x_d, &
290 n, host_to_device, sync = .false.)
291 call device_memcpy(this%advs%x, this%advs%x_d, &
292 n, host_to_device, sync = .false.)
295 call this%gs_Xh%op(this%s, gs_op_add)
296 call this%gs_Xh%op(this%slag%lf(1), gs_op_add)
297 call this%gs_Xh%op(this%slag%lf(2), gs_op_add)
299 end subroutine scalar_pnpn_restart
301 subroutine scalar_pnpn_free(this)
302 class(scalar_pnpn_t),
intent(inout) :: this
305 call this%scheme_free()
307 call this%bc_res%free()
308 call this%bclst_ds%free()
309 call this%proj_s%free()
311 call this%s_res%free()
315 call this%abx1%free()
316 call this%abx2%free()
318 call this%advs%free()
320 if (
allocated(this%adv))
then
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 scalar_pnpn_free
347 subroutine scalar_pnpn_step(this, time, ext_bdf, dt_controller, &
349 class(scalar_pnpn_t),
intent(inout) :: this
350 type(time_state_t),
intent(in) :: time
351 type(time_scheme_controller_t),
intent(in) :: ext_bdf
352 type(time_step_controller_t),
intent(in) :: dt_controller
353 type(ksp_monitor_t),
intent(inout) :: ksp_results
357 if (this%freeze)
return
359 n = this%dm_Xh%size()
361 call profiler_start_region(trim(this%name), 2)
362 associate(u => this%u, v => this%v, w => this%w, s => this%s, &
363 cp => this%cp, rho => this%rho, lambda_tot => this%lambda_tot, &
365 s_res => this%s_res, &
366 ax => this%Ax, f_xh => this%f_Xh, xh => this%Xh, &
367 c_xh => this%c_Xh, dm_xh => this%dm_Xh, gs_xh => this%gs_Xh, &
368 slag => this%slag, oifs => this%oifs, &
369 projection_dim => this%projection_dim, &
370 msh => this%msh, res => this%res, makeoifs => this%makeoifs, &
371 makeext => this%makeext, makebdf => this%makebdf, &
372 t => time%t, tstep => time%tstep, dt => time%dt)
375 call print_debug(this)
377 call this%source_term%compute(time)
380 call this%bcs%apply_scalar(this%f_Xh%x, dm_xh%size(), time, .false.)
384 call this%adv%compute_scalar(u, v, w, s, this%advs, &
385 xh, this%c_Xh, dm_xh%size())
387 call makeext%compute_scalar(this%abx1, this%abx2, f_xh%x, &
388 rho%x(1,1,1,1), ext_bdf%advection_coeffs, n)
390 call makeoifs%compute_scalar(this%advs%x, f_xh%x, rho%x(1,1,1,1), dt,&
394 call this%adv%compute_scalar(u, v, w, s, f_xh, &
395 xh, this%c_Xh, dm_xh%size())
401 call makeext%compute_scalar(this%abx1, this%abx2, f_xh%x, &
402 rho%x(1,1,1,1), ext_bdf%advection_coeffs, n)
405 call makebdf%compute_scalar(slag, f_xh%x, s, c_xh%B, rho%x(1,1,1,1), &
406 dt, ext_bdf%diffusion_coeffs, ext_bdf%ndiff, n)
412 call this%apply_strong_bcs(time)
415 call this%update_material_properties(time)
418 call profiler_start_region(trim(this%name) //
'_residual', 20)
419 call res%compute(ax, s, s_res, f_xh, c_xh, msh, xh, lambda_tot, &
420 rho%x(1,1,1,1)*cp%x(1,1,1,1), ext_bdf%diffusion_coeffs(1), dt, &
423 call gs_xh%op(s_res, gs_op_add)
426 call this%bclst_ds%apply_scalar(s_res%x, dm_xh%size())
428 call profiler_end_region(trim(this%name) //
'_residual', 20)
430 call this%proj_s%pre_solving(s_res%x, tstep, c_xh, n, dt_controller)
432 call this%pc%update()
433 call profiler_start_region(trim(this%name) //
'_solve', 21)
434 ksp_results = this%ksp%solve(ax, ds, s_res%x, n, &
435 c_xh, this%bclst_ds, gs_xh)
436 ksp_results%name = trim(this%name)
437 call profiler_end_region(trim(this%name) //
'_solve', 21)
439 call this%proj_s%post_solving(ds%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%x_d, ds%x_d, 1.0_rp, n)
446 call add2s2(s%x, ds%x, 1.0_rp, n)
450 call profiler_end_region(trim(this%name), 2)
451 end subroutine scalar_pnpn_step
453 subroutine print_debug(this)
454 class(scalar_pnpn_t),
intent(inout) :: this
455 character(len=LOG_SIZE) :: log_buf
458 n = this%dm_Xh%size()
460 write(log_buf,
'(A, A, E15.7, A, E15.7, A, E15.7)')
'Scalar debug', &
461 ' l2norm s', glsc2(this%s%x, this%s%x, n), &
462 ' slag1', glsc2(this%slag%lf(1)%x, this%slag%lf(1)%x, n), &
463 ' slag2', glsc2(this%slag%lf(2)%x, this%slag%lf(2)%x, n)
464 call neko_log%message(log_buf, lvl = neko_log_debug)
465 write(log_buf,
'(A, A, E15.7, A, E15.7)')
'Scalar debug2', &
466 ' l2norm abx1', glsc2(this%abx1%x, this%abx1%x, n), &
467 ' abx2', glsc2(this%abx2%x, this%abx2%x, n)
468 call neko_log%message(log_buf, lvl = neko_log_debug)
469 end subroutine print_debug
473 subroutine scalar_pnpn_setup_bcs_(this, user)
474 class(scalar_pnpn_t),
intent(inout) :: this
475 type(user_t),
target,
intent(in) :: user
476 integer :: i, j, n_bcs, zone_size, global_zone_size, ierr
477 type(json_core) :: core
478 type(json_value),
pointer :: bc_object
479 type(json_file) :: bc_subdict
480 class(bc_t),
pointer :: bc_i
483 logical,
allocatable :: marked_zones(:)
484 integer,
allocatable :: zone_indices(:)
486 if (this%params%valid_path(
'boundary_conditions'))
then
487 call this%params%info(
'boundary_conditions', &
489 call this%params%get_core(core)
490 call this%params%get(
'boundary_conditions', bc_object, found)
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 exist in the mesh. Check scalar boundary condition ", &
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 "scalar and make sure that each zone index appears only ",&
526 "in a single boundary condition."
529 marked_zones(zone_indices(j)) = .true.
535 call 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 scalar boundary condition assigned to zone ", i
550 do i = 1,
size(this%msh%labeled_zones)
551 if (this%msh%labeled_zones(i)%size .gt. 0)
then
552 write(error_unit,
'(A, A, A)')
"*** ERROR ***: ", &
553 "No boundary_conditions entry in the case file for scalar ", &
564 end subroutine scalar_pnpn_setup_bcs_
568 subroutine scalar_scheme_apply_strong_bcs(this, time)
569 class(scalar_pnpn_t),
intent(inout) :: this
570 type(time_state_t),
intent(in) :: time
573 class(bc_t),
pointer :: bc_i
577 call this%bcs%apply(this%s, time = time, strong = .true.)
583 call this%gs_Xh%op(this%s, gs_op_min, glb_cmd_event)
584 call device_event_sync(glb_cmd_event)
588 call this%bcs%apply(this%s, time = time, strong = .true.)
591 call this%gs_Xh%op(this%s, gs_op_max, glb_cmd_event)
592 call device_event_sync(glb_cmd_event)
595 do i = 1, this%bcs%size()
596 bc_i => this%bcs%get(i)
597 bc_i%updated = .false.
601 end subroutine scalar_scheme_apply_strong_bcs
Copy data between host and device (or device and device)
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Retrieves a parameter by name or throws an error.
Subroutines to add advection terms to the RHS of a transport equation.
Defines a Matrix-vector product.
Defines a boundary condition.
type(mpi_comm), public neko_comm
MPI communicator.
subroutine, public device_add2s2(a_d, b_d, c1, n, strm)
Vector addition with scalar multiplication (multiplication on first argument)
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
Device abstraction, common interface for various accelerators.
subroutine, public device_event_sync(event)
Synchronize an event.
integer, parameter, public host_to_device
type(c_ptr), bind(C), public glb_cmd_event
Event for the global command queue.
Dirichlet condition applied in the facet normal direction.
Contains the field_serties_t type.
Utilities for retrieving parameters from the case files.
Implements the base abstract type for Krylov solvers plus helper types.
integer, parameter, public neko_log_debug
Debug log level.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
subroutine, public col2(a, b, n)
Vector multiplication .
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
subroutine, public profiler_end_region(name, region_id)
End the most recently started profiler region.
Project x onto X, the space of old solutions and back again.
Defines a registry for storing solution fields.
type(registry_t), target, public neko_registry
Global field registry.
Routines to generate the right-hand sides for the convection-diffusion equation. Employs the EXT/BDF ...
Contains the scalar_pnpn_t type.
subroutine scalar_pnpn_step(this, time, ext_bdf, dt_controller, ksp_results)
subroutine scalar_pnpn_setup_bcs_(this, user)
Initialize boundary conditions.
subroutine scalar_scheme_apply_strong_bcs(this, time)
Apply strong boundary conditions.
subroutine scalar_pnpn_restart(this, chkp)
subroutine scalar_pnpn_init(this, msh, coef, gs, params, numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
Boundary condition factory. Both constructs and initializes the object.
subroutine scalar_pnpn_free(this)
Defines the residual for the scalar transport equation.
Contains the scalar_scheme_t type.
Compound scheme for the advection and diffusion operators in a transport equation.
Base class for time integration schemes.
Module with things related to the simulation time.
Implements type time_step_controller.
Interfaces for user interaction with NEKO.
Defines a zero-valued Dirichlet boundary condition.
Base abstract type for computing the advection operator.
Base type for a matrix-vector product providing .
Base type for a boundary condition.
A list of allocatable `bc_t`. Follows the standard interface of lists.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Dirichlet condition in facet normal direction.
Stores a series (sequence) of fields, logically connected to a base field, and arranged according to ...
Type for storing initial and final residuals in a Krylov solver.
Abstract type to add contributions to F from lagged BD terms.
Abstract type to sum up contributions to kth order extrapolation scheme.
Abstract type to add contributions of kth order OIFS scheme.
Abstract type to compute scalar residual.
Base type for a scalar advection-diffusion solver.
Implements the logic to compute the time coefficients for the advection and diffusion operators in a ...
A struct that contains all info about the time, expand as needed.
Provides a tool to set time step dt.
A type collecting all the overridable user routines and flag to suppress type injection from custom m...
Zero-valued Dirichlet boundary condition. Used for no-slip walls, but also for various auxillary cond...