Neko 1.99.2
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, &
61 c_d, d, gdim, n) &
62 bind(c, name = 'hip_opcolv3c')
63 use, intrinsic :: iso_c_binding
64 import c_rp
65 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
66 real(c_rp) :: d
67 integer(c_int) :: gdim, n
68 end subroutine hip_opcolv3c
69 end interface
70
71 interface
72 subroutine hip_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, &
73 c, gdim, n) &
74 bind(c, name = 'hip_opadd2cm')
75 use, intrinsic :: iso_c_binding
76 import c_rp
77 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d
78 real(c_rp) :: c
79 integer(c_int) :: gdim, n
80 end subroutine hip_opadd2cm
81 end interface
82
83 interface
84 subroutine hip_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, &
85 c_d, gdim, n) &
86 bind(c, name = 'hip_opadd2col')
87 use, intrinsic :: iso_c_binding
88 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
89 integer(c_int) :: gdim, n
90 end subroutine hip_opadd2col
91 end interface
92#elif HAVE_CUDA
93 interface
94 subroutine cuda_opchsign(a1_d, a2_d, a3_d, gdim, n) &
95 bind(c, name = 'cuda_opchsign')
96 use, intrinsic :: iso_c_binding
97 type(c_ptr), value :: a1_d, a2_d, a3_d
98 integer(c_int) :: gdim, n
99 end subroutine cuda_opchsign
100 end interface
101
102 interface
103 subroutine cuda_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n) &
104 bind(c, name = 'cuda_opcolv')
105 use, intrinsic :: iso_c_binding
106 type(c_ptr), value :: a1_d, a2_d, a3_d, c_d
107 integer(c_int) :: gdim, n
108 end subroutine cuda_opcolv
109 end interface
110
111 interface
112 subroutine cuda_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, &
113 c_d, d, gdim, n) &
114 bind(c, name = 'cuda_opcolv3c')
115 use, intrinsic :: iso_c_binding
116 import c_rp
117 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
118 real(c_rp) :: d
119 integer(c_int) :: gdim, n
120 end subroutine cuda_opcolv3c
121 end interface
122
123 interface
124 subroutine cuda_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, &
125 c, gdim, n) &
126 bind(c, name = 'cuda_opadd2cm')
127 use, intrinsic :: iso_c_binding
128 import c_rp
129 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d
130 real(c_rp) :: c
131 integer(c_int) :: gdim, n
132 end subroutine cuda_opadd2cm
133 end interface
134
135 interface
136 subroutine cuda_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, &
137 c_d, gdim, n) &
138 bind(c, name = 'cuda_opadd2col')
139 use, intrinsic :: iso_c_binding
140 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
141 integer(c_int) :: gdim, n
142 end subroutine cuda_opadd2col
143 end interface
144#elif HAVE_OPENCL
145 interface
146 subroutine opencl_opchsign(a1_d, a2_d, a3_d, gdim, n) &
147 bind(c, name = 'opencl_opchsign')
148 use, intrinsic :: iso_c_binding
149 type(c_ptr), value :: a1_d, a2_d, a3_d
150 integer(c_int) :: gdim, n
151 end subroutine opencl_opchsign
152 end interface
153
154 interface
155 subroutine opencl_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n) &
156 bind(c, name = 'opencl_opcolv')
157 use, intrinsic :: iso_c_binding
158 type(c_ptr), value :: a1_d, a2_d, a3_d, c_d
159 integer(c_int) :: gdim, n
160 end subroutine opencl_opcolv
161 end interface
162
163 interface
164 subroutine opencl_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, &
165 c_d, d, gdim, n) &
166 bind(c, name = 'opencl_opcolv3c')
167 use, intrinsic :: iso_c_binding
168 import c_rp
169 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
170 real(c_rp) :: d
171 integer(c_int) :: gdim, n
172 end subroutine opencl_opcolv3c
173 end interface
174
175 interface
176 subroutine opencl_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, &
177 c, gdim, n) &
178 bind(c, name = 'opencl_opadd2cm')
179 use, intrinsic :: iso_c_binding
180 import c_rp
181 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d
182 real(c_rp) :: c
183 integer(c_int) :: gdim, n
184 end subroutine opencl_opadd2cm
185 end interface
186
187 interface
188 subroutine opencl_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, &
189 c_d, gdim, n) &
190 bind(c, name = 'opencl_opadd2col')
191 use, intrinsic :: iso_c_binding
192 type(c_ptr), value :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
193 integer(c_int) :: gdim, n
194 end subroutine opencl_opadd2col
195 end interface
196#endif
197
200
201contains
202
204 subroutine device_opchsign(a1_d, a2_d, a3_d, gdim, n)
205 type(c_ptr) :: a1_d, a2_d, a3_d
206 integer :: n, gdim
207#ifdef HAVE_HIP
208 call hip_opchsign(a1_d, a2_d, a3_d, gdim, n)
209#elif HAVE_CUDA
210 call cuda_opchsign(a1_d, a2_d, a3_d, gdim, n)
211#elif HAVE_OPENCL
212 call opencl_opchsign(a1_d, a2_d, a3_d, gdim, n)
213#else
214 call neko_error('No device backend configured')
215#endif
216 end subroutine device_opchsign
217
219 subroutine device_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n)
220 type(c_ptr) :: a1_d, a2_d, a3_d, c_d
221 integer :: n, gdim
222#ifdef HAVE_HIP
223 call hip_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n)
224#elif HAVE_CUDA
225 call cuda_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n)
226#elif HAVE_OPENCL
227 call opencl_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n)
228#else
229 call neko_error('No device backend configured')
230#endif
231 end subroutine device_opcolv
232
234 subroutine device_opcolv3c(a1_d, a2_d, a3_d, &
235 b1_d, b2_d, b3_d, c_d, d, n, gdim)
236 type(c_ptr) :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
237 real(kind=rp) :: d
238 integer :: n, gdim
239#ifdef HAVE_HIP
240 call hip_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, gdim, n)
241#elif HAVE_CUDA
242 call cuda_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, gdim, n)
243#elif HAVE_OPENCL
244 call opencl_opcolv3c(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, d, gdim, n)
245#else
246 call neko_error('No device backend configured')
247#endif
248 end subroutine device_opcolv3c
249
251 subroutine device_opadd2cm (a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, n, gdim)
252 type(c_ptr) :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d
253 real(kind=rp) :: c
254 integer :: n, gdim
255#ifdef HAVE_HIP
256 call hip_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, gdim, n)
257#elif HAVE_CUDA
258 call cuda_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, gdim, n)
259#elif HAVE_OPENCL
260 call opencl_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, gdim, n)
261#else
262 call neko_error('No device backend configured')
263#endif
264 end subroutine device_opadd2cm
265
267 subroutine device_opadd2col (a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, n, gdim)
268 type(c_ptr) :: a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d
269 integer :: n, gdim
270#ifdef HAVE_HIP
271 call hip_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, gdim, n)
272#elif HAVE_CUDA
273 call cuda_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, gdim, n)
274#elif HAVE_OPENCL
275 call opencl_opadd2col(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c_d, gdim, n)
276#else
277 call neko_error('No device backend configured')
278#endif
279 end subroutine device_opadd2col
280
281end 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:103
void opencl_opadd2cm(void *a1, void *a2, void *a3, void *b1, void *b2, void *b3, real *c, int *gdim, int *n)
Definition mathops.c:136
void opencl_opcolv(void *a1, void *a2, void *a3, void *c, int *gdim, int *n)
Definition mathops.c:76
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:168
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