Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
scalar_stats_simcomp.f90
Go to the documentation of this file.
1! Copyright (c) 2025-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
48 use logger, only : log_size, neko_log
51 use comm, only : neko_comm
52 use mpi_f08, only : mpi_wtime, mpi_barrier
53 implicit none
54 private
55
68 type(scalar_stats_output_t) :: stats_output
70 real(kind=rp) :: start_time
71 real(kind=rp) :: time
72 logical :: default_fname = .true.
73
74 contains
76 procedure, pass(this) :: init => scalar_stats_simcomp_init_from_json
78 procedure, pass(this) :: init_from_components => &
81 procedure, pass(this) :: free => scalar_stats_simcomp_free
83 procedure, pass(this) :: compute_ => scalar_stats_simcomp_compute
85 procedure, pass(this) :: output_ => scalar_stats_simcomp_compute
87 procedure, pass(this) :: restart_ => scalar_stats_simcomp_restart
89
90contains
91
95 subroutine scalar_stats_simcomp_init_from_json(this, json, case)
96 class(scalar_stats_simcomp_t), target, intent(inout) :: this
97 type(json_file), intent(inout) :: json
98 class(case_t), intent(inout), target :: case
99 character(len=:), allocatable :: filename
100 character(len=NEKO_VARNAME_LEN), allocatable :: fields(:)
101 character(len=:), allocatable :: hom_dir
102 character(len=:), allocatable :: stat_set
103 character(len=:), allocatable :: sname
104 character(len=:), allocatable :: name
105 real(kind=rp) :: start_time
106 type(field_t), pointer :: s, u, v, w, p
107 type(coef_t), pointer :: coef
108 logical :: sname_provided
109
110 sname_provided = json%valid_path('field')
111
112 call json_get_or_default(json, 'field', &
113 sname, 's')
114 if (sname_provided) then
115 call json_get_or_default(json, "name", &
116 name, "scalar_stats_" // trim(sname))
117 else
118 call json_get_or_default(json, "name", &
119 name, "scalar_stats")
120 end if
121 call this%init_base(json, case)
122 call json_get_or_default(json, 'avg_direction', &
123 hom_dir, 'none')
124 call json_get_or_lookup_or_default(json, 'start_time', &
125 start_time, 0.0_rp)
126 call json_get_or_default(json, 'set_of_stats', &
127 stat_set, 'full')
128
129 s => neko_registry%get_field_by_name(sname)
130 u => neko_registry%get_field("u")
131 v => neko_registry%get_field("v")
132 w => neko_registry%get_field("w")
133 p => neko_registry%get_field("p")
134 coef => case%fluid%c_Xh
135 this%name = name
136
137 if (json%valid_path("output_filename")) then
138 call json_get(json, "output_filename", filename)
139 call scalar_stats_simcomp_init_from_components(this, name, s, u, v, w, &
140 p, coef, start_time, hom_dir, stat_set, filename)
141 else if (sname_provided) then
142 call scalar_stats_simcomp_init_from_components(this, name, s, u, v, w, &
143 p, coef, start_time, hom_dir, stat_set, "scalar_stats_" // &
144 trim(sname) // "0")
145 else
146 call scalar_stats_simcomp_init_from_components(this, name, s, u, v, w, &
147 p, coef, start_time, hom_dir, stat_set)
148 end if
149
151
163 subroutine scalar_stats_simcomp_init_from_components(this, name, s, u, v, w, &
164 p, coef, start_time, hom_dir, stat_set, fname)
165 class(scalar_stats_simcomp_t), target, intent(inout) :: this
166 character(len=*), intent(in) :: name
167 character(len=*), intent(in) :: hom_dir
168 character(len=*), intent(in) :: stat_set
169 real(kind=rp), intent(in) :: start_time
170 type(field_t), intent(in), target :: s, u, v, w, p
171 type(coef_t), intent(in), target :: coef
172 character(len=*), intent(in), optional :: fname
173 character(len=NEKO_FNAME_LEN) :: stats_fname
174 character(len=LOG_SIZE) :: log_buf
175 character(len=5) :: prefix
176
177 this%name = name
178 call neko_log%section('Scalar stats')
179 write(log_buf, '(A,A)') 'Scalar field: ', trim(s%name)
180 call neko_log%message(log_buf)
181 write(log_buf, '(A,E15.7)') 'Start time: ', start_time
182 call neko_log%message(log_buf)
183 write(log_buf, '(A,A)') 'Set of statistics: ', trim(stat_set)
184 call neko_log%message(log_buf)
185 write(log_buf, '(A,A)') 'Averaging in direction: ', trim(hom_dir)
186 call neko_log%message(log_buf)
187
188
189 call this%stats%init(coef, s, u, v, w, p, stat_set, name)
190
191 this%start_time = start_time
192 this%time = start_time
193 if (present(fname)) then
194 this%default_fname = .false.
195 stats_fname = fname
196 else
197 stats_fname = "scalar_stats0"
198 this%default_fname = .true.
199 end if
200
201 call this%stats_output%init(this%stats, this%start_time, &
202 hom_dir = hom_dir, name = stats_fname, &
203 path = this%case%output_directory)
204
205 call this%case%output_controller%add(this%stats_output, &
206 this%output_controller%control_value, &
207 this%output_controller%control_mode)
208
209 call neko_log%end_section()
210
212
215 class(scalar_stats_simcomp_t), intent(inout) :: this
216 call this%free_base()
217 call this%stats%free()
218 end subroutine scalar_stats_simcomp_free
219
220 subroutine scalar_stats_simcomp_restart(this, time)
221 class(scalar_stats_simcomp_t), intent(inout) :: this
222 type(time_state_t), intent(in) :: time
223 character(len=NEKO_FNAME_LEN) :: fname
224 character(len=5) :: prefix, suffix
225 integer :: last_slash_pos
226 real(kind=rp) :: t
227 t = time%t
228 if (t .gt. this%time) this%time = t
229 if (this%default_fname) then
230 fname = this%stats_output%file_%get_base_fname()
231 write (prefix, '(I5)') this%stats_output%file_%get_counter()
232 call filename_suffix(fname, suffix)
233 last_slash_pos = &
235 if (last_slash_pos .ne. 0) then
236 fname = &
237 trim(fname(1:last_slash_pos))// &
238 "scalar_stats"//trim(adjustl(prefix))//"."//suffix
239 else
240 fname = "scalar_stats"// &
241 trim(adjustl(prefix))//"."//suffix
242 end if
243 call this%stats_output%init_base(fname)
244 end if
245 end subroutine scalar_stats_simcomp_restart
246
249 subroutine scalar_stats_simcomp_compute(this, time)
250 class(scalar_stats_simcomp_t), intent(inout) :: this
251 type(time_state_t), intent(in) :: time
252 real(kind=rp) :: delta_t, t
253 real(kind=rp) :: sample_start_time, sample_time
254 character(len=LOG_SIZE) :: log_buf
255 integer :: ierr
256
257 if (time%start_time .gt. this%start_time) then
258 write(log_buf, '(A)') 'Simulation start time is later than the ' &
259 // 'scalar stats start time.'
260 call neko_log%warning(log_buf)
261 write(log_buf, '(A,E15.7)') 'Simulation start time:', time%start_time
262 call neko_log%warning(log_buf)
263 write(log_buf, '(A,E15.7)') 'Scalar stats start time:', this%start_time
264 call neko_log%warning(log_buf)
265 write(log_buf, '(A)') 'Assigning the statistics start time to ' &
266 // 'the simulation start time.'
267 call neko_log%warning(log_buf)
268 this%start_time = time%start_time
269 this%time = time%start_time
270 end if
271
272 t = time%t
273
274 if (t .ge. this%start_time) then
275 delta_t = t - this%time !This is only a real number
276
277 call mpi_barrier(neko_comm, ierr)
278
279 sample_start_time = mpi_wtime()
280
281 call this%stats%update(delta_t)
282 call mpi_barrier(neko_comm, ierr)
283 this%time = t
284
285 sample_time = mpi_wtime() - sample_start_time
286
287 call neko_log%section('Scalar stats')
288 write(log_buf, '(A,E15.7)') 'Sampling at time:', t
289 call neko_log%message(log_buf)
290 write(log_buf, '(A33,E15.7)') 'Simulationtime since last sample:', &
291 delta_t
292 call neko_log%message(log_buf)
293 write(log_buf, '(A,E15.7)') 'Sampling time (s):', sample_time
294 call neko_log%message(log_buf)
295 call neko_log%end_section()
296 end if
297
298 end subroutine scalar_stats_simcomp_compute
299
300end module scalar_stats_simcomp
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_stats_ouput_t.
Implements the scalar_stats_simcomp_t type.
subroutine scalar_stats_simcomp_free(this)
Destructor.
subroutine scalar_stats_simcomp_init_from_json(this, json, case)
Constructor from json.
subroutine scalar_stats_simcomp_restart(this, time)
subroutine scalar_stats_simcomp_init_from_components(this, name, s, u, v, w, p, coef, start_time, hom_dir, stat_set, fname)
Actual constructor.
subroutine scalar_stats_simcomp_compute(this, time)
scalar_stats, called depending on compute_control and compute_value
Computes various statistics for the scalar fields. We use the Reynolds decomposition for a field u = ...
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:41
subroutine, public filename_suffix(fname, suffix)
Extract a filename's suffix.
Definition utils.f90:118
integer, parameter, public neko_varname_len
Definition utils.f90:42
pure integer function, public filename_tslash_pos(fname)
Find position (in the string) of a filename's trailing slash.
Definition utils.f90:75
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
Defines an output for the scalar statistics computed using the scalar_stats_t object.
A simulation component that computes the scalar statistics for the skewness, kurtosis,...
Base abstract class for simulation components.
A struct that contains all info about the time, expand as needed.