40 use,
intrinsic :: iso_c_binding, only : c_ptr
64 uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2, ab1, ab2, ab3, nab, n)&
65 bind(c, name=
'rhs_maker_sumab_hip')
66 use,
intrinsic :: iso_c_binding
68 type(c_ptr),
value :: u_d, v_d, w_d, uu_d, vv_d, ww_d
69 type(c_ptr),
value :: uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2
70 real(c_rp) :: ab1, ab2, ab3
71 integer(c_int) :: nab, n
77 abx2_d, aby2_d, abz2_d, &
78 bfx_d, bfy_d, bfz_d, &
79 rho, ab1, ab2, ab3, n) &
80 bind(c, name=
'rhs_maker_ext_hip')
81 use,
intrinsic :: iso_c_binding
83 type(c_ptr),
value :: abx1_d, aby1_d, abz1_d
84 type(c_ptr),
value :: abx2_d, aby2_d, abz2_d
85 type(c_ptr),
value :: bfx_d, bfy_d, bfz_d
86 real(c_rp) :: rho, ab1, ab2, ab3
93 ext1, ext2, ext3, n) &
94 bind(c, name=
'scalar_rhs_maker_ext_hip')
95 use,
intrinsic :: iso_c_binding
97 type(c_ptr),
value :: fs_lag_d, fs_laglag_d, fs_d
98 real(c_rp) :: rho, ext1, ext2, ext3
105 wlag1_d, wlag2_d, bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d, &
106 rho, dt, bd2, bd3, bd4, nbd, n) bind(c, name='rhs_maker_bdf_hip')
107 use,
intrinsic :: iso_c_binding
109 type(c_ptr),
value :: ulag1_d, ulag2_d, vlag1_d
110 type(c_ptr),
value :: vlag2_d, wlag1_d, wlag2_d
111 type(c_ptr),
value :: bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d
112 reaL(c_rp) :: rho, dt, bd2, bd3, bd4
113 integer(c_int) :: nbd, n
119 rho, dt, bd2, bd3, bd4, nbd, n) &
120 bind(c, name=
'scalar_rhs_maker_bdf_hip')
121 use,
intrinsic :: iso_c_binding
123 type(c_ptr),
value :: s_lag_d, s_laglag_d
124 type(c_ptr),
value :: fs_d, s_d, B_d
125 reaL(c_rp) :: rho, dt, bd2, bd3, bd4
126 integer(c_int) :: nbd, n
132 uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2, ab1, ab2, ab3, nab, n)&
133 bind(c, name=
'rhs_maker_sumab_cuda')
134 use,
intrinsic :: iso_c_binding
136 type(c_ptr),
value :: u_d, v_d, w_d, uu_d, vv_d, ww_d
137 type(c_ptr),
value :: uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2
138 real(c_rp) :: ab1, ab2, ab3
139 integer(c_int) :: nab, n
145 abx2_d, aby2_d, abz2_d, &
146 bfx_d, bfy_d, bfz_d, &
147 rho, ab1, ab2, ab3, n) &
148 bind(c, name=
'rhs_maker_ext_cuda')
149 use,
intrinsic :: iso_c_binding
151 type(c_ptr),
value :: abx1_d, aby1_d, abz1_d
152 type(c_ptr),
value :: abx2_d, aby2_d, abz2_d
153 type(c_ptr),
value :: bfx_d, bfy_d, bfz_d
154 real(c_rp) :: rho, ab1, ab2, ab3
161 ext1, ext2, ext3, n) &
162 bind(c, name=
'scalar_rhs_maker_ext_cuda')
163 use,
intrinsic :: iso_c_binding
165 type(c_ptr),
value :: fs_lag_d, fs_laglag_d, fs_d
166 real(c_rp) :: rho, ext1, ext2, ext3
173 wlag1_d, wlag2_d, bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d, &
174 rho, dt, bd2, bd3, bd4, nbd, n) bind(c, name='rhs_maker_bdf_cuda')
175 use,
intrinsic :: iso_c_binding
177 type(c_ptr),
value :: ulag1_d, ulag2_d, vlag1_d
178 type(c_ptr),
value :: vlag2_d, wlag1_d, wlag2_d
179 type(c_ptr),
value :: bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d
180 reaL(c_rp) :: rho, dt, bd2, bd3, bd4
181 integer(c_int) :: nbd, n
187 rho, dt, bd2, bd3, bd4, nbd, n) &
188 bind(c, name=
'scalar_rhs_maker_bdf_cuda')
189 use,
intrinsic :: iso_c_binding
191 type(c_ptr),
value :: s_lag_d, s_laglag_d
192 type(c_ptr),
value :: fs_d, s_d, B_d
193 reaL(c_rp) :: rho, dt, bd2, bd3, bd4
194 integer(c_int) :: nbd, n
200 uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2, ab1, ab2, ab3, nab, n)&
201 bind(c, name=
'rhs_maker_sumab_opencl')
202 use,
intrinsic :: iso_c_binding
204 type(c_ptr),
value :: u_d, v_d, w_d, uu_d, vv_d, ww_d
205 type(c_ptr),
value :: uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2
206 real(c_rp) :: ab1, ab2, ab3
207 integer(c_int) :: nab, n
213 abx2_d, aby2_d, abz2_d, &
214 bfx_d, bfy_d, bfz_d, &
215 rho, ab1, ab2, ab3, n) &
216 bind(c, name=
'rhs_maker_ext_opencl')
217 use,
intrinsic :: iso_c_binding
219 type(c_ptr),
value :: abx1_d, aby1_d, abz1_d
220 type(c_ptr),
value :: abx2_d, aby2_d, abz2_d
221 type(c_ptr),
value :: bfx_d, bfy_d, bfz_d
222 real(c_rp) :: rho, ab1, ab2, ab3
229 ext1, ext2, ext3, n) &
230 bind(c, name=
'scalar_rhs_maker_ext_opencl')
231 use,
intrinsic :: iso_c_binding
233 type(c_ptr),
value :: fs_lag_d, fs_laglag_d, fs_d
234 real(c_rp) :: rho, ext1, ext2, ext3
241 wlag1_d, wlag2_d, bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d, &
242 rho, dt, bd2, bd3, bd4, nbd, n) bind(c, name='rhs_maker_bdf_opencl')
243 use,
intrinsic :: iso_c_binding
245 type(c_ptr),
value :: ulag1_d, ulag2_d, vlag1_d
246 type(c_ptr),
value :: vlag2_d, wlag1_d, wlag2_d
247 type(c_ptr),
value :: bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d
248 reaL(c_rp) :: rho, dt, bd2, bd3, bd4
249 integer(c_int) :: nbd, n
255 rho, dt, bd2, bd3, bd4, nbd, n) &
256 bind(c, name=
'scalar_rhs_maker_bdf_opencl')
257 use,
intrinsic :: iso_c_binding
259 type(c_ptr),
value :: s_lag_d, s_laglag_d
260 type(c_ptr),
value :: fs_d, s_d, B_d
261 reaL(c_rp) :: rho, dt, bd2, bd3, bd4
262 integer(c_int) :: nbd, n
269 subroutine rhs_maker_sumab_device(u, v, w, uu, vv, ww, uulag, vvlag, wwlag, ab, nab)
270 type(
field_t),
intent(inout) :: u,v, w
271 type(
field_t),
intent(inout) :: uu, vv, ww
273 real(kind=
rp),
dimension(3),
intent(in) :: ab
274 integer,
intent(in) :: nab
278 uulag%lf(1)%x_d, uulag%lf(2)%x_d, vvlag%lf(1)%x_d, vvlag%lf(2)%x_d, &
279 wwlag%lf(1)%x_d, wwlag%lf(2)%x_d, ab(1), ab(2), ab(3), nab, &
283 uulag%lf(1)%x_d, uulag%lf(2)%x_d, vvlag%lf(1)%x_d, vvlag%lf(2)%x_d, &
284 wwlag%lf(1)%x_d, wwlag%lf(2)%x_d, ab(1), ab(2), ab(3), nab, &
288 uulag%lf(1)%x_d, uulag%lf(2)%x_d, vvlag%lf(1)%x_d, vvlag%lf(2)%x_d, &
289 wwlag%lf(1)%x_d, wwlag%lf(2)%x_d, ab(1), ab(2), ab(3), nab, &
296 fx_laglag, fy_laglag, fz_laglag, fx, fy, fz, &
298 type(
field_t),
intent(inout) :: fx_lag, fy_lag, fz_lag
299 type(
field_t),
intent(inout) :: fx_laglag, fy_laglag, fz_laglag
300 real(kind=
rp),
intent(inout) :: rho, ext_coeffs(4)
301 integer,
intent(in) :: n
302 real(kind=
rp),
intent(inout) :: fx(n), fy(n), fz(n)
303 type(c_ptr) :: fx_d, fy_d, fz_d
311 fx_laglag%x_d, fy_laglag%x_d, fz_laglag%x_d, &
312 fx_d, fy_d, fz_d, rho, &
313 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
316 fx_laglag%x_d, fy_laglag%x_d, fz_laglag%x_d, &
317 fx_d, fy_d, fz_d, rho, &
318 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
321 fx_laglag%x_d, fy_laglag%x_d, fz_laglag%x_d, &
322 fx_d, fy_d, fz_d, rho, &
323 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
330 type(
field_t),
intent(inout) :: fs_lag
331 type(
field_t),
intent(inout) :: fs_laglag
332 real(kind=
rp),
intent(inout) :: rho, ext_coeffs(4)
333 integer,
intent(in) :: n
334 real(kind=
rp),
intent(inout) :: fs(n)
341 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
344 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
347 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
353 u, v, w, B, rho, dt, bd, nbd, n)
354 integer,
intent(in) :: n, nbd
355 type(
field_t),
intent(in) :: u, v, w
357 real(kind=
rp),
intent(inout) :: bfx(n), bfy(n), bfz(n)
358 real(kind=
rp),
intent(in) :: b(n)
359 real(kind=
rp),
intent(in) :: dt, rho, bd(4)
360 type(c_ptr) :: bfx_d, bfy_d, bfz_d, b_d
369 vlag%lf(1)%x_d, vlag%lf(2)%x_d, &
370 wlag%lf(1)%x_d, wlag%lf(2)%x_d, &
371 bfx_d, bfy_d, bfz_d, u%x_d, v%x_d, w%x_d, &
372 b_d, rho, dt, bd(2), bd(3), bd(4), nbd, n)
375 vlag%lf(1)%x_d, vlag%lf(2)%x_d, &
376 wlag%lf(1)%x_d, wlag%lf(2)%x_d, &
377 bfx_d, bfy_d, bfz_d, u%x_d, v%x_d, w%x_d, &
378 b_d, rho, dt, bd(2), bd(3), bd(4), nbd, n)
381 vlag%lf(1)%x_d, vlag%lf(2)%x_d, &
382 wlag%lf(1)%x_d, wlag%lf(2)%x_d, &
383 bfx_d, bfy_d, bfz_d, u%x_d, v%x_d, w%x_d, &
384 b_d, rho, dt, bd(2), bd(3), bd(4), nbd, n)
391 integer,
intent(in) :: n, nbd
394 real(kind=
rp),
intent(inout) :: fs(n)
395 real(kind=
rp),
intent(in) :: b(n)
396 real(kind=
rp),
intent(in) :: dt, rho, bd(4)
397 type(c_ptr) :: fs_d, b_d
404 fs_d, s%x_d, b_d, rho, dt, bd(2), bd(3), bd(4), &
408 fs_d, s%x_d, b_d, rho, dt, bd(2), bd(3), bd(4), &
412 fs_d, s%x_d, b_d, rho, dt, bd(2), bd(3), bd(4), &
Return the device pointer for an associated Fortran array.
Device abstraction, common interface for various accelerators.
integer, parameter, public c_rp
integer, parameter, public rp
Global precision used in computations.
subroutine rhs_maker_bdf_device(ulag, vlag, wlag, bfx, bfy, bfz, u, v, w, B, rho, dt, bd, nbd, n)
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 scalar_rhs_maker_ext_device(fs_lag, fs_laglag, fs, rho, ext_coeffs, n)
subroutine rhs_maker_ext_device(fx_lag, fy_lag, fz_lag, fx_laglag, fy_laglag, fz_laglag, fx, fy, fz, rho, ext_coeffs, n)
Routines to generate the right-hand sides for the convection-diffusion equation. Employs the EXT/BDF ...
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)
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)
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)
void scalar_rhs_maker_ext_opencl(void *fs_lag, void *fs_laglag, void *fs, real *rho, real *ext1, real *ext2, real *ext3, int *n)
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)
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)
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)
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)
void scalar_rhs_maker_ext_cuda(void *fs_lag, void *fs_laglag, void *fs, real *rho, real *ext1, real *ext2, real *ext3, int *n)
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)
Abstract type to add contributions to F from lagged BD terms.
Abstract type to sum up contributions to kth order extrapolation scheme.
Abstract type to compute extrapolated velocity field for the pressure equation.