Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
fluid_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
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
64 type(fluid_sgs_stats_output_t) :: stats_output
66 real(kind=rp) :: start_time
67 real(kind=rp) :: time
68 logical :: default_fname = .true.
69
70 contains
72 procedure, pass(this) :: init => fluid_sgs_stats_simcomp_init_from_json
74 procedure, pass(this) :: init_from_components => &
77 procedure, pass(this) :: free => fluid_sgs_stats_simcomp_free
79 procedure, pass(this) :: compute_ => fluid_sgs_stats_simcomp_compute
81 procedure, pass(this) :: output_ => fluid_sgs_stats_simcomp_compute
83 procedure, pass(this) :: restart_ => fluid_sgs_stats_simcomp_restart
85
86contains
87
91 subroutine fluid_sgs_stats_simcomp_init_from_json(this, json, case)
92 class(fluid_sgs_stats_simcomp_t), target, intent(inout) :: this
93 type(json_file), intent(inout) :: json
94 class(case_t), intent(inout), target :: case
95 character(len=:), allocatable :: filename
96 character(len=NEKO_VARNAME_LEN), allocatable :: fields(:)
97 character(len=:), allocatable :: hom_dir
98 character(len=:), allocatable :: nut_field
99 character(len=:), allocatable :: name
100 real(kind=rp) :: start_time
101 type(field_t), pointer :: u, v, w
102 type(coef_t), pointer :: coef
103
104 call json_get_or_default(json, "name", name, "fluid_sgs_stats")
105 call this%init_base(json, case)
106 call json_get_or_default(json, 'avg_direction', &
107 hom_dir, 'none')
108 call json_get_or_lookup_or_default(json, 'start_time', &
109 start_time, 0.0_rp)
110 call json_get_or_default(json, 'nut_field', &
111 nut_field, 'nut')
112
113
114 u => neko_registry%get_field("u")
115 v => neko_registry%get_field("v")
116 w => neko_registry%get_field("w")
117 coef => case%fluid%c_Xh
118 this%name = name
119
120 if (json%valid_path("output_filename")) then
121 call json_get(json, "output_filename", filename)
122 call fluid_sgs_stats_simcomp_init_from_components(this, u, v, w, coef, &
123 start_time, hom_dir, nut_field, filename)
124 else
125 call fluid_sgs_stats_simcomp_init_from_components(this, u, v, w, coef, &
126 start_time, hom_dir, nut_field)
127 end if
128
130
140 subroutine fluid_sgs_stats_simcomp_init_from_components(this, u, v, w, coef, &
141 start_time, hom_dir, nut_field, fname)
142 class(fluid_sgs_stats_simcomp_t), target, intent(inout) :: this
143 character(len=*), intent(in) :: hom_dir
144 character(len=*), intent(in) :: nut_field
145 real(kind=rp), intent(in) :: start_time
146 type(field_t), intent(in), target :: u, v, w
147 type(coef_t), intent(in), target :: coef
148 character(len=*), intent(in), optional :: fname
149 character(len=NEKO_FNAME_LEN) :: stats_fname
150 character(len=LOG_SIZE) :: log_buf
151 character(len=5) :: prefix
152
153 call neko_log%section('Fluid stats')
154 write(log_buf, '(A,E15.7)') 'Start time: ', start_time
155 call neko_log%message(log_buf)
156 write(log_buf, '(A,A)') 'Nut field name: ', trim(nut_field)
157 call neko_log%message(log_buf)
158 write(log_buf, '(A,A)') 'Averaging in direction: ', trim(hom_dir)
159 call neko_log%message(log_buf)
160
161
162 call this%stats%init(coef, u, v, w, nut_field)
163
164 this%start_time = start_time
165 this%time = start_time
166 if (present(fname)) then
167 this%default_fname = .false.
168 stats_fname = fname
169 else
170 stats_fname = "fluid_sgs_stats0"
171 this%default_fname = .true.
172 end if
173
174 call this%stats_output%init(this%stats, this%start_time, &
175 hom_dir = hom_dir, name = stats_fname, &
176 path = this%case%output_directory)
177
178 call this%case%output_controller%add(this%stats_output, &
179 this%output_controller%control_value, &
180 this%output_controller%control_mode)
181
182 call neko_log%end_section()
183
185
188 class(fluid_sgs_stats_simcomp_t), intent(inout) :: this
189 call this%free_base()
190 call this%stats%free()
191 call this%stats_output%free()
192 end subroutine fluid_sgs_stats_simcomp_free
193
194 subroutine fluid_sgs_stats_simcomp_restart(this, time)
195 class(fluid_sgs_stats_simcomp_t), intent(inout) :: this
196 type(time_state_t), intent(in) :: time
197 character(len=NEKO_FNAME_LEN) :: fname
198 character(len=5) :: prefix, suffix
199 integer :: last_slash_pos
200 real(kind=rp) :: t
201 t = time%t
202 if (t .gt. this%time) this%time = t
203 if (this%default_fname) then
204 fname = this%stats_output%file_%get_base_fname()
205 write (prefix, '(I5)') this%stats_output%file_%get_counter()
206 call filename_suffix(fname, suffix)
207 last_slash_pos = &
209 if (last_slash_pos .ne. 0) then
210 fname = &
211 trim(fname(1:last_slash_pos))// &
212 "fluid_sgs_stats"//trim(adjustl(prefix))//"."//suffix
213 else
214 fname = "fluid_sgs_stats"// &
215 trim(adjustl(prefix))//"."//suffix
216 end if
217 call this%stats_output%init_base(fname)
218 end if
220
223 subroutine fluid_sgs_stats_simcomp_compute(this, time)
224 class(fluid_sgs_stats_simcomp_t), intent(inout) :: this
225 type(time_state_t), intent(in) :: time
226 real(kind=rp) :: delta_t, t
227 real(kind=rp) :: sample_start_time, sample_time
228 character(len=LOG_SIZE) :: log_buf
229 integer :: ierr
230
231 if (time%start_time .gt. this%start_time) then
232 write(log_buf, '(A)') 'Simulation start time is later than the ' &
233 // 'fluid stats start time.'
234 call neko_log%warning(log_buf)
235 write(log_buf, '(A,E15.7)') 'Simulation start time:', time%start_time
236 call neko_log%warning(log_buf)
237 write(log_buf, '(A,E15.7)') 'Fluid stats start time:', this%start_time
238 call neko_log%warning(log_buf)
239 write(log_buf, '(A)') 'Assigning the statistics start time to ' &
240 // 'the simulation start time.'
241 call neko_log%warning(log_buf)
242 this%start_time = time%start_time
243 this%time = time%start_time
244 end if
245
246 t = time%t
247
248 if (t .ge. this%start_time) then
249 delta_t = t - this%time !This is only a real number
250
251 call mpi_barrier(neko_comm, ierr)
252
253 sample_start_time = mpi_wtime()
254
255 call this%stats%update(delta_t)
256 call mpi_barrier(neko_comm, ierr)
257 this%time = t
258
259 sample_time = mpi_wtime() - sample_start_time
260
261 call neko_log%section('Fluid SGS stats')
262 write(log_buf, '(A,E15.7)') 'Sampling at time:', t
263 call neko_log%message(log_buf)
264 write(log_buf, '(A33,E15.7)') 'Simulation time since last sample:', &
265 delta_t
266 call neko_log%message(log_buf)
267 write(log_buf, '(A,E15.7)') 'Sampling time (s):', sample_time
268 call neko_log%message(log_buf)
269 call neko_log%end_section()
270 end if
271
273
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
Implements fluid_sgs_stats_ouput_t.
Implements the fluid_sgs_stats_simcomp_t type.
subroutine fluid_sgs_stats_simcomp_init_from_components(this, u, v, w, coef, start_time, hom_dir, nut_field, fname)
Actual constructor.
subroutine fluid_sgs_stats_simcomp_compute(this, time)
fluid_sgs_stats, called depending on compute_control and compute_value
subroutine fluid_sgs_stats_simcomp_restart(this, time)
subroutine fluid_sgs_stats_simcomp_free(this)
Destructor.
subroutine fluid_sgs_stats_simcomp_init_from_json(this, json, case)
Constructor from json.
Computes the subgrid-scale contributions for Reynolds stresses. We use the Reynolds decomposition for...
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
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 sgs statistics for fluid computed using the fluid_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.