41 use,
intrinsic :: iso_c_binding, only : c_ptr
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
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
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
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
100 ext1, ext2, ext3, n) &
101 bind(c, name =
'scalar_rhs_maker_ext_hip')
102 use,
intrinsic :: iso_c_binding
104 type(c_ptr),
value :: fs_lag_d, fs_laglag_d, fs_d
105 real(c_rp) :: rho, ext1, ext2, ext3
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
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
126 rho, dt, bd2, bd3, bd4, nbd, n) &
127 bind(c, name =
'scalar_rhs_maker_bdf_hip')
128 use,
intrinsic :: iso_c_binding
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
139 bf_y_d, bf_z_d, rho, dt, n) bind(c, name = 'rhs_maker_oifs_hip')
140 use,
intrinsic :: iso_c_binding
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
151 bind(c, name =
'scalar_rhs_maker_oifs_hip')
152 use,
intrinsic :: iso_c_binding
154 type(c_ptr),
value :: bf_s_d
155 type(c_ptr),
value :: phi_s_d
156 reaL(c_rp) :: rho, dt
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
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
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
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
192 ext1, ext2, ext3, n) &
193 bind(c, name =
'scalar_rhs_maker_ext_cuda')
194 use,
intrinsic :: iso_c_binding
196 type(c_ptr),
value :: fs_lag_d, fs_laglag_d, fs_d
197 real(c_rp) :: rho, ext1, ext2, ext3
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
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
218 rho, dt, bd2, bd3, bd4, nbd, n) &
219 bind(c, name =
'scalar_rhs_maker_bdf_cuda')
220 use,
intrinsic :: iso_c_binding
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
231 bf_y_d, bf_z_d, rho, dt, n) bind(c, name = 'rhs_maker_oifs_cuda')
232 use,
intrinsic :: iso_c_binding
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
243 bind(c, name =
'scalar_rhs_maker_oifs_cuda')
244 use,
intrinsic :: iso_c_binding
246 type(c_ptr),
value :: bf_s_d
247 type(c_ptr),
value :: phi_s_d
248 reaL(c_rp) :: rho, dt
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
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
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
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
284 ext1, ext2, ext3, n) &
285 bind(c, name =
'scalar_rhs_maker_ext_opencl')
286 use,
intrinsic :: iso_c_binding
288 type(c_ptr),
value :: fs_lag_d, fs_laglag_d, fs_d
289 real(c_rp) :: rho, ext1, ext2, ext3
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
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
310 rho, dt, bd2, bd3, bd4, nbd, n) &
311 bind(c, name =
'scalar_rhs_maker_bdf_opencl')
312 use,
intrinsic :: iso_c_binding
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
323 bf_y_d, bf_z_d, rho, dt, n) bind(c, name = 'rhs_maker_oifs_opencl')
324 use,
intrinsic :: iso_c_binding
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
335 bind(c, name =
'scalar_rhs_maker_oifs_opencl')
336 use,
intrinsic :: iso_c_binding
338 type(c_ptr),
value :: bf_s_d
339 type(c_ptr),
value :: phi_s_d
340 reaL(c_rp) :: rho, dt
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
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
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
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
371 end subroutine rhs_maker_ext_metal
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
380 type(c_ptr),
value :: fs_lag_d, fs_laglag_d, fs_d
381 real(c_rp) :: rho, ext1, ext2, ext3
383 end subroutine scalar_rhs_maker_ext_metal
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
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
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
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
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
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
422 end subroutine rhs_maker_oifs_metal
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
430 type(c_ptr),
value :: bf_s_d
431 type(c_ptr),
value :: phi_s_d
432 reaL(c_rp) :: rho, dt
434 end subroutine scalar_rhs_maker_oifs_metal
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
444 real(kind=
rp),
dimension(3),
intent(in) :: ab
445 integer,
intent(in) :: nab
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, &
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, &
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, &
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, &
472 fx_laglag, fy_laglag, fz_laglag, fx, fy, fz, &
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
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)
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)
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)
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)
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)
522 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
525 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
528 ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
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)
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
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
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)
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)
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)
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)
582 integer,
intent(in) :: n, nbd
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
595 fs_d, s%x_d, b_d, rho, dt, bd(2), bd(3), bd(4), &
599 fs_d, s%x_d, b_d, rho, dt, bd(2), bd(3), bd(4), &
603 fs_d, s%x_d, b_d, rho, dt, bd(2), bd(3), bd(4), &
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), &
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)
620 type(c_ptr) :: phi_x_d, phi_y_d, phi_z_d, bf_x_d, bf_y_d, bf_z_d
631 bf_x_d, bf_y_d, bf_z_d, rho, dt, n)
634 bf_x_d, bf_y_d, bf_z_d, rho, dt, n)
637 bf_x_d, bf_y_d, bf_z_d, rho, dt, n)
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)
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
662 call scalar_rhs_maker_oifs_metal(phi_s_d, bf_s_d, rho, dt, n)
Return the device pointer for an associated Fortran array.
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
Device abstraction, common interface for various accelerators.
Contains the field_serties_t type.
integer, parameter, public c_rp
integer, parameter, public rp
Global precision used in computations.
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 ...
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 scalar_rhs_maker_oifs_opencl(void *phi_s, void *bf_s, real *rho, real *dt, 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 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)
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_oifs_cuda(void *phi_s, void *bf_s, real *rho, real *dt, 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_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)
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)
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.
Abstract type to sum up contributions to kth order extrapolation scheme.
Abstract type to add contributions of kth order OIFS scheme.
Abstract type to compute extrapolated velocity field for the pressure equation.