45 use json_module,
only: json_file
65 class(
ksp_t),
allocatable :: shared_ksp
87 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
89 integer,
intent(in) :: n_scalars
90 type(
mesh_t),
target,
intent(in) :: msh
91 type(
coef_t),
target,
intent(in) :: coef
92 type(
gs_t),
target,
intent(inout) :: gs
93 type(json_file),
target,
intent(inout) :: params
94 type(json_file),
target,
intent(inout) :: numerics_params
95 type(
user_t),
target,
intent(in) :: user
98 TYPE(
field_t),
TARGET,
INTENT(IN) :: rho
99 type(
chkp_t),
target,
intent(inout) :: chkp
100 type(json_file) :: json_subdict
102 character(len=:),
allocatable :: field_name
103 character(len=:),
allocatable :: field_names(:)
104 character(len=256) :: error_msg, buffer
111 allocate(
character(len=256) :: field_names(n_scalars))
121 if (len_trim(field_name) == 0)
then
122 if (n_scalars == 1)
then
125 write(
buffer,
'(A,I0)')
's_', i
130 field_names(i) = trim(field_name)
133 if (n_scalars > 1)
then
136 if (trim(field_names(i)) == trim(field_names(j)))
then
137 call neko_error(
"Duplicate scalar field name detected: "// &
138 trim(field_names(i)) // &
139 ". Please provide unique names for each scalar field.")
151 call json_subdict%add(
'name', trim(field_names(i)))
153 call this%scalar_fields(i)%init(msh, coef, gs, json_subdict, &
158 if (n_scalars > 1)
then
159 call this%register_lags_with_checkpoint(chkp)
162 select type(scalar => this%scalar_fields(1))
164 call chkp%add_scalar(scalar%s, scalar%slag, scalar%abx1, scalar%abx2)
170 user, chkp, ulag, vlag, wlag, time_scheme, rho)
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
176 type(json_file),
target,
intent(inout) :: numerics_params
177 type(
user_t),
target,
intent(in) :: user
178 type(
chkp_t),
target,
intent(inout) :: chkp
181 TYPE(
field_t),
TARGET,
INTENT(IN) :: rho
187 if (.not. params%valid_path(
'name'))
then
188 call params%add(
'name',
's')
192 call this%scalar_fields(1)%init(msh, coef, gs, params, numerics_params, &
196 select type(scalar => this%scalar_fields(1))
198 call chkp%add_scalar(scalar%s, scalar%slag, scalar%abx1, scalar%abx2)
209 type(
ksp_monitor_t),
dimension(size(this%scalar_fields)) :: ksp_results
210 logical :: all_frozen
215 do i = 1,
size(this%scalar_fields)
216 all_frozen = all_frozen .and. this%scalar_fields(i)%freeze
217 call this%scalar_fields(i)%step(time, ext_bdf, dt_controller, &
221 if (.not. all_frozen)
then
222 call ksp_results(i)%print_header()
225 do i = 1,
size(this%scalar_fields)
226 if (this%scalar_fields(i)%freeze) cycle
234 type(
chkp_t),
intent(inout) :: chkp
235 integer :: i, n_scalars
237 n_scalars =
size(this%scalar_fields)
238 do i = 1,
size(this%scalar_fields)
239 call this%scalar_fields(i)%restart(chkp)
248 do i = 1,
size(this%scalar_fields)
249 call this%scalar_fields(i)%slag%set(this%scalar_fields(i)%s)
250 call this%scalar_fields(i)%validate()
260 if (
allocated(this%scalar_fields))
then
261 do i = 1,
size(this%scalar_fields)
262 call this%scalar_fields(i)%free()
264 deallocate(this%scalar_fields)
267 if (
allocated(this%shared_ksp))
then
268 call this%shared_ksp%free()
269 deallocate(this%shared_ksp)
276 type(
chkp_t),
intent(inout) :: chkp
277 integer :: i, n_scalars
279 n_scalars =
size(this%scalar_fields)
282 allocate(chkp%scalar_abx1(n_scalars))
283 allocate(chkp%scalar_abx2(n_scalars))
287 call chkp%scalar_lags%append(this%scalar_fields(i)%slag)
290 select type(scalar_field => this%scalar_fields(i))
300 type(
chkp_t),
intent(inout) :: chkp
301 integer,
intent(in) :: index
304 chkp%scalar_abx1(index)%ptr => scalar_field%abx1
305 chkp%scalar_abx2(index)%ptr => scalar_field%abx2
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.
Generic buffer that is extended with buffers of varying rank.
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_verbose
Verbose log level.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
integer, parameter, public rp
Global precision used in computations.
Auxiliary routines for fluid solvers.
subroutine, public scalar_step_info(time, ksp_results, strict_convergence, allow_stabilization)
Prints for prs, velx, vely, velz the following: Number of iterations, start residual,...
Contains the scalar_pnpn_t type.
Contains the scalar_scheme_t type.
Contains the scalars_t type that manages multiple scalar fields.
subroutine scalars_validate(this)
Check if the configuration is valid.
subroutine associate_scalar_abx_fields(chkp, index, scalar_field)
Helper subroutine to associate ABX field pointers with proper TARGET attribute.
subroutine register_lags_with_checkpoint(this, chkp)
Register scalar lag fields with checkpoint.
subroutine scalars_init(this, n_scalars, msh, coef, gs, params, numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
Initialize the scalars container.
subroutine scalars_init_single(this, msh, coef, gs, params, numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
subroutine scalars_free(this)
Clean up all resources.
subroutine scalars_step(this, time, ext_bdf, dt_controller)
Perform a time step for all scalar fields.
subroutine scalars_restart(this, chkp)
Restart from checkpoint data.
Defines a function space.
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.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
field_list_t, To be able to group fields together
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.
Base abstract type for a canonical Krylov method, solving .
Base type for a scalar advection-diffusion solver.
Type to manage multiple scalar transport equations.
The function space for the SEM solution fields.
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...