Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
device_math_ext.f90
Go to the documentation of this file.
1
34!
35module device_math_ext
36 use utils, only: neko_error
37 use num_types, only: rp, c_rp
38 use, intrinsic :: iso_c_binding, only: c_ptr, c_int
39 implicit none
40
41#if HAVE_HIP
42
43 interface
44 subroutine hip_copy_mask(a_d, b_d, size, mask_d, mask_size) &
45 bind(c, name = 'hip_copy_mask')
46 import c_rp, c_int, c_ptr
47 type(c_ptr), value :: a_d
48 type(c_ptr), value :: b_d
49 integer(c_int) :: size
50 type(c_ptr), value :: mask_d
51 integer(c_int) :: mask_size
52 end subroutine hip_copy_mask
53 end interface
54 interface
55 subroutine hip_cadd_mask(a_d, c, size, mask_d, mask_size) &
56 bind(c, name = 'hip_cadd_mask')
57 import c_rp, c_int, c_ptr
58 type(c_ptr), value :: a_d
59 real(c_rp) :: c
60 integer(c_int) :: size
61 type(c_ptr), value :: mask_d
62 integer(c_int) :: mask_size
63 end subroutine hip_cadd_mask
64 end interface
65 interface
66 subroutine hip_invcol1_mask(a_d, size, mask_d, mask_size) &
67 bind(c, name = 'hip_invcol1_mask')
68 import c_rp, c_int, c_ptr
69 type(c_ptr), value :: a_d
70 integer(c_int) :: size
71 type(c_ptr), value :: mask_d
72 integer(c_int) :: mask_size
73 end subroutine hip_invcol1_mask
74 end interface
75 interface
76 subroutine hip_col2_mask(a_d, b_d, size, mask_d, mask_size) &
77 bind(c, name = 'hip_col2_mask')
78 import c_rp, c_int, c_ptr
79 type(c_ptr), value :: a_d
80 type(c_ptr), value :: b_d
81 integer(c_int) :: size
82 type(c_ptr), value :: mask_d
83 integer(c_int) :: mask_size
84 end subroutine hip_col2_mask
85 end interface
86 interface
87 subroutine hip_col3_mask(a_d, b_d, c_d, size, mask_d, mask_size) &
88 bind(c, name = 'hip_col3_mask')
89 import c_rp, c_int, c_ptr
90 type(c_ptr), value :: a_d
91 type(c_ptr), value :: b_d
92 type(c_ptr), value :: c_d
93 integer(c_int) :: size
94 type(c_ptr), value :: mask_d
95 integer(c_int) :: mask_size
96 end subroutine hip_col3_mask
97 end interface
98 interface
99 subroutine hip_sub3_mask(a_d, b_d, c_d, size, mask_d, mask_size) &
100 bind(c, name = 'hip_sub3_mask')
101 import c_rp, c_int, c_ptr
102 type(c_ptr), value :: a_d
103 type(c_ptr), value :: b_d
104 type(c_ptr), value :: c_d
105 integer(c_int) :: size
106 type(c_ptr), value :: mask_d
107 integer(c_int) :: mask_size
108 end subroutine hip_sub3_mask
109 end interface
110
111#elif HAVE_CUDA
112
113 interface
114 subroutine cuda_copy_mask(a_d, b_d, size, mask_d, mask_size) &
115 bind(c, name = 'cuda_copy_mask')
116 import c_rp, c_int, c_ptr
117 type(c_ptr), value :: a_d
118 type(c_ptr), value :: b_d
119 integer(c_int) :: size
120 type(c_ptr), value :: mask_d
121 integer(c_int) :: mask_size
122 end subroutine cuda_copy_mask
123 end interface
124 interface
125 subroutine cuda_cadd_mask(a_d, c, size, mask_d, mask_size) &
126 bind(c, name = 'cuda_cadd_mask')
127 import c_rp, c_int, c_ptr
128 type(c_ptr), value :: a_d
129 real(c_rp) :: c
130 integer(c_int) :: size
131 type(c_ptr), value :: mask_d
132 integer(c_int) :: mask_size
133 end subroutine cuda_cadd_mask
134 end interface
135 interface
136 subroutine cuda_invcol1_mask(a_d, size, mask_d, mask_size) &
137 bind(c, name = 'cuda_invcol1_mask')
138 import c_rp, c_int, c_ptr
139 type(c_ptr), value :: a_d
140 integer(c_int) :: size
141 type(c_ptr), value :: mask_d
142 integer(c_int) :: mask_size
143 end subroutine cuda_invcol1_mask
144 end interface
145 interface
146 subroutine cuda_col2_mask(a_d, b_d, size, mask_d, mask_size) &
147 bind(c, name = 'cuda_col2_mask')
148 import c_rp, c_int, c_ptr
149 type(c_ptr), value :: a_d
150 type(c_ptr), value :: b_d
151 integer(c_int) :: size
152 type(c_ptr), value :: mask_d
153 integer(c_int) :: mask_size
154 end subroutine cuda_col2_mask
155 end interface
156 interface
157 subroutine cuda_col3_mask(a_d, b_d, c_d, size, mask_d, mask_size) &
158 bind(c, name = 'cuda_col3_mask')
159 import c_rp, c_int, c_ptr
160 type(c_ptr), value :: a_d
161 type(c_ptr), value :: b_d
162 type(c_ptr), value :: c_d
163 integer(c_int) :: size
164 type(c_ptr), value :: mask_d
165 integer(c_int) :: mask_size
166 end subroutine cuda_col3_mask
167 end interface
168 interface
169 subroutine cuda_sub3_mask(a_d, b_d, c_d, size, mask_d, mask_size) &
170 bind(c, name = 'cuda_sub3_mask')
171 import c_rp, c_int, c_ptr
172 type(c_ptr), value :: a_d
173 type(c_ptr), value :: b_d
174 type(c_ptr), value :: c_d
175 integer(c_int) :: size
176 type(c_ptr), value :: mask_d
177 integer(c_int) :: mask_size
178 end subroutine cuda_sub3_mask
179 end interface
180
181#elif HAVE_OPENCL
182
183#endif
184
185contains
186
187 subroutine device_copy_mask(a_d, b_d, size, mask_d, mask_size)
188 type(c_ptr) :: a_d
189 type(c_ptr) :: b_d
190 integer :: size
191 type(c_ptr) :: mask_d
192 integer :: mask_size
193#if HAVE_HIP
194 call hip_copy_mask(a_d, b_d, size, mask_d, mask_size)
195#elif HAVE_CUDA
196 call cuda_copy_mask(a_d, b_d, size, mask_d, mask_size)
197#else
198 call neko_error('No device backend configured for device_copy_mask')
199#endif
200 end subroutine device_copy_mask
201
202 subroutine device_cadd_mask(a_d, c, size, mask_d, mask_size)
203 type(c_ptr) :: a_d
204 real(kind=rp), intent(in) :: c
205 integer :: size
206 type(c_ptr) :: mask_d
207 integer :: mask_size
208#if HAVE_HIP
209 call hip_cadd_mask(a_d, c, size, mask_d, mask_size)
210#elif HAVE_CUDA
211 call cuda_cadd_mask(a_d, c, size, mask_d, mask_size)
212#else
213 call neko_error('No device backend configured for device_cadd_mask')
214#endif
215 end subroutine device_cadd_mask
216
217 subroutine device_invcol1_mask(a_d, size, mask_d, mask_size)
218 type(c_ptr) :: a_d
219 integer :: size
220 type(c_ptr) :: mask_d
221 integer :: mask_size
222#if HAVE_HIP
223 call hip_invcol1_mask(a_d, size, mask_d, mask_size)
224#elif HAVE_CUDA
225 call cuda_invcol1_mask(a_d, size, mask_d, mask_size)
226#else
227 call neko_error('No device backend configured for device_invcol1_mask')
228#endif
229 end subroutine device_invcol1_mask
230
231 subroutine device_col2_mask(a_d, b_d, size, mask_d, mask_size)
232 type(c_ptr) :: a_d
233 type(c_ptr) :: b_d
234 integer :: size
235 type(c_ptr) :: mask_d
236 integer :: mask_size
237#if HAVE_HIP
238 call hip_col2_mask(a_d, b_d, size, mask_d, mask_size)
239#elif HAVE_CUDA
240 call cuda_col2_mask(a_d, b_d, size, mask_d, mask_size)
241#else
242 call neko_error('No device backend configured for device_col2_mask')
243#endif
244 end subroutine device_col2_mask
245
246 subroutine device_col3_mask(a_d, b_d, c_d, size, mask_d, mask_size)
247 type(c_ptr) :: a_d
248 type(c_ptr) :: b_d
249 type(c_ptr) :: c_d
250 integer :: size
251 type(c_ptr) :: mask_d
252 integer :: mask_size
253#if HAVE_HIP
254 call hip_col3_mask(a_d, b_d, c_d, size, mask_d, mask_size)
255#elif HAVE_CUDA
256 call cuda_col3_mask(a_d, b_d, c_d, size, mask_d, mask_size)
257#else
258 call neko_error('No device backend configured for device_col3_mask')
259#endif
260 end subroutine device_col3_mask
261
262 subroutine device_sub3_mask(a_d, b_d, c_d, size, mask_d, mask_size)
263 type(c_ptr) :: a_d
264 type(c_ptr) :: b_d
265 type(c_ptr) :: c_d
266 integer :: size
267 type(c_ptr) :: mask_d
268 integer :: mask_size
269#if HAVE_HIP
270 call hip_sub3_mask(a_d, b_d, c_d, size, mask_d, mask_size)
271#elif HAVE_CUDA
272 call cuda_sub3_mask(a_d, b_d, c_d, size, mask_d, mask_size)
273#else
274 call neko_error('No device backend configured for device_sub3_mask')
275#endif
276 end subroutine device_sub3_mask
277
278end module device_math_ext