Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
device_mathops.F90
Go to the documentation of this file.
1! Copyright (c) 2021-2023, 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, c_rp
35 use utils, only : neko_error
36 use, intrinsic :: iso_c_binding, only : c_int, c_ptr
37 implicit none
38 private
39
40#ifdef HAVE_HIP
41 interface
42 subroutine hip_opchsign(a1_d, a2_d, a3_d, gdim, n) &
43 bind(c, name='hip_opchsign')
44 use, intrinsic :: iso_c_binding
45 type(c_ptr), value :: a1_d, a2_d, a3_d
46 integer(c_int) :: gdim, n
47 end subroutine hip_opchsign
48 end interface
49
50 interface
51 subroutine hip_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n) &
52 bind(c, name='hip_opcolv')
53 use, intrinsic :: iso_c_binding
54 type(c_ptr), value :: a1_d, a2_d, a3_d, c_d
55 integer(c_int) :: gdim, n
56 end subroutine hip_opcolv
57 end interface
58
59 interface
60 subroutine hip_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, gdim, n) &
61 bind(c, name='hip_opcolv3c')
62 use, intrinsic :: iso_c_binding
63 import c_rp
64 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
65 real(c_rp) :: d
66 integer(c_int) :: gdim, n
67 end subroutine hip_opcolv3c
68 end interface
69
70 interface
71 subroutine hip_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, gdim, n) &
72 bind(c, name='hip_opadd2cm')
73 use, intrinsic :: iso_c_binding
74 import c_rp
75 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d
76 real(c_rp) :: c
77 integer(c_int) :: gdim, n
78 end subroutine hip_opadd2cm
79 end interface
80
81 interface
82 subroutine hip_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, gdim, n) &
83 bind(c, name='hip_opadd2col')
84 use, intrinsic :: iso_c_binding
85 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
86 integer(c_int) :: gdim, n
87 end subroutine hip_opadd2col
88 end interface
89#elif HAVE_CUDA
90 interface
91 subroutine cuda_opchsign(a1_d, a2_d, a3_d, gdim, n) &
92 bind(c, name='cuda_opchsign')
93 use, intrinsic :: iso_c_binding
94 type(c_ptr), value :: a1_d, a2_d, a3_d
95 integer(c_int) :: gdim, n
96 end subroutine cuda_opchsign
97 end interface
98
99 interface
100 subroutine cuda_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n) &
101 bind(c, name='cuda_opcolv')
102 use, intrinsic :: iso_c_binding
103 type(c_ptr), value :: a1_d, a2_d, a3_d, c_d
104 integer(c_int) :: gdim, n
105 end subroutine cuda_opcolv
106 end interface
107
108 interface
109 subroutine cuda_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, gdim, n) &
110 bind(c, name='cuda_opcolv3c')
111 use, intrinsic :: iso_c_binding
112 import c_rp
113 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
114 real(c_rp) :: d
115 integer(c_int) :: gdim, n
116 end subroutine cuda_opcolv3c
117 end interface
118
119 interface
120 subroutine cuda_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, gdim, n) &
121 bind(c, name='cuda_opadd2cm')
122 use, intrinsic :: iso_c_binding
123 import c_rp
124 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d
125 real(c_rp) :: c
126 integer(c_int) :: gdim, n
127 end subroutine cuda_opadd2cm
128 end interface
129
130 interface
131 subroutine cuda_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, gdim, n) &
132 bind(c, name='cuda_opadd2col')
133 use, intrinsic :: iso_c_binding
134 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
135 integer(c_int) :: gdim, n
136 end subroutine cuda_opadd2col
137 end interface
138#elif HAVE_OPENCL
139 interface
140 subroutine opencl_opchsign(a1_d, a2_d, a3_d, gdim, n) &
141 bind(c, name='opencl_opchsign')
142 use, intrinsic :: iso_c_binding
143 type(c_ptr), value :: a1_d, a2_d, a3_d
144 integer(c_int) :: gdim, n
145 end subroutine opencl_opchsign
146 end interface
147
148 interface
149 subroutine opencl_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n) &
150 bind(c, name='opencl_opcolv')
151 use, intrinsic :: iso_c_binding
152 type(c_ptr), value :: a1_d, a2_d, a3_d, c_d
153 integer(c_int) :: gdim, n
154 end subroutine opencl_opcolv
155 end interface
156
157 interface
158 subroutine opencl_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, gdim, n) &
159 bind(c, name='opencl_opcolv3c')
160 use, intrinsic :: iso_c_binding
161 import c_rp
162 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
163 real(c_rp) :: d
164 integer(c_int) :: gdim, n
165 end subroutine opencl_opcolv3c
166 end interface
167
168 interface
169 subroutine opencl_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, gdim, n) &
170 bind(c, name='opencl_opadd2cm')
171 use, intrinsic :: iso_c_binding
172 import c_rp
173 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d
174 real(c_rp) :: c
175 integer(c_int) :: gdim, n
176 end subroutine opencl_opadd2cm
177 end interface
178
179 interface
180 subroutine opencl_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, gdim, n) &
181 bind(c, name='opencl_opadd2col')
182 use, intrinsic :: iso_c_binding
183 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
184 integer(c_int) :: gdim, n
185 end subroutine opencl_opadd2col
186 end interface
187#endif
188
191
192contains
193
195 subroutine device_opchsign(a1_d, a2_d, a3_d, gdim, n)
196 type(c_ptr) :: a1_d, a2_d, a3_d
197 integer :: n, gdim
198#ifdef HAVE_HIP
199 call hip_opchsign(a1_d, a2_d, a3_d, gdim, n)
200#elif HAVE_CUDA
201 call cuda_opchsign(a1_d, a2_d, a3_d, gdim, n)
202#elif HAVE_OPENCL
203 call opencl_opchsign(a1_d, a2_d, a3_d, gdim, n)
204#else
205 call neko_error('No device backend configured')
206#endif
207 end subroutine device_opchsign
208
210 subroutine device_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n)
211 type(c_ptr) :: a1_d, a2_d, a3_d, c_d
212 integer :: n, gdim
213#ifdef HAVE_HIP
214 call hip_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n)
215#elif HAVE_CUDA
216 call cuda_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n)
217#elif HAVE_OPENCL
218 call opencl_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n)
219#else
220 call neko_error('No device backend configured')
221#endif
222 end subroutine device_opcolv
223
225 subroutine device_opcolv3c(a1_d, a2_d, a3_d, &
226 b1_d, b2_d, b3_d, c_d, d, n, gdim)
227 type(c_ptr) :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
228 real(kind=rp) :: d
229 integer :: n, gdim
230#ifdef HAVE_HIP
231 call hip_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, gdim, n)
232#elif HAVE_CUDA
233 call cuda_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, gdim, n)
234#elif HAVE_OPENCL
235 call opencl_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, gdim, n)
236#else
237 call neko_error('No device backend configured')
238#endif
239 end subroutine device_opcolv3c
240
242 subroutine device_opadd2cm (a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, n, gdim)
243 type(c_ptr) :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d
244 real(kind=rp) :: c
245 integer :: n, gdim
246#ifdef HAVE_HIP
247 call hip_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, gdim, n)
248#elif HAVE_CUDA
249 call cuda_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, gdim, n)
250#elif HAVE_OPENCL
251 call opencl_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, gdim, n)
252#else
253 call neko_error('No device backend configured')
254#endif
255 end subroutine device_opadd2cm
256
258 subroutine device_opadd2col (a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, n, gdim)
259 type(c_ptr) :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
260 integer :: n, gdim
261#ifdef HAVE_HIP
262 call hip_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, gdim, n)
263#elif HAVE_CUDA
264 call cuda_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, gdim, n)
265#elif HAVE_OPENCL
266 call opencl_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, gdim, n)
267#else
268 call neko_error('No device backend configured')
269#endif
270 end subroutine device_opadd2col
271
272end module device_mathops
void opencl_opcolv3c(void *a1, void *a2, void *a3, void *b1, void *b2, void *b3, void *c, real *d, int *gdim, int *n)
Definition mathops.c:101
void opencl_opadd2cm(void *a1, void *a2, void *a3, void *b1, void *b2, void *b3, real *c, int *gdim, int *n)
Definition mathops.c:133
void opencl_opcolv(void *a1, void *a2, void *a3, void *c, int *gdim, int *n)
Definition mathops.c:75
void opencl_opchsign(void *a1, void *a2, void *a3, int *gdim, int *n)
Definition mathops.c:50
void opencl_opadd2col(void *a1, void *a2, void *a3, void *b1, void *b2, void *b3, void *c, int *gdim, int *n)
Definition mathops.c:164
void cuda_opchsign(void *a1, void *a2, void *a3, int *gdim, int *n)
Definition mathops.cu:42
void cuda_opcolv3c(void *a1, void *a2, void *a3, void *b1, void *b2, void *b3, void *c, real *d, int *gdim, int *n)
Definition mathops.cu:69
void cuda_opadd2cm(void *a1, void *a2, void *a3, void *b1, void *b2, void *b3, real *c, int *gdim, int *n)
Definition mathops.cu:85
void cuda_opcolv(void *a1, void *a2, void *a3, void *c, int *gdim, int *n)
Definition mathops.cu:55
void cuda_opadd2col(void *a1, void *a2, void *a3, void *b1, void *b2, void *b3, void *c, int *gdim, int *n)
Definition mathops.cu:101
subroutine, public device_opchsign(a1_d, a2_d, a3_d, gdim, n)
subroutine, public device_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, n, gdim)
subroutine, public device_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, n, gdim)
subroutine, public device_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, n, gdim)
subroutine, public device_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n)
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