Neko 1.99.2
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 integer :: i, n
149
150 n = this%dof%size()
151
152 call this%compute_residual(tstep, dt, time%dtlag)
153 call this%compute_viscosity(tstep)
154
155 end subroutine entropy_viscosity_compute
156
157 subroutine entropy_viscosity_compute_residual(this, tstep, dt, dt_lag)
158 class(entropy_viscosity_t), intent(inout) :: this
159 integer, intent(in) :: tstep
160 real(kind=rp), intent(in) :: dt
161 real(kind=rp), intent(in) :: dt_lag(10)
162 integer :: i, n
163 type(field_t), pointer :: us_field, vs_field, ws_field, div_field
164 integer :: temp_indices(4)
165 real(kind=rp) :: bdf_coeffs(4)
166 type(bdf_time_scheme_t) :: bdf_scheme
167 real(kind=rp) :: dt_local(10)
168
169 if (tstep .le. 3) then
170 return
171 end if
172
173 n = this%dof%size()
174 call field_cfill(this%entropy_residual, 0.0_rp, n)
175
176 bdf_coeffs = 0.0_rp
177 dt_local = dt_lag
178
179 call bdf_scheme%compute_coeffs(bdf_coeffs, dt_local, 3)
180
181 if (neko_bcknd_device .eq. 1) then
183 this%entropy_residual%x_d, &
184 this%S%x_d, this%S_lag%lf(1)%x_d, &
185 this%S_lag%lf(2)%x_d, this%S_lag%lf(3)%x_d, &
186 bdf_coeffs, dt, n)
187 else
189 this%entropy_residual%x, &
190 this%S%x, this%S_lag%lf(1)%x, &
191 this%S_lag%lf(2)%x, this%S_lag%lf(3)%x, &
192 bdf_coeffs, dt, n)
193 end if
194
195 call neko_scratch_registry%request_field(us_field, temp_indices(1), .false.)
196 call neko_scratch_registry%request_field(vs_field, temp_indices(2), .false.)
197 call neko_scratch_registry%request_field(ws_field, temp_indices(3), .false.)
198 call neko_scratch_registry%request_field(div_field, temp_indices(4), .false.)
199
200 if (neko_bcknd_device .eq. 1) then
201 call device_col3(us_field%x_d, this%u%x_d, this%S%x_d, n)
202 call device_col3(vs_field%x_d, this%v%x_d, this%S%x_d, n)
203 call device_col3(ws_field%x_d, this%w%x_d, this%S%x_d, n)
204 else
205 do i = 1, n
206 us_field%x(i,1,1,1) = this%u%x(i,1,1,1) * this%S%x(i,1,1,1)
207 vs_field%x(i,1,1,1) = this%v%x(i,1,1,1) * this%S%x(i,1,1,1)
208 ws_field%x(i,1,1,1) = this%w%x(i,1,1,1) * this%S%x(i,1,1,1)
209 end do
210 end if
211
212 call div(div_field%x, us_field%x, vs_field%x, ws_field%x, this%coef)
213
214 if (neko_bcknd_device .eq. 1) then
215 call device_memcpy(this%entropy_residual%x, this%entropy_residual%x_d, &
216 n, device_to_host, sync = .false.)
217 call device_memcpy(div_field%x, div_field%x_d, n, device_to_host, &
218 sync = .true.)
219 end if
220
221 do i = 1, n
222 this%entropy_residual%x(i,1,1,1) = abs(this%entropy_residual%x(i,1,1,1) &
223 + div_field%x(i,1,1,1))
224 end do
225
226 if (neko_bcknd_device .eq. 1) then
227 call device_memcpy(this%entropy_residual%x, this%entropy_residual%x_d, &
228 n, host_to_device, sync = .false.)
229 end if
230
231 call neko_scratch_registry%relinquish_field(temp_indices)
232
234
236 class(entropy_viscosity_t), intent(inout) :: this
237 integer, intent(in) :: tstep
238 integer :: i, n, temp_indices(1)
239 real(kind=rp) :: s_mean, n_s
240 type(field_t), pointer :: temp_field
241
242 n = this%dof%size()
243
244 if (tstep .le. 3) then
245 call field_cfill(this%reg_coeff, 0.0_rp, n)
246 return
247 end if
248
249 call neko_scratch_registry%request_field(temp_field, temp_indices(1), .false.)
250
251 call field_cfill(temp_field, 1.0_rp, n)
252 s_mean = field_glsum(this%S, n) / field_glsum(temp_field, n)
253
254 call field_copy(temp_field, this%S, n)
255 call field_cadd(temp_field, -s_mean, n)
256
257 if (neko_bcknd_device .eq. 1) then
258 call device_absval(temp_field%x_d, n)
259 call device_memcpy(temp_field%x, temp_field%x_d, n, device_to_host, &
260 sync = .true.)
261 else
262 call absval(temp_field%x, n)
263 end if
264
265 ! Normalization factor n_S = |S-S_mean|_inf
266 n_s = glmax(temp_field%x, n)
267
268 call neko_scratch_registry%relinquish_field(temp_indices)
269
270 if (n_s < 1.0e-12_rp) then
271 n_s = 1.0e-12_rp
272 end if
273
274 ! entropy viscosity = c_avisc_entropy * h^2 * entropy_residual / n_S
275 if (neko_bcknd_device .eq. 1) then
277 this%reg_coeff%x_d, this%entropy_residual%x_d, &
278 this%h%x_d, this%c_avisc_entropy, n_s, n)
279 else
281 this%reg_coeff%x, this%entropy_residual%x, &
282 this%h%x, this%c_avisc_entropy, n_s, n)
283 end if
284
285 ! effective viscosity = min(entropy viscosity, low-order viscosity)
286 if (neko_bcknd_device .eq. 1) then
288 this%reg_coeff%x_d, this%h%x_d, this%max_wave_speed%x_d, &
289 this%c_avisc_low, n)
290 else
292 this%reg_coeff%x, this%h%x, this%max_wave_speed%x, &
293 this%c_avisc_low, n)
294 end if
295
296 call this%apply_element_max()
297
298 call this%smooth_viscosity()
299
301
305 class(entropy_viscosity_t), intent(inout) :: this
306 integer :: n
307 type(field_t), pointer :: temp_field, mult_field
308 integer :: temp_indices(2)
309
310 n = this%dof%size()
311
312 call neko_scratch_registry%request_field(temp_field, temp_indices(1), .false.)
313 call neko_scratch_registry%request_field(mult_field, temp_indices(2), .false.)
314
315 call field_copy(temp_field, this%reg_coeff, n)
316 call this%gs%op(temp_field, gs_op_add)
317
318 call field_cfill(mult_field, 1.0_rp, n)
319 call this%gs%op(mult_field, gs_op_add)
320
321 if (neko_bcknd_device .eq. 1) then
323 this%reg_coeff%x_d, temp_field%x_d, mult_field%x_d, n)
324 else
326 this%reg_coeff%x, temp_field%x, mult_field%x, n)
327 end if
328
329 call neko_scratch_registry%relinquish_field(temp_indices)
330
332
334 class(entropy_viscosity_t), intent(inout) :: this
335 integer :: lx
336
337 lx = this%Xh%lx
338
339 if (neko_bcknd_device .eq. 1) then
341 this%reg_coeff%x_d, lx, this%msh%nelv)
342 else
344 this%reg_coeff%x, lx, this%msh%nelv)
345 end if
346
348
349 subroutine entropy_viscosity_set_fields(this, S, u, v, w, h, max_wave_speed, &
350 msh, Xh, gs)
351 class(entropy_viscosity_t), intent(inout) :: this
352 type(field_t), target, intent(inout) :: s
353 type(field_t), target, intent(in) :: u, v, w, h, max_wave_speed
354 type(mesh_t), target, intent(in) :: msh
355 type(space_t), target, intent(in) :: xh
356 type(gs_t), target, intent(in) :: gs
357
358 this%S => s
359 this%u => u
360 this%v => v
361 this%w => w
362 this%h => h
363 this%max_wave_speed => max_wave_speed
364 this%msh => msh
365 this%Xh => xh
366 this%gs => gs
367
368 call this%S_lag%init(s, 3)
369
370 end subroutine entropy_viscosity_set_fields
371
373 class(entropy_viscosity_t), intent(inout) :: this
374
375 call this%S_lag%update()
376
377 end subroutine entropy_viscosity_update_lag
378
380 pure function entropy_viscosity_low_order(this, i) result(visc)
381 class(entropy_viscosity_t), intent(in) :: this
382 integer, intent(in) :: i
383 real(kind=rp) :: visc
384
385 visc = this%c_avisc_low * this%h%x(i,1,1,1) * this%max_wave_speed%x(i,1,1,1)
386
387 end function entropy_viscosity_low_order
388
389end module entropy_viscosity
Copy data between host and device (or device and device)
Definition device.F90:71
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
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:47
integer, parameter, public device_to_host
Definition device.F90:47
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_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:1372
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
Definition math.f90:516
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
subroutine, public div(res, ux, uy, uz, coef)
Compute the divergence of a vector field.
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:56
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.