Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
device_vreman_nut.F90
Go to the documentation of this file.
1! Copyright (c) 2024, 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, intrinsic :: iso_c_binding, only: c_ptr, c_int
35 use num_types, only: rp, c_rp
36 use utils, only: neko_error
38 use mpi_f08, only: mpi_sum, mpi_in_place, mpi_allreduce
39
40 implicit none
41 private
42
43#ifdef HAVE_HIP
44 interface
45 subroutine hip_vreman_nut_compute(a11_d, a12_d, a13_d, &
46 a21_d, a22_d, a23_d, &
47 a31_d, a32_d, a33_d, &
48 delta_d, nut_d, mult_d, c, eps, n) &
49 bind(c, name = 'hip_vreman_nut_compute')
50 use, intrinsic :: iso_c_binding, only: c_ptr, c_int
51 import c_rp
52 type(c_ptr), value :: a11_d, a12_d, a13_d, &
53 a21_d, a22_d, a23_d, &
54 a31_d, a32_d, a33_d, &
55 delta_d, nut_d, mult_d
56 integer(c_int) :: n
57 real(c_rp) :: c, eps
58 end subroutine hip_vreman_nut_compute
59
60
61 subroutine hip_vreman_nut_compute_buoy(a11_d, a12_d, a13_d, &
62 a21_d, a22_d, a23_d, &
63 a31_d, a32_d, a33_d, &
64 delta_d, nut_d, mult_d, c, eps, n, &
65 dTdx_d, dTdy_d, dTdz_d, g, ri_c, ref_temp) &
66 bind(c, name = 'hip_vreman_nut_compute_buoy')
67 use, intrinsic :: iso_c_binding, only: c_ptr, c_int
68 import c_rp
69 type(c_ptr), value :: a11_d, a12_d, a13_d, &
70 a21_d, a22_d, a23_d, &
71 a31_d, a32_d, a33_d, &
72 delta_d, nut_d, mult_d
73 integer(c_int) :: n
74 real(c_rp) :: c, eps
75 type(c_ptr), value :: dTdx_d, dTdy_d, dTdz_d
76 real(c_rp) :: g(3)
77 real(c_rp) :: ri_c, ref_temp
78 end subroutine hip_vreman_nut_compute_buoy
79
80 end interface
81#elif HAVE_CUDA
82 interface
83 subroutine cuda_vreman_nut_compute(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, mult_d, c, eps, n) &
87 bind(c, name = 'cuda_vreman_nut_compute')
88 use, intrinsic :: iso_c_binding, only: c_ptr, c_int
89 import c_rp
90 type(c_ptr), value :: a11_d, a12_d, a13_d, &
91 a21_d, a22_d, a23_d, &
92 a31_d, a32_d, a33_d, &
93 delta_d, nut_d, mult_d
94 integer(c_int) :: n
95 real(c_rp) :: c, eps
96 end subroutine cuda_vreman_nut_compute
97
98
99 subroutine cuda_vreman_nut_compute_buoy(a11_d, a12_d, a13_d, &
100 a21_d, a22_d, a23_d, &
101 a31_d, a32_d, a33_d, &
102 delta_d, nut_d, mult_d, c, eps, n, &
103 dTdx_d, dTdy_d, dTdz_d, g, ri_c, ref_temp) &
104 bind(c, name = 'cuda_vreman_nut_compute_buoy')
105 use, intrinsic :: iso_c_binding, only: c_ptr, c_int
106 import c_rp
107 type(c_ptr), value :: a11_d, a12_d, a13_d, &
108 a21_d, a22_d, a23_d, &
109 a31_d, a32_d, a33_d, &
110 delta_d, nut_d, mult_d
111 integer(c_int) :: n
112 real(c_rp) :: c, eps
113 type(c_ptr), value :: dTdx_d, dTdy_d, dTdz_d
114 real(c_rp) :: g(3)
115 real(c_rp) :: ri_c, ref_temp
116 end subroutine cuda_vreman_nut_compute_buoy
117
118 end interface
119#elif HAVE_OPENCL
120#endif
121
124
125contains
126
128 subroutine device_vreman_nut_compute(a11_d, a12_d, a13_d, &
129 a21_d, a22_d, a23_d, &
130 a31_d, a32_d, a33_d, &
131 delta_d, nut_d, mult_d, c, eps, n)
132 type(c_ptr) :: a11_d, a12_d, a13_d, &
133 a21_d, a22_d, a23_d, &
134 a31_d, a32_d, a33_d, &
135 delta_d, nut_d, mult_d
136 integer :: n
137 real(kind=rp) :: c, eps
138#if HAVE_HIP
139 call hip_vreman_nut_compute(a11_d, a12_d, a13_d, &
140 a21_d, a22_d, a23_d, &
141 a31_d, a32_d, a33_d, &
142 delta_d, nut_d, mult_d, c, eps, n)
143#elif HAVE_CUDA
144 call cuda_vreman_nut_compute(a11_d, a12_d, a13_d, &
145 a21_d, a22_d, a23_d, &
146 a31_d, a32_d, a33_d, &
147 delta_d, nut_d, mult_d, c, eps, n)
148#elif HAVE_OPENCL
149 call neko_error('opencl backend is not supported for device_vreman_nut')
150#else
151 call neko_error('no device backend configured')
152#endif
153 end subroutine device_vreman_nut_compute
154
155
156 subroutine device_vreman_nut_compute_buoy(a11_d, a12_d, a13_d, &
157 a21_d, a22_d, a23_d, &
158 a31_d, a32_d, a33_d, &
159 delta_d, nut_d, mult_d, c, eps, n, &
160 dTdx_d, dTdy_d, dTdz_d, g, ri_c, ref_temp)
161 type(c_ptr) :: a11_d, a12_d, a13_d, &
162 a21_d, a22_d, a23_d, &
163 a31_d, a32_d, a33_d, &
164 delta_d, nut_d, mult_d
165 type(c_ptr) :: dtdx_d, dtdy_d, dtdz_d
166 integer :: n
167 real(kind=rp) :: c, eps
168 real(kind=rp) :: g(3), ri_c, ref_temp
169#if HAVE_HIP
170 call hip_vreman_nut_compute_buoy(a11_d, a12_d, a13_d, &
171 a21_d, a22_d, a23_d, &
172 a31_d, a32_d, a33_d, &
173 delta_d, nut_d, mult_d, c, eps, n, &
174 dtdx_d, dtdy_d, dtdz_d, g, ri_c, ref_temp)
175#elif HAVE_CUDA
176 call cuda_vreman_nut_compute_buoy(a11_d, a12_d, a13_d, &
177 a21_d, a22_d, a23_d, &
178 a31_d, a32_d, a33_d, &
179 delta_d, nut_d, mult_d, c, eps, n, &
180 dtdx_d, dtdy_d, dtdz_d, g, ri_c, ref_temp)
181#elif HAVE_OPENCL
182 call neko_error('opencl backend is not supported for device_vreman_nut (buoyancy)')
183#else
184 call neko_error('no device backend configured (vreman buoyancy)')
185#endif
186 end subroutine device_vreman_nut_compute_buoy
187
188
189end module device_vreman_nut
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
subroutine, public device_vreman_nut_compute(a11_d, a12_d, a13_d, a21_d, a22_d, a23_d, a31_d, a32_d, a33_d, delta_d, nut_d, mult_d, c, eps, n)
Compute the eddy viscosity field for the Sigma model indevice.
subroutine, public device_vreman_nut_compute_buoy(a11_d, a12_d, a13_d, a21_d, a22_d, a23_d, a31_d, a32_d, a33_d, delta_d, nut_d, mult_d, c, eps, n, dtdx_d, dtdy_d, dtdz_d, g, ri_c, ref_temp)
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
void cuda_vreman_nut_compute(void *a11, void *a12, void *a13, void *a21, void *a22, void *a23, void *a31, void *a32, void *a33, void *delta, void *nut, void *mult, real *c, real *eps, int *n)
Definition vreman_nut.cu:50
void cuda_vreman_nut_compute_buoy(void *a11, void *a12, void *a13, void *a21, void *a22, void *a23, void *a31, void *a32, void *a33, void *delta, void *nut, void *mult, real *c, real *eps, int *n, void *dTdx, void *dTdy, void *dTdz, real *g, real *ri_c, real *ref_temp)
Definition vreman_nut.cu:70