Neko  0.8.99
A portable framework for high-order spectral element flow simulations
fluid_abbdf_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 !
33 module rhs_maker_device
34  use rhs_maker
35  use device
36  use utils
37  use, intrinsic :: iso_c_binding
38  implicit none
39  private
40 
41  type, public, extends(rhs_maker_sumab_t) :: rhs_maker_sumab_device_t
42  contains
43  procedure, nopass :: compute_fluid => rhs_maker_sumab_device
45 
46  type, public, extends(rhs_maker_ext_t) :: rhs_maker_ext_device_t
47  contains
48  procedure, nopass :: compute_fluid => rhs_maker_ext_device
49  end type rhs_maker_ext_device_t
50 
51  type, public, extends(rhs_maker_bdf_t) :: rhs_maker_bdf_device_t
52  contains
53  procedure, nopass :: compute_fluid => rhs_maker_bdf_device
54  end type rhs_maker_bdf_device_t
55 
56 #ifdef HAVE_HIP
57  interface
58  subroutine rhs_maker_sumab_hip(u_d, v_d, w_d, uu_d, vv_d, ww_d, &
59  uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2, ab1, ab2, ab3, nab, n)&
60  bind(c, name='rhs_maker_sumab_hip')
61  use, intrinsic :: iso_c_binding
62  import c_rp
63  type(c_ptr), value :: u_d, v_d, w_d, uu_d, vv_d, ww_d
64  type(c_ptr), value :: uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2
65  real(c_rp) :: ab1, ab2, ab3
66  integer(c_int) :: nab, n
67  end subroutine rhs_maker_sumab_hip
68  end interface
69 
70  interface
71  subroutine rhs_maker_ext_hip(abx1_d, aby1_d, abz1_d, &
72  abx2_d, aby2_d, abz2_d, &
73  bfx_d, bfy_d, bfz_d, &
74  rho, ab1, ab2, ab3, n) &
75  bind(c, name='rhs_maker_ext_hip')
76  use, intrinsic :: iso_c_binding
77  import c_rp
78  type(c_ptr), value :: abx1_d, aby1_d, abz1_d
79  type(c_ptr), value :: abx2_d, aby2_d, abz2_d
80  type(c_ptr), value :: bfx_d, bfy_d, bfz_d
81  real(c_rp) :: rho, ab1, ab2, ab3
82  integer(c_int) :: n
83  end subroutine rhs_maker_ext_hip
84  end interface
85 
86  interface
87  subroutine rhs_maker_bdf_hip(ulag1_d, ulag2_d, vlag1_d, vlag2_d, &
88  wlag1_d, wlag2_d, bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d, &
89  rho, dt, bd2, bd3, bd4, nbd, n) bind(c, name='rhs_maker_bdf_hip')
90  use, intrinsic :: iso_c_binding
91  import c_rp
92  type(c_ptr), value :: ulag1_d, ulag2_d, vlag1_d
93  type(c_ptr), value :: vlag2_d, wlag1_d, wlag2_d
94  type(c_ptr), value :: bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d
95  reaL(c_rp) :: rho, dt, bd2, bd3, bd4
96  integer(c_int) :: nbd, n
97  end subroutine rhs_maker_bdf_hip
98  end interface
99 #elif HAVE_CUDA
100  interface
101  subroutine rhs_maker_sumab_cuda(u_d, v_d, w_d, uu_d, vv_d, ww_d, &
102  uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2, ab1, ab2, ab3, nab, n)&
103  bind(c, name='rhs_maker_sumab_cuda')
104  use, intrinsic :: iso_c_binding
105  import c_rp
106  type(c_ptr), value :: u_d, v_d, w_d, uu_d, vv_d, ww_d
107  type(c_ptr), value :: uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2
108  real(c_rp) :: ab1, ab2, ab3
109  integer(c_int) :: nab, n
110  end subroutine rhs_maker_sumab_cuda
111  end interface
112 
113  interface
114  subroutine rhs_maker_ext_cuda(abx1_d, aby1_d, abz1_d, &
115  abx2_d, aby2_d, abz2_d, &
116  bfx_d, bfy_d, bfz_d, &
117  rho, ab1, ab2, ab3, n) &
118  bind(c, name='rhs_maker_ext_cuda')
119  use, intrinsic :: iso_c_binding
120  import c_rp
121  type(c_ptr), value :: abx1_d, aby1_d, abz1_d
122  type(c_ptr), value :: abx2_d, aby2_d, abz2_d
123  type(c_ptr), value :: bfx_d, bfy_d, bfz_d
124  real(c_rp) :: rho, ab1, ab2, ab3
125  integer(c_int) :: n
126  end subroutine rhs_maker_ext_cuda
127  end interface
128 
129  interface
130  subroutine rhs_maker_bdf_cuda(ulag1_d, ulag2_d, vlag1_d, vlag2_d, &
131  wlag1_d, wlag2_d, bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d, &
132  rho, dt, bd2, bd3, bd4, nbd, n) bind(c, name='rhs_maker_bdf_cuda')
133  use, intrinsic :: iso_c_binding
134  import c_rp
135  type(c_ptr), value :: ulag1_d, ulag2_d, vlag1_d
136  type(c_ptr), value :: vlag2_d, wlag1_d, wlag2_d
137  type(c_ptr), value :: bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d
138  reaL(c_rp) :: rho, dt, bd2, bd3, bd4
139  integer(c_int) :: nbd, n
140  end subroutine rhs_maker_bdf_cuda
141  end interface
142 #elif HAVE_OPENCL
143  interface
144  subroutine rhs_maker_sumab_opencl(u_d, v_d, w_d, uu_d, vv_d, ww_d, &
145  uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2, ab1, ab2, ab3, nab, n)&
146  bind(c, name='rhs_maker_sumab_opencl')
147  use, intrinsic :: iso_c_binding
148  import c_rp
149  type(c_ptr), value :: u_d, v_d, w_d, uu_d, vv_d, ww_d
150  type(c_ptr), value :: uulag1, uulag2, vvlag1, vvlag2, wwlag1, wwlag2
151  real(c_rp) :: ab1, ab2, ab3
152  integer(c_int) :: nab, n
153  end subroutine rhs_maker_sumab_opencl
154  end interface
155 
156  interface
157  subroutine rhs_maker_ext_opencl(abx1_d, aby1_d, abz1_d, &
158  abx2_d, aby2_d, abz2_d, &
159  bfx_d, bfy_d, bfz_d, &
160  rho, ab1, ab2, ab3, n) &
161  bind(c, name='rhs_maker_ext_opencl')
162  use, intrinsic :: iso_c_binding
163  import c_rp
164  type(c_ptr), value :: abx1_d, aby1_d, abz1_d
165  type(c_ptr), value :: abx2_d, aby2_d, abz2_d
166  type(c_ptr), value :: bfx_d, bfy_d, bfz_d
167  real(c_rp) :: rho, ab1, ab2, ab3
168  integer(c_int) :: n
169  end subroutine rhs_maker_ext_opencl
170  end interface
171 
172  interface
173  subroutine rhs_maker_bdf_opencl(ulag1_d, ulag2_d, vlag1_d, vlag2_d, &
174  wlag1_d, wlag2_d, bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d, &
175  rho, dt, bd2, bd3, bd4, nbd, n) bind(c, name='rhs_maker_bdf_opencl')
176  use, intrinsic :: iso_c_binding
177  import c_rp
178  type(c_ptr), value :: ulag1_d, ulag2_d, vlag1_d
179  type(c_ptr), value :: vlag2_d, wlag1_d, wlag2_d
180  type(c_ptr), value :: bfx_d, bfy_d, bfz_d, u_d, v_d, w_d, B_d
181  reaL(c_rp) :: rho, dt, bd2, bd3, bd4
182  integer(c_int) :: nbd, n
183  end subroutine rhs_maker_bdf_opencl
184  end interface
185 #endif
186 
187 contains
188 
189  subroutine rhs_maker_sumab_device(u, v, w, uu, vv, ww, uulag, vvlag, wwlag, ab, nab)
190  type(field_t), intent(inout) :: u,v, w
191  type(field_t), intent(inout) :: uu, vv, ww
192  type(field_series_t), intent(inout) :: uulag, vvlag, wwlag
193  real(kind=rp), dimension(3), intent(in) :: ab
194  integer, intent(in) :: nab
195 
196 #ifdef HAVE_HIP
197  call rhs_maker_sumab_hip(u%x_d, v%x_d, w%x_d, uu%x_d, vv%x_d, ww%x_d, &
198  uulag%lf(1)%x_d, uulag%lf(2)%x_d, vvlag%lf(1)%x_d, vvlag%lf(2)%x_d, &
199  wwlag%lf(1)%x_d, wwlag%lf(2)%x_d, ab(1), ab(2), ab(3), nab, &
200  uu%dof%size())
201 #elif HAVE_CUDA
202  call rhs_maker_sumab_cuda(u%x_d, v%x_d, w%x_d, uu%x_d, vv%x_d, ww%x_d, &
203  uulag%lf(1)%x_d, uulag%lf(2)%x_d, vvlag%lf(1)%x_d, vvlag%lf(2)%x_d, &
204  wwlag%lf(1)%x_d, wwlag%lf(2)%x_d, ab(1), ab(2), ab(3), nab, &
205  uu%dof%size())
206 #elif HAVE_OPENCL
207  call rhs_maker_sumab_opencl(u%x_d, v%x_d, w%x_d, uu%x_d, vv%x_d, ww%x_d, &
208  uulag%lf(1)%x_d, uulag%lf(2)%x_d, vvlag%lf(1)%x_d, vvlag%lf(2)%x_d, &
209  wwlag%lf(1)%x_d, wwlag%lf(2)%x_d, ab(1), ab(2), ab(3), nab, &
210  uu%dof%size())
211 #endif
212 
213  end subroutine rhs_maker_sumab_device
214 
215  subroutine rhs_maker_ext_device(temp1, temp2, temp3, fx_lag, fy_lag, fz_lag, &
216  fx_laglag, fy_laglag, fz_laglag, fx, fy, fz, &
217  rho, ext_coeffs, n)
218  type(field_t), intent(inout) :: temp1, temp2, temp3
219  type(field_t), intent(inout) :: fx_lag, fy_lag, fz_lag
220  type(field_t), intent(inout) :: fx_laglag, fy_laglag, fz_laglag
221  real(kind=rp), intent(inout) :: rho, ext_coeffs(4)
222  integer, intent(in) :: n
223  real(kind=rp), intent(inout) :: fx(n), fy(n), fz(n)
224  type(c_ptr) :: fx_d, fy_d, fz_d
225 
226  fx_d = device_get_ptr(fx)
227  fy_d = device_get_ptr(fy)
228  fz_d = device_get_ptr(fz)
229 
230 #ifdef HAVE_HIP
231  call rhs_maker_ext_hip(fx_lag%x_d, fy_lag%x_d, fz_lag%x_d, &
232  fx_laglag%x_d, fy_laglag%x_d, fz_laglag%x_d, &
233  fx_d, fy_d, fz_d, rho, &
234  ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
235 #elif HAVE_CUDA
236  call rhs_maker_ext_cuda(fx_lag%x_d, fy_lag%x_d, fz_lag%x_d, &
237  fx_laglag%x_d, fy_laglag%x_d, fz_laglag%x_d, &
238  fx_d, fy_d, fz_d, rho, &
239  ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
240 #elif HAVE_OPENCL
241  call rhs_maker_ext_opencl(fx_lag%x_d, fy_lag%x_d, fz_lag%x_d, &
242  fx_laglag%x_d, fy_laglag%x_d, fz_laglag%x_d, &
243  fx_d, fy_d, fz_d, rho, &
244  ext_coeffs(1), ext_coeffs(2), ext_coeffs(3), n)
245 #endif
246 
247  end subroutine rhs_maker_ext_device
248 
249  subroutine rhs_maker_bdf_device(ta1, ta2, ta3, tb1, tb2, tb3, &
250  ulag, vlag, wlag, bfx, bfy, bfz, &
251  u, v, w, B, rho, dt, bd, nbd, n)
252  integer, intent(in) :: n, nbd
253  type(field_t), intent(inout) :: ta1, ta2, ta3
254  type(field_t), intent(in) :: u, v, w
255  type(field_t), intent(inout) :: tb1, tb2, tb3
256  type(field_series_t), intent(in) :: ulag, vlag, wlag
257  real(kind=rp), intent(inout) :: bfx(n), bfy(n), bfz(n)
258  real(kind=rp), intent(in) :: b(n)
259  real(kind=rp), intent(in) :: dt, rho, bd(10)
260  type(c_ptr) :: bfx_d, bfy_d, bfz_d, b_d
261 
262  bfx_d = device_get_ptr(bfx)
263  bfy_d = device_get_ptr(bfy)
264  bfz_d = device_get_ptr(bfz)
265  b_d = device_get_ptr(b)
266 
267 #ifdef HAVE_HIP
268  call rhs_maker_bdf_hip(ulag%lf(1)%x_d, ulag%lf(2)%x_d, &
269  vlag%lf(1)%x_d, vlag%lf(2)%x_d, &
270  wlag%lf(1)%x_d, wlag%lf(2)%x_d, &
271  bfx_d, bfy_d, bfz_d, u%x_d, v%x_d, w%x_d, &
272  b_d, rho, dt, bd(2), bd(3), bd(4), nbd, n)
273 #elif HAVE_CUDA
274  call rhs_maker_bdf_cuda(ulag%lf(1)%x_d, ulag%lf(2)%x_d, &
275  vlag%lf(1)%x_d, vlag%lf(2)%x_d, &
276  wlag%lf(1)%x_d, wlag%lf(2)%x_d, &
277  bfx_d, bfy_d, bfz_d, u%x_d, v%x_d, w%x_d, &
278  b_d, rho, dt, bd(2), bd(3), bd(4), nbd, n)
279 #elif HAVE_OPENCL
280  call rhs_maker_bdf_opencl(ulag%lf(1)%x_d, ulag%lf(2)%x_d, &
281  vlag%lf(1)%x_d, vlag%lf(2)%x_d, &
282  wlag%lf(1)%x_d, wlag%lf(2)%x_d, &
283  bfx_d, bfy_d, bfz_d, u%x_d, v%x_d, w%x_d, &
284  b_d, rho, dt, bd(2), bd(3), bd(4), nbd, n)
285 #endif
286 
287  end subroutine rhs_maker_bdf_device
288 
289 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
subroutine rhs_maker_bdf_device(ulag, vlag, wlag, bfx, bfy, bfz, u, v, w, 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_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 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 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 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 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