Neko 1.99.4
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
entropy_viscosity.f90
Go to the documentation of this file.
1! Copyright (c) 2025, 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!
34 use num_types, only : rp
36 use json_module, only : json_file
38 use field, only : field_t
41 use coefs, only : coef_t
42 use dofmap, only : dofmap_t
43 use time_state, only : time_state_t
45 use operators, only : div
46 use math, only : glmax, absval
48 use mesh, only : mesh_t
49 use space, only : space_t
50 use gather_scatter, only : gs_t
51 use gs_ops, only : gs_op_add
65 implicit none
66 private
67
68 type, public, extends(regularization_t) :: entropy_viscosity_t
69 real(kind=rp) :: c_avisc_entropy
70 real(kind=rp) :: c_avisc_low
71 type(field_t) :: entropy_residual
72 type(field_series_t) :: s_lag
73 type(field_t), pointer :: s => null()
74 type(field_t), pointer :: u => null()
75 type(field_t), pointer :: v => null()
76 type(field_t), pointer :: w => null()
77 type(field_t), pointer :: h => null()
78 type(field_t), pointer :: max_wave_speed => null()
79 type(mesh_t), pointer :: msh => null()
80 type(space_t), pointer :: xh => null()
81 type(gs_t), pointer :: gs => null()
82 contains
83 procedure, pass(this) :: init => entropy_viscosity_init
84 procedure, pass(this) :: free => entropy_viscosity_free
85 procedure, pass(this) :: compute => entropy_viscosity_compute
86 procedure, pass(this) :: update_lag => entropy_viscosity_update_lag
87 procedure, pass(this), private :: compute_residual => entropy_viscosity_compute_residual
88 procedure, pass(this), private :: compute_viscosity => entropy_viscosity_compute_viscosity
89 procedure, pass(this), private :: smooth_viscosity => entropy_viscosity_smooth_viscosity
90 procedure, pass(this), private :: apply_element_max => entropy_viscosity_apply_element_max
91 procedure, pass(this), private :: low_order_viscosity => entropy_viscosity_low_order
92 end type entropy_viscosity_t
93
95
96contains
97
98 subroutine entropy_viscosity_init(this, json, coef, dof, reg_coeff)
99 class(entropy_viscosity_t), intent(inout) :: this
100 type(json_file), intent(inout) :: json
101 type(coef_t), intent(in), target :: coef
102 type(dofmap_t), intent(in), target :: dof
103 type(field_t), intent(in), target :: reg_coeff
104
105 call this%init_base(json, coef, dof, reg_coeff)
106
107 call json_get_or_default(json, 'c_avisc_low', this%c_avisc_low, 1.0_rp)
108 call json_get_or_default(json, 'c_avisc_entropy', this%c_avisc_entropy, 1.0_rp)
109
110 call this%entropy_residual%init(dof, 'entropy_residual')
111
112 nullify(this%S)
113 nullify(this%u)
114 nullify(this%v)
115 nullify(this%w)
116 nullify(this%h)
117 nullify(this%max_wave_speed)
118 nullify(this%msh)
119 nullify(this%Xh)
120 nullify(this%gs)
121
122 end subroutine entropy_viscosity_init
123
124 subroutine entropy_viscosity_free(this)
125 class(entropy_viscosity_t), intent(inout) :: this
126
127 call this%free_base()
128 call this%entropy_residual%free()
129 call this%S_lag%free()
130
131 nullify(this%S)
132 nullify(this%u)
133 nullify(this%v)
134 nullify(this%w)
135 nullify(this%h)
136 nullify(this%max_wave_speed)
137 nullify(this%msh)
138 nullify(this%Xh)
139 nullify(this%gs)
140
141 end subroutine entropy_viscosity_free
142
143 subroutine entropy_viscosity_compute(this, time, tstep, dt)
144 class(entropy_viscosity_t), intent(inout) :: this
145 type(time_state_t), intent(in) :: time
146 integer, intent(in) :: tstep
147 real(kind=rp), intent(in) :: dt
148
149 call this%compute_residual(tstep, dt, time%dtlag)
150 call this%compute_viscosity(tstep)
151
152 end subroutine entropy_viscosity_compute
153
154 subroutine entropy_viscosity_compute_residual(this, tstep, dt, dt_lag)
155 class(entropy_viscosity_t), intent(inout) :: this
156 integer, intent(in) :: tstep
157 real(kind=rp), intent(in) :: dt
158 real(kind=rp), intent(in) :: dt_lag(10)
159 integer :: n
160 type(field_t), pointer :: us_field, vs_field, ws_field, div_field
161 integer :: temp_indices(4)
162 real(kind=rp) :: bdf_coeffs(4)
163 type(bdf_time_scheme_t) :: bdf_scheme
164 real(kind=rp) :: dt_local(10)
165
166 if (tstep .le. 3) then
167 return
168 end if
169
170 n = this%dof%size()
171 call field_cfill(this%entropy_residual, 0.0_rp, n)
172
173 bdf_coeffs = 0.0_rp
174 dt_local = dt_lag
175
176 call bdf_scheme%compute_coeffs(bdf_coeffs, dt_local, 3)
177
178 if (neko_bcknd_device .eq. 1) then
180 this%entropy_residual%x_d, &
181 this%S%x_d, this%S_lag%lf(1)%x_d, &
182 this%S_lag%lf(2)%x_d, this%S_lag%lf(3)%x_d, &
183 bdf_coeffs, dt, n)
184 else
186 this%entropy_residual%x, &
187 this%S%x, this%S_lag%lf(1)%x, &
188 this%S_lag%lf(2)%x, this%S_lag%lf(3)%x, &
189 bdf_coeffs, dt, n)
190 end if
191
192 call neko_scratch_registry%request_field(us_field, temp_indices(1), .false.)
193 call neko_scratch_registry%request_field(vs_field, temp_indices(2), .false.)
194 call neko_scratch_registry%request_field(ws_field, temp_indices(3), .false.)
195 call neko_scratch_registry%request_field(div_field, temp_indices(4), .false.)
196
197 if (neko_bcknd_device .eq. 1) then
198 call device_col3(us_field%x_d, this%u%x_d, this%S%x_d, n)
199 call device_col3(vs_field%x_d, this%v%x_d, this%S%x_d, n)
200 call device_col3(ws_field%x_d, this%w%x_d, this%S%x_d, n)
201 else
202 call entropy_viscosity_col3_vector_cpu(us_field%x, vs_field%x, &
203 ws_field%x, this%u%x, this%v%x, this%w%x, this%S%x, n)
204 end if
205
206 call div(div_field%x, us_field%x, vs_field%x, ws_field%x, this%coef)
207
208 if (neko_bcknd_device .eq. 1) then
209 call device_memcpy(this%entropy_residual%x, this%entropy_residual%x_d, &
210 n, device_to_host, sync = .false.)
211 call device_memcpy(div_field%x, div_field%x_d, n, device_to_host, &
212 sync = .true.)
213 end if
214
215 call entropy_viscosity_abs_add_cpu(this%entropy_residual%x, div_field%x, n)
216
217 if (neko_bcknd_device .eq. 1) then
218 call device_memcpy(this%entropy_residual%x, this%entropy_residual%x_d, &
219 n, host_to_device, sync = .false.)
220 end if
221
222 call neko_scratch_registry%relinquish_field(temp_indices)
223
225
227 class(entropy_viscosity_t), intent(inout) :: this
228 integer, intent(in) :: tstep
229 integer :: n, temp_indices(1)
230 real(kind=rp) :: s_mean, n_s
231 type(field_t), pointer :: temp_field
232
233 n = this%dof%size()
234
235 if (tstep .le. 3) then
236 call field_cfill(this%reg_coeff, 0.0_rp, n)
237 return
238 end if
239
240 call neko_scratch_registry%request_field(temp_field, temp_indices(1), .false.)
241
242 call field_cfill(temp_field, 1.0_rp, n)
243 s_mean = field_glsum(this%S, n) / field_glsum(temp_field, n)
244
245 call field_copy(temp_field, this%S, n)
246 call field_cadd(temp_field, -s_mean, n)
247
248 if (neko_bcknd_device .eq. 1) then
249 call device_absval(temp_field%x_d, n)
250 call device_memcpy(temp_field%x, temp_field%x_d, n, device_to_host, &
251 sync = .true.)
252 else
253 call absval(temp_field%x, n)
254 end if
255
256 ! Normalization factor n_S = |S-S_mean|_inf
257 n_s = glmax(temp_field%x, n)
258
259 call neko_scratch_registry%relinquish_field(temp_indices)
260
261 if (n_s < 1.0e-12_rp) then
262 n_s = 1.0e-12_rp
263 end if
264
265 ! entropy viscosity = c_avisc_entropy * h^2 * entropy_residual / n_S
266 if (neko_bcknd_device .eq. 1) then
268 this%reg_coeff%x_d, this%entropy_residual%x_d, &
269 this%h%x_d, this%c_avisc_entropy, n_s, n)
270 else
272 this%reg_coeff%x, this%entropy_residual%x, &
273 this%h%x, this%c_avisc_entropy, n_s, n)
274 end if
275
276 ! effective viscosity = min(entropy viscosity, low-order viscosity)
277 if (neko_bcknd_device .eq. 1) then
279 this%reg_coeff%x_d, this%h%x_d, this%max_wave_speed%x_d, &
280 this%c_avisc_low, n)
281 else
283 this%reg_coeff%x, this%h%x, this%max_wave_speed%x, &
284 this%c_avisc_low, n)
285 end if
286
287 call this%apply_element_max()
288
289 call this%smooth_viscosity()
290
292
296 class(entropy_viscosity_t), intent(inout) :: this
297 integer :: n
298 type(field_t), pointer :: temp_field, mult_field
299 integer :: temp_indices(2)
300
301 n = this%dof%size()
302
303 call neko_scratch_registry%request_field(temp_field, temp_indices(1), .false.)
304 call neko_scratch_registry%request_field(mult_field, temp_indices(2), .false.)
305
306 call field_copy(temp_field, this%reg_coeff, n)
307 call this%gs%op(temp_field, gs_op_add)
308
309 call field_cfill(mult_field, 1.0_rp, n)
310 call this%gs%op(mult_field, gs_op_add)
311
312 if (neko_bcknd_device .eq. 1) then
314 this%reg_coeff%x_d, temp_field%x_d, mult_field%x_d, n)
315 else
317 this%reg_coeff%x, temp_field%x, mult_field%x, n)
318 end if
319
320 call neko_scratch_registry%relinquish_field(temp_indices)
321
323
325 class(entropy_viscosity_t), intent(inout) :: this
326 integer :: lx
327
328 lx = this%Xh%lx
329
330 if (neko_bcknd_device .eq. 1) then
332 this%reg_coeff%x_d, lx, this%msh%nelv)
333 else
335 this%reg_coeff%x, lx, this%msh%nelv)
336 end if
337
339
340 subroutine entropy_viscosity_set_fields(this, S, u, v, w, h, max_wave_speed, &
341 msh, Xh, gs)
342 class(entropy_viscosity_t), intent(inout) :: this
343 type(field_t), target, intent(inout) :: s
344 type(field_t), target, intent(in) :: u, v, w, h, max_wave_speed
345 type(mesh_t), target, intent(in) :: msh
346 type(space_t), target, intent(in) :: xh
347 type(gs_t), target, intent(in) :: gs
348
349 this%S => s
350 this%u => u
351 this%v => v
352 this%w => w
353 this%h => h
354 this%max_wave_speed => max_wave_speed
355 this%msh => msh
356 this%Xh => xh
357 this%gs => gs
358
359 call this%S_lag%init(s, 3)
360
361 end subroutine entropy_viscosity_set_fields
362
364 class(entropy_viscosity_t), intent(inout) :: this
365
366 call this%S_lag%update()
367
368 end subroutine entropy_viscosity_update_lag
369
370 subroutine entropy_viscosity_col3_vector_cpu(us, vs, ws, u, v, w, S, n)
371 integer, intent(in) :: n
372 real(kind=rp), intent(out) :: us(n), vs(n), ws(n)
373 real(kind=rp), intent(in) :: u(n), v(n), w(n), s(n)
374 integer :: i
375
376 !OCL NORECURRENCE, NOVREC, NOALIAS
377 !DIR$ CONCURRENT
378 !GCC$ ivdep
379 !$omp parallel do simd
380 do i = 1, n
381 us(i) = u(i) * s(i)
382 vs(i) = v(i) * s(i)
383 ws(i) = w(i) * s(i)
384 end do
385 !$omp end parallel do simd
387
388 subroutine entropy_viscosity_abs_add_cpu(entropy_residual, div_field, n)
389 integer, intent(in) :: n
390 real(kind=rp), intent(inout) :: entropy_residual(n)
391 real(kind=rp), intent(in) :: div_field(n)
392 integer :: i
393
394 !OCL NORECURRENCE, NOVREC, NOALIAS
395 !DIR$ CONCURRENT
396 !GCC$ ivdep
397 !$omp parallel do simd
398 do i = 1, n
399 entropy_residual(i) = abs(entropy_residual(i) + div_field(i))
400 end do
401 !$omp end parallel do simd
402 end subroutine entropy_viscosity_abs_add_cpu
403
405 pure function entropy_viscosity_low_order(this, i) result(visc)
406 class(entropy_viscosity_t), intent(in) :: this
407 integer, intent(in) :: i
408 real(kind=rp) :: visc
409
410 visc = this%c_avisc_low * this%h%x(i,1,1,1) * this%max_wave_speed%x(i,1,1,1)
411
412 end function entropy_viscosity_low_order
413
414end module entropy_viscosity
Copy data between host and device (or device and device)
Definition device.F90:72
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Compute the divergence of a vector field.
Definition operators.f90:85
Backward-differencing scheme for time integration.
Coefficients.
Definition coef.f90:34
real(kind=rp) function, public device_glsum(a_d, n, strm)
Sum a vector of length n.
subroutine, public device_absval(a_d, n, strm)
subroutine, public device_col3(a_d, b_d, c_d, n, strm)
Vector multiplication with 3 vectors .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:48
integer, parameter, public device_to_host
Definition device.F90:48
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
CPU backend for entropy viscosity regularization.
subroutine, public entropy_viscosity_compute_viscosity_cpu(reg_coeff, entropy_residual, h, c_avisc_entropy, n_s, n)
Compute viscosity from entropy residual on CPU.
subroutine, public entropy_viscosity_compute_residual_cpu(entropy_residual, s, s_lag1, s_lag2, s_lag3, bdf_coeffs, dt, n)
Compute entropy residual on CPU.
subroutine, public entropy_viscosity_smooth_divide_cpu(reg_coeff, temp_field, mult_field, n)
Divide by multiplicity for smoothing on CPU.
subroutine, public entropy_viscosity_clamp_to_low_order_cpu(reg_coeff, h, max_wave_speed, c_avisc_low, n)
Clamp regularization coefficient to low-order viscosity on CPU.
subroutine, public entropy_viscosity_apply_element_max_cpu(reg_coeff, lx, nelv)
Apply element-wise maximum on CPU.
Device backend for entropy viscosity regularization.
subroutine, public entropy_viscosity_smooth_divide_device(reg_coeff_d, temp_field_d, mult_field_d, n)
Divide by multiplicity for smoothing on device.
subroutine, public entropy_viscosity_clamp_to_low_order_device(reg_coeff_d, h_d, max_wave_speed_d, c_avisc_low, n)
Clamp regularization coefficient to low-order viscosity on device.
subroutine, public entropy_viscosity_apply_element_max_device(reg_coeff_d, lx, nelv)
Apply element-wise maximum on device.
subroutine, public entropy_viscosity_compute_residual_device(entropy_residual_d, s_d, s_lag1_d, s_lag2_d, s_lag3_d, bdf_coeffs, dt, n)
Compute entropy residual on device.
subroutine, public entropy_viscosity_compute_viscosity_device(reg_coeff_d, entropy_residual_d, h_d, c_avisc_entropy, n_s, n)
Compute viscosity from entropy residual on device.
subroutine entropy_viscosity_apply_element_max(this)
subroutine entropy_viscosity_abs_add_cpu(entropy_residual, div_field, n)
subroutine entropy_viscosity_col3_vector_cpu(us, vs, ws, u, v, w, s, n)
subroutine entropy_viscosity_update_lag(this)
subroutine entropy_viscosity_smooth_viscosity(this)
Cross-element smoothing via gather-scatter averaging. Averages viscosity values at shared nodes betwe...
subroutine entropy_viscosity_compute_viscosity(this, tstep)
subroutine entropy_viscosity_compute_residual(this, tstep, dt, dt_lag)
subroutine entropy_viscosity_compute(this, time, tstep, dt)
subroutine entropy_viscosity_init(this, json, coef, dof, reg_coeff)
subroutine entropy_viscosity_free(this)
pure real(kind=rp) function entropy_viscosity_low_order(this, i)
Compute low-order viscosity at point i: c_avisc_low * h * max_wave_speed.
subroutine, public entropy_viscosity_set_fields(this, s, u, v, w, h, max_wave_speed, msh, xh, gs)
subroutine, public field_cadd(a, s, n)
Add a scalar to vector .
subroutine, public field_cfill(a, c, n)
Set all elements to a constant c .
real(kind=rp) function, public field_glsum(a, n)
subroutine, public field_copy(a, b, n)
Copy a vector .
Contains the field_serties_t type.
Defines a field.
Definition field.f90:34
Gather-scatter.
Defines Gather-scatter operations.
Definition gs_ops.f90:34
integer, parameter, public gs_op_add
Definition gs_ops.f90:36
Utilities for retrieving parameters from the case files.
Definition math.f90:60
subroutine, public absval(a, n)
Take the absolute value of an array.
Definition math.f90:1638
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
Definition math.f90:648
Defines a mesh.
Definition mesh.f90:34
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
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 function space.
Definition space.f90:34
Module with things related to the simulation time.
Implicit backward-differencing scheme for time integration.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
Stores a series (sequence) of fields, logically connected to a base field, and arranged according to ...
The function space for the SEM solution fields.
Definition space.f90:63
A struct that contains all info about the time, expand as needed.