Neko  0.9.0
A portable framework for high-order spectral element flow simulations
device_dirichlet.F90
Go to the documentation of this file.
1 ! Copyright (c) 2021-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 num_types
35  use utils
36  use, intrinsic :: iso_c_binding, only : c_ptr, c_int
37  implicit none
38  private
39 
40 #ifdef HAVE_HIP
41  interface
42  subroutine hip_dirichlet_apply_scalar(msk, x, g, m) &
43  bind(c, name='hip_dirichlet_apply_scalar')
44  use, intrinsic :: iso_c_binding
45  import c_rp
46  implicit none
47  real(c_rp) :: g
48  integer(c_int) :: m
49  type(c_ptr), value :: msk, x
50  end subroutine hip_dirichlet_apply_scalar
51  end interface
52 
53  interface
54  subroutine hip_dirichlet_apply_vector(msk, x, y, z, g, m) &
55  bind(c, name='hip_dirichlet_apply_vector')
56  use, intrinsic :: iso_c_binding
57  import c_rp
58  implicit none
59  real(c_rp) :: g
60  integer(c_int) :: m
61  type(c_ptr), value :: msk, x, y, z
62  end subroutine hip_dirichlet_apply_vector
63  end interface
64 #elif HAVE_CUDA
65  interface
66  subroutine cuda_dirichlet_apply_scalar(msk, x, g, m) &
67  bind(c, name='cuda_dirichlet_apply_scalar')
68  use, intrinsic :: iso_c_binding
69  import c_rp
70  implicit none
71  real(c_rp) :: g
72  integer(c_int) :: m
73  type(c_ptr), value :: msk, x
74  end subroutine cuda_dirichlet_apply_scalar
75  end interface
76 
77  interface
78  subroutine cuda_dirichlet_apply_vector(msk, x, y, z, g, m) &
79  bind(c, name='cuda_dirichlet_apply_vector')
80  use, intrinsic :: iso_c_binding
81  import c_rp
82  implicit none
83  real(c_rp) :: g
84  integer(c_int) :: m
85  type(c_ptr), value :: msk, x, y, z
86  end subroutine cuda_dirichlet_apply_vector
87  end interface
88 #elif HAVE_OPENCL
89  interface
90  subroutine opencl_dirichlet_apply_scalar(msk, x, g, m) &
91  bind(c, name='opencl_dirichlet_apply_scalar')
92  use, intrinsic :: iso_c_binding
93  import c_rp
94  implicit none
95  real(c_rp) :: g
96  integer(c_int) :: m
97  type(c_ptr), value :: msk, x
98  end subroutine opencl_dirichlet_apply_scalar
99  end interface
100 
101  interface
102  subroutine opencl_dirichlet_apply_vector(msk, x, y, z, g, m) &
103  bind(c, name='opencl_dirichlet_apply_vector')
104  use, intrinsic :: iso_c_binding
105  import c_rp
106  implicit none
107  real(c_rp) :: g
108  integer(c_int) :: m
109  type(c_ptr), value :: msk, x, y, z
110  end subroutine opencl_dirichlet_apply_vector
111  end interface
112 #endif
113 
115 
116 contains
117 
118  subroutine device_dirichlet_apply_scalar(msk, x, g, m)
119  integer, intent(in) :: m
120  type(c_ptr) :: msk, x
121  real(kind=rp), intent(in) :: g
122 
123 #ifdef HAVE_HIP
124  call hip_dirichlet_apply_scalar(msk, x, g, m)
125 #elif HAVE_CUDA
126  call cuda_dirichlet_apply_scalar(msk, x, g, m)
127 #elif HAVE_OPENCL
128  call opencl_dirichlet_apply_scalar(msk, x, g, m)
129 #else
130  call neko_error('No device backend configured')
131 #endif
132 
133  end subroutine device_dirichlet_apply_scalar
134 
135  subroutine device_dirichlet_apply_vector(msk, x, y, z, g, m)
136  integer, intent(in) :: m
137  type(c_ptr) :: msk, x, y, z
138  real(kind=rp), intent(in) :: g
139 
140 #ifdef HAVE_HIP
141  call hip_dirichlet_apply_vector(msk, x, y, z, g, m)
142 #elif HAVE_CUDA
143  call cuda_dirichlet_apply_vector(msk, x, y, z, g, m)
144 #elif HAVE_OPENCL
145  call opencl_dirichlet_apply_vector(msk, x, y, z, g, m)
146 #else
147  call neko_error('No device backend configured')
148 #endif
149 
150  end subroutine device_dirichlet_apply_vector
151 
152 end module device_dirichlet
void opencl_dirichlet_apply_vector(void *msk, void *x, void *y, void *z, real *g, int *m)
Definition: dirichlet.c:80
void opencl_dirichlet_apply_scalar(void *msk, void *x, real *g, int *m)
Definition: dirichlet.c:52
void cuda_dirichlet_apply_vector(void *msk, void *x, void *y, void *z, real *g, int *m)
Definition: dirichlet.cu:60
void cuda_dirichlet_apply_scalar(void *msk, void *x, real *g, int *m)
Definition: dirichlet.cu:44
subroutine, public device_dirichlet_apply_scalar(msk, x, g, m)
subroutine, public device_dirichlet_apply_vector(msk, x, y, z, g, m)
integer, parameter, public c_rp
Definition: num_types.f90:13
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Utilities.
Definition: utils.f90:35