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
49 use comm, only : neko_comm
50 use mpi_f08, only : mpi_wtime, mpi_barrier
51 implicit none
52 private
53
66 type(scalar_stats_output_t) :: stats_output
68 real(kind=rp) :: start_time
69 real(kind=rp) :: time
70 logical :: default_fname = .true.
71
72 contains
74 procedure, pass(this) :: init => scalar_stats_simcomp_init_from_json
76 procedure, pass(this) :: init_from_components => &
79 procedure, pass(this) :: free => scalar_stats_simcomp_free
81 procedure, pass(this) :: compute_ => scalar_stats_simcomp_compute
83 procedure, pass(this) :: output_ => scalar_stats_simcomp_compute
85 procedure, pass(this) :: restart_ => scalar_stats_simcomp_restart
87
88contains
89
93 subroutine scalar_stats_simcomp_init_from_json(this, json, case)
94 class(scalar_stats_simcomp_t), target, intent(inout) :: this
95 type(json_file), intent(inout) :: json
96 class(case_t), intent(inout), target :: case
97 character(len=:), allocatable :: filename
98 character(len=20), allocatable :: fields(:)
99 character(len=:), allocatable :: hom_dir
100 character(len=:), allocatable :: stat_set
101 character(len=:), allocatable :: sname
102 character(len=:), allocatable :: name
103 real(kind=rp) :: start_time
104 type(field_t), pointer :: s, u, v, w, p
105 type(coef_t), pointer :: coef
106
107 call json_get_or_default(json, "name", name, "scalar_stats")
108 call this%init_base(json, case)
109 call json_get_or_default(json, 'avg_direction', &
110 hom_dir, 'none')
111 call json_get_or_default(json, 'start_time', &
112 start_time, 0.0_rp)
113 call json_get_or_default(json, 'set_of_stats', &
114 stat_set, 'full')
115 call json_get_or_default(json, 'field', &
116 sname, 's')
117
118 s => neko_registry%get_field_by_name(sname)
119 u => neko_registry%get_field("u")
120 v => neko_registry%get_field("v")
121 w => neko_registry%get_field("w")
122 p => neko_registry%get_field("p")
123 coef => case%fluid%c_Xh
124 this%name = name
125
126 if (json%valid_path("output_filename")) then
127 call json_get(json, "output_filename", filename)
128 call scalar_stats_simcomp_init_from_components(this, name, s, u, v, w, &
129 p, coef, start_time, hom_dir, stat_set, filename)
130 else
131 call scalar_stats_simcomp_init_from_components(this, name, s, u, v, w, &
132 p, coef, start_time, hom_dir, stat_set)
133 end if
134
136
148 subroutine scalar_stats_simcomp_init_from_components(this, name, s, u, v, w, &
149 p, coef, start_time, hom_dir, stat_set, fname)
150 class(scalar_stats_simcomp_t), target, intent(inout) :: this
151 character(len=*), intent(in) :: name
152 character(len=*), intent(in) :: hom_dir
153 character(len=*), intent(in) :: stat_set
154 real(kind=rp), intent(in) :: start_time
155 type(field_t), intent(in), target :: s, u, v, w, p
156 type(coef_t), intent(in), target :: coef
157 character(len=*), intent(in), optional :: fname
158 character(len=NEKO_FNAME_LEN) :: stats_fname
159 character(len=LOG_SIZE) :: log_buf
160 character(len=5) :: prefix
161
162 this%name = name
163 call neko_log%section('Scalar stats')
164 write(log_buf, '(A,A)') 'Scalar field: ', trim(s%name)
165 call neko_log%message(log_buf)
166 write(log_buf, '(A,E15.7)') 'Start time: ', start_time
167 call neko_log%message(log_buf)
168 write(log_buf, '(A,A)') 'Set of statistics: ', trim(stat_set)
169 call neko_log%message(log_buf)
170 write(log_buf, '(A,A)') 'Averaging in direction: ', trim(hom_dir)
171 call neko_log%message(log_buf)
172
173
174 call this%stats%init(coef, s, u, v, w, p, stat_set, name)
175
176 this%start_time = start_time
177 this%time = start_time
178 if (present(fname)) then
179 this%default_fname = .false.
180 stats_fname = fname
181 else
182 stats_fname = "scalar_stats0"
183 this%default_fname = .true.
184 end if
185
186 call this%stats_output%init(this%stats, this%start_time, &
187 hom_dir = hom_dir, name = stats_fname, &
188 path = this%case%output_directory)
189
190 call this%case%output_controller%add(this%stats_output, &
191 this%output_controller%control_value, &
192 this%output_controller%control_mode)
193
194 call neko_log%end_section()
195
197
200 class(scalar_stats_simcomp_t), intent(inout) :: this
201 call this%free_base()
202 call this%stats%free()
203 end subroutine scalar_stats_simcomp_free
204
205 subroutine scalar_stats_simcomp_restart(this, time)
206 class(scalar_stats_simcomp_t), intent(inout) :: this
207 type(time_state_t), intent(in) :: time
208 character(len=NEKO_FNAME_LEN) :: fname
209 character(len=5) :: prefix, suffix
210 integer :: last_slash_pos
211 real(kind=rp) :: t
212 t = time%t
213 if (t .gt. this%time) this%time = t
214 if (this%default_fname) then
215 fname = this%stats_output%file_%get_base_fname()
216 write (prefix, '(I5)') this%stats_output%file_%get_counter()
217 call filename_suffix(fname, suffix)
218 last_slash_pos = &
220 if (last_slash_pos .ne. 0) then
221 fname = &
222 trim(fname(1:last_slash_pos))// &
223 "scalar_stats"//trim(adjustl(prefix))//"."//suffix
224 else
225 fname = "scalar_stats"// &
226 trim(adjustl(prefix))//"."//suffix
227 end if
228 call this%stats_output%init_base(fname)
229 end if
230 end subroutine scalar_stats_simcomp_restart
231
234 subroutine scalar_stats_simcomp_compute(this, time)
235 class(scalar_stats_simcomp_t), intent(inout) :: this
236 type(time_state_t), intent(in) :: time
237 real(kind=rp) :: delta_t, t
238 real(kind=rp) :: sample_start_time, sample_time
239 character(len=LOG_SIZE) :: log_buf
240 integer :: ierr
241
242 if (time%start_time .gt. this%start_time) then
243 write(log_buf, '(A)') 'Simulation start time is later than the ' &
244 // 'scalar stats start time.'
245 call neko_log%warning(log_buf)
246 write(log_buf, '(A,E15.7)') 'Simulation start time:', time%start_time
247 call neko_log%warning(log_buf)
248 write(log_buf, '(A,E15.7)') 'Scalar stats start time:', this%start_time
249 call neko_log%warning(log_buf)
250 write(log_buf, '(A)') 'Assigning the statistics start time to ' &
251 // 'the simulation start time.'
252 call neko_log%warning(log_buf)
253 this%start_time = time%start_time
254 this%time = time%start_time
255 end if
256
257 t = time%t
258
259 if (t .ge. this%start_time) then
260 delta_t = t - this%time !This is only a real number
261
262 call mpi_barrier(neko_comm, ierr)
263
264 sample_start_time = mpi_wtime()
265
266 call this%stats%update(delta_t)
267 call mpi_barrier(neko_comm, ierr)
268 this%time = t
269
270 sample_time = mpi_wtime() - sample_start_time
271
272 call neko_log%section('Scalar stats')
273 write(log_buf, '(A,E15.7)') 'Sampling at time:', t
274 call neko_log%message(log_buf)
275 write(log_buf, '(A33,E15.7)') 'Simulationtime since last sample:', &
276 delta_t
277 call neko_log%message(log_buf)
278 write(log_buf, '(A,E15.7)') 'Sampling time (s):', sample_time
279 call neko_log%message(log_buf)
280 call neko_log%end_section()
281 end if
282
283 end subroutine scalar_stats_simcomp_compute
284
285end 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:58
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.