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