Neko 1.99.2
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
47 use field, only: field_t
48 use field_list, only: field_list_t
50 use checkpoint, only: chkp_t
51 use krylov, only: ksp_t, ksp_monitor_t
53 use user_intf, only: user_t
54 use utils, only: neko_error
55 use coefs, only : coef_t
56 use time_state, only : time_state_t
57 implicit none
58 private
59
61 type, public :: scalars_t
63 class(scalar_scheme_t), allocatable :: scalar_fields(:)
65 class(ksp_t), allocatable :: shared_ksp
66 contains
69 procedure, private :: scalars_init
70 procedure, private :: scalars_init_single
72 procedure :: step => scalars_step
74 procedure :: restart => scalars_restart
76 procedure :: validate => scalars_validate
78 procedure :: free => scalars_free
80 procedure, private :: register_lags_with_checkpoint
81 end type scalars_t
82
83contains
84
86 subroutine scalars_init(this, n_scalars, msh, coef, gs, params, &
87 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
88 class(scalars_t), intent(inout) :: this
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
96 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
97 type(time_scheme_controller_t), target, intent(in) :: time_scheme
98 TYPE(field_t), TARGET, INTENT(IN) :: rho
99 type(chkp_t), target, intent(inout) :: chkp
100 type(json_file) :: json_subdict
101 integer :: i, j
102 character(len=:), allocatable :: field_name
103 character(len=:), allocatable :: field_names(:)
104 character(len=256) :: error_msg, buffer
105
106 ! Allocate the scalar fields
107 ! If there are more scalar_scheme_t types, add a factory function here
108 allocate(scalar_pnpn_t::this%scalar_fields(n_scalars))
109
110 ! Collect and validate field names for all scalars
111 allocate(character(len=256) :: field_names(n_scalars))
112
113 do i = 1, n_scalars
114 ! Extract element i from the "scalars" array
115 call json_extract_item(params, "", i, json_subdict)
116
117 ! Try to get name from JSON, generate one if not found or empty
118 call json_get_or_default(json_subdict, 'name', field_name, '')
119
120 ! If name is empty or not provided, generate a default one
121 if (len_trim(field_name) == 0) then
122 if (n_scalars == 1) then
123 field_name = 's' ! Single scalar gets default name 's'
124 else
125 write(buffer, '(A,I0)') 's_', i
126 field_name = trim(buffer)
127 end if
128 end if
129
130 field_names(i) = trim(field_name)
131
132 ! If there's a duplicate, append a number until unique
133 if (n_scalars > 1) then
134 j = 1
135 do while (j < i)
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.")
140 else
141 j = j + 1
142 end if
143 end do
144 end if
145 end do
146
147 do i = 1, n_scalars
148 call json_extract_item(params, "", i, json_subdict)
149
150 ! Use the processed field names for all scalars
151 call json_subdict%add('name', trim(field_names(i)))
152
153 call this%scalar_fields(i)%init(msh, coef, gs, json_subdict, &
154 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
155 end do
156
157 ! Register all scalar lag fields with checkpoint using scalable approach
158 if (n_scalars > 1) then
159 call this%register_lags_with_checkpoint(chkp)
160 else
161 ! For single scalar, use legacy interface
162 select type(scalar => this%scalar_fields(1))
163 type is (scalar_pnpn_t)
164 call chkp%add_scalar(scalar%s, scalar%slag, scalar%abx1, scalar%abx2)
165 end select
166 end if
167 end subroutine scalars_init
168
169 subroutine scalars_init_single(this, msh, coef, gs, params, numerics_params, &
170 user, chkp, ulag, vlag, wlag, time_scheme, rho)
171 class(scalars_t), 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
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
179 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
180 type(time_scheme_controller_t), target, intent(in) :: time_scheme
181 TYPE(field_t), TARGET, INTENT(IN) :: rho
182
183 ! Allocate a single scalar field
184 allocate(scalar_pnpn_t::this%scalar_fields(1))
185
186 ! Set the scalar name to "s"
187 if (.not. params%valid_path('name')) then
188 call params%add('name', 's')
189 end if
190
191 ! Initialize it directly with the params
192 call this%scalar_fields(1)%init(msh, coef, gs, params, numerics_params, &
193 user, chkp, ulag, vlag, wlag, time_scheme, rho)
194
195 ! Register single scalar with checkpoint
196 select type(scalar => this%scalar_fields(1))
197 type is (scalar_pnpn_t)
198 call chkp%add_scalar(scalar%s, scalar%slag, scalar%abx1, scalar%abx2)
199 end select
200 end subroutine scalars_init_single
201
203 subroutine scalars_step(this, time, ext_bdf, dt_controller)
204 class(scalars_t), intent(inout) :: this
205 type(time_state_t), intent(in) :: time
206 type(time_scheme_controller_t), intent(inout) :: ext_bdf
207 type(time_step_controller_t), intent(inout) :: dt_controller
208 integer :: i
209 type(ksp_monitor_t), dimension(size(this%scalar_fields)) :: ksp_results
210 logical :: all_frozen
211
212 all_frozen = .true.
213
214 ! Iterate through all scalar fields
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, &
218 ksp_results(i))
219 end do
220
221 if (.not. all_frozen) then
222 call ksp_results(i)%print_header()
223 end if
224
225 do i = 1, size(this%scalar_fields)
226 if (this%scalar_fields(i)%freeze) cycle
227 call scalar_step_info(time, ksp_results(i))
228 end do
229 end subroutine scalars_step
230
232 subroutine scalars_restart(this, chkp)
233 class(scalars_t), intent(inout) :: this
234 type(chkp_t), intent(inout) :: chkp
235 integer :: i, n_scalars
236
237 n_scalars = size(this%scalar_fields)
238 do i = 1, size(this%scalar_fields)
239 call this%scalar_fields(i)%restart(chkp)
240 end do
241 end subroutine scalars_restart
242
244 subroutine scalars_validate(this)
245 class(scalars_t), intent(inout) :: this
246 integer :: i
247 ! Iterate through all scalar fields
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()
251 end do
252 end subroutine scalars_validate
253
255 subroutine scalars_free(this)
256 class(scalars_t), intent(inout) :: this
257 integer :: i
258
259 ! Iterate through all scalar fields
260 if (allocated(this%scalar_fields)) then
261 do i = 1, size(this%scalar_fields)
262 call this%scalar_fields(i)%free()
263 end do
264 deallocate(this%scalar_fields)
265 end if
266
267 if (allocated(this%shared_ksp)) then
268 call this%shared_ksp%free()
269 deallocate(this%shared_ksp)
270 end if
271 end subroutine scalars_free
272
274 subroutine register_lags_with_checkpoint(this, chkp)
275 class(scalars_t), intent(inout) :: this
276 type(chkp_t), intent(inout) :: chkp
277 integer :: i, n_scalars
278
279 n_scalars = size(this%scalar_fields)
280
281 ! Allocate ABX field arrays
282 allocate(chkp%scalar_abx1(n_scalars))
283 allocate(chkp%scalar_abx2(n_scalars))
284
285 ! Add all scalar lag fields to the checkpoint list and populate ABX fields
286 do i = 1, n_scalars
287 call chkp%scalar_lags%append(this%scalar_fields(i)%slag)
288
289 ! Cast to scalar_pnpn_t to access ABX fields
290 select type(scalar_field => this%scalar_fields(i))
291 type is(scalar_pnpn_t)
292 call associate_scalar_abx_fields(chkp, i, scalar_field)
293 end select
294 end do
295
296 end subroutine register_lags_with_checkpoint
297
299 subroutine associate_scalar_abx_fields(chkp, index, scalar_field)
300 type(chkp_t), intent(inout) :: chkp
301 integer, intent(in) :: index
302 type(scalar_pnpn_t), target, intent(in) :: scalar_field
303
304 chkp%scalar_abx1(index)%ptr => scalar_field%abx1
305 chkp%scalar_abx2(index)%ptr => scalar_field%abx2
306 end subroutine associate_scalar_abx_fields
307
308end 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.
Generic buffer that is extended with buffers of varying rank.
Definition buffer.F90:34
Defines a checkpoint.
Coefficients.
Definition coef.f90:34
Contains the field_serties_t type.
Defines a field.
Definition field.f90:34
Gather-scatter.
Utilities for retrieving parameters from the case files.
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:82
type(log_t), public neko_log
Global log stream.
Definition log.f90:76
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.
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.
Definition scalars.f90:35
subroutine scalars_validate(this)
Check if the configuration is valid.
Definition scalars.f90:245
subroutine associate_scalar_abx_fields(chkp, index, scalar_field)
Helper subroutine to associate ABX field pointers with proper TARGET attribute.
Definition scalars.f90:300
subroutine register_lags_with_checkpoint(this, chkp)
Register scalar lag fields with checkpoint.
Definition scalars.f90:275
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:88
subroutine scalars_init_single(this, msh, coef, gs, params, numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
Definition scalars.f90:171
subroutine scalars_free(this)
Clean up all resources.
Definition scalars.f90:256
subroutine scalars_step(this, time, ext_bdf, dt_controller)
Perform a time step for all scalar fields.
Definition scalars.f90:204
subroutine scalars_restart(this, chkp)
Restart from checkpoint data.
Definition scalars.f90:233
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:56
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:61
The function space for the SEM solution fields.
Definition space.f90:63
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...