Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
device_deardorff_nut.F90
Go to the documentation of this file.
1! Copyright (c) 2025-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!
34
36 use, intrinsic :: iso_c_binding, only: c_ptr, c_int
37 use num_types, only: rp, c_rp
38 use utils, only: neko_error
40 use mpi_f08, only: mpi_sum, mpi_in_place, mpi_allreduce
41
42 implicit none
43 private
44
45#ifdef HAVE_HIP
46 interface
47 subroutine hip_deardorff_nut_compute(TKE_d, &
48 dTdx_d, dTdy_d, dTdz_d, &
49 a11_d, a12_d, a13_d, &
50 a21_d, a22_d, a23_d, &
51 a31_d, a32_d, a33_d, &
52 delta_d, nut_d, temperature_alphat, TKE_alphat, TKE_source, &
53 c_k, T0, g1, g2, g3, &
54 eps, n) bind(C,name="hip_deardorff_nut_compute")
55 use, intrinsic :: iso_c_binding, only: c_ptr, c_int
56 import c_rp
57 type(c_ptr), value :: TKE_d, &
58 dTdx_d, dTdy_d, dTdz_d, &
59 a11_d, a12_d, a13_d, &
60 a21_d, a22_d, a23_d, &
61 a31_d, a32_d, a33_d, &
62 delta_d, nut_d, temperature_alphat, &
63 TKE_alphat, TKE_source
64 integer(c_int) :: n
65 real(c_rp) :: c_k, T0, g1, g2, g3, eps
66 end subroutine hip_deardorff_nut_compute
67 end interface
68#elif HAVE_CUDA
69 interface
70 subroutine cuda_deardorff_nut_compute(TKE_d, &
71 dTdx_d, dTdy_d, dTdz_d, &
72 a11_d, a12_d, a13_d, &
73 a21_d, a22_d, a23_d, &
74 a31_d, a32_d, a33_d, &
75 delta_d, nut_d, temperature_alphat, TKE_alphat, TKE_source, &
76 c_k, T0, g1, g2, g3, &
77 eps, n) &
78 bind(c, name = 'cuda_deardorff_nut_compute')
79 use, intrinsic :: iso_c_binding, only: c_ptr, c_int
80 import c_rp
81 type(c_ptr), value :: TKE_d, &
82 dTdx_d, dTdy_d, dTdz_d, &
83 a11_d, a12_d, a13_d, &
84 a21_d, a22_d, a23_d, &
85 a31_d, a32_d, a33_d, &
86 delta_d, nut_d, temperature_alphat, &
87 TKE_alphat, TKE_source
88 integer(c_int) :: n
89 real(c_rp) :: c_k, T0, g1, g2, g3, eps
90 end subroutine cuda_deardorff_nut_compute
91 end interface
92#elif HAVE_OPENCL
93#endif
94
96
97contains
98
124 dTdx_d, dTdy_d, dTdz_d, a11_d, a12_d, a13_d, &
125 a21_d, a22_d, a23_d, &
126 a31_d, a32_d, a33_d, &
127 delta_d, &
128 nut_d, temperature_alphat, TKE_alphat, TKE_source, &
129 c_k, T0, g, eps, n)
130 type(c_ptr) :: tke_d, dtdx_d, dtdy_d, dtdz_d, &
131 a11_d, a12_d, a13_d, &
132 a21_d, a22_d, a23_d, &
133 a31_d, a32_d, a33_d, &
134 delta_d, nut_d, temperature_alphat, tke_alphat, tke_source
135 integer :: n
136 real(kind=rp) :: c_k, t0, g(3), eps
137#if HAVE_HIP
138 call hip_deardorff_nut_compute(tke_d, &
139 dtdx_d, dtdy_d, dtdz_d, &
140 a11_d, a12_d, a13_d, &
141 a21_d, a22_d, a23_d, &
142 a31_d, a32_d, a33_d, &
143 delta_d, nut_d, temperature_alphat, tke_alphat, tke_source, &
144 c_k, t0, g(1), g(2), g(3), eps, n)
145#elif HAVE_CUDA
146 call cuda_deardorff_nut_compute(tke_d, &
147 dtdx_d, dtdy_d, dtdz_d, &
148 a11_d, a12_d, a13_d, &
149 a21_d, a22_d, a23_d, &
150 a31_d, a32_d, a33_d, &
151 delta_d, nut_d, temperature_alphat, tke_alphat, tke_source, &
152 c_k, t0, g(1), g(2), g(3), eps, n)
153#elif HAVE_OPENCL
154 call neko_error('opencl backend is not supported for device_deardorff_nut')
155#else
156 call neko_error('no device backend configured')
157#endif
158 end subroutine device_deardorff_nut_compute
159
160
161end module device_deardorff_nut
void cuda_deardorff_nut_compute(void *TKE, void *dTdx, void *dTdy, void *dTdz, void *a11, void *a12, void *a13, void *a21, void *a22, void *a23, void *a31, void *a32, void *a33, void *delta, void *nut, void *temperature_alphat, void *TKE_alphat, void *TKE_source, real *c_k, real *T0, real *g1, real *g2, real *g3, real *eps, int *n)
Definition comm.F90:1
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
Definition comm.F90:51
integer, public pe_size
MPI size of communicator.
Definition comm.F90:59
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
Device kernel wrapper for computing Deardorff SGS quantities.
subroutine, public device_deardorff_nut_compute(tke_d, dtdx_d, dtdy_d, dtdz_d, a11_d, a12_d, a13_d, a21_d, a22_d, a23_d, a31_d, a32_d, a33_d, delta_d, nut_d, temperature_alphat, tke_alphat, tke_source, c_k, t0, g, eps, n)
Compute Deardorff SGS quantities on the device backend.
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