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
49 use comm, only : neko_comm
50 use mpi_f08, only : mpi_wtime, mpi_barrier
51 implicit none
52 private
53
62 type(scalar_sgs_stats_output_t) :: stats_output
64 real(kind=rp) :: start_time
65 real(kind=rp) :: time
66 logical :: default_fname = .true.
67
68 contains
69 generic :: init_from_components => &
70 init_from_components_alphat, &
71 init_from_components_nut
73 procedure, pass(this) :: init => scalar_sgs_stats_simcomp_init_from_json
75 procedure, pass(this) :: init_from_components_alphat => &
78 procedure, pass(this) :: init_from_components_nut => &
81 procedure, pass(this) :: free => scalar_sgs_stats_simcomp_free
83 procedure, pass(this) :: compute_ => scalar_sgs_stats_simcomp_compute
85 procedure, pass(this) :: output_ => scalar_sgs_stats_simcomp_compute
87 procedure, pass(this) :: restart_ => scalar_sgs_stats_simcomp_restart
89
90contains
91
95 subroutine scalar_sgs_stats_simcomp_init_from_json(this, json, case)
96 class(scalar_sgs_stats_simcomp_t), target, intent(inout) :: this
97 type(json_file), intent(inout) :: json
98 class(case_t), intent(inout), target :: case
99 type(json_file) :: json_subdict
100 character(len=:), allocatable :: filename
101 character(len=20), allocatable :: fields(:)
102 character(len=:), allocatable :: hom_dir
103 character(len=:), allocatable :: sname
104 character(len=:), allocatable :: name
105 real(kind=rp) :: start_time
106 type(field_t), pointer :: s
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_default(json, 'start_time', &
117 start_time, 0.0_rp)
118 call json_get_or_default(json, 'field', &
119 sname, '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(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 s => neko_registry%get_field(sname)
131 coef => case%fluid%c_Xh
132 this%name = name
133
134 if (json%valid_path("output_filename")) then
135 call json_get(json, "output_filename", filename)
136 if (nut_dependency) then
137 call this%init_from_components(s, coef, &
138 start_time, hom_dir, nut_field, pr_turb, filename)
139 else
140 call this%init_from_components(s, coef, &
141 start_time, hom_dir, alphat_field, filename)
142 end if
143 else
144 if (nut_dependency) then
145 call this%init_from_components(s, coef, &
146 start_time, hom_dir, nut_field, pr_turb)
147 else
148 call this%init_from_components(s, coef, &
149 start_time, hom_dir, alphat_field)
150 end if
151 end if
152
154
163 s, coef, start_time, hom_dir, alphat_field, fname)
164 class(scalar_sgs_stats_simcomp_t), target, intent(inout) :: this
165 character(len=*), intent(in) :: hom_dir
166 real(kind=rp), intent(in) :: start_time
167 type(field_t), intent(in), target :: s
168 type(coef_t), intent(in), target :: coef
169 character(len=*), intent(in) :: alphat_field
170 character(len=*), intent(in), optional :: fname
171 character(len=NEKO_FNAME_LEN) :: stats_fname
172 character(len=LOG_SIZE) :: log_buf
173 character(len=5) :: prefix
174
175 call neko_log%section('scalar SGS stats')
176 write(log_buf, '(A,A)') 'Scalar field: ', trim(s%name)
177 call neko_log%message(log_buf)
178 write(log_buf, '(A,A)') 'Eddy diffusivity field: ', trim(alphat_field)
179 call neko_log%message(log_buf)
180 write(log_buf, '(A,E15.7)') 'Start time: ', start_time
181 call neko_log%message(log_buf)
182 write(log_buf, '(A,A)') 'Averaging in direction: ', trim(hom_dir)
183 call neko_log%message(log_buf)
184
185
186 call this%stats%init(coef, s, alphat_field)
187
188 this%start_time = start_time
189 this%time = start_time
190 if (present(fname)) then
191 this%default_fname = .false.
192 stats_fname = fname
193 else
194 stats_fname = "scalar_sgs_stats0"
195 this%default_fname = .true.
196 end if
197
198 call this%stats_output%init(this%stats, this%start_time, &
199 hom_dir = hom_dir, name = stats_fname, &
200 path = this%case%output_directory)
201
202 call this%case%output_controller%add(this%stats_output, &
203 this%output_controller%control_value, &
204 this%output_controller%control_mode)
205
206 call neko_log%end_section()
207
209
220 s, coef, start_time, hom_dir, nut_field, pr_turb, fname)
221 class(scalar_sgs_stats_simcomp_t), target, intent(inout) :: this
222 character(len=*), intent(in) :: hom_dir
223 real(kind=rp), intent(in) :: start_time
224 type(field_t), intent(in), target :: s
225 type(coef_t), intent(in), target :: coef
226 character(len=*), intent(in) :: nut_field
227 real(kind=rp), intent(in) :: pr_turb
228 character(len=*), intent(in), optional :: fname
229 character(len=NEKO_FNAME_LEN) :: stats_fname
230 character(len=LOG_SIZE) :: log_buf
231 character(len=5) :: prefix
232
233 call neko_log%section('scalar stats')
234 write(log_buf, '(A,A)') 'Scalar field: ', trim(s%name)
235 call neko_log%message(log_buf)
236 write(log_buf, '(A,A)') 'Eddy viscosity field: ', trim(nut_field)
237 call neko_log%message(log_buf)
238 write(log_buf, '(A,E15.7)') 'Turbulent Prandtl number: ', pr_turb
239 call neko_log%message(log_buf)
240 write(log_buf, '(A,E15.7)') 'Start time: ', start_time
241 call neko_log%message(log_buf)
242 write(log_buf, '(A,A)') 'Averaging in direction: ', trim(hom_dir)
243 call neko_log%message(log_buf)
244
245
246 call this%stats%init(coef, s, nut_field, pr_turb)
247
248 this%start_time = start_time
249 this%time = start_time
250 if (present(fname)) then
251 this%default_fname = .false.
252 stats_fname = fname
253 else
254 stats_fname = "scalar_sgs_stats0"
255 this%default_fname = .true.
256 end if
257
258 call this%stats_output%init(this%stats, this%start_time, &
259 hom_dir = hom_dir, name = stats_fname, &
260 path = this%case%output_directory)
261
262 call this%case%output_controller%add(this%stats_output, &
263 this%output_controller%control_value, &
264 this%output_controller%control_mode)
265
266 call neko_log%end_section()
267
269
272 class(scalar_sgs_stats_simcomp_t), intent(inout) :: this
273 call this%free_base()
274 call this%stats%free()
275 call this%stats_output%free()
276 end subroutine scalar_sgs_stats_simcomp_free
277
278 subroutine scalar_sgs_stats_simcomp_restart(this, time)
279 class(scalar_sgs_stats_simcomp_t), intent(inout) :: this
280 type(time_state_t), intent(in) :: time
281 character(len=NEKO_FNAME_LEN) :: fname
282 character(len=5) :: prefix, suffix
283 integer :: last_slash_pos
284 real(kind=rp) :: t
285 t = time%t
286 if (t .gt. this%time) this%time = t
287 if (this%default_fname) then
288 fname = this%stats_output%file_%get_base_fname()
289 write (prefix, '(I5)') this%stats_output%file_%get_counter()
290 call filename_suffix(fname, suffix)
291 last_slash_pos = &
293 if (last_slash_pos .ne. 0) then
294 fname = &
295 trim(fname(1:last_slash_pos))// &
296 "scalar_sgs_stats"//trim(adjustl(prefix))//"."//suffix
297 else
298 fname = "scalar_sgs_stats"// &
299 trim(adjustl(prefix))//"."//suffix
300 end if
301 call this%stats_output%init_base(fname)
302 end if
304
307 subroutine scalar_sgs_stats_simcomp_compute(this, time)
308 class(scalar_sgs_stats_simcomp_t), intent(inout) :: this
309 type(time_state_t), intent(in) :: time
310 real(kind=rp) :: delta_t, t
311 real(kind=rp) :: sample_start_time, sample_time
312 character(len=LOG_SIZE) :: log_buf
313 integer :: ierr
314
315 if (time%start_time .gt. this%start_time) then
316 write(log_buf, '(A)') 'Simulation start time is later than the ' &
317 // 'scalar stats start time.'
318 call neko_log%warning(log_buf)
319 write(log_buf, '(A,E15.7)') 'Simulation start time:', time%start_time
320 call neko_log%warning(log_buf)
321 write(log_buf, '(A,E15.7)') 'scalar stats start time:', this%start_time
322 call neko_log%warning(log_buf)
323 write(log_buf, '(A)') 'Assigning the statistics start time to ' &
324 // 'the simulation start time.'
325 call neko_log%warning(log_buf)
326 this%start_time = time%start_time
327 this%time = time%start_time
328 end if
329
330 t = time%t
331
332 if (t .ge. this%start_time) then
333 delta_t = t - this%time !This is only a real number
334
335 call mpi_barrier(neko_comm, ierr)
336
337 sample_start_time = mpi_wtime()
338
339 call this%stats%update(delta_t)
340 call mpi_barrier(neko_comm, ierr)
341 this%time = t
342
343 sample_time = mpi_wtime() - sample_start_time
344
345 call neko_log%section('scalar stats')
346 write(log_buf, '(A,E15.7)') 'Sampling at time:', t
347 call neko_log%message(log_buf)
348 write(log_buf, '(A33,E15.7)') 'Simulation time since last sample:', &
349 delta_t
350 call neko_log%message(log_buf)
351 write(log_buf, '(A,E15.7)') 'Sampling time (s):', sample_time
352 call neko_log%message(log_buf)
353 call neko_log%end_section()
354 end if
355
357
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_alphat(this, s, coef, start_time, hom_dir, alphat_field, fname)
Actual constructor using directly the alphat field.
subroutine scalar_sgs_stats_simcomp_init_from_components_nut(this, s, 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_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:58
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.