Neko  0.8.99
A portable framework for high-order spectral element flow simulations
fdm_device.F90
Go to the documentation of this file.
1 ! Copyright (c) 2023, 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 fdm_device
34  use num_types
35  use utils
36  use device, only : device_get_ptr, glb_cmd_queue
37  use, intrinsic :: iso_c_binding, only : c_ptr, c_int
38  implicit none
39  private
40 
41 #ifdef HAVE_HIP
42  interface
43  subroutine hip_fdm_do_fast(e_d, r_d, s_d, d_d, nl, nelv, stream) &
44  bind(c, name='hip_fdm_do_fast')
45  use, intrinsic :: iso_c_binding
46  type(c_ptr), value :: e_d, r_d, s_d, d_d, stream
47  integer(c_int) :: nl, nelv
48  end subroutine hip_fdm_do_fast
49  end interface
50 #elif HAVE_CUDA
51  interface
52  subroutine cuda_fdm_do_fast(e_d, r_d, s_d, d_d, nl, nelv, stream) &
53  bind(c, name='cuda_fdm_do_fast')
54  use, intrinsic :: iso_c_binding
55  type(c_ptr), value :: e_d, r_d, s_d, d_d, stream
56  integer(c_int) :: nl, nelv
57  end subroutine cuda_fdm_do_fast
58  end interface
59 #elif HAVE_OPENCL
60  interface
61  subroutine opencl_fdm_do_fast(e_d, r_d, s_d, d_d, nl, nelv, stream) &
62  bind(c, name='opencl_fdm_do_fast')
63  use, intrinsic :: iso_c_binding
64  type(c_ptr), value :: e_d, r_d, s_d, d_d, stream
65  integer(c_int) :: nl, nelv
66  end subroutine opencl_fdm_do_fast
67  end interface
68 #endif
69 
70  public :: fdm_do_fast_device
71 
72 contains
73 
74  subroutine fdm_do_fast_device(e, r, s, d, nl, ldim, nelv, stream)
75  integer, intent(in) :: nl, nelv, ldim
76  real(kind=rp), intent(inout) :: e(nl**ldim, nelv)
77  real(kind=rp), intent(inout) :: r(nl**ldim, nelv)
78  real(kind=rp), intent(inout) :: s(nl*nl,2,ldim, nelv)
79  real(kind=rp), intent(inout) :: d(nl**ldim, nelv)
80  type(c_ptr) :: e_d, r_d, s_d, d_d
81  type(c_ptr), optional :: stream
82 
83  e_d = device_get_ptr(e)
84  r_d = device_get_ptr(r)
85  s_d = device_get_ptr(s)
86  d_d = device_get_ptr(d)
87  if (.not. present(stream)) stream = glb_cmd_queue
88  if (ldim .ne. 3) call neko_error('fdm dim not supported')
89 
90 #ifdef HAVE_HIP
91  call hip_fdm_do_fast(e_d, r_d, s_d, d_d, nl, nelv, stream)
92 #elif HAVE_CUDA
93  call cuda_fdm_do_fast(e_d, r_d, s_d, d_d, nl, nelv, stream)
94 #elif HAVE_OPENCL
95  call opencl_fdm_do_fast(e_d, r_d, s_d, d_d, nl, nelv, stream)
96 #else
97  call neko_error('No device backend configured')
98 #endif
99  end subroutine fdm_do_fast_device
100 
101 end module fdm_device
void opencl_fdm_do_fast(void *e, void *r, void *s, void *d, int *nl, int *nel)
Definition: fdm.c:49
void cuda_fdm_do_fast(void *e, void *r, void *s, void *d, int *nl, int *nel, cudaStream_t stream)
Definition: fdm.cu:42
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, public fdm_do_fast_device(e, r, s, d, nl, ldim, nelv, stream)
Definition: fdm_device.F90:75
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Utilities.
Definition: utils.f90:35