37 use json_module,
only : json_file
51 use mpi_f08,
only : mpi_wtime, mpi_barrier
65 real(kind=
rp) :: start_time
67 logical :: default_fname = .true.
73 procedure, pass(this) :: init_from_components => &
92 type(json_file),
intent(inout) :: json
93 class(
case_t),
intent(inout),
target :: case
94 character(len=:),
allocatable :: filename
95 character(len=20),
allocatable :: fields(:)
96 character(len=:),
allocatable :: hom_dir
97 character(len=:),
allocatable :: nut_field
98 character(len=:),
allocatable :: name
99 real(kind=
rp) :: start_time
100 type(
field_t),
pointer :: u, v, w
101 type(
coef_t),
pointer :: coef
104 call this%init_base(json,
case)
116 coef =>
case%fluid%c_Xh
119 if (json%valid_path(
"output_filename"))
then
120 call json_get(json,
"output_filename", filename)
122 start_time, hom_dir, nut_field, filename)
125 start_time, hom_dir, nut_field)
140 start_time, hom_dir, nut_field, fname)
142 character(len=*),
intent(in) :: hom_dir
143 character(len=*),
intent(in) :: nut_field
144 real(kind=
rp),
intent(in) :: start_time
145 type(
field_t),
intent(in),
target :: u, v, w
146 type(
coef_t),
intent(in),
target :: coef
147 character(len=*),
intent(in),
optional :: fname
148 character(len=NEKO_FNAME_LEN) :: stats_fname
149 character(len=LOG_SIZE) :: log_buf
150 character(len=5) :: prefix
152 call neko_log%section(
'Fluid stats')
153 write(log_buf,
'(A,E15.7)')
'Start time: ', start_time
155 write(log_buf,
'(A,A)')
'Nut field name: ', trim(nut_field)
157 write(log_buf,
'(A,A)')
'Averaging in direction: ', trim(hom_dir)
161 call this%stats%init(coef, u, v, w, nut_field)
163 this%start_time = start_time
164 this%time = start_time
165 if (
present(fname))
then
166 this%default_fname = .false.
169 stats_fname =
"fluid_sgs_stats0"
170 this%default_fname = .true.
173 call this%stats_output%init(this%stats, this%start_time, &
174 hom_dir = hom_dir, name = stats_fname, &
175 path = this%case%output_directory)
177 call this%case%output_controller%add(this%stats_output, &
178 this%output_controller%control_value, &
179 this%output_controller%control_mode)
188 call this%free_base()
189 call this%stats%free()
190 call this%stats_output%free()
196 character(len=NEKO_FNAME_LEN) :: fname
197 character(len=5) :: prefix, suffix
198 integer :: last_slash_pos
201 if (t .gt. this%time) this%time = t
202 if (this%default_fname)
then
203 fname = this%stats_output%file_%get_base_fname()
204 write (prefix,
'(I5)') this%stats_output%file_%get_counter()
208 if (last_slash_pos .ne. 0)
then
210 trim(fname(1:last_slash_pos))// &
211 "fluid_sgs_stats"//trim(adjustl(prefix))//
"."//suffix
213 fname =
"fluid_sgs_stats"// &
214 trim(adjustl(prefix))//
"."//suffix
216 call this%stats_output%init_base(fname)
225 real(kind=
rp) :: delta_t, t
226 real(kind=
rp) :: sample_start_time, sample_time
227 character(len=LOG_SIZE) :: log_buf
230 if (time%start_time .gt. this%start_time)
then
231 write(log_buf,
'(A)')
'Simulation start time is later than the ' &
232 //
'fluid stats start time.'
234 write(log_buf,
'(A,E15.7)')
'Simulation start time:', time%start_time
236 write(log_buf,
'(A,E15.7)')
'Fluid stats start time:', this%start_time
238 write(log_buf,
'(A)')
'Assigning the statistics start time to ' &
239 //
'the simulation start time.'
241 this%start_time = time%start_time
242 this%time = time%start_time
247 if (t .ge. this%start_time)
then
248 delta_t = t - this%time
252 sample_start_time = mpi_wtime()
254 call this%stats%update(delta_t)
258 sample_time = mpi_wtime() - sample_start_time
260 call neko_log%section(
'Fluid SGS stats')
261 write(log_buf,
'(A,E15.7)')
'Sampling at time:', t
263 write(log_buf,
'(A33,E15.7)')
'Simulation time since last sample:', &
266 write(log_buf,
'(A,E15.7)')
'Sampling time (s):', sample_time
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.
type(mpi_comm), public neko_comm
MPI communicator.
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.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
integer, parameter, public dp
integer, parameter, public sp
integer, parameter, public rp
Global precision used in computations.
Defines a registry for storing solution fields.
type(registry_t), target, public neko_registry
Global field registry.
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.
integer, parameter, public neko_fname_len
subroutine, public filename_suffix(fname, suffix)
Extract a filename's suffix.
pure integer function, public filename_tslash_pos(fname)
Find position (in the string) of a filename's trailing slash.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
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.