Neko 1.99.4
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
entropy_viscosity_device.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!
35 use, intrinsic :: iso_c_binding, only: c_ptr, c_int
36 use num_types, only: rp, c_rp
37 use utils, only: neko_error
38 implicit none
39 private
40
41#ifdef HAVE_HIP
42 interface
43 subroutine hip_entropy_visc_compute_residual(entropy_residual_d, &
44 S_d, S_lag1_d, S_lag2_d, S_lag3_d, &
45 bdf1, bdf2, bdf3, bdf4, dt, n) &
46 bind(c, name = 'hip_entropy_visc_compute_residual')
47 use, intrinsic :: iso_c_binding
48 import c_rp
49 type(c_ptr), value :: entropy_residual_d
50 type(c_ptr), value :: S_d, S_lag1_d, S_lag2_d, S_lag3_d
51 real(c_rp) :: bdf1, bdf2, bdf3, bdf4, dt
52 integer(c_int) :: n
54 end interface
55
56 interface
57 subroutine hip_entropy_visc_compute_viscosity(reg_coeff_d, &
58 entropy_residual_d, h_d, c_avisc_entropy, n_S, n) &
59 bind(c, name = 'hip_entropy_visc_compute_viscosity')
60 use, intrinsic :: iso_c_binding
61 import c_rp
62 type(c_ptr), value :: reg_coeff_d, entropy_residual_d, h_d
63 real(c_rp) :: c_avisc_entropy, n_S
64 integer(c_int) :: n
66 end interface
67
68 interface
69 subroutine hip_entropy_visc_apply_element_max(reg_coeff_d, lx, nelv) &
70 bind(c, name = 'hip_entropy_visc_apply_element_max')
71 use, intrinsic :: iso_c_binding
72 type(c_ptr), value :: reg_coeff_d
73 integer(c_int) :: lx, nelv
75 end interface
76
77 interface
78 subroutine hip_entropy_visc_clamp_to_low_order(reg_coeff_d, &
79 h_d, max_wave_speed_d, c_avisc_low, n) &
80 bind(c, name = 'hip_entropy_visc_clamp_to_low_order')
81 use, intrinsic :: iso_c_binding
82 import c_rp
83 type(c_ptr), value :: reg_coeff_d, h_d, max_wave_speed_d
84 real(c_rp) :: c_avisc_low
85 integer(c_int) :: n
87 end interface
88
89 interface
90 subroutine hip_entropy_visc_smooth_divide(reg_coeff_d, &
91 temp_field_d, mult_field_d, n) &
92 bind(c, name = 'hip_entropy_visc_smooth_divide')
93 use, intrinsic :: iso_c_binding
94 type(c_ptr), value :: reg_coeff_d, temp_field_d, mult_field_d
95 integer(c_int) :: n
97 end interface
98
99#elif HAVE_CUDA
100 interface
101 subroutine cuda_entropy_visc_compute_residual(entropy_residual_d, &
102 S_d, S_lag1_d, S_lag2_d, S_lag3_d, &
103 bdf1, bdf2, bdf3, bdf4, dt, n) &
104 bind(c, name = 'cuda_entropy_visc_compute_residual')
105 use, intrinsic :: iso_c_binding
106 import c_rp
107 type(c_ptr), value :: entropy_residual_d
108 type(c_ptr), value :: S_d, S_lag1_d, S_lag2_d, S_lag3_d
109 real(c_rp) :: bdf1, bdf2, bdf3, bdf4, dt
110 integer(c_int) :: n
112 end interface
113
114 interface
115 subroutine cuda_entropy_visc_compute_viscosity(reg_coeff_d, &
116 entropy_residual_d, h_d, c_avisc_entropy, n_S, n) &
117 bind(c, name = 'cuda_entropy_visc_compute_viscosity')
118 use, intrinsic :: iso_c_binding
119 import c_rp
120 type(c_ptr), value :: reg_coeff_d, entropy_residual_d, h_d
121 real(c_rp) :: c_avisc_entropy, n_S
122 integer(c_int) :: n
124 end interface
125
126 interface
127 subroutine cuda_entropy_visc_apply_element_max(reg_coeff_d, lx, nelv) &
128 bind(c, name = 'cuda_entropy_visc_apply_element_max')
129 use, intrinsic :: iso_c_binding
130 type(c_ptr), value :: reg_coeff_d
131 integer(c_int) :: lx, nelv
133 end interface
134
135 interface
136 subroutine cuda_entropy_visc_clamp_to_low_order(reg_coeff_d, &
137 h_d, max_wave_speed_d, c_avisc_low, n) &
138 bind(c, name = 'cuda_entropy_visc_clamp_to_low_order')
139 use, intrinsic :: iso_c_binding
140 import c_rp
141 type(c_ptr), value :: reg_coeff_d, h_d, max_wave_speed_d
142 real(c_rp) :: c_avisc_low
143 integer(c_int) :: n
145 end interface
146
147 interface
148 subroutine cuda_entropy_visc_smooth_divide(reg_coeff_d, &
149 temp_field_d, mult_field_d, n) &
150 bind(c, name = 'cuda_entropy_visc_smooth_divide')
151 use, intrinsic :: iso_c_binding
152 type(c_ptr), value :: reg_coeff_d, temp_field_d, mult_field_d
153 integer(c_int) :: n
155 end interface
156
157#elif HAVE_OPENCL
158 interface
159 subroutine opencl_entropy_visc_compute_residual(entropy_residual_d, &
160 S_d, S_lag1_d, S_lag2_d, S_lag3_d, &
161 bdf1, bdf2, bdf3, bdf4, dt, n) &
162 bind(c, name = 'opencl_entropy_visc_compute_residual')
163 use, intrinsic :: iso_c_binding
164 import c_rp
165 type(c_ptr), value :: entropy_residual_d
166 type(c_ptr), value :: S_d, S_lag1_d, S_lag2_d, S_lag3_d
167 real(c_rp), value :: bdf1, bdf2, bdf3, bdf4, dt
168 integer(c_int), value :: n
170 end interface
171
172 interface
173 subroutine opencl_entropy_visc_compute_viscosity(reg_coeff_d, &
174 entropy_residual_d, h_d, c_avisc_entropy, n_S, n) &
175 bind(c, name = 'opencl_entropy_visc_compute_viscosity')
176 use, intrinsic :: iso_c_binding
177 import c_rp
178 type(c_ptr), value :: reg_coeff_d, entropy_residual_d, h_d
179 real(c_rp), value :: c_avisc_entropy, n_S
180 integer(c_int), value :: n
182 end interface
183
184 interface
185 subroutine opencl_entropy_visc_apply_element_max(reg_coeff_d, lx, nelv) &
186 bind(c, name = 'opencl_entropy_visc_apply_element_max')
187 use, intrinsic :: iso_c_binding
188 type(c_ptr), value :: reg_coeff_d
189 integer(c_int), value :: lx, nelv
191 end interface
192
193 interface
194 subroutine opencl_entropy_visc_clamp_to_low_order(reg_coeff_d, &
195 h_d, max_wave_speed_d, c_avisc_low, n) &
196 bind(c, name = 'opencl_entropy_visc_clamp_to_low_order')
197 use, intrinsic :: iso_c_binding
198 import c_rp
199 type(c_ptr), value :: reg_coeff_d, h_d, max_wave_speed_d
200 real(c_rp), value :: c_avisc_low
201 integer(c_int), value :: n
203 end interface
204
205 interface
206 subroutine opencl_entropy_visc_smooth_divide(reg_coeff_d, &
207 temp_field_d, mult_field_d, n) &
208 bind(c, name = 'opencl_entropy_visc_smooth_divide')
209 use, intrinsic :: iso_c_binding
210 type(c_ptr), value :: reg_coeff_d, temp_field_d, mult_field_d
211 integer(c_int), value :: n
213 end interface
214#elif HAVE_METAL
215 interface
216 subroutine metal_entropy_visc_compute_residual(entropy_residual_d, &
217 S_d, S_lag1_d, S_lag2_d, S_lag3_d, &
218 bdf1, bdf2, bdf3, bdf4, dt, n) &
219 bind(c, name = 'metal_entropy_visc_compute_residual')
220 use, intrinsic :: iso_c_binding
221 import c_rp
222 type(c_ptr), value :: entropy_residual_d
223 type(c_ptr), value :: S_d, S_lag1_d, S_lag2_d, S_lag3_d
224 real(c_rp), value :: bdf1, bdf2, bdf3, bdf4, dt
225 integer(c_int), value :: n
226 end subroutine metal_entropy_visc_compute_residual
227 end interface
228
229 interface
230 subroutine metal_entropy_visc_compute_viscosity(reg_coeff_d, &
231 entropy_residual_d, h_d, c_avisc_entropy, n_S, n) &
232 bind(c, name = 'metal_entropy_visc_compute_viscosity')
233 use, intrinsic :: iso_c_binding
234 import c_rp
235 type(c_ptr), value :: reg_coeff_d, entropy_residual_d, h_d
236 real(c_rp), value :: c_avisc_entropy, n_S
237 integer(c_int), value :: n
238 end subroutine metal_entropy_visc_compute_viscosity
239 end interface
240
241 interface
242 subroutine metal_entropy_visc_apply_element_max(reg_coeff_d, lx, nelv) &
243 bind(c, name = 'metal_entropy_visc_apply_element_max')
244 use, intrinsic :: iso_c_binding
245 type(c_ptr), value :: reg_coeff_d
246 integer(c_int), value :: lx, nelv
247 end subroutine metal_entropy_visc_apply_element_max
248 end interface
249
250 interface
251 subroutine metal_entropy_visc_clamp_to_low_order(reg_coeff_d, &
252 h_d, max_wave_speed_d, c_avisc_low, n) &
253 bind(c, name = 'metal_entropy_visc_clamp_to_low_order')
254 use, intrinsic :: iso_c_binding
255 import c_rp
256 type(c_ptr), value :: reg_coeff_d, h_d, max_wave_speed_d
257 real(c_rp), value :: c_avisc_low
258 integer(c_int), value :: n
259 end subroutine metal_entropy_visc_clamp_to_low_order
260 end interface
261
262 interface
263 subroutine metal_entropy_visc_smooth_divide(reg_coeff_d, &
264 temp_field_d, mult_field_d, n) &
265 bind(c, name = 'metal_entropy_visc_smooth_divide')
266 use, intrinsic :: iso_c_binding
267 type(c_ptr), value :: reg_coeff_d, temp_field_d, mult_field_d
268 integer(c_int), value :: n
269 end subroutine metal_entropy_visc_smooth_divide
270 end interface
271#endif
272
278
279contains
280
282 subroutine entropy_viscosity_compute_residual_device(entropy_residual_d, &
283 S_d, S_lag1_d, S_lag2_d, S_lag3_d, bdf_coeffs, dt, n)
284 type(c_ptr), intent(in) :: entropy_residual_d
285 type(c_ptr), intent(in) :: s_d, s_lag1_d, s_lag2_d, s_lag3_d
286 real(kind=rp), intent(in) :: bdf_coeffs(4)
287 real(kind=rp), intent(in) :: dt
288 integer, intent(in) :: n
289
290#ifdef HAVE_HIP
291 call hip_entropy_visc_compute_residual(entropy_residual_d, &
292 s_d, s_lag1_d, s_lag2_d, s_lag3_d, &
293 bdf_coeffs(1), bdf_coeffs(2), bdf_coeffs(3), bdf_coeffs(4), dt, n)
294#elif HAVE_CUDA
295 call cuda_entropy_visc_compute_residual(entropy_residual_d, &
296 s_d, s_lag1_d, s_lag2_d, s_lag3_d, &
297 bdf_coeffs(1), bdf_coeffs(2), bdf_coeffs(3), bdf_coeffs(4), dt, n)
298#elif HAVE_OPENCL
299 call opencl_entropy_visc_compute_residual(entropy_residual_d, &
300 s_d, s_lag1_d, s_lag2_d, s_lag3_d, &
301 bdf_coeffs(1), bdf_coeffs(2), bdf_coeffs(3), bdf_coeffs(4), dt, n)
302#elif HAVE_METAL
303 call metal_entropy_visc_compute_residual(entropy_residual_d, &
304 s_d, s_lag1_d, s_lag2_d, s_lag3_d, &
305 bdf_coeffs(1), bdf_coeffs(2), bdf_coeffs(3), bdf_coeffs(4), dt, n)
306#else
307 call neko_error('No device backend configured')
308#endif
310
313 entropy_residual_d, h_d, c_avisc_entropy, n_S, n)
314 type(c_ptr), intent(in) :: reg_coeff_d, entropy_residual_d, h_d
315 real(kind=rp), intent(in) :: c_avisc_entropy, n_s
316 integer, intent(in) :: n
317
318#ifdef HAVE_HIP
319 call hip_entropy_visc_compute_viscosity(reg_coeff_d, &
320 entropy_residual_d, h_d, c_avisc_entropy, n_s, n)
321#elif HAVE_CUDA
322 call cuda_entropy_visc_compute_viscosity(reg_coeff_d, &
323 entropy_residual_d, h_d, c_avisc_entropy, n_s, n)
324#elif HAVE_OPENCL
325 call opencl_entropy_visc_compute_viscosity(reg_coeff_d, &
326 entropy_residual_d, h_d, c_avisc_entropy, n_s, n)
327#elif HAVE_METAL
328 call metal_entropy_visc_compute_viscosity(reg_coeff_d, &
329 entropy_residual_d, h_d, c_avisc_entropy, n_s, n)
330#else
331 call neko_error('No device backend configured')
332#endif
334
336 subroutine entropy_viscosity_apply_element_max_device(reg_coeff_d, lx, nelv)
337 type(c_ptr), intent(in) :: reg_coeff_d
338 integer, intent(in) :: lx, nelv
339
340#ifdef HAVE_HIP
341 call hip_entropy_visc_apply_element_max(reg_coeff_d, lx, nelv)
342#elif HAVE_CUDA
343 call cuda_entropy_visc_apply_element_max(reg_coeff_d, lx, nelv)
344#elif HAVE_OPENCL
345 call opencl_entropy_visc_apply_element_max(reg_coeff_d, lx, nelv)
346#elif HAVE_METAL
347 call metal_entropy_visc_apply_element_max(reg_coeff_d, lx, nelv)
348#else
349 call neko_error('No device backend configured')
350#endif
352
355 h_d, max_wave_speed_d, c_avisc_low, n)
356 type(c_ptr), intent(in) :: reg_coeff_d, h_d, max_wave_speed_d
357 real(kind=rp), intent(in) :: c_avisc_low
358 integer, intent(in) :: n
359
360#ifdef HAVE_HIP
361 call hip_entropy_visc_clamp_to_low_order(reg_coeff_d, &
362 h_d, max_wave_speed_d, c_avisc_low, n)
363#elif HAVE_CUDA
364 call cuda_entropy_visc_clamp_to_low_order(reg_coeff_d, &
365 h_d, max_wave_speed_d, c_avisc_low, n)
366#elif HAVE_OPENCL
368 h_d, max_wave_speed_d, c_avisc_low, n)
369#elif HAVE_METAL
370 call metal_entropy_visc_clamp_to_low_order(reg_coeff_d, &
371 h_d, max_wave_speed_d, c_avisc_low, n)
372#else
373 call neko_error('No device backend configured')
374#endif
376
379 temp_field_d, mult_field_d, n)
380 type(c_ptr), intent(in) :: reg_coeff_d, temp_field_d, mult_field_d
381 integer, intent(in) :: n
382
383#ifdef HAVE_HIP
384 call hip_entropy_visc_smooth_divide(reg_coeff_d, &
385 temp_field_d, mult_field_d, n)
386#elif HAVE_CUDA
387 call cuda_entropy_visc_smooth_divide(reg_coeff_d, &
388 temp_field_d, mult_field_d, n)
389#elif HAVE_OPENCL
390 call opencl_entropy_visc_smooth_divide(reg_coeff_d, &
391 temp_field_d, mult_field_d, n)
392#elif HAVE_METAL
393 call metal_entropy_visc_smooth_divide(reg_coeff_d, &
394 temp_field_d, mult_field_d, n)
395#else
396 call neko_error('No device backend configured')
397#endif
399
401
void opencl_entropy_visc_compute_viscosity(void *reg_coeff, void *entropy_residual, void *h, real c_avisc_entropy, real n_S, int n)
void opencl_entropy_visc_compute_residual(void *entropy_residual, void *S, void *S_lag1, void *S_lag2, void *S_lag3, real bdf1, real bdf2, real bdf3, real bdf4, real dt, int n)
void opencl_entropy_visc_smooth_divide(void *reg_coeff, void *temp_field, void *mult_field, int n)
void opencl_entropy_visc_clamp_to_low_order(void *reg_coeff, void *h, void *max_wave_speed, real c_avisc_low, int n)
void opencl_entropy_visc_apply_element_max(void *reg_coeff, int lx, int nelv)
void cuda_entropy_visc_apply_element_max(void *reg_coeff, int *lx, int *nelv)
void cuda_entropy_visc_smooth_divide(void *reg_coeff, void *temp_field, void *mult_field, int *n)
void cuda_entropy_visc_compute_residual(void *entropy_residual, void *S, void *S_lag1, void *S_lag2, void *S_lag3, real *bdf1, real *bdf2, real *bdf3, real *bdf4, real *dt, int *n)
void cuda_entropy_visc_clamp_to_low_order(void *reg_coeff, void *h, void *max_wave_speed, real *c_avisc_low, int *n)
void cuda_entropy_visc_compute_viscosity(void *reg_coeff, void *entropy_residual, void *h, real *c_avisc_entropy, real *n_S, int *n)
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.
integer, parameter, public c_rp
Definition num_types.f90:13
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Utilities.
Definition utils.f90:35