Neko 1.99.1
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, 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
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 real(kind=rp) :: start_time
103 type(field_t), pointer :: s, u, v, w, p
104 type(coef_t), pointer :: coef
105
106 call this%init_base(json, case)
107 call json_get_or_default(json, 'avg_direction', &
108 hom_dir, 'none')
109 call json_get_or_default(json, 'start_time', &
110 start_time, 0.0_rp)
111 call json_get_or_default(json, 'set_of_stats', &
112 stat_set, 'full')
113 call json_get_or_default(json, 'field', &
114 sname, 's')
115
116 s => neko_field_registry%get_field_by_name(sname)
117 u => neko_field_registry%get_field("u")
118 v => neko_field_registry%get_field("v")
119 w => neko_field_registry%get_field("w")
120 p => neko_field_registry%get_field("p")
121 coef => case%fluid%c_Xh
122
123 if (json%valid_path("output_filename")) then
124 call json_get(json, "output_filename", filename)
125 call scalar_stats_simcomp_init_from_components(this, s, u, v, w, p, coef, &
126 start_time, hom_dir, stat_set, filename)
127 else
128 call scalar_stats_simcomp_init_from_components(this, s, u, v, w, p, coef, &
129 start_time, hom_dir, stat_set)
130 end if
131
133
143 subroutine scalar_stats_simcomp_init_from_components(this, s, u, v, w, p, coef, &
144 start_time, hom_dir, stat_set, fname)
145 class(scalar_stats_simcomp_t), target, intent(inout) :: this
146 character(len=*), intent(in) :: hom_dir
147 character(len=*), intent(in) :: stat_set
148 real(kind=rp), intent(in) :: start_time
149 type(field_t), intent(in), target :: s, u, v, w, p
150 type(coef_t), intent(in), target :: coef
151 character(len=*), intent(in), optional :: fname
152 character(len=NEKO_FNAME_LEN) :: stats_fname
153 character(len=LOG_SIZE) :: log_buf
154 character(len=5) :: prefix
155
156 call neko_log%section('Scalar stats')
157 write(log_buf, '(A,A)') 'Scalar field: ', trim(s%name)
158 call neko_log%message(log_buf)
159 write(log_buf, '(A,E15.7)') 'Start time: ', start_time
160 call neko_log%message(log_buf)
161 write(log_buf, '(A,A)') 'Set of statistics: ', trim(stat_set)
162 call neko_log%message(log_buf)
163 write(log_buf, '(A,A)') 'Averaging in direction: ', trim(hom_dir)
164 call neko_log%message(log_buf)
165
166
167 call this%stats%init(coef, s, u, v, w, p, stat_set)
168
169 this%start_time = start_time
170 this%time = start_time
171 if (present(fname)) then
172 this%default_fname = .false.
173 stats_fname = fname
174 else
175 stats_fname = "scalar_stats0"
176 this%default_fname = .true.
177 end if
178
179 call this%stats_output%init(this%stats, this%start_time, &
180 hom_dir = hom_dir,name = stats_fname, &
181 path = this%case%output_directory)
182
183 call this%case%output_controller%add(this%stats_output, &
184 this%output_controller%control_value, &
185 this%output_controller%control_mode)
186
187 call neko_log%end_section()
188
190
193 class(scalar_stats_simcomp_t), intent(inout) :: this
194 call this%free_base()
195 call this%stats%free()
196 end subroutine scalar_stats_simcomp_free
197
198 subroutine scalar_stats_simcomp_restart(this, time)
199 class(scalar_stats_simcomp_t), intent(inout) :: this
200 type(time_state_t), intent(in) :: time
201 character(len=NEKO_FNAME_LEN) :: fname
202 character(len=5) :: prefix,suffix
203 integer :: last_slash_pos
204 real(kind=rp) :: t
205 t = time%t
206 if (t .gt. this%time) this%time = t
207 if (this%default_fname) then
208 fname = this%stats_output%file_%get_base_fname()
209 write (prefix, '(I5)') this%stats_output%file_%get_counter()
210 call filename_suffix(fname,suffix)
211 last_slash_pos = &
213 if (last_slash_pos .ne. 0) then
214 fname = &
215 trim(fname(1:last_slash_pos))// &
216 "scalar_stats"//trim(adjustl(prefix))//"."//suffix
217 else
218 fname = "scalar_stats"// &
219 trim(adjustl(prefix))//"."//suffix
220 end if
221 call this%stats_output%init_base(fname)
222 end if
223 end subroutine scalar_stats_simcomp_restart
224
228 subroutine scalar_stats_simcomp_compute(this, time)
229 class(scalar_stats_simcomp_t), intent(inout) :: this
230 type(time_state_t), intent(in) :: time
231 real(kind=rp) :: delta_t, t
232 real(kind=rp) :: sample_start_time, sample_time
233 character(len=LOG_SIZE) :: log_buf
234 integer :: ierr
235
236 if (time%start_time .gt. this%start_time) then
237 write(log_buf, '(A)') 'Simulation start time is later than the ' &
238 // 'scalar stats start time.'
239 call neko_log%warning(log_buf)
240 write(log_buf, '(A,E15.7)') 'Simulation start time:', time%start_time
241 call neko_log%warning(log_buf)
242 write(log_buf, '(A,E15.7)') 'Scalar stats start time:', this%start_time
243 call neko_log%warning(log_buf)
244 write(log_buf, '(A)') 'Assigning the statistics start time to ' &
245 // 'the simulation start time.'
246 call neko_log%warning(log_buf)
247 this%start_time = time%start_time
248 this%time = time%start_time
249 end if
250
251 t = time%t
252
253 if (t .ge. this%start_time) then
254 delta_t = t - this%time !This is only a real number
255
256 call mpi_barrier(neko_comm, ierr)
257
258 sample_start_time = mpi_wtime()
259
260 call this%stats%update(delta_t)
261 call mpi_barrier(neko_comm, ierr)
262 this%time = t
263
264 sample_time = mpi_wtime() - sample_start_time
265
266 call neko_log%section('Scalar stats')
267 write(log_buf, '(A,E15.7)') 'Sampling at time:', t
268 call neko_log%message(log_buf)
269 write(log_buf, '(A33,E15.7)') 'Simulationtime since last sample:', &
270 delta_t
271 call neko_log%message(log_buf)
272 write(log_buf, '(A,E15.7)') 'Sampling time (s):', sample_time
273 call neko_log%message(log_buf)
274 call neko_log%end_section()
275 end if
276
277 end subroutine scalar_stats_simcomp_compute
278
279end 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 registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
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:70
integer, parameter, public log_size
Definition log.f90:46
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
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_compute(this, time)
scalar_stats, called depending on compute_control and compute_value
subroutine scalar_stats_simcomp_init_from_components(this, s, u, v, w, p, coef, start_time, hom_dir, stat_set, fname)
Actual constructor.
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:55
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.