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
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
67 type(scalar_stats_output_t) :: stats_output
69 real(kind=rp) :: start_time
70 real(kind=rp) :: time
71 logical :: default_fname = .true.
72
73 contains
75 procedure, pass(this) :: init => scalar_stats_simcomp_init_from_json
77 procedure, pass(this) :: init_from_components => &
80 procedure, pass(this) :: free => scalar_stats_simcomp_free
82 procedure, pass(this) :: compute_ => scalar_stats_simcomp_compute
84 procedure, pass(this) :: output_ => scalar_stats_simcomp_compute
86 procedure, pass(this) :: restart_ => scalar_stats_simcomp_restart
88
89contains
90
94 subroutine scalar_stats_simcomp_init_from_json(this, json, case)
95 class(scalar_stats_simcomp_t), target, intent(inout) :: this
96 type(json_file), intent(inout) :: json
97 class(case_t), intent(inout), target :: case
98 character(len=:), allocatable :: filename
99 character(len=20), allocatable :: fields(:)
100 character(len=:), allocatable :: hom_dir
101 character(len=:), allocatable :: stat_set
102 character(len=:), allocatable :: sname
103 character(len=:), allocatable :: name
104 real(kind=rp) :: start_time
105 type(field_t), pointer :: s, u, v, w, p
106 type(coef_t), pointer :: coef
107 logical :: sname_provided
108
109 sname_provided = json%valid_path('field')
110
111 call json_get_or_default(json, 'field', &
112 sname, 's')
113 if (sname_provided) then
114 call json_get_or_default(json, "name", &
115 name, "scalar_stats_" // trim(sname))
116 else
117 call json_get_or_default(json, "name", &
118 name, "scalar_stats")
119 endif
120 call this%init_base(json, case)
121 call json_get_or_default(json, 'avg_direction', &
122 hom_dir, 'none')
123 call json_get_or_lookup_or_default(json, 'start_time', &
124 start_time, 0.0_rp)
125 call json_get_or_default(json, 'set_of_stats', &
126 stat_set, 'full')
127
128 s => neko_registry%get_field_by_name(sname)
129 u => neko_registry%get_field("u")
130 v => neko_registry%get_field("v")
131 w => neko_registry%get_field("w")
132 p => neko_registry%get_field("p")
133 coef => case%fluid%c_Xh
134 this%name = name
135
136 if (json%valid_path("output_filename")) then
137 call json_get(json, "output_filename", filename)
138 call scalar_stats_simcomp_init_from_components(this, name, s, u, v, w, &
139 p, coef, start_time, hom_dir, stat_set, filename)
140 else if (sname_provided) then
141 call scalar_stats_simcomp_init_from_components(this, name, s, u, v, w, &
142 p, coef, start_time, hom_dir, stat_set, "scalar_stats_" // trim(sname) // "0")
143 else
144 call scalar_stats_simcomp_init_from_components(this, name, s, u, v, w, &
145 p, coef, start_time, hom_dir, stat_set)
146 end if
147
149
161 subroutine scalar_stats_simcomp_init_from_components(this, name, s, u, v, w, &
162 p, coef, start_time, hom_dir, stat_set, fname)
163 class(scalar_stats_simcomp_t), target, intent(inout) :: this
164 character(len=*), intent(in) :: name
165 character(len=*), intent(in) :: hom_dir
166 character(len=*), intent(in) :: stat_set
167 real(kind=rp), intent(in) :: start_time
168 type(field_t), intent(in), target :: s, u, v, w, p
169 type(coef_t), intent(in), target :: coef
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 this%name = name
176 call neko_log%section('Scalar stats')
177 write(log_buf, '(A,A)') 'Scalar field: ', trim(s%name)
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)') 'Set of statistics: ', trim(stat_set)
182 call neko_log%message(log_buf)
183 write(log_buf, '(A,A)') 'Averaging in direction: ', trim(hom_dir)
184 call neko_log%message(log_buf)
185
186
187 call this%stats%init(coef, s, u, v, w, p, stat_set, name)
188
189 this%start_time = start_time
190 this%time = start_time
191 if (present(fname)) then
192 this%default_fname = .false.
193 stats_fname = fname
194 else
195 stats_fname = "scalar_stats0"
196 this%default_fname = .true.
197 end if
198
199 call this%stats_output%init(this%stats, this%start_time, &
200 hom_dir = hom_dir, name = stats_fname, &
201 path = this%case%output_directory)
202
203 call this%case%output_controller%add(this%stats_output, &
204 this%output_controller%control_value, &
205 this%output_controller%control_mode)
206
207 call neko_log%end_section()
208
210
213 class(scalar_stats_simcomp_t), intent(inout) :: this
214 call this%free_base()
215 call this%stats%free()
216 end subroutine scalar_stats_simcomp_free
217
218 subroutine scalar_stats_simcomp_restart(this, time)
219 class(scalar_stats_simcomp_t), intent(inout) :: this
220 type(time_state_t), intent(in) :: time
221 character(len=NEKO_FNAME_LEN) :: fname
222 character(len=5) :: prefix, suffix
223 integer :: last_slash_pos
224 real(kind=rp) :: t
225 t = time%t
226 if (t .gt. this%time) this%time = t
227 if (this%default_fname) then
228 fname = this%stats_output%file_%get_base_fname()
229 write (prefix, '(I5)') this%stats_output%file_%get_counter()
230 call filename_suffix(fname, suffix)
231 last_slash_pos = &
233 if (last_slash_pos .ne. 0) then
234 fname = &
235 trim(fname(1:last_slash_pos))// &
236 "scalar_stats"//trim(adjustl(prefix))//"."//suffix
237 else
238 fname = "scalar_stats"// &
239 trim(adjustl(prefix))//"."//suffix
240 end if
241 call this%stats_output%init_base(fname)
242 end if
243 end subroutine scalar_stats_simcomp_restart
244
247 subroutine scalar_stats_simcomp_compute(this, time)
248 class(scalar_stats_simcomp_t), intent(inout) :: this
249 type(time_state_t), intent(in) :: time
250 real(kind=rp) :: delta_t, t
251 real(kind=rp) :: sample_start_time, sample_time
252 character(len=LOG_SIZE) :: log_buf
253 integer :: ierr
254
255 if (time%start_time .gt. this%start_time) then
256 write(log_buf, '(A)') 'Simulation start time is later than the ' &
257 // 'scalar stats start time.'
258 call neko_log%warning(log_buf)
259 write(log_buf, '(A,E15.7)') 'Simulation start time:', time%start_time
260 call neko_log%warning(log_buf)
261 write(log_buf, '(A,E15.7)') 'Scalar stats start time:', this%start_time
262 call neko_log%warning(log_buf)
263 write(log_buf, '(A)') 'Assigning the statistics start time to ' &
264 // 'the simulation start time.'
265 call neko_log%warning(log_buf)
266 this%start_time = time%start_time
267 this%time = time%start_time
268 end if
269
270 t = time%t
271
272 if (t .ge. this%start_time) then
273 delta_t = t - this%time !This is only a real number
274
275 call mpi_barrier(neko_comm, ierr)
276
277 sample_start_time = mpi_wtime()
278
279 call this%stats%update(delta_t)
280 call mpi_barrier(neko_comm, ierr)
281 this%time = t
282
283 sample_time = mpi_wtime() - sample_start_time
284
285 call neko_log%section('Scalar stats')
286 write(log_buf, '(A,E15.7)') 'Sampling at time:', t
287 call neko_log%message(log_buf)
288 write(log_buf, '(A33,E15.7)') 'Simulationtime since last sample:', &
289 delta_t
290 call neko_log%message(log_buf)
291 write(log_buf, '(A,E15.7)') 'Sampling time (s):', sample_time
292 call neko_log%message(log_buf)
293 call neko_log%end_section()
294 end if
295
296 end subroutine scalar_stats_simcomp_compute
297
298end 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: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 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.