Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
rhs_maker_device.F90
Go to the documentation of this file.
1! Copyright (c) 2022, 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!
36 use device, only : device_get_ptr
37 use device_math, only : device_copy
39 use field, only : field_t
40 use num_types, only : rp, c_rp
41 use, intrinsic :: iso_c_binding, only : c_ptr
42 implicit none
43 private
44
46 contains
47 procedure, nopass :: compute_fluid => rhs_maker_sumab_device
49
50 type, public, extends(rhs_maker_ext_t) :: rhs_maker_ext_device_t
51 contains
52 procedure, nopass :: compute_fluid => rhs_maker_ext_device
53 procedure, nopass :: compute_scalar => scalar_rhs_maker_ext_device
55
56 type, public, extends(rhs_maker_bdf_t) :: rhs_maker_bdf_device_t
57 contains
58 procedure, nopass :: compute_fluid => rhs_maker_bdf_device
59 procedure, nopass :: compute_scalar => scalar_rhs_maker_bdf_device
61
62 type, public, extends(rhs_maker_oifs_t) :: rhs_maker_oifs_device_t
63 contains
64 procedure, nopass :: compute_fluid => rhs_maker_oifs_device
65 procedure, nopass :: compute_scalar => scalar_rhs_maker_oifs_device
67
68#ifdef HAVE_HIP
69 interface
70 subroutine rhs_maker_sumab_hip(u_d, v_d, w_d, uu_d, vv_d, ww_d, &
71 uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2, ab1, ab2, ab3, nab, n)&
72 bind(c, name = 'rhs_maker_sumab_hip')
73 use, intrinsic :: iso_c_binding
74 import c_rp
75 type(c_ptr), value :: u_d, v_d, w_d, uu_d, vv_d, ww_d
76 type(c_ptr), value :: uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2
77 real(c_rp) :: ab1, ab2, ab3
78 integer(c_int) :: nab, n
79 end subroutine rhs_maker_sumab_hip
80 end interface
81
82 interface
83 subroutine rhs_maker_ext_hip(abx1_d, aby1_d, abz1_d, &
84 abx2_d, aby2_d, abz2_d, &
85 bfx_d, bfy_d, bfz_d, &
86 rho, ab1, ab2, ab3, n) &
87 bind(c, name = 'rhs_maker_ext_hip')
88 use, intrinsic :: iso_c_binding
89 import c_rp
90 type(c_ptr), value :: abx1_d, aby1_d, abz1_d
91 type(c_ptr), value :: abx2_d, aby2_d, abz2_d
92 type(c_ptr), value :: bfx_d, bfy_d, bfz_d
93 real(c_rp) :: rho, ab1, ab2, ab3
94 integer(c_int) :: n
95 end subroutine rhs_maker_ext_hip
96 end interface
97
98 interface
99 subroutine scalar_rhs_maker_ext_hip(fs_lag_d, fs_laglag_d, fs_d, rho, &
100 ext1, ext2, ext3, n) &
101 bind(c, name = 'scalar_rhs_maker_ext_hip')
102 use, intrinsic :: iso_c_binding
103 import c_rp
104 type(c_ptr), value :: fs_lag_d, fs_laglag_d, fs_d
105 real(c_rp) :: rho, ext1, ext2, ext3
106 integer(c_int) :: n
107 end subroutine scalar_rhs_maker_ext_hip
108 end interface
109
110 interface
111 subroutine rhs_maker_bdf_hip(ulag1_d, ulag2_d, vlag1_d, vlag2_d, &
112 wlag1_d, wlag2_d, bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d, &
113 rho, dt, bd2, bd3, bd4, nbd, n) bind(c, name = 'rhs_maker_bdf_hip')
114 use, intrinsic :: iso_c_binding
115 import c_rp
116 type(c_ptr), value :: ulag1_d, ulag2_d, vlag1_d
117 type(c_ptr), value :: vlag2_d, wlag1_d, wlag2_d
118 type(c_ptr), value :: bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d
119 reaL(c_rp) :: rho, dt, bd2, bd3, bd4
120 integer(c_int) :: nbd, n
121 end subroutine rhs_maker_bdf_hip
122 end interface
123
124 interface
125 subroutine scalar_rhs_maker_bdf_hip(s_lag_d, s_laglag_d, fs_d, s_d, B_d, &
126 rho, dt, bd2, bd3, bd4, nbd, n) &
127 bind(c, name = 'scalar_rhs_maker_bdf_hip')
128 use, intrinsic :: iso_c_binding
129 import c_rp
130 type(c_ptr), value :: s_lag_d, s_laglag_d
131 type(c_ptr), value :: fs_d, s_d, B_d
132 reaL(c_rp) :: rho, dt, bd2, bd3, bd4
133 integer(c_int) :: nbd, n
134 end subroutine scalar_rhs_maker_bdf_hip
135 end interface
136
137 interface
138 subroutine rhs_maker_oifs_hip(phi_x_d, phi_y_d, phi_z_d, bf_x_d, &
139 bf_y_d, bf_z_d, rho, dt, n) bind(c, name = 'rhs_maker_oifs_hip')
140 use, intrinsic :: iso_c_binding
141 import c_rp
142 type(c_ptr), value :: bf_x_d, bf_y_d, bf_z_d
143 type(c_ptr), value :: phi_x_d, phi_y_d, phi_z_d
144 reaL(c_rp) :: rho, dt
145 integer(c_int) :: n
146 end subroutine rhs_maker_oifs_hip
147 end interface
148
149 interface
150 subroutine scalar_rhs_maker_oifs_hip(phi_s_d, bf_s_d, rho, dt, n) &
151 bind(c, name = 'scalar_rhs_maker_oifs_hip')
152 use, intrinsic :: iso_c_binding
153 import c_rp
154 type(c_ptr), value :: bf_s_d
155 type(c_ptr), value :: phi_s_d
156 reaL(c_rp) :: rho, dt
157 integer(c_int) :: n
158 end subroutine scalar_rhs_maker_oifs_hip
159 end interface
160#elif HAVE_CUDA
161 interface
162 subroutine rhs_maker_sumab_cuda(u_d, v_d, w_d, uu_d, vv_d, ww_d, &
163 uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2, ab1, ab2, ab3, nab, n)&
164 bind(c, name = 'rhs_maker_sumab_cuda')
165 use, intrinsic :: iso_c_binding
166 import c_rp
167 type(c_ptr), value :: u_d, v_d, w_d, uu_d, vv_d, ww_d
168 type(c_ptr), value :: uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2
169 real(c_rp) :: ab1, ab2, ab3
170 integer(c_int) :: nab, n
171 end subroutine rhs_maker_sumab_cuda
172 end interface
173
174 interface
175 subroutine rhs_maker_ext_cuda(abx1_d, aby1_d, abz1_d, &
176 abx2_d, aby2_d, abz2_d, &
177 bfx_d, bfy_d, bfz_d, &
178 rho, ab1, ab2, ab3, n) &
179 bind(c, name = 'rhs_maker_ext_cuda')
180 use, intrinsic :: iso_c_binding
181 import c_rp
182 type(c_ptr), value :: abx1_d, aby1_d, abz1_d
183 type(c_ptr), value :: abx2_d, aby2_d, abz2_d
184 type(c_ptr), value :: bfx_d, bfy_d, bfz_d
185 real(c_rp) :: rho, ab1, ab2, ab3
186 integer(c_int) :: n
187 end subroutine rhs_maker_ext_cuda
188 end interface
189
190 interface
191 subroutine scalar_rhs_maker_ext_cuda(fs_lag_d, fs_laglag_d, fs_d, rho, &
192 ext1, ext2, ext3, n) &
193 bind(c, name = 'scalar_rhs_maker_ext_cuda')
194 use, intrinsic :: iso_c_binding
195 import c_rp
196 type(c_ptr), value :: fs_lag_d, fs_laglag_d, fs_d
197 real(c_rp) :: rho, ext1, ext2, ext3
198 integer(c_int) :: n
199 end subroutine scalar_rhs_maker_ext_cuda
200 end interface
201
202 interface
203 subroutine rhs_maker_bdf_cuda(ulag1_d, ulag2_d, vlag1_d, vlag2_d, &
204 wlag1_d, wlag2_d, bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d, &
205 rho, dt, bd2, bd3, bd4, nbd, n) bind(c, name = 'rhs_maker_bdf_cuda')
206 use, intrinsic :: iso_c_binding
207 import c_rp
208 type(c_ptr), value :: ulag1_d, ulag2_d, vlag1_d
209 type(c_ptr), value :: vlag2_d, wlag1_d, wlag2_d
210 type(c_ptr), value :: bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d
211 reaL(c_rp) :: rho, dt, bd2, bd3, bd4
212 integer(c_int) :: nbd, n
213 end subroutine rhs_maker_bdf_cuda
214 end interface
215
216 interface
217 subroutine scalar_rhs_maker_bdf_cuda(s_lag_d, s_laglag_d, fs_d, s_d, B_d, &
218 rho, dt, bd2, bd3, bd4, nbd, n) &
219 bind(c, name = 'scalar_rhs_maker_bdf_cuda')
220 use, intrinsic :: iso_c_binding
221 import c_rp
222 type(c_ptr), value :: s_lag_d, s_laglag_d
223 type(c_ptr), value :: fs_d, s_d, B_d
224 reaL(c_rp) :: rho, dt, bd2, bd3, bd4
225 integer(c_int) :: nbd, n
226 end subroutine scalar_rhs_maker_bdf_cuda
227 end interface
228
229 interface
230 subroutine rhs_maker_oifs_cuda(phi_x_d, phi_y_d, phi_z_d, bf_x_d, &
231 bf_y_d, bf_z_d, rho, dt, n) bind(c, name = 'rhs_maker_oifs_cuda')
232 use, intrinsic :: iso_c_binding
233 import c_rp
234 type(c_ptr), value :: bf_x_d, bf_y_d, bf_z_d
235 type(c_ptr), value :: phi_x_d, phi_y_d, phi_z_d
236 reaL(c_rp) :: rho, dt
237 integer(c_int) :: n
238 end subroutine rhs_maker_oifs_cuda
239 end interface
240
241 interface
242 subroutine scalar_rhs_maker_oifs_cuda(phi_s_d, bf_s_d, rho, dt, n) &
243 bind(c, name = 'scalar_rhs_maker_oifs_cuda')
244 use, intrinsic :: iso_c_binding
245 import c_rp
246 type(c_ptr), value :: bf_s_d
247 type(c_ptr), value :: phi_s_d
248 reaL(c_rp) :: rho, dt
249 integer(c_int) :: n
250 end subroutine scalar_rhs_maker_oifs_cuda
251 end interface
252#elif HAVE_OPENCL
253 interface
254 subroutine rhs_maker_sumab_opencl(u_d, v_d, w_d, uu_d, vv_d, ww_d, &
255 uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2, ab1, ab2, ab3, nab, n)&
256 bind(c, name = 'rhs_maker_sumab_opencl')
257 use, intrinsic :: iso_c_binding
258 import c_rp
259 type(c_ptr), value :: u_d, v_d, w_d, uu_d, vv_d, ww_d
260 type(c_ptr), value :: uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2
261 real(c_rp) :: ab1, ab2, ab3
262 integer(c_int) :: nab, n
263 end subroutine rhs_maker_sumab_opencl
264 end interface
265
266 interface
267 subroutine rhs_maker_ext_opencl(abx1_d, aby1_d, abz1_d, &
268 abx2_d, aby2_d, abz2_d, &
269 bfx_d, bfy_d, bfz_d, &
270 rho, ab1, ab2, ab3, n) &
271 bind(c, name = 'rhs_maker_ext_opencl')
272 use, intrinsic :: iso_c_binding
273 import c_rp
274 type(c_ptr), value :: abx1_d, aby1_d, abz1_d
275 type(c_ptr), value :: abx2_d, aby2_d, abz2_d
276 type(c_ptr), value :: bfx_d, bfy_d, bfz_d
277 real(c_rp) :: rho, ab1, ab2, ab3
278 integer(c_int) :: n
279 end subroutine rhs_maker_ext_opencl
280 end interface
281
282 interface
283 subroutine scalar_rhs_maker_ext_opencl(fs_lag_d, fs_laglag_d, fs_d, rho, &
284 ext1, ext2, ext3, n) &
285 bind(c, name = 'scalar_rhs_maker_ext_opencl')
286 use, intrinsic :: iso_c_binding
287 import c_rp
288 type(c_ptr), value :: fs_lag_d, fs_laglag_d, fs_d
289 real(c_rp) :: rho, ext1, ext2, ext3
290 integer(c_int) :: n
291 end subroutine scalar_rhs_maker_ext_opencl
292 end interface
293
294 interface
295 subroutine rhs_maker_bdf_opencl(ulag1_d, ulag2_d, vlag1_d, vlag2_d, &
296 wlag1_d, wlag2_d, bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d, &
297 rho, dt, bd2, bd3, bd4, nbd, n) bind(c, name = 'rhs_maker_bdf_opencl')
298 use, intrinsic :: iso_c_binding
299 import c_rp
300 type(c_ptr), value :: ulag1_d, ulag2_d, vlag1_d
301 type(c_ptr), value :: vlag2_d, wlag1_d, wlag2_d
302 type(c_ptr), value :: bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d
303 reaL(c_rp) :: rho, dt, bd2, bd3, bd4
304 integer(c_int) :: nbd, n
305 end subroutine rhs_maker_bdf_opencl
306 end interface
307
308 interface
309 subroutine scalar_rhs_maker_bdf_opencl(s_lag_d, s_laglag_d, fs_d, s_d, B_d, &
310 rho, dt, bd2, bd3, bd4, nbd, n) &
311 bind(c, name = 'scalar_rhs_maker_bdf_opencl')
312 use, intrinsic :: iso_c_binding
313 import c_rp
314 type(c_ptr), value :: s_lag_d, s_laglag_d
315 type(c_ptr), value :: fs_d, s_d, B_d
316 reaL(c_rp) :: rho, dt, bd2, bd3, bd4
317 integer(c_int) :: nbd, n
318 end subroutine scalar_rhs_maker_bdf_opencl
319 end interface
320
321 interface
322 subroutine rhs_maker_oifs_opencl(phi_x_d, phi_y_d, phi_z_d, bf_x_d, &
323 bf_y_d, bf_z_d, rho, dt, n) bind(c, name = 'rhs_maker_oifs_opencl')
324 use, intrinsic :: iso_c_binding
325 import c_rp
326 type(c_ptr), value :: bf_x_d, bf_y_d, bf_z_d
327 type(c_ptr), value :: phi_x_d, phi_y_d, phi_z_d
328 reaL(c_rp) :: rho, dt
329 integer(c_int) :: n
330 end subroutine rhs_maker_oifs_opencl
331 end interface
332
333 interface
334 subroutine scalar_rhs_maker_oifs_opencl(phi_s_d, bf_s_d, rho, dt, n) &
335 bind(c, name = 'scalar_rhs_maker_oifs_opencl')
336 use, intrinsic :: iso_c_binding
337 import c_rp
338 type(c_ptr), value :: bf_s_d
339 type(c_ptr), value :: phi_s_d
340 reaL(c_rp) :: rho, dt
341 integer(c_int) :: n
342 end subroutine scalar_rhs_maker_oifs_opencl
343 end interface
344#elif HAVE_METAL
345 interface
346 subroutine rhs_maker_sumab_metal(u_d, v_d, w_d, uu_d, vv_d, ww_d, &
347 uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2, ab1, ab2, ab3, nab, n)&
348 bind(c, name = 'rhs_maker_sumab_metal')
349 use, intrinsic :: iso_c_binding
350 import c_rp
351 type(c_ptr), value :: u_d, v_d, w_d, uu_d, vv_d, ww_d
352 type(c_ptr), value :: uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2
353 real(c_rp) :: ab1, ab2, ab3
354 integer(c_int) :: nab, n
355 end subroutine rhs_maker_sumab_metal
356 end interface
357
358 interface
359 subroutine rhs_maker_ext_metal(abx1_d, aby1_d, abz1_d, &
360 abx2_d, aby2_d, abz2_d, &
361 bfx_d, bfy_d, bfz_d, &
362 rho, ab1, ab2, ab3, n) &
363 bind(c, name = 'rhs_maker_ext_metal')
364 use, intrinsic :: iso_c_binding
365 import c_rp
366 type(c_ptr), value :: abx1_d, aby1_d, abz1_d
367 type(c_ptr), value :: abx2_d, aby2_d, abz2_d
368 type(c_ptr), value :: bfx_d, bfy_d, bfz_d
369 real(c_rp) :: rho, ab1, ab2, ab3
370 integer(c_int) :: n
371 end subroutine rhs_maker_ext_metal
372 end interface
373
374 interface
375 subroutine scalar_rhs_maker_ext_metal(fs_lag_d, fs_laglag_d, fs_d, rho, &
376 ext1, ext2, ext3, n) &
377 bind(c, name = 'scalar_rhs_maker_ext_metal')
378 use, intrinsic :: iso_c_binding
379 import c_rp
380 type(c_ptr), value :: fs_lag_d, fs_laglag_d, fs_d
381 real(c_rp) :: rho, ext1, ext2, ext3
382 integer(c_int) :: n
383 end subroutine scalar_rhs_maker_ext_metal
384 end interface
385
386 interface
387 subroutine rhs_maker_bdf_metal(ulag1_d, ulag2_d, vlag1_d, vlag2_d, &
388 wlag1_d, wlag2_d, bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d, &
389 rho, dt, bd2, bd3, bd4, nbd, n) bind(c, name = 'rhs_maker_bdf_metal')
390 use, intrinsic :: iso_c_binding
391 import c_rp
392 type(c_ptr), value :: ulag1_d, ulag2_d, vlag1_d
393 type(c_ptr), value :: vlag2_d, wlag1_d, wlag2_d
394 type(c_ptr), value :: bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d
395 reaL(c_rp) :: rho, dt, bd2, bd3, bd4
396 integer(c_int) :: nbd, n
397 end subroutine rhs_maker_bdf_metal
398 end interface
399
400 interface
401 subroutine scalar_rhs_maker_bdf_metal(s_lag_d, s_laglag_d, fs_d, s_d, B_d, &
402 rho, dt, bd2, bd3, bd4, nbd, n) &
403 bind(c, name = 'scalar_rhs_maker_bdf_metal')
404 use, intrinsic :: iso_c_binding
405 import c_rp
406 type(c_ptr), value :: s_lag_d, s_laglag_d
407 type(c_ptr), value :: fs_d, s_d, B_d
408 reaL(c_rp) :: rho, dt, bd2, bd3, bd4
409 integer(c_int) :: nbd, n
410 end subroutine scalar_rhs_maker_bdf_metal
411 end interface
412
413 interface
414 subroutine rhs_maker_oifs_metal(phi_x_d, phi_y_d, phi_z_d, bf_x_d, &
415 bf_y_d, bf_z_d, rho, dt, n) bind(c, name = 'rhs_maker_oifs_metal')
416 use, intrinsic :: iso_c_binding
417 import c_rp
418 type(c_ptr), value :: bf_x_d, bf_y_d, bf_z_d
419 type(c_ptr), value :: phi_x_d, phi_y_d, phi_z_d
420 reaL(c_rp) :: rho, dt
421 integer(c_int) :: n
422 end subroutine rhs_maker_oifs_metal
423 end interface
424
425 interface
426 subroutine scalar_rhs_maker_oifs_metal(phi_s_d, bf_s_d, rho, dt, n) &
427 bind(c, name = 'scalar_rhs_maker_oifs_metal')
428 use, intrinsic :: iso_c_binding
429 import c_rp
430 type(c_ptr), value :: bf_s_d
431 type(c_ptr), value :: phi_s_d
432 reaL(c_rp) :: rho, dt
433 integer(c_int) :: n
434 end subroutine scalar_rhs_maker_oifs_metal
435 end interface
436#endif
437
438contains
439
440 subroutine rhs_maker_sumab_device(u, v, w, uu, vv, ww, uulag, vvlag, wwlag, ab, nab)
441 type(field_t), intent(inout) :: u,v, w
442 type(field_t), intent(inout) :: uu, vv, ww
443 type(field_series_t), intent(inout) :: uulag, vvlag, wwlag
444 real(kind=rp), dimension(3), intent(in) :: ab
445 integer, intent(in) :: nab
446
447#ifdef HAVE_HIP
448 call rhs_maker_sumab_hip(u%x_d, v%x_d, w%x_d, uu%x_d, vv%x_d, ww%x_d, &
449 uulag%lf(1)%x_d, uulag%lf(2)%x_d, vvlag%lf(1)%x_d, vvlag%lf(2)%x_d, &
450 wwlag%lf(1)%x_d, wwlag%lf(2)%x_d, ab(1), ab(2), ab(3), nab, &
451 uu%dof%size())
452#elif HAVE_CUDA
453 call rhs_maker_sumab_cuda(u%x_d, v%x_d, w%x_d, uu%x_d, vv%x_d, ww%x_d, &
454 uulag%lf(1)%x_d, uulag%lf(2)%x_d, vvlag%lf(1)%x_d, vvlag%lf(2)%x_d, &
455 wwlag%lf(1)%x_d, wwlag%lf(2)%x_d, ab(1), ab(2), ab(3), nab, &
456 uu%dof%size())
457#elif HAVE_OPENCL
458 call rhs_maker_sumab_opencl(u%x_d, v%x_d, w%x_d, uu%x_d, vv%x_d, ww%x_d, &
459 uulag%lf(1)%x_d, uulag%lf(2)%x_d, vvlag%lf(1)%x_d, vvlag%lf(2)%x_d, &
460 wwlag%lf(1)%x_d, wwlag%lf(2)%x_d, ab(1), ab(2), ab(3), nab, &
461 uu%dof%size())
462#elif HAVE_METAL
463 call rhs_maker_sumab_metal(u%x_d, v%x_d, w%x_d, uu%x_d, vv%x_d, ww%x_d, &
464 uulag%lf(1)%x_d, uulag%lf(2)%x_d, vvlag%lf(1)%x_d, vvlag%lf(2)%x_d, &
465 wwlag%lf(1)%x_d, wwlag%lf(2)%x_d, ab(1), ab(2), ab(3), nab, &
466 uu%dof%size())
467#endif
468
469 end subroutine rhs_maker_sumab_device
470
471 subroutine rhs_maker_ext_device(fx_lag, fy_lag, fz_lag, &
472 fx_laglag, fy_laglag, fz_laglag, fx, fy, fz, &
473 rho, ext_coeffs, n)
474 type(field_t), intent(inout) :: fx_lag, fy_lag, fz_lag
475 type(field_t), intent(inout) :: fx_laglag, fy_laglag, fz_laglag
476 real(kind=rp), intent(in) :: rho, ext_coeffs(4)
477 integer, intent(in) :: n
478 real(kind=rp), intent(inout) :: fx(n), fy(n), fz(n)
479 type(c_ptr) :: fx_d, fy_d, fz_d
480
481 fx_d = device_get_ptr(fx)
482 fy_d = device_get_ptr(fy)
483 fz_d = device_get_ptr(fz)
484
485#ifdef HAVE_HIP
486 call rhs_maker_ext_hip(fx_lag%x_d, fy_lag%x_d, fz_lag%x_d, &
487 fx_laglag%x_d, fy_laglag%x_d, fz_laglag%x_d, &
488 fx_d, fy_d, fz_d, rho, &
489 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
490#elif HAVE_CUDA
491 call rhs_maker_ext_cuda(fx_lag%x_d, fy_lag%x_d, fz_lag%x_d, &
492 fx_laglag%x_d, fy_laglag%x_d, fz_laglag%x_d, &
493 fx_d, fy_d, fz_d, rho, &
494 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
495#elif HAVE_OPENCL
496 call rhs_maker_ext_opencl(fx_lag%x_d, fy_lag%x_d, fz_lag%x_d, &
497 fx_laglag%x_d, fy_laglag%x_d, fz_laglag%x_d, &
498 fx_d, fy_d, fz_d, rho, &
499 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
500#elif HAVE_METAL
501 call rhs_maker_ext_metal(fx_lag%x_d, fy_lag%x_d, fz_lag%x_d, &
502 fx_laglag%x_d, fy_laglag%x_d, fz_laglag%x_d, &
503 fx_d, fy_d, fz_d, rho, &
504 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
505#endif
506
507 end subroutine rhs_maker_ext_device
508
509 subroutine scalar_rhs_maker_ext_device(fs_lag, fs_laglag, fs, &
510 rho, ext_coeffs, n)
511 type(field_t), intent(inout) :: fs_lag
512 type(field_t), intent(inout) :: fs_laglag
513 real(kind=rp), intent(in) :: rho, ext_coeffs(4)
514 integer, intent(in) :: n
515 real(kind=rp), intent(inout) :: fs(n)
516 type(c_ptr) :: fs_d
517
518 fs_d = device_get_ptr(fs)
519
520#ifdef HAVE_HIP
521 call scalar_rhs_maker_ext_hip(fs_lag%x_d, fs_laglag%x_d, fs_d, rho, &
522 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
523#elif HAVE_CUDA
524 call scalar_rhs_maker_ext_cuda(fs_lag%x_d, fs_laglag%x_d, fs_d, rho, &
525 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
526#elif HAVE_OPENCL
527 call scalar_rhs_maker_ext_opencl(fs_lag%x_d, fs_laglag%x_d, fs_d, rho, &
528 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
529#elif HAVE_METAL
530 call scalar_rhs_maker_ext_metal(fs_lag%x_d, fs_laglag%x_d, fs_d, rho, &
531 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
532#endif
533
534 end subroutine scalar_rhs_maker_ext_device
535
536 subroutine rhs_maker_bdf_device(ulag, vlag, wlag, bfx, bfy, bfz, &
537 u, v, w, B, rho, dt, bd, nbd, n, Blag, Blaglag)
538 integer, intent(in) :: n, nbd
539 type(field_t), intent(in) :: u, v, w
540 type(field_series_t), intent(in) :: ulag, vlag, wlag
541 real(kind=rp), intent(in) :: blag(n), blaglag(n)
542 real(kind=rp), intent(inout) :: bfx(n), bfy(n), bfz(n)
543 real(kind=rp), intent(in) :: b(n)
544 real(kind=rp), intent(in) :: dt, rho, bd(4)
545 type(c_ptr) :: bfx_d, bfy_d, bfz_d, b_d
546
547 bfx_d = device_get_ptr(bfx)
548 bfy_d = device_get_ptr(bfy)
549 bfz_d = device_get_ptr(bfz)
550 b_d = device_get_ptr(b)
551
552#ifdef HAVE_HIP
553 call rhs_maker_bdf_hip(ulag%lf(1)%x_d, ulag%lf(2)%x_d, &
554 vlag%lf(1)%x_d, vlag%lf(2)%x_d, &
555 wlag%lf(1)%x_d, wlag%lf(2)%x_d, &
556 bfx_d, bfy_d, bfz_d, u%x_d, v%x_d, w%x_d, &
557 b_d, rho, dt, bd(2), bd(3), bd(4), nbd, n)
558#elif HAVE_CUDA
559 call rhs_maker_bdf_cuda(ulag%lf(1)%x_d, ulag%lf(2)%x_d, &
560 vlag%lf(1)%x_d, vlag%lf(2)%x_d, &
561 wlag%lf(1)%x_d, wlag%lf(2)%x_d, &
562 bfx_d, bfy_d, bfz_d, u%x_d, v%x_d, w%x_d, &
563 b_d, rho, dt, bd(2), bd(3), bd(4), nbd, n)
564#elif HAVE_OPENCL
565 call rhs_maker_bdf_opencl(ulag%lf(1)%x_d, ulag%lf(2)%x_d, &
566 vlag%lf(1)%x_d, vlag%lf(2)%x_d, &
567 wlag%lf(1)%x_d, wlag%lf(2)%x_d, &
568 bfx_d, bfy_d, bfz_d, u%x_d, v%x_d, w%x_d, &
569 b_d, rho, dt, bd(2), bd(3), bd(4), nbd, n)
570#elif HAVE_METAL
571 call rhs_maker_bdf_metal(ulag%lf(1)%x_d, ulag%lf(2)%x_d, &
572 vlag%lf(1)%x_d, vlag%lf(2)%x_d, &
573 wlag%lf(1)%x_d, wlag%lf(2)%x_d, &
574 bfx_d, bfy_d, bfz_d, u%x_d, v%x_d, w%x_d, &
575 b_d, rho, dt, bd(2), bd(3), bd(4), nbd, n)
576#endif
577
578 end subroutine rhs_maker_bdf_device
579
580 subroutine scalar_rhs_maker_bdf_device(s_lag, fs, s, B, rho, dt, &
581 bd, nbd, n)
582 integer, intent(in) :: n, nbd
583 type(field_t), intent(in) :: s
584 type(field_series_t), intent(in) :: s_lag
585 real(kind=rp), intent(inout) :: fs(n)
586 real(kind=rp), intent(in) :: b(n)
587 real(kind=rp), intent(in) :: dt, rho, bd(4)
588 type(c_ptr) :: fs_d, b_d
589
590 fs_d = device_get_ptr(fs)
591 b_d = device_get_ptr(b)
592
593#ifdef HAVE_HIP
594 call scalar_rhs_maker_bdf_hip(s_lag%lf(1)%x_d, s_lag%lf(2)%x_d, &
595 fs_d, s%x_d, b_d, rho, dt, bd(2), bd(3), bd(4), &
596 nbd, n)
597#elif HAVE_CUDA
598 call scalar_rhs_maker_bdf_cuda(s_lag%lf(1)%x_d, s_lag%lf(2)%x_d, &
599 fs_d, s%x_d, b_d, rho, dt, bd(2), bd(3), bd(4), &
600 nbd, n)
601#elif HAVE_OPENCL
602 call scalar_rhs_maker_bdf_opencl(s_lag%lf(1)%x_d, s_lag%lf(2)%x_d, &
603 fs_d, s%x_d, b_d, rho, dt, bd(2), bd(3), bd(4), &
604 nbd, n)
605#elif HAVE_METAL
606 call scalar_rhs_maker_bdf_metal(s_lag%lf(1)%x_d, s_lag%lf(2)%x_d, &
607 fs_d, s%x_d, b_d, rho, dt, bd(2), bd(3), bd(4), &
608 nbd, n)
609#endif
610
611 end subroutine scalar_rhs_maker_bdf_device
612
613 subroutine rhs_maker_oifs_device(phi_x, phi_y, phi_z, bf_x, bf_y, bf_z, &
614 rho, dt, n)
615 real(kind=rp), intent(in) :: rho, dt
616 integer, intent(in) :: n
617 real(kind=rp), intent(inout) :: bf_x(n), bf_y(n), bf_z(n)
618 real(kind=rp), intent(inout) :: phi_x(n), phi_y(n), phi_z(n)
619
620 type(c_ptr) :: phi_x_d, phi_y_d, phi_z_d, bf_x_d, bf_y_d, bf_z_d
621
622 phi_x_d = device_get_ptr(phi_x)
623 phi_y_d = device_get_ptr(phi_y)
624 phi_z_d = device_get_ptr(phi_z)
625 bf_x_d = device_get_ptr(bf_x)
626 bf_y_d = device_get_ptr(bf_y)
627 bf_z_d = device_get_ptr(bf_z)
628
629#ifdef HAVE_HIP
630 call rhs_maker_oifs_hip(phi_x_d, phi_y_d, phi_z_d, &
631 bf_x_d, bf_y_d, bf_z_d, rho, dt, n)
632#elif HAVE_CUDA
633 call rhs_maker_oifs_cuda(phi_x_d, phi_y_d, phi_z_d, &
634 bf_x_d, bf_y_d, bf_z_d, rho, dt, n)
635#elif HAVE_OPENCL
636 call rhs_maker_oifs_opencl(phi_x_d, phi_y_d, phi_z_d, &
637 bf_x_d, bf_y_d, bf_z_d, rho, dt, n)
638#elif HAVE_METAL
639 call rhs_maker_oifs_metal(phi_x_d, phi_y_d, phi_z_d, &
640 bf_x_d, bf_y_d, bf_z_d, rho, dt, n)
641#endif
642
643 end subroutine rhs_maker_oifs_device
644
645 subroutine scalar_rhs_maker_oifs_device(phi_s,bf_s, rho, dt, n)
646 real(kind=rp), intent(in) :: rho, dt
647 integer, intent(in) :: n
648 real(kind=rp), intent(inout) :: bf_s(n)
649 real(kind=rp), intent(inout) :: phi_s(n)
650 type(c_ptr) :: phi_s_d, bf_s_d
651
652 phi_s_d = device_get_ptr(phi_s)
653 bf_s_d = device_get_ptr(bf_s)
654
655#ifdef HAVE_HIP
656 call scalar_rhs_maker_oifs_hip(phi_s_d, bf_s_d, rho, dt, n)
657#elif HAVE_CUDA
658 call scalar_rhs_maker_oifs_cuda(phi_s_d, bf_s_d, rho, dt, n)
659#elif HAVE_OPENCL
660 call scalar_rhs_maker_oifs_opencl(phi_s_d, bf_s_d, rho, dt, n)
661#elif HAVE_METAL
662 call scalar_rhs_maker_oifs_metal(phi_s_d, bf_s_d, rho, dt, n)
663#endif
664
665 end subroutine scalar_rhs_maker_oifs_device
666
667end module rhs_maker_device
Return the device pointer for an associated Fortran array.
Definition device.F90:108
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
Contains the field_serties_t type.
Defines a field.
Definition field.f90:34
integer, parameter, public c_rp
Definition num_types.f90:13
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
subroutine scalar_rhs_maker_bdf_device(s_lag, fs, s, b, rho, dt, bd, nbd, n)
subroutine rhs_maker_sumab_device(u, v, w, uu, vv, ww, uulag, vvlag, wwlag, ab, nab)
subroutine rhs_maker_oifs_device(phi_x, phi_y, phi_z, bf_x, bf_y, bf_z, rho, dt, n)
subroutine scalar_rhs_maker_ext_device(fs_lag, fs_laglag, fs, rho, ext_coeffs, n)
subroutine rhs_maker_bdf_device(ulag, vlag, wlag, bfx, bfy, bfz, u, v, w, b, rho, dt, bd, nbd, n, blag, blaglag)
subroutine rhs_maker_ext_device(fx_lag, fy_lag, fz_lag, fx_laglag, fy_laglag, fz_laglag, fx, fy, fz, rho, ext_coeffs, n)
subroutine scalar_rhs_maker_oifs_device(phi_s, bf_s, rho, dt, n)
Routines to generate the right-hand sides for the convection-diffusion equation. Employs the EXT/BDF ...
Definition rhs_maker.f90:38
void scalar_rhs_maker_bdf_opencl(void *s_lag, void *s_laglag, void *fs, void *s, void *B, real *rho, real *dt, real *bd2, real *bd3, real *bd4, int *nbd, int *n)
Definition rhs_maker.c:205
void rhs_maker_sumab_opencl(void *u, void *v, void *w, void *uu, void *vv, void *ww, void *ulag1, void *ulag2, void *vlag1, void *vlag2, void *wlag1, void *wlag2, real *ext1, real *ext2, real *ext3, int *nab, int *n)
Definition rhs_maker.c:49
void scalar_rhs_maker_oifs_opencl(void *phi_s, void *bf_s, real *rho, real *dt, int *n)
Definition rhs_maker.c:273
void rhs_maker_ext_opencl(void *abx1, void *aby1, void *abz1, void *abx2, void *aby2, void *abz2, void *bfx, void *bfy, void *bfz, real *rho, real *ext1, real *ext2, real *ext3, int *n)
Definition rhs_maker.c:90
void rhs_maker_oifs_opencl(void *phi_x, void *phi_y, void *phi_z, void *bf_x, void *bf_y, void *bf_z, real *rho, real *dt, int *n)
Definition rhs_maker.c:241
void scalar_rhs_maker_ext_opencl(void *fs_lag, void *fs_laglag, void *fs, real *rho, real *ext1, real *ext2, real *ext3, int *n)
Definition rhs_maker.c:128
void rhs_maker_bdf_opencl(void *ulag1, void *ulag2, void *vlag1, void *vlag2, void *wlag1, void *wlag2, void *bfx, void *bfy, void *bfz, void *u, void *v, void *w, void *B, real *rho, real *dt, real *bd2, real *bd3, real *bd4, int *nbd, int *n)
Definition rhs_maker.c:159
void scalar_rhs_maker_oifs_cuda(void *phi_s, void *bf_s, real *rho, real *dt, int *n)
Definition rhs_maker.cu:157
void scalar_rhs_maker_bdf_cuda(void *s_lag, void *s_laglag, void *fs, void *s, void *B, real *rho, real *dt, real *bd2, real *bd3, real *bd4, int *nbd, int *n)
Definition rhs_maker.cu:120
void rhs_maker_oifs_cuda(void *phi_x, void *phi_y, void *phi_z, void *bf_x, void *bf_y, void *bf_z, real *rho, real *dt, int *n)
Definition rhs_maker.cu:140
void rhs_maker_ext_cuda(void *abx1, void *aby1, void *abz1, void *abx2, void *aby2, void *abz2, void *bfx, void *bfy, void *bfz, real *rho, real *ab1, real *ab2, real *ab3, int *n)
Definition rhs_maker.cu:63
void rhs_maker_bdf_cuda(void *ulag1, void *ulag2, void *vlag1, void *vlag2, void *wlag1, void *wlag2, void *bfx, void *bfy, void *bfz, void *u, void *v, void *w, void *B, real *rho, real *dt, real *bd2, real *bd3, real *bd4, int *nbd, int *n)
Definition rhs_maker.cu:98
void scalar_rhs_maker_ext_cuda(void *fs_lag, void *fs_laglag, void *fs, real *rho, real *ext1, real *ext2, real *ext3, int *n)
Definition rhs_maker.cu:82
void rhs_maker_sumab_cuda(void *u, void *v, void *w, void *uu, void *vv, void *ww, void *ulag1, void *ulag2, void *vlag1, void *vlag2, void *wlag1, void *wlag2, real *ab1, real *ab2, real *ab3, int *nab, int *n)
Definition rhs_maker.cu:44
Stores a series (sequence) of fields, logically connected to a base field, and arranged according to ...
Abstract type to add contributions to F from lagged BD terms.
Definition rhs_maker.f90:59
Abstract type to sum up contributions to kth order extrapolation scheme.
Definition rhs_maker.f90:52
Abstract type to add contributions of kth order OIFS scheme.
Definition rhs_maker.f90:66
Abstract type to compute extrapolated velocity field for the pressure equation.
Definition rhs_maker.f90:46