Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
fluid_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 : strain_rate
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) :: fluid_sgs_stats_t
51 type(field_t) :: stats_work
52 type(field_t), pointer :: s11_work
53 type(field_t), pointer :: s22_work
54 type(field_t), pointer :: s33_work
55 type(field_t), pointer :: s12_work
56 type(field_t), pointer :: s13_work
57 type(field_t), pointer :: s23_work
58
60 type(field_t), pointer :: nut
61 type(field_t), pointer :: u
62 type(field_t), pointer :: v
63 type(field_t), pointer :: w
64
65 type(mean_field_t) :: nut_mean
66
67 type(mean_field_t) :: uu_sgs
68 type(mean_field_t) :: vv_sgs
69 type(mean_field_t) :: ww_sgs
70 type(mean_field_t) :: uv_sgs
71 type(mean_field_t) :: uw_sgs
72 type(mean_field_t) :: vw_sgs
73
75 type(coef_t), pointer :: coef
76
78 integer :: n_stats = 7
79
82 type(field_list_t) :: stat_fields
83 contains
85 procedure, pass(this) :: init => fluid_sgs_stats_init
87 procedure, pass(this) :: free => fluid_sgs_stats_free
89 procedure, pass(this) :: update => fluid_sgs_stats_update
91 procedure, pass(this) :: reset => fluid_sgs_stats_reset
92 end type fluid_sgs_stats_t
93
94contains
95
103 subroutine fluid_sgs_stats_init(this, coef, u, v, w, nut_field)
104 class(fluid_sgs_stats_t), intent(inout), target:: this
105 type(coef_t), target :: coef
106 type(field_t), target, intent(in) :: u, v, w
107 character(*), intent(in), optional :: nut_field
108
109 call this%free()
110 this%coef => coef
111
112 this%u => u
113 this%v => v
114 this%w => w
115
116 if (present(nut_field)) then
117 this%nut => neko_registry%get_field_by_name(trim(nut_field))
118 else
119 this%nut => neko_registry%get_field_by_name('nut')
120 end if
121
122
123 ! Initialize work fields
124 call this%stats_work%init(this%u%dof, 'stats')
125
126 ! Initialize mean fields
127 call this%nut_mean%init(this%nut)
128
129 call this%uu_sgs%init(this%stats_work, 'uu_sgs')
130 call this%vv_sgs%init(this%stats_work, 'vv_sgs')
131 call this%ww_sgs%init(this%stats_work, 'ww_sgs')
132 call this%uv_sgs%init(this%stats_work, 'uv_sgs')
133 call this%uw_sgs%init(this%stats_work, 'uw_sgs')
134 call this%vw_sgs%init(this%stats_work, 'vw_sgs')
135
136 allocate(this%stat_fields%items(this%n_stats))
137
138 call this%stat_fields%assign_to_field(1, this%nut_mean%mf)
139
140 call this%stat_fields%assign_to_field(2, this%uu_sgs%mf)
141 call this%stat_fields%assign_to_field(3, this%vv_sgs%mf)
142 call this%stat_fields%assign_to_field(4, this%ww_sgs%mf)
143 call this%stat_fields%assign_to_field(5, this%uv_sgs%mf)
144 call this%stat_fields%assign_to_field(6, this%uw_sgs%mf)
145 call this%stat_fields%assign_to_field(7, this%vw_sgs%mf)
146
147 end subroutine fluid_sgs_stats_init
148
151 subroutine fluid_sgs_stats_update(this, k)
152 class(fluid_sgs_stats_t), intent(inout) :: this
153 real(kind=rp), intent(in) :: k
154 integer :: n
155 integer :: temp_indices(6)
156
157 associate(stats_work => this%stats_work)
158 n = stats_work%dof%size()
159
160 call neko_scratch_registry%request_field(this%s11_work, &
161 temp_indices(1), .false.)
162 call neko_scratch_registry%request_field(this%s22_work, &
163 temp_indices(2), .false.)
164 call neko_scratch_registry%request_field(this%s33_work, &
165 temp_indices(3), .false.)
166 call neko_scratch_registry%request_field(this%s12_work, &
167 temp_indices(4), .false.)
168 call neko_scratch_registry%request_field(this%s13_work, &
169 temp_indices(5), .false.)
170 call neko_scratch_registry%request_field(this%s23_work, &
171 temp_indices(6), .false.)
172
173 call this%nut_mean%update(k)
174
175 call strain_rate(this%s11_work%x, &
176 this%s22_work%x, &
177 this%s33_work%x, &
178 this%s12_work%x, &
179 this%s13_work%x, &
180 this%s23_work%x, this%u, this%v, this%w, this%coef)
181
182 ! form the double sij tensor
183 call field_cmult(this%s11_work, 2.0_rp)
184 call field_cmult(this%s22_work, 2.0_rp)
185 call field_cmult(this%s33_work, 2.0_rp)
186 call field_cmult(this%s12_work, 2.0_rp)
187 call field_cmult(this%s13_work, 2.0_rp)
188 call field_cmult(this%s23_work, 2.0_rp)
189
190 call field_col3(stats_work, this%nut, this%s11_work)
191 call this%uu_sgs%update(k)
192 call field_col3(stats_work, this%nut, this%s22_work)
193 call this%vv_sgs%update(k)
194 call field_col3(stats_work, this%nut, this%s33_work)
195 call this%ww_sgs%update(k)
196 call field_col3(stats_work, this%nut, this%s12_work)
197 call this%uv_sgs%update(k)
198 call field_col3(stats_work, this%nut, this%s13_work)
199 call this%uw_sgs%update(k)
200 call field_col3(stats_work, this%nut, this%s23_work)
201 call this%vw_sgs%update(k)
202
203 end associate
204
205 call neko_scratch_registry%relinquish_field(temp_indices)
206 end subroutine fluid_sgs_stats_update
207
208
210 subroutine fluid_sgs_stats_free(this)
211 class(fluid_sgs_stats_t), intent(inout) :: this
212
213 call this%stats_work%free()
214
215 call this%nut_mean%free()
216
217 call this%uu_sgs%free()
218 call this%vv_sgs%free()
219 call this%ww_sgs%free()
220 call this%uv_sgs%free()
221 call this%uw_sgs%free()
222 call this%vw_sgs%free()
223
224 nullify(this%coef)
225 nullify(this%u)
226 nullify(this%v)
227 nullify(this%w)
228 nullify(this%nut)
229
230 call this%stat_fields%free()
231
232 end subroutine fluid_sgs_stats_free
233
235 subroutine fluid_sgs_stats_reset(this)
236 class(fluid_sgs_stats_t), intent(inout), target:: this
237
238 call this%nut_mean%reset()
239
240 call this%uu_sgs%reset()
241 call this%vv_sgs%reset()
242 call this%ww_sgs%reset()
243 call this%uv_sgs%reset()
244 call this%uw_sgs%reset()
245 call this%vw_sgs%reset()
246
247 end subroutine fluid_sgs_stats_reset
248
249end module fluid_sgs_stats
Coefficients.
Definition coef.f90:34
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
Computes the subgrid-scale contributions for Reynolds stresses. We use the Reynolds decomposition for...
subroutine fluid_sgs_stats_update(this, k)
Updates all fields with a new sample.
subroutine fluid_sgs_stats_free(this)
Destructor.
subroutine fluid_sgs_stats_init(this, coef, u, v, w, nut_field)
Constructor. Initialize the fields associated with fluid_sgs_stats.
subroutine fluid_sgs_stats_reset(this)
Resets all the computed means values and sampling times to zero.
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 strain_rate(s11, s22, s33, s12, s13, s23, u, v, w, coef)
Compute the strain rate tensor, i.e 0.5 * du_i/dx_j + du_j/dx_i.
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
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.