37 use json_module,
only : json_file
50 use mpi_f08,
only : mpi_wtime, mpi_barrier
64 real(kind=
rp) :: start_time
66 logical :: default_fname = .true.
69 generic :: init_from_components => &
70 init_from_components_alphat, &
71 init_from_components_nut
75 procedure, pass(this) :: init_from_components_alphat => &
78 procedure, pass(this) :: init_from_components_nut => &
97 type(json_file),
intent(inout) :: json
98 class(
case_t),
intent(inout),
target :: case
99 type(json_file) :: json_subdict
100 character(len=:),
allocatable :: filename
101 character(len=20),
allocatable :: fields(:)
102 character(len=:),
allocatable :: hom_dir
103 character(len=:),
allocatable :: sname
104 character(len=:),
allocatable :: name
105 real(kind=
rp) :: start_time
107 type(
coef_t),
pointer :: coef
108 character(len=:),
allocatable :: alphat_field, nut_field
109 real(kind=
rp) :: pr_turb
110 logical :: nut_dependency
113 call this%init_base(json,
case)
121 call json_get(json,
'alphat', json_subdict)
122 call json_get(json_subdict,
'nut_dependency', nut_dependency)
123 if (nut_dependency)
then
124 call json_get(json_subdict,
'Pr_t', pr_turb)
125 call json_get(json_subdict,
'nut_field', nut_field)
127 call json_get(json_subdict,
'alphat_field', alphat_field)
131 coef =>
case%fluid%c_Xh
134 if (json%valid_path(
"output_filename"))
then
135 call json_get(json,
"output_filename", filename)
136 if (nut_dependency)
then
137 call this%init_from_components(s, coef, &
138 start_time, hom_dir, nut_field, pr_turb, filename)
140 call this%init_from_components(s, coef, &
141 start_time, hom_dir, alphat_field, filename)
144 if (nut_dependency)
then
145 call this%init_from_components(s, coef, &
146 start_time, hom_dir, nut_field, pr_turb)
148 call this%init_from_components(s, coef, &
149 start_time, hom_dir, alphat_field)
163 s, coef, start_time, hom_dir, alphat_field, fname)
165 character(len=*),
intent(in) :: hom_dir
166 real(kind=
rp),
intent(in) :: start_time
167 type(
field_t),
intent(in),
target :: s
168 type(
coef_t),
intent(in),
target :: coef
169 character(len=*),
intent(in) :: alphat_field
170 character(len=*),
intent(in),
optional :: fname
171 character(len=NEKO_FNAME_LEN) :: stats_fname
172 character(len=LOG_SIZE) :: log_buf
173 character(len=5) :: prefix
175 call neko_log%section(
'scalar SGS stats')
176 write(log_buf,
'(A,A)')
'Scalar field: ', trim(s%name)
178 write(log_buf,
'(A,A)')
'Eddy diffusivity field: ', trim(alphat_field)
180 write(log_buf,
'(A,E15.7)')
'Start time: ', start_time
182 write(log_buf,
'(A,A)')
'Averaging in direction: ', trim(hom_dir)
186 call this%stats%init(coef, s, alphat_field)
188 this%start_time = start_time
189 this%time = start_time
190 if (
present(fname))
then
191 this%default_fname = .false.
194 stats_fname =
"scalar_sgs_stats0"
195 this%default_fname = .true.
198 call this%stats_output%init(this%stats, this%start_time, &
199 hom_dir = hom_dir, name = stats_fname, &
200 path = this%case%output_directory)
202 call this%case%output_controller%add(this%stats_output, &
203 this%output_controller%control_value, &
204 this%output_controller%control_mode)
220 s, coef, start_time, hom_dir, nut_field, pr_turb, fname)
222 character(len=*),
intent(in) :: hom_dir
223 real(kind=
rp),
intent(in) :: start_time
224 type(
field_t),
intent(in),
target :: s
225 type(
coef_t),
intent(in),
target :: coef
226 character(len=*),
intent(in) :: nut_field
227 real(kind=
rp),
intent(in) :: pr_turb
228 character(len=*),
intent(in),
optional :: fname
229 character(len=NEKO_FNAME_LEN) :: stats_fname
230 character(len=LOG_SIZE) :: log_buf
231 character(len=5) :: prefix
233 call neko_log%section(
'scalar stats')
234 write(log_buf,
'(A,A)')
'Scalar field: ', trim(s%name)
236 write(log_buf,
'(A,A)')
'Eddy viscosity field: ', trim(nut_field)
238 write(log_buf,
'(A,E15.7)')
'Turbulent Prandtl number: ', pr_turb
240 write(log_buf,
'(A,E15.7)')
'Start time: ', start_time
242 write(log_buf,
'(A,A)')
'Averaging in direction: ', trim(hom_dir)
246 call this%stats%init(coef, s, nut_field, pr_turb)
248 this%start_time = start_time
249 this%time = start_time
250 if (
present(fname))
then
251 this%default_fname = .false.
254 stats_fname =
"scalar_sgs_stats0"
255 this%default_fname = .true.
258 call this%stats_output%init(this%stats, this%start_time, &
259 hom_dir = hom_dir, name = stats_fname, &
260 path = this%case%output_directory)
262 call this%case%output_controller%add(this%stats_output, &
263 this%output_controller%control_value, &
264 this%output_controller%control_mode)
273 call this%free_base()
274 call this%stats%free()
275 call this%stats_output%free()
281 character(len=NEKO_FNAME_LEN) :: fname
282 character(len=5) :: prefix, suffix
283 integer :: last_slash_pos
286 if (t .gt. this%time) this%time = t
287 if (this%default_fname)
then
288 fname = this%stats_output%file_%get_base_fname()
289 write (prefix,
'(I5)') this%stats_output%file_%get_counter()
293 if (last_slash_pos .ne. 0)
then
295 trim(fname(1:last_slash_pos))// &
296 "scalar_sgs_stats"//trim(adjustl(prefix))//
"."//suffix
298 fname =
"scalar_sgs_stats"// &
299 trim(adjustl(prefix))//
"."//suffix
301 call this%stats_output%init_base(fname)
310 real(kind=
rp) :: delta_t, t
311 real(kind=
rp) :: sample_start_time, sample_time
312 character(len=LOG_SIZE) :: log_buf
315 if (time%start_time .gt. this%start_time)
then
316 write(log_buf,
'(A)')
'Simulation start time is later than the ' &
317 //
'scalar stats start time.'
319 write(log_buf,
'(A,E15.7)')
'Simulation start time:', time%start_time
321 write(log_buf,
'(A,E15.7)')
'scalar stats start time:', this%start_time
323 write(log_buf,
'(A)')
'Assigning the statistics start time to ' &
324 //
'the simulation start time.'
326 this%start_time = time%start_time
327 this%time = time%start_time
332 if (t .ge. this%start_time)
then
333 delta_t = t - this%time
337 sample_start_time = mpi_wtime()
339 call this%stats%update(delta_t)
343 sample_time = mpi_wtime() - sample_start_time
345 call neko_log%section(
'scalar stats')
346 write(log_buf,
'(A,E15.7)')
'Sampling at time:', t
348 write(log_buf,
'(A33,E15.7)')
'Simulation time since last sample:', &
351 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.
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.
Implements scalar_sgs_stats_ouput_t.
Implements the scalar_sgs_stats_simcomp_t type.
subroutine scalar_sgs_stats_simcomp_init_from_components_alphat(this, s, coef, start_time, hom_dir, alphat_field, fname)
Actual constructor using directly the alphat field.
subroutine scalar_sgs_stats_simcomp_init_from_components_nut(this, s, coef, start_time, hom_dir, nut_field, pr_turb, fname)
Actual constructor using directly the nut field and the turbulent Prandtl number.
subroutine scalar_sgs_stats_simcomp_free(this)
Destructor.
subroutine scalar_sgs_stats_simcomp_compute(this, time)
scalar_sgs_stats, called depending on compute_control and compute_value
subroutine scalar_sgs_stats_simcomp_init_from_json(this, json, case)
Constructor from json.
subroutine scalar_sgs_stats_simcomp_restart(this, time)
Computes the subgrid-scale contributions for the scalar fluxes. We use the Reynolds decomposition for...
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 scalar computed using the scalar_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.