Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
scalars.f90
Go to the documentation of this file.
1! Copyright (c) 2022-2025, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
34
35module scalars
36 use num_types, only: rp
37 use scalar_pnpn, only: scalar_pnpn_t
40 use mesh, only: mesh_t
41 use space, only: space_t
42 use gather_scatter, only: gs_t
45 use json_module, only: json_file
48 use field, only: field_t
49 use field_list, only: field_list_t
52 use checkpoint, only: chkp_t
53 use krylov, only: ksp_t, ksp_monitor_t
55 use user_intf, only: user_t
56 use utils, only: neko_error
57 use coefs, only : coef_t
58 use time_state, only : time_state_t
59 implicit none
60 private
61
63 type, public :: scalars_t
65 class(scalar_scheme_t), allocatable :: scalar_fields(:)
67 class(ksp_t), allocatable :: shared_ksp
68 contains
71 procedure, private :: scalars_init
72 procedure, private :: scalars_init_single
74 procedure :: step => scalars_step
76 procedure :: restart => scalars_restart
78 procedure :: validate => scalars_validate
80 procedure :: free => scalars_free
82 procedure, private :: register_lags_with_checkpoint
83 end type scalars_t
84
85contains
86
88 subroutine scalars_init(this, n_scalars, msh, coef, gs, params, &
89 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
90 class(scalars_t), intent(inout) :: this
91 integer, intent(in) :: n_scalars
92 type(mesh_t), target, intent(in) :: msh
93 type(coef_t), target, intent(in) :: coef
94 type(gs_t), target, intent(inout) :: gs
95 type(json_file), target, intent(inout) :: params
96 type(json_file), target, intent(inout) :: numerics_params
97 type(user_t), target, intent(in) :: user
98 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
99 type(time_scheme_controller_t), target, intent(in) :: time_scheme
100 TYPE(field_t), TARGET, INTENT(IN) :: rho
101 type(chkp_t), target, intent(inout) :: chkp
102 type(json_file) :: json_subdict
103 integer :: i, j
104 character(len=:), allocatable :: field_name
105 character(len=:), allocatable :: field_names(:)
106 character(len=256) :: error_msg
107
108 ! Allocate the scalar fields
109 ! If there are more scalar_scheme_t types, add a factory function here
110 allocate(scalar_pnpn_t::this%scalar_fields(n_scalars))
111
112 ! Collect and validate field names for all scalars
113 allocate(character(len=256) :: field_names(n_scalars))
114
115 do i = 1, n_scalars
116 ! Extract element i from the "scalars" array
117 call json_extract_item(params, "", i, json_subdict)
118
119 ! Try to get name from JSON, generate one if not found or empty
120 if (json_subdict%valid_path('name')) then
121 call json_get(json_subdict, 'name', field_name)
122 else
123 field_name = ''
124 end if
125
126 ! If name is empty or not provided, generate a default one
127 if (len_trim(field_name) == 0) then
128 if (n_scalars == 1) then
129 field_name = 's' ! Single scalar gets default name 's'
130 else
131 write(field_name, '(A,I0)') 's_', i
132 end if
133 end if
134
135 field_names(i) = trim(field_name)
136
137 ! If there's a duplicate, append a number until unique
138 if (n_scalars > 1) then
139 j = 1
140 do while (j < i)
141 if (trim(field_names(i)) == trim(field_names(j))) then
142 write(field_name, '(A,I0)') trim(field_names(i))//'_', j
143 field_names(i) = trim(field_name)
144 j = 1 ! Start over to check if new name is unique
145 else
146 j = j + 1
147 end if
148 end do
149 end if
150 end do
151
152 do i = 1, n_scalars
153 call json_extract_item(params, "", i, json_subdict)
154
155 ! Use the processed field names for all scalars
156 call json_subdict%add('name', trim(field_names(i)))
157
158 call this%scalar_fields(i)%init(msh, coef, gs, json_subdict, &
159 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
160 end do
161
162 ! Register all scalar lag fields with checkpoint using scalable approach
163 if (n_scalars > 1) then
164 call this%register_lags_with_checkpoint(chkp)
165 else
166 ! For single scalar, use legacy interface
167 call chkp%add_scalar(this%scalar_fields(1)%s, this%scalar_fields(1)%slag)
168 end if
169 end subroutine scalars_init
170
171 subroutine scalars_init_single(this, msh, coef, gs, params, numerics_params, &
172 user, chkp, ulag, vlag, wlag, time_scheme, rho)
173 class(scalars_t), 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
181 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
182 type(time_scheme_controller_t), target, intent(in) :: time_scheme
183 TYPE(field_t), TARGET, INTENT(IN) :: rho
184
185 ! Allocate a single scalar field
186 allocate(scalar_pnpn_t::this%scalar_fields(1))
187
188 ! Set the scalar name to "s"
189 if (.not. params%valid_path('name')) then
190 call params%add('name', 's')
191 end if
192
193 ! Initialize it directly with the params
194 call this%scalar_fields(1)%init(msh, coef, gs, params, numerics_params, &
195 user, chkp, ulag, vlag, wlag, time_scheme, rho)
196
197 ! Register single scalar with checkpoint
198 call chkp%add_scalar(this%scalar_fields(1)%s, this%scalar_fields(1)%slag)
199 end subroutine scalars_init_single
200
202 subroutine scalars_step(this, time, ext_bdf, dt_controller)
203 class(scalars_t), intent(inout) :: this
204 type(time_state_t), intent(in) :: time
205 type(time_scheme_controller_t), intent(inout) :: ext_bdf
206 type(time_step_controller_t), intent(inout) :: dt_controller
207 integer :: i
208 type(ksp_monitor_t), dimension(size(this%scalar_fields)) :: ksp_results
209
210 ! Iterate through all scalar fields
211 do i = 1, size(this%scalar_fields)
212 call this%scalar_fields(i)%step(time, ext_bdf, dt_controller, &
213 ksp_results(i))
214 end do
215
216 call scalar_step_info(time, ksp_results)
217 end subroutine scalars_step
218
220 subroutine scalars_restart(this, chkp)
221 class(scalars_t), intent(inout) :: this
222 type(chkp_t), intent(inout) :: chkp
223 integer :: i, n_scalars
224
225 n_scalars = size(this%scalar_fields)
226 do i = 1, size(this%scalar_fields)
227 call this%scalar_fields(i)%restart(chkp)
228 end do
229 end subroutine scalars_restart
230
232 subroutine scalars_validate(this)
233 class(scalars_t), intent(inout) :: this
234 integer :: i
235 ! Iterate through all scalar fields
236 do i = 1, size(this%scalar_fields)
237 call this%scalar_fields(i)%slag%set(this%scalar_fields(i)%s)
238 call this%scalar_fields(i)%validate()
239 end do
240 end subroutine scalars_validate
241
243 subroutine scalars_free(this)
244 class(scalars_t), intent(inout) :: this
245 integer :: i
246
247 ! Iterate through all scalar fields
248 if (allocated(this%scalar_fields)) then
249 do i = 1, size(this%scalar_fields)
250 call this%scalar_fields(i)%free()
251 end do
252 deallocate(this%scalar_fields)
253 end if
254 end subroutine scalars_free
255
257 subroutine register_lags_with_checkpoint(this, chkp)
258 class(scalars_t), intent(inout) :: this
259 type(chkp_t), intent(inout) :: chkp
260 integer :: i, n_scalars
261
262 n_scalars = size(this%scalar_fields)
263
264 ! Allocate ABX field arrays
265 allocate(chkp%scalar_abx1(n_scalars))
266 allocate(chkp%scalar_abx2(n_scalars))
267
268 ! Add all scalar lag fields to the checkpoint list and populate ABX fields
269 do i = 1, n_scalars
270 call chkp%scalar_lags%append(this%scalar_fields(i)%slag)
271
272 ! Cast to scalar_pnpn_t to access ABX fields
273 select type(scalar_field => this%scalar_fields(i))
274 type is(scalar_pnpn_t)
275 call associate_scalar_abx_fields(chkp, i, scalar_field)
276 end select
277 end do
278
279 end subroutine register_lags_with_checkpoint
280
282 subroutine associate_scalar_abx_fields(chkp, index, scalar_field)
283 type(chkp_t), intent(inout) :: chkp
284 integer, intent(in) :: index
285 type(scalar_pnpn_t), target, intent(in) :: scalar_field
286
287 chkp%scalar_abx1(index)%ptr => scalar_field%abx1
288 chkp%scalar_abx2(index)%ptr => scalar_field%abx2
289 end subroutine associate_scalar_abx_fields
290
291end module scalars
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.
Defines a checkpoint.
Coefficients.
Definition coef.f90:34
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Contains the field_serties_t type.
Defines a field.
Definition field.f90:34
Gather-scatter.
Utilities for retrieving parameters from the case files.
subroutine, public json_extract_object(json, name, object)
Extract object as a separate JSON dictionary.
Implements the base abstract type for Krylov solvers plus helper types.
Definition krylov.f90:34
Logging routines.
Definition log.f90:34
integer, parameter, public neko_log_verbose
Verbose log level.
Definition log.f90:76
type(log_t), public neko_log
Global log stream.
Definition log.f90:70
integer, parameter, public log_size
Definition log.f90:46
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Auxiliary routines for fluid solvers.
Definition scalar_aux.f90:2
subroutine scalar_step_info(time, ksp_results, strict_convergence)
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.
Definition scalars.f90:35
subroutine scalars_validate(this)
Check if the configuration is valid.
Definition scalars.f90:233
subroutine associate_scalar_abx_fields(chkp, index, scalar_field)
Helper subroutine to associate ABX field pointers with proper TARGET attribute.
Definition scalars.f90:283
subroutine register_lags_with_checkpoint(this, chkp)
Register scalar lag fields with checkpoint.
Definition scalars.f90:258
subroutine scalars_init(this, n_scalars, msh, coef, gs, params, numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
Initialize the scalars container.
Definition scalars.f90:90
subroutine scalars_init_single(this, msh, coef, gs, params, numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
Definition scalars.f90:173
subroutine scalars_free(this)
Clean up all resources.
Definition scalars.f90:244
subroutine scalars_step(this, time, ext_bdf, dt_controller)
Perform a time step for all scalar fields.
Definition scalars.f90:203
subroutine scalars_restart(this, chkp)
Restart from checkpoint data.
Definition scalars.f90:221
Defines a function space.
Definition space.f90:34
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.
Definition user_intf.f90:34
Utilities.
Definition utils.f90:35
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
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.
Definition krylov.f90:56
Base abstract type for a canonical Krylov method, solving .
Definition krylov.f90:73
Base type for a scalar advection-diffusion solver.
Type to manage multiple scalar transport equations.
Definition scalars.f90:63
The function space for the SEM solution fields.
Definition space.f90:62
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.
A type collecting all the overridable user routines and flag to suppress type injection from custom m...