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
47 use field, only: field_t
48 use field_list, only: field_list_t
51 use checkpoint, only: chkp_t
52 use krylov, only: ksp_t, ksp_monitor_t
54 use user_intf, only: user_t
55 use utils, only: neko_error
56 use coefs, only : coef_t
57 use time_state, only : time_state_t
58 implicit none
59 private
60
62 type, public :: scalars_t
64 class(scalar_scheme_t), allocatable :: scalar_fields(:)
66 class(ksp_t), allocatable :: shared_ksp
67 contains
70 procedure, private :: scalars_init
71 procedure, private :: scalars_init_single
73 procedure :: step => scalars_step
75 procedure :: restart => scalars_restart
77 procedure :: validate => scalars_validate
79 procedure :: free => scalars_free
81 procedure, private :: register_lags_with_checkpoint
82 end type scalars_t
83
84contains
85
87 subroutine scalars_init(this, n_scalars, msh, coef, gs, params, &
88 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
89 class(scalars_t), intent(inout) :: this
90 integer, intent(in) :: n_scalars
91 type(mesh_t), target, intent(in) :: msh
92 type(coef_t), target, intent(in) :: coef
93 type(gs_t), target, intent(inout) :: gs
94 type(json_file), target, intent(inout) :: params
95 type(json_file), target, intent(inout) :: numerics_params
96 type(user_t), target, intent(in) :: user
97 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
98 type(time_scheme_controller_t), target, intent(in) :: time_scheme
99 TYPE(field_t), TARGET, INTENT(IN) :: rho
100 type(chkp_t), target, intent(inout) :: chkp
101 type(json_file) :: json_subdict
102 integer :: i, j
103 character(len=:), allocatable :: field_name
104 character(len=:), allocatable :: field_names(:)
105 character(len=256) :: error_msg
106
107 ! Allocate the scalar fields
108 ! If there are more scalar_scheme_t types, add a factory function here
109 allocate(scalar_pnpn_t::this%scalar_fields(n_scalars))
110
111 ! Collect and validate field names for all scalars
112 allocate(character(len=256) :: field_names(n_scalars))
113
114 do i = 1, n_scalars
115 ! Extract element i from the "scalars" array
116 call json_extract_item(params, "", i, json_subdict)
117
118 ! Try to get name from JSON, generate one if not found or empty
119 if (json_subdict%valid_path('name')) then
120 call json_get(json_subdict, 'name', field_name)
121 else
122 field_name = ''
123 end if
124
125 ! If name is empty or not provided, generate a default one
126 if (len_trim(field_name) == 0) then
127 if (n_scalars == 1) then
128 field_name = 's' ! Single scalar gets default name 's'
129 else
130 write(field_name, '(A,I0)') 's_', i
131 end if
132 end if
133
134 field_names(i) = trim(field_name)
135
136 ! If there's a duplicate, append a number until unique
137 if (n_scalars > 1) then
138 j = 1
139 do while (j < i)
140 if (trim(field_names(i)) == trim(field_names(j))) then
141 write(field_name, '(A,I0)') trim(field_names(i))//'_', j
142 field_names(i) = trim(field_name)
143 j = 1 ! Start over to check if new name is unique
144 else
145 j = j + 1
146 end if
147 end do
148 end if
149 end do
150
151 do i = 1, n_scalars
152 call json_extract_item(params, "", i, json_subdict)
153
154 ! Use the processed field names for all scalars
155 call json_subdict%add('name', trim(field_names(i)))
156
157 call this%scalar_fields(i)%init(msh, coef, gs, json_subdict, &
158 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
159 end do
160
161 ! Register all scalar lag fields with checkpoint using scalable approach
162 if (n_scalars > 1) then
163 call this%register_lags_with_checkpoint(chkp)
164 else
165 ! For single scalar, use legacy interface
166 select type(scalar => this%scalar_fields(1))
167 type is (scalar_pnpn_t)
168 call chkp%add_scalar(scalar%s, scalar%slag, scalar%abx1, scalar%abx2)
169 end select
170 end if
171 end subroutine scalars_init
172
173 subroutine scalars_init_single(this, msh, coef, gs, params, numerics_params, &
174 user, chkp, ulag, vlag, wlag, time_scheme, rho)
175 class(scalars_t), intent(inout) :: this
176 type(mesh_t), target, intent(in) :: msh
177 type(coef_t), target, intent(in) :: coef
178 type(gs_t), target, intent(inout) :: gs
179 type(json_file), target, intent(inout) :: params
180 type(json_file), target, intent(inout) :: numerics_params
181 type(user_t), target, intent(in) :: user
182 type(chkp_t), target, intent(inout) :: chkp
183 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
184 type(time_scheme_controller_t), target, intent(in) :: time_scheme
185 TYPE(field_t), TARGET, INTENT(IN) :: rho
186
187 ! Allocate a single scalar field
188 allocate(scalar_pnpn_t::this%scalar_fields(1))
189
190 ! Set the scalar name to "s"
191 if (.not. params%valid_path('name')) then
192 call params%add('name', 's')
193 end if
194
195 ! Initialize it directly with the params
196 call this%scalar_fields(1)%init(msh, coef, gs, params, numerics_params, &
197 user, chkp, ulag, vlag, wlag, time_scheme, rho)
198
199 ! Register single scalar with checkpoint
200 select type(scalar => this%scalar_fields(1))
201 type is (scalar_pnpn_t)
202 call chkp%add_scalar(scalar%s, scalar%slag, scalar%abx1, scalar%abx2)
203 end select
204 end subroutine scalars_init_single
205
207 subroutine scalars_step(this, time, ext_bdf, dt_controller)
208 class(scalars_t), intent(inout) :: this
209 type(time_state_t), intent(in) :: time
210 type(time_scheme_controller_t), intent(inout) :: ext_bdf
211 type(time_step_controller_t), intent(inout) :: dt_controller
212 integer :: i
213 type(ksp_monitor_t), dimension(size(this%scalar_fields)) :: ksp_results
214
215 ! Iterate through all scalar fields
216 do i = 1, size(this%scalar_fields)
217 call this%scalar_fields(i)%step(time, ext_bdf, dt_controller, &
218 ksp_results(i))
219 end do
220
221 call scalar_step_info(time, ksp_results)
222 end subroutine scalars_step
223
225 subroutine scalars_restart(this, chkp)
226 class(scalars_t), intent(inout) :: this
227 type(chkp_t), intent(inout) :: chkp
228 integer :: i, n_scalars
229
230 n_scalars = size(this%scalar_fields)
231 do i = 1, size(this%scalar_fields)
232 call this%scalar_fields(i)%restart(chkp)
233 end do
234 end subroutine scalars_restart
235
237 subroutine scalars_validate(this)
238 class(scalars_t), intent(inout) :: this
239 integer :: i
240 ! Iterate through all scalar fields
241 do i = 1, size(this%scalar_fields)
242 call this%scalar_fields(i)%slag%set(this%scalar_fields(i)%s)
243 call this%scalar_fields(i)%validate()
244 end do
245 end subroutine scalars_validate
246
248 subroutine scalars_free(this)
249 class(scalars_t), intent(inout) :: this
250 integer :: i
251
252 ! Iterate through all scalar fields
253 if (allocated(this%scalar_fields)) then
254 do i = 1, size(this%scalar_fields)
255 call this%scalar_fields(i)%free()
256 end do
257 deallocate(this%scalar_fields)
258 end if
259 end subroutine scalars_free
260
262 subroutine register_lags_with_checkpoint(this, chkp)
263 class(scalars_t), intent(inout) :: this
264 type(chkp_t), intent(inout) :: chkp
265 integer :: i, n_scalars
266
267 n_scalars = size(this%scalar_fields)
268
269 ! Allocate ABX field arrays
270 allocate(chkp%scalar_abx1(n_scalars))
271 allocate(chkp%scalar_abx2(n_scalars))
272
273 ! Add all scalar lag fields to the checkpoint list and populate ABX fields
274 do i = 1, n_scalars
275 call chkp%scalar_lags%append(this%scalar_fields(i)%slag)
276
277 ! Cast to scalar_pnpn_t to access ABX fields
278 select type(scalar_field => this%scalar_fields(i))
279 type is(scalar_pnpn_t)
280 call associate_scalar_abx_fields(chkp, i, scalar_field)
281 end select
282 end do
283
284 end subroutine register_lags_with_checkpoint
285
287 subroutine associate_scalar_abx_fields(chkp, index, scalar_field)
288 type(chkp_t), intent(inout) :: chkp
289 integer, intent(in) :: index
290 type(scalar_pnpn_t), target, intent(in) :: scalar_field
291
292 chkp%scalar_abx1(index)%ptr => scalar_field%abx1
293 chkp%scalar_abx2(index)%ptr => scalar_field%abx2
294 end subroutine associate_scalar_abx_fields
295
296end 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.
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:238
subroutine associate_scalar_abx_fields(chkp, index, scalar_field)
Helper subroutine to associate ABX field pointers with proper TARGET attribute.
Definition scalars.f90:288
subroutine register_lags_with_checkpoint(this, chkp)
Register scalar lag fields with checkpoint.
Definition scalars.f90:263
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:89
subroutine scalars_init_single(this, msh, coef, gs, params, numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
Definition scalars.f90:175
subroutine scalars_free(this)
Clean up all resources.
Definition scalars.f90:249
subroutine scalars_step(this, time, ext_bdf, dt_controller)
Perform a time step for all scalar fields.
Definition scalars.f90:208
subroutine scalars_restart(this, chkp)
Restart from checkpoint data.
Definition scalars.f90:226
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:62
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...