Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
scalar_sgs_stats_simcomp.f90
Go to the documentation of this file.
1! Copyright (c) 2026, 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!
33!
36 use num_types, only : rp, dp, sp
37 use json_module, only : json_file
39 use registry, only : neko_registry
40 use time_state, only : time_state_t
41 use field, only : field_t
44 use case, only : case_t
45 use coefs, only : coef_t
47 use logger, only : log_size, neko_log
50 use comm, only : neko_comm
51 use mpi_f08, only : mpi_wtime, mpi_barrier
52 implicit none
53 private
54
63 type(scalar_sgs_stats_output_t) :: stats_output
65 real(kind=rp) :: start_time
66 real(kind=rp) :: time
67 logical :: default_fname = .true.
68
69 contains
70 generic :: init_from_components => &
71 init_from_components_alphat, &
72 init_from_components_nut
74 procedure, pass(this) :: init => scalar_sgs_stats_simcomp_init_from_json
76 procedure, pass(this) :: init_from_components_alphat => &
79 procedure, pass(this) :: init_from_components_nut => &
82 procedure, pass(this) :: free => scalar_sgs_stats_simcomp_free
84 procedure, pass(this) :: compute_ => scalar_sgs_stats_simcomp_compute
86 procedure, pass(this) :: output_ => scalar_sgs_stats_simcomp_compute
88 procedure, pass(this) :: restart_ => scalar_sgs_stats_simcomp_restart
90
91contains
92
96 subroutine scalar_sgs_stats_simcomp_init_from_json(this, json, case)
97 class(scalar_sgs_stats_simcomp_t), target, intent(inout) :: this
98 type(json_file), intent(inout) :: json
99 class(case_t), intent(inout), target :: case
100 type(json_file) :: json_subdict
101 character(len=:), allocatable :: filename
102 character(len=20), allocatable :: fields(:)
103 character(len=:), allocatable :: hom_dir
104 character(len=:), allocatable :: s_name
105 character(len=:), allocatable :: name
106 real(kind=rp) :: start_time
107 type(coef_t), pointer :: coef
108 character(len=:), allocatable :: alphat_field, nut_field
109 real(kind=rp) :: pr_turb
110 logical :: nut_dependency
111
112 call json_get_or_default(json, "name", name, "scalar_sgs_stats")
113 call this%init_base(json, case)
114 call json_get_or_default(json, 'avg_direction', &
115 hom_dir, 'none')
116 call json_get_or_lookup_or_default(json, 'start_time', &
117 start_time, 0.0_rp)
118 call json_get_or_default(json, 'field', &
119 s_name, 's')
120
121 call json_get(json, 'alphat', json_subdict)
122 call json_get(json_subdict, 'nut_dependency', nut_dependency)
123 if (nut_dependency) then
124 call json_get_or_lookup(json_subdict, 'Pr_t', pr_turb)
125 call json_get(json_subdict, 'nut_field', nut_field)
126 else
127 call json_get(json_subdict, 'alphat_field', alphat_field)
128 end if
129
130 coef => case%fluid%c_Xh
131 this%name = name
132
133 if (json%valid_path("output_filename")) then
134 call json_get(json, "output_filename", filename)
135 if (nut_dependency) then
136 call this%init_from_components(s_name, coef, &
137 start_time, hom_dir, nut_field, pr_turb, filename)
138 else
139 call this%init_from_components(s_name, coef, &
140 start_time, hom_dir, alphat_field, filename)
141 end if
142 else
143 if (nut_dependency) then
144 call this%init_from_components(s_name, coef, &
145 start_time, hom_dir, nut_field, pr_turb)
146 else
147 call this%init_from_components(s_name, coef, &
148 start_time, hom_dir, alphat_field)
149 end if
150 end if
151
153
162 s_name, coef, start_time, hom_dir, alphat_field, fname)
163 class(scalar_sgs_stats_simcomp_t), target, intent(inout) :: this
164 character(len=*), intent(in) :: s_name
165 character(len=*), intent(in) :: hom_dir
166 real(kind=rp), intent(in) :: start_time
167 type(coef_t), intent(in), target :: coef
168 character(len=*), intent(in) :: alphat_field
169 character(len=*), intent(in), optional :: fname
170 character(len=NEKO_FNAME_LEN) :: stats_fname
171 character(len=LOG_SIZE) :: log_buf
172 character(len=5) :: prefix
173
174 call neko_log%section('scalar SGS stats')
175 write(log_buf, '(A,A)') 'Scalar field: ', trim(s_name)
176 call neko_log%message(log_buf)
177 write(log_buf, '(A,A)') 'Eddy diffusivity field: ', trim(alphat_field)
178 call neko_log%message(log_buf)
179 write(log_buf, '(A,E15.7)') 'Start time: ', start_time
180 call neko_log%message(log_buf)
181 write(log_buf, '(A,A)') 'Averaging in direction: ', trim(hom_dir)
182 call neko_log%message(log_buf)
183
184
185 call this%stats%init(coef, s_name, alphat_field)
186
187 this%start_time = start_time
188 this%time = start_time
189 if (present(fname)) then
190 this%default_fname = .false.
191 stats_fname = fname
192 else
193 stats_fname = "scalar_sgs_stats0"
194 this%default_fname = .true.
195 end if
196
197 call this%stats_output%init(this%stats, this%start_time, &
198 hom_dir = hom_dir, name = stats_fname, &
199 path = this%case%output_directory)
200
201 call this%case%output_controller%add(this%stats_output, &
202 this%output_controller%control_value, &
203 this%output_controller%control_mode)
204
205 call neko_log%end_section()
206
208
219 s_name, coef, start_time, hom_dir, nut_field, pr_turb, fname)
220 class(scalar_sgs_stats_simcomp_t), target, intent(inout) :: this
221 character(len=*), intent(in) :: s_name
222 character(len=*), intent(in) :: hom_dir
223 real(kind=rp), intent(in) :: start_time
224 type(coef_t), intent(in), target :: coef
225 character(len=*), intent(in) :: nut_field
226 real(kind=rp), intent(in) :: pr_turb
227 character(len=*), intent(in), optional :: fname
228 character(len=NEKO_FNAME_LEN) :: stats_fname
229 character(len=LOG_SIZE) :: log_buf
230 character(len=5) :: prefix
231
232 call neko_log%section('scalar stats')
233 write(log_buf, '(A,A)') 'Scalar field: ', trim(s_name)
234 call neko_log%message(log_buf)
235 write(log_buf, '(A,A)') 'Eddy viscosity field: ', trim(nut_field)
236 call neko_log%message(log_buf)
237 write(log_buf, '(A,E15.7)') 'Turbulent Prandtl number: ', pr_turb
238 call neko_log%message(log_buf)
239 write(log_buf, '(A,E15.7)') 'Start time: ', start_time
240 call neko_log%message(log_buf)
241 write(log_buf, '(A,A)') 'Averaging in direction: ', trim(hom_dir)
242 call neko_log%message(log_buf)
243
244
245 call this%stats%init(coef, s_name, nut_field, pr_turb)
246
247 this%start_time = start_time
248 this%time = start_time
249 if (present(fname)) then
250 this%default_fname = .false.
251 stats_fname = fname
252 else
253 stats_fname = "scalar_sgs_stats0"
254 this%default_fname = .true.
255 end if
256
257 call this%stats_output%init(this%stats, this%start_time, &
258 hom_dir = hom_dir, name = stats_fname, &
259 path = this%case%output_directory)
260
261 call this%case%output_controller%add(this%stats_output, &
262 this%output_controller%control_value, &
263 this%output_controller%control_mode)
264
265 call neko_log%end_section()
266
268
271 class(scalar_sgs_stats_simcomp_t), intent(inout) :: this
272 call this%free_base()
273 call this%stats%free()
274 call this%stats_output%free()
275 end subroutine scalar_sgs_stats_simcomp_free
276
277 subroutine scalar_sgs_stats_simcomp_restart(this, time)
278 class(scalar_sgs_stats_simcomp_t), intent(inout) :: this
279 type(time_state_t), intent(in) :: time
280 character(len=NEKO_FNAME_LEN) :: fname
281 character(len=5) :: prefix, suffix
282 integer :: last_slash_pos
283 real(kind=rp) :: t
284 t = time%t
285 if (t .gt. this%time) this%time = t
286 if (this%default_fname) then
287 fname = this%stats_output%file_%get_base_fname()
288 write (prefix, '(I5)') this%stats_output%file_%get_counter()
289 call filename_suffix(fname, suffix)
290 last_slash_pos = &
292 if (last_slash_pos .ne. 0) then
293 fname = &
294 trim(fname(1:last_slash_pos))// &
295 "scalar_sgs_stats"//trim(adjustl(prefix))//"."//suffix
296 else
297 fname = "scalar_sgs_stats"// &
298 trim(adjustl(prefix))//"."//suffix
299 end if
300 call this%stats_output%init_base(fname)
301 end if
303
306 subroutine scalar_sgs_stats_simcomp_compute(this, time)
307 class(scalar_sgs_stats_simcomp_t), intent(inout) :: this
308 type(time_state_t), intent(in) :: time
309 real(kind=rp) :: delta_t, t
310 real(kind=rp) :: sample_start_time, sample_time
311 character(len=LOG_SIZE) :: log_buf
312 integer :: ierr
313
314 if (time%start_time .gt. this%start_time) then
315 write(log_buf, '(A)') 'Simulation start time is later than the ' &
316 // 'scalar stats start time.'
317 call neko_log%warning(log_buf)
318 write(log_buf, '(A,E15.7)') 'Simulation start time:', time%start_time
319 call neko_log%warning(log_buf)
320 write(log_buf, '(A,E15.7)') 'scalar stats start time:', this%start_time
321 call neko_log%warning(log_buf)
322 write(log_buf, '(A)') 'Assigning the statistics start time to ' &
323 // 'the simulation start time.'
324 call neko_log%warning(log_buf)
325 this%start_time = time%start_time
326 this%time = time%start_time
327 end if
328
329 t = time%t
330
331 if (t .ge. this%start_time) then
332 delta_t = t - this%time !This is only a real number
333
334 call mpi_barrier(neko_comm, ierr)
335
336 sample_start_time = mpi_wtime()
337
338 call this%stats%update(delta_t)
339 call mpi_barrier(neko_comm, ierr)
340 this%time = t
341
342 sample_time = mpi_wtime() - sample_start_time
343
344 call neko_log%section('scalar stats')
345 write(log_buf, '(A,E15.7)') 'Sampling at time:', t
346 call neko_log%message(log_buf)
347 write(log_buf, '(A33,E15.7)') 'Simulation time since last sample:', &
348 delta_t
349 call neko_log%message(log_buf)
350 write(log_buf, '(A,E15.7)') 'Sampling time (s):', sample_time
351 call neko_log%message(log_buf)
352 call neko_log%end_section()
353 end if
354
356
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 simulation case.
Definition case.f90:34
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
Defines a field.
Definition field.f90:34
Utilities for retrieving parameters from the case files.
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public sp
Definition num_types.f90:8
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:144
Implements scalar_sgs_stats_ouput_t.
Implements the scalar_sgs_stats_simcomp_t type.
subroutine scalar_sgs_stats_simcomp_init_from_components_nut(this, s_name, coef, start_time, hom_dir, nut_field, pr_turb, fname)
Actual constructor using directly the nut field and the turbulent Prandtl number.
subroutine scalar_sgs_stats_simcomp_free(this)
Destructor.
subroutine scalar_sgs_stats_simcomp_init_from_components_alphat(this, s_name, coef, start_time, hom_dir, alphat_field, fname)
Actual constructor using directly the alphat field.
subroutine scalar_sgs_stats_simcomp_compute(this, time)
scalar_sgs_stats, called depending on compute_control and compute_value
subroutine scalar_sgs_stats_simcomp_init_from_json(this, json, case)
Constructor from json.
subroutine scalar_sgs_stats_simcomp_restart(this, time)
Computes the subgrid-scale contributions for the scalar fluxes. We use the Reynolds decomposition for...
Simulation components are objects that encapsulate functionality that can be fit to a particular comp...
subroutine restart_(this, time)
Dummy restart function.
subroutine compute_(this, time)
Dummy compute function.
Defines a container for all statistics.
Module with things related to the simulation time.
Utilities.
Definition utils.f90:35
integer, parameter, public neko_fname_len
Definition utils.f90:40
subroutine, public filename_suffix(fname, suffix)
Extract a filename's suffix.
Definition utils.f90:108
pure integer function, public filename_tslash_pos(fname)
Find position (in the string) of a filename's trailing slash.
Definition utils.f90:65
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:62
Defines an output for the sgs statistics for scalar computed using the scalar_sgs_stats_t object.
A simulation component that computes the subgrid-scale contributions to the Reynolds stresses in LES.
Base abstract class for simulation components.
A struct that contains all info about the time, expand as needed.