Neko  0.8.99
A portable framework for high-order spectral element flow simulations
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 !
34  use rhs_maker
35  use device
36  use utils
37  use field_series, only : field_series_t
38  use field, only : field_t
39  use num_types, only : rp, c_rp
40  use, intrinsic :: iso_c_binding, only : c_ptr
41  implicit none
42  private
43 
44  type, public, extends(rhs_maker_sumab_t) :: rhs_maker_sumab_device_t
45  contains
46  procedure, nopass :: compute_fluid => rhs_maker_sumab_device
48 
49  type, public, extends(rhs_maker_ext_t) :: rhs_maker_ext_device_t
50  contains
51  procedure, nopass :: compute_fluid => rhs_maker_ext_device
52  procedure, nopass :: compute_scalar => scalar_rhs_maker_ext_device
53  end type rhs_maker_ext_device_t
54 
55  type, public, extends(rhs_maker_bdf_t) :: rhs_maker_bdf_device_t
56  contains
57  procedure, nopass :: compute_fluid => rhs_maker_bdf_device
58  procedure, nopass :: compute_scalar => scalar_rhs_maker_bdf_device
59  end type rhs_maker_bdf_device_t
60 
61 #ifdef HAVE_HIP
62  interface
63  subroutine rhs_maker_sumab_hip(u_d, v_d, w_d, uu_d, vv_d, ww_d, &
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
67  import c_rp
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
72  end subroutine rhs_maker_sumab_hip
73  end interface
74 
75  interface
76  subroutine rhs_maker_ext_hip(abx1_d, aby1_d, abz1_d, &
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
82  import c_rp
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
87  integer(c_int) :: n
88  end subroutine rhs_maker_ext_hip
89  end interface
90 
91  interface
92  subroutine scalar_rhs_maker_ext_hip(fs_lag_d, fs_laglag_d, fs_d, rho, &
93  ext1, ext2, ext3, n) &
94  bind(c, name='scalar_rhs_maker_ext_hip')
95  use, intrinsic :: iso_c_binding
96  import c_rp
97  type(c_ptr), value :: fs_lag_d, fs_laglag_d, fs_d
98  real(c_rp) :: rho, ext1, ext2, ext3
99  integer(c_int) :: n
100  end subroutine scalar_rhs_maker_ext_hip
101  end interface
102 
103  interface
104  subroutine rhs_maker_bdf_hip(ulag1_d, ulag2_d, vlag1_d, vlag2_d, &
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
108  import c_rp
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
114  end subroutine rhs_maker_bdf_hip
115  end interface
116 
117  interface
118  subroutine scalar_rhs_maker_bdf_hip(s_lag_d, s_laglag_d, fs_d, s_d, B_d, &
119  rho, dt, bd2, bd3, bd4, nbd, n) &
120  bind(c, name='scalar_rhs_maker_bdf_hip')
121  use, intrinsic :: iso_c_binding
122  import c_rp
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
127  end subroutine scalar_rhs_maker_bdf_hip
128  end interface
129 #elif HAVE_CUDA
130  interface
131  subroutine rhs_maker_sumab_cuda(u_d, v_d, w_d, uu_d, vv_d, ww_d, &
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
135  import c_rp
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
140  end subroutine rhs_maker_sumab_cuda
141  end interface
142 
143  interface
144  subroutine rhs_maker_ext_cuda(abx1_d, aby1_d, abz1_d, &
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
150  import c_rp
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
155  integer(c_int) :: n
156  end subroutine rhs_maker_ext_cuda
157  end interface
158 
159  interface
160  subroutine scalar_rhs_maker_ext_cuda(fs_lag_d, fs_laglag_d, fs_d, rho, &
161  ext1, ext2, ext3, n) &
162  bind(c, name='scalar_rhs_maker_ext_cuda')
163  use, intrinsic :: iso_c_binding
164  import c_rp
165  type(c_ptr), value :: fs_lag_d, fs_laglag_d, fs_d
166  real(c_rp) :: rho, ext1, ext2, ext3
167  integer(c_int) :: n
168  end subroutine scalar_rhs_maker_ext_cuda
169  end interface
170 
171  interface
172  subroutine rhs_maker_bdf_cuda(ulag1_d, ulag2_d, vlag1_d, vlag2_d, &
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
176  import c_rp
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
182  end subroutine rhs_maker_bdf_cuda
183  end interface
184 
185  interface
186  subroutine scalar_rhs_maker_bdf_cuda(s_lag_d, s_laglag_d, fs_d, s_d, B_d, &
187  rho, dt, bd2, bd3, bd4, nbd, n) &
188  bind(c, name='scalar_rhs_maker_bdf_cuda')
189  use, intrinsic :: iso_c_binding
190  import c_rp
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
195  end subroutine scalar_rhs_maker_bdf_cuda
196  end interface
197 #elif HAVE_OPENCL
198  interface
199  subroutine rhs_maker_sumab_opencl(u_d, v_d, w_d, uu_d, vv_d, ww_d, &
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
203  import c_rp
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
208  end subroutine rhs_maker_sumab_opencl
209  end interface
210 
211  interface
212  subroutine rhs_maker_ext_opencl(abx1_d, aby1_d, abz1_d, &
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
218  import c_rp
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
223  integer(c_int) :: n
224  end subroutine rhs_maker_ext_opencl
225  end interface
226 
227  interface
228  subroutine scalar_rhs_maker_ext_opencl(fs_lag_d, fs_laglag_d, fs_d, rho, &
229  ext1, ext2, ext3, n) &
230  bind(c, name='scalar_rhs_maker_ext_opencl')
231  use, intrinsic :: iso_c_binding
232  import c_rp
233  type(c_ptr), value :: fs_lag_d, fs_laglag_d, fs_d
234  real(c_rp) :: rho, ext1, ext2, ext3
235  integer(c_int) :: n
236  end subroutine scalar_rhs_maker_ext_opencl
237  end interface
238 
239  interface
240  subroutine rhs_maker_bdf_opencl(ulag1_d, ulag2_d, vlag1_d, vlag2_d, &
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
244  import c_rp
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
250  end subroutine rhs_maker_bdf_opencl
251  end interface
252 
253  interface
254  subroutine scalar_rhs_maker_bdf_opencl(s_lag_d, s_laglag_d, fs_d, s_d, B_d, &
255  rho, dt, bd2, bd3, bd4, nbd, n) &
256  bind(c, name='scalar_rhs_maker_bdf_opencl')
257  use, intrinsic :: iso_c_binding
258  import c_rp
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
263  end subroutine scalar_rhs_maker_bdf_opencl
264  end interface
265 #endif
266 
267 contains
268 
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
272  type(field_series_t), intent(inout) :: uulag, vvlag, wwlag
273  real(kind=rp), dimension(3), intent(in) :: ab
274  integer, intent(in) :: nab
275 
276 #ifdef HAVE_HIP
277  call rhs_maker_sumab_hip(u%x_d, v%x_d, w%x_d, uu%x_d, vv%x_d, ww%x_d, &
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, &
280  uu%dof%size())
281 #elif HAVE_CUDA
282  call rhs_maker_sumab_cuda(u%x_d, v%x_d, w%x_d, uu%x_d, vv%x_d, ww%x_d, &
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, &
285  uu%dof%size())
286 #elif HAVE_OPENCL
287  call rhs_maker_sumab_opencl(u%x_d, v%x_d, w%x_d, uu%x_d, vv%x_d, ww%x_d, &
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, &
290  uu%dof%size())
291 #endif
292 
293  end subroutine rhs_maker_sumab_device
294 
295  subroutine rhs_maker_ext_device(fx_lag, fy_lag, fz_lag, &
296  fx_laglag, fy_laglag, fz_laglag, fx, fy, fz, &
297  rho, ext_coeffs, n)
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
304 
305  fx_d = device_get_ptr(fx)
306  fy_d = device_get_ptr(fy)
307  fz_d = device_get_ptr(fz)
308 
309 #ifdef HAVE_HIP
310  call rhs_maker_ext_hip(fx_lag%x_d, fy_lag%x_d, fz_lag%x_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)
314 #elif HAVE_CUDA
315  call rhs_maker_ext_cuda(fx_lag%x_d, fy_lag%x_d, fz_lag%x_d, &
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)
319 #elif HAVE_OPENCL
320  call rhs_maker_ext_opencl(fx_lag%x_d, fy_lag%x_d, fz_lag%x_d, &
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)
324 #endif
325 
326  end subroutine rhs_maker_ext_device
327 
328  subroutine scalar_rhs_maker_ext_device(fs_lag, fs_laglag, fs, &
329  rho, ext_coeffs, 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)
335  type(c_ptr) :: fs_d
336 
337  fs_d = device_get_ptr(fs)
338 
339 #ifdef HAVE_HIP
340  call scalar_rhs_maker_ext_hip(fs_lag%x_d, fs_laglag%x_d, fs_d, rho, &
341  ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
342 #elif HAVE_CUDA
343  call scalar_rhs_maker_ext_cuda(fs_lag%x_d, fs_laglag%x_d, fs_d, rho, &
344  ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
345 #elif HAVE_OPENCL
346  call scalar_rhs_maker_ext_opencl(fs_lag%x_d, fs_laglag%x_d, fs_d, rho, &
347  ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
348 #endif
349 
350  end subroutine scalar_rhs_maker_ext_device
351 
352  subroutine rhs_maker_bdf_device(ulag, vlag, wlag, bfx, bfy, bfz, &
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
356  type(field_series_t), intent(in) :: ulag, vlag, wlag
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
361 
362  bfx_d = device_get_ptr(bfx)
363  bfy_d = device_get_ptr(bfy)
364  bfz_d = device_get_ptr(bfz)
365  b_d = device_get_ptr(b)
366 
367 #ifdef HAVE_HIP
368  call rhs_maker_bdf_hip(ulag%lf(1)%x_d, ulag%lf(2)%x_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)
373 #elif HAVE_CUDA
374  call rhs_maker_bdf_cuda(ulag%lf(1)%x_d, ulag%lf(2)%x_d, &
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)
379 #elif HAVE_OPENCL
380  call rhs_maker_bdf_opencl(ulag%lf(1)%x_d, ulag%lf(2)%x_d, &
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)
385 #endif
386 
387  end subroutine rhs_maker_bdf_device
388 
389  subroutine scalar_rhs_maker_bdf_device(s_lag, fs, s, B, rho, dt, &
390  bd, nbd, n)
391  integer, intent(in) :: n, nbd
392  type(field_t), intent(in) :: s
393  type(field_series_t), intent(in) :: s_lag
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
398 
399  fs_d = device_get_ptr(fs)
400  b_d = device_get_ptr(b)
401 
402 #ifdef HAVE_HIP
403  call scalar_rhs_maker_bdf_hip(s_lag%lf(1)%x_d, s_lag%lf(2)%x_d, &
404  fs_d, s%x_d, b_d, rho, dt, bd(2), bd(3), bd(4), &
405  nbd, n)
406 #elif HAVE_CUDA
407  call scalar_rhs_maker_bdf_cuda(s_lag%lf(1)%x_d, s_lag%lf(2)%x_d, &
408  fs_d, s%x_d, b_d, rho, dt, bd(2), bd(3), bd(4), &
409  nbd, n)
410 #elif HAVE_OPENCL
411  call scalar_rhs_maker_bdf_opencl(s_lag%lf(1)%x_d, s_lag%lf(2)%x_d, &
412  fs_d, s%x_d, b_d, rho, dt, bd(2), bd(3), bd(4), &
413  nbd, n)
414 #endif
415 
416  end subroutine scalar_rhs_maker_bdf_device
417 
418 end module rhs_maker_device
Return the device pointer for an associated Fortran array.
Definition: device.F90:81
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
Stores a series fields.
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 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 ...
Definition: rhs_maker.f90:38
Utilities.
Definition: utils.f90:35
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:201
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 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:89
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:126
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:156
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:119
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:62
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:97
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:81
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:43
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 compute extrapolated velocity field for the pressure equation.
Definition: rhs_maker.f90:46