Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
scalar_sgs_stats.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!
36 use mean_field, only : mean_field_t
37 use num_types, only : rp
39 use operators, only : grad
40 use coefs, only : coef_t
41 use field, only : field_t
42 use field_list, only : field_list_t
43 use stats_quant, only : stats_quant_t
44 use registry, only : neko_registry
46 implicit none
47 private
48
49 type, public, extends(stats_quant_t) :: scalar_sgs_stats_t
51 type(field_t) :: stats_work
52
54 type(field_t), pointer :: alphat => null()
55 type(field_t), pointer :: nut => null()
56 logical :: nut_dependency = .false.
57 real(kind=rp) :: pr_turb
58 type(field_t), pointer :: s => null()
59
60 type(mean_field_t) :: alphat_mean
61
62 type(mean_field_t) :: alphatdsdx
63 type(mean_field_t) :: alphatdsdy
64 type(mean_field_t) :: alphatdsdz
65
67 type(field_t), pointer :: dsdx_work
68 type(field_t), pointer :: dsdy_work
69 type(field_t), pointer :: dsdz_work
70
72 type(coef_t), pointer :: coef => null()
74 integer :: n_stats = 4
77 type(field_list_t) :: stat_fields
78 contains
79 generic :: init => init_alphat, init_nut
81 procedure, pass(this) :: init_alphat => scalar_sgs_stats_init_alphat
82 procedure, pass(this) :: init_nut => scalar_sgs_stats_init_nut
84 procedure, pass(this) :: free => scalar_sgs_stats_free
86 procedure, pass(this) :: update => scalar_sgs_stats_update
88 procedure, pass(this) :: reset => scalar_sgs_stats_reset
89 end type scalar_sgs_stats_t
90
91contains
92
98 subroutine scalar_sgs_stats_init_alphat(this, coef, s, alphat_field)
99 class(scalar_sgs_stats_t), intent(inout), target:: this
100 type(coef_t), target, optional :: coef
101 type(field_t), target, intent(in) :: s
102 character(len=*), intent(in) :: alphat_field
103
104 call this%free()
105 this%coef => coef
106
107 this%s => s
108 this%alphat => neko_registry%get_field(alphat_field)
109 this%nut_dependency = .false.
110
111 ! Initialize work fields
112 call this%stats_work%init(this%s%dof, 'stats')
113
114 ! Initialize mean fields
115 call this%alphat_mean%init(this%alphat)
116
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')
120
121 allocate(this%stat_fields%items(this%n_stats))
122
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)
127
128 end subroutine scalar_sgs_stats_init_alphat
129
136 subroutine scalar_sgs_stats_init_nut(this, coef, s, nut_field, pr_turb)
137 class(scalar_sgs_stats_t), intent(inout), target:: this
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
142
143 call this%free()
144 this%coef => coef
145
146 this%s => s
147 this%nut => neko_registry%get_field(nut_field)
148 this%pr_turb = pr_turb
149 this%nut_dependency = .true.
150
151 allocate(this%alphat)
152 call this%alphat%init(this%nut%dof, 'alphat_temp')
153
154 ! Initialize work fields
155 call this%stats_work%init(this%s%dof, 'stats')
156
157 ! Initialize mean fields
158 call this%alphat_mean%init(this%alphat)
159
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')
163
164 allocate(this%stat_fields%items(this%n_stats))
165
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)
170
171 end subroutine scalar_sgs_stats_init_nut
172
175 subroutine scalar_sgs_stats_update(this, k)
176 class(scalar_sgs_stats_t), intent(inout) :: this
177 real(kind=rp), intent(in) :: k
178 integer :: n
179 integer :: temp_indices(3)
180
181 associate(stats_work => this%stats_work)
182 n = stats_work%dof%size()
183
184 call neko_scratch_registry%request_field(this%dsdx_work, &
185 temp_indices(1), .false.)
186 call neko_scratch_registry%request_field(this%dsdy_work, &
187 temp_indices(2), .false.)
188 call neko_scratch_registry%request_field(this%dsdz_work, &
189 temp_indices(3), .false.)
190
191 if (this%nut_dependency) then
192 call field_cmult2(this%alphat, this%nut, 1.0_rp / this%pr_turb)
193 end if
194 call this%alphat_mean%update(k)
195
196 call grad(this%dsdx_work%x, &
197 this%dsdy_work%x, &
198 this%dsdz_work%x, &
199 this%s%x, this%coef)
200
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)
207
208 end associate
209
210 call neko_scratch_registry%relinquish_field(temp_indices)
211 end subroutine scalar_sgs_stats_update
212
213
215 subroutine scalar_sgs_stats_free(this)
216 class(scalar_sgs_stats_t), intent(inout) :: this
217
218 call this%stats_work%free()
219
220 call this%alphat_mean%free()
221
222 call this%alphatdsdx%free()
223 call this%alphatdsdy%free()
224 call this%alphatdsdz%free()
225
226 if (this%nut_dependency .and. associated(this%alphat)) then
227 call this%alphat%free()
228 deallocate(this%alphat)
229 end if
230
231 nullify(this%coef)
232 nullify(this%s)
233 nullify(this%alphat)
234
235 call this%stat_fields%free()
236
237 end subroutine scalar_sgs_stats_free
238
240 subroutine scalar_sgs_stats_reset(this)
241 class(scalar_sgs_stats_t), intent(inout), target:: this
242
243 call this%alphat_mean%reset()
244
245 call this%alphatdsdx%reset()
246 call this%alphatdsdy%reset()
247 call this%alphatdsdz%reset()
248
249 end subroutine scalar_sgs_stats_reset
250
251end module scalar_sgs_stats
Coefficients.
Definition coef.f90:34
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 .
Defines a field.
Definition field.f90:34
Implements mean_field_t.
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
subroutine, public grad(ux, uy, uz, u, coef)
Compute the gradient of a scalar field.
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
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,...
Definition coef.f90:58
field_list_t, To be able to group fields together
Computes the temporal mean of a field.
Abstract type defining a statistical quantity.