54 type(
field_t),
pointer :: alphat => null()
56 logical :: nut_dependency = .false.
57 real(kind=
rp) :: pr_turb
72 type(
coef_t),
pointer :: coef => null()
74 integer :: n_stats = 4
79 generic :: init => init_alphat, init_nut
100 type(
coef_t),
target,
optional :: coef
101 type(
field_t),
target,
intent(in) :: s
102 character(len=*),
intent(in) :: alphat_field
109 this%nut_dependency = .false.
112 call this%stats_work%init(this%s%dof,
'stats')
115 call this%alphat_mean%init(this%alphat)
117 call this%alphatdsdx%init(this%stats_work,
'alphatdsdx')
118 call this%alphatdsdy%init(this%stats_work,
'alphatdsdy')
119 call this%alphatdsdz%init(this%stats_work,
'alphatdsdz')
121 allocate(this%stat_fields%items(this%n_stats))
123 call this%stat_fields%assign_to_field(1, this%alphat_mean%mf)
124 call this%stat_fields%assign_to_field(2, this%alphatdsdx%mf)
125 call this%stat_fields%assign_to_field(3, this%alphatdsdy%mf)
126 call this%stat_fields%assign_to_field(4, this%alphatdsdz%mf)
138 type(
coef_t),
target,
optional :: coef
139 type(
field_t),
target,
intent(in) :: s
140 character(len=*),
intent(in) :: nut_field
141 real(kind=
rp),
intent(in) :: pr_turb
148 this%pr_turb = pr_turb
149 this%nut_dependency = .true.
151 allocate(this%alphat)
152 call this%alphat%init(this%nut%dof,
'alphat_temp')
155 call this%stats_work%init(this%s%dof,
'stats')
158 call this%alphat_mean%init(this%alphat)
160 call this%alphatdsdx%init(this%stats_work,
'alphatdsdx')
161 call this%alphatdsdy%init(this%stats_work,
'alphatdsdy')
162 call this%alphatdsdz%init(this%stats_work,
'alphatdsdz')
164 allocate(this%stat_fields%items(this%n_stats))
166 call this%stat_fields%assign_to_field(1, this%alphat_mean%mf)
167 call this%stat_fields%assign_to_field(2, this%alphatdsdx%mf)
168 call this%stat_fields%assign_to_field(3, this%alphatdsdy%mf)
169 call this%stat_fields%assign_to_field(4, this%alphatdsdz%mf)
177 real(kind=
rp),
intent(in) :: k
179 integer :: temp_indices(3)
181 associate(stats_work => this%stats_work)
182 n = stats_work%dof%size()
185 temp_indices(1), .false.)
187 temp_indices(2), .false.)
189 temp_indices(3), .false.)
191 if (this%nut_dependency)
then
192 call field_cmult2(this%alphat, this%nut, 1.0_rp / this%pr_turb)
194 call this%alphat_mean%update(k)
196 call grad(this%dsdx_work%x, &
201 call field_col3(stats_work, this%alphat, this%dsdx_work)
202 call this%alphatdsdx%update(k)
203 call field_col3(stats_work, this%alphat, this%dsdy_work)
204 call this%alphatdsdy%update(k)
205 call field_col3(stats_work, this%alphat, this%dsdz_work)
206 call this%alphatdsdz%update(k)
218 call this%stats_work%free()
220 call this%alphat_mean%free()
222 call this%alphatdsdx%free()
223 call this%alphatdsdy%free()
224 call this%alphatdsdz%free()
226 if (this%nut_dependency .and.
associated(this%alphat))
then
227 call this%alphat%free()
228 deallocate(this%alphat)
235 call this%stat_fields%free()
243 call this%alphat_mean%reset()
245 call this%alphatdsdx%reset()
246 call this%alphatdsdy%reset()
247 call this%alphatdsdz%reset()
subroutine, public field_cmult2(a, b, c, n)
Multiplication by constant c .
subroutine, public field_col3(a, b, c, n)
Vector multiplication with 3 vectors .
subroutine, public field_cmult(a, c, n)
Multiplication by constant c .
integer, parameter, public rp
Global precision used in computations.
subroutine, public grad(ux, uy, uz, u, coef)
Compute the gradient of a scalar field.
Defines a registry for storing solution fields.
type(registry_t), target, public neko_registry
Global field registry.
Computes the subgrid-scale contributions for the scalar fluxes. We use the Reynolds decomposition for...
subroutine scalar_sgs_stats_free(this)
Destructor.
subroutine scalar_sgs_stats_reset(this)
Resets all the computed means values and sampling times to zero.
subroutine scalar_sgs_stats_init_nut(this, coef, s, nut_field, pr_turb)
Constructor. Initialize the fields associated with scalar_sgs_stats. This version uses nut field and ...
subroutine scalar_sgs_stats_init_alphat(this, coef, s, alphat_field)
Constructor. Initialize the fields associated with scalar_sgs_stats. This version uses alphat directl...
subroutine scalar_sgs_stats_update(this, k)
Updates all fields with a new sample.
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Defines a statistical quantity.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
field_list_t, To be able to group fields together
Computes the temporal mean of a field.
Abstract type defining a statistical quantity.