Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
device_local_interpolation.F90
Go to the documentation of this file.
1! Copyright (c) 2022-2025, 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, only : rp, c_rp
35 use utils, only : neko_error
36 use, intrinsic :: iso_c_binding, only : c_ptr, c_int
37 implicit none
38 private
39
41
42#ifdef HAVE_HIP
43 interface
44 subroutine hip_find_rst_legendre(rst, pt_x, pt_y, pt_z, &
45 x_hat, y_hat, z_hat, &
46 resx, resy, resz, &
47 lx, el_ids, n_pt, tol, conv_pts) &
48 bind(c, name='hip_find_rst_legendre')
49 use, intrinsic :: iso_c_binding
50 use num_types
51 implicit none
52 type(c_ptr), value :: rst
53 type(c_ptr), value :: pt_x, pt_y, pt_z
54 type(c_ptr), value :: x_hat, y_hat, z_hat
55 type(c_ptr), value :: resx, resy, resz
56 type(c_ptr), value :: el_ids, conv_pts
57 integer(c_int) :: lx, n_pt
58 real(c_rp) :: tol
59 end subroutine hip_find_rst_legendre
60
61 end interface
62#elif HAVE_CUDA
63 interface
64 subroutine cuda_find_rst_legendre(rst, pt_x, pt_y, pt_z, &
65 x_hat, y_hat, z_hat, &
66 resx, resy, resz, &
67 lx, el_ids, n_pt, tol, conv_pts) &
68 bind(c, name='cuda_find_rst_legendre')
69 use, intrinsic :: iso_c_binding
70 use num_types
71 implicit none
72 type(c_ptr), value :: rst
73 type(c_ptr), value :: pt_x, pt_y, pt_z
74 type(c_ptr), value :: x_hat, y_hat, z_hat
75 type(c_ptr), value :: resx, resy, resz
76 type(c_ptr), value :: el_ids, conv_pts
77 integer(c_int) :: lx, n_pt
78 real(c_rp) :: tol
79 end subroutine cuda_find_rst_legendre
80 end interface
81#elif HAVE_OPENCL
82 interface
83 subroutine opencl_find_rst_legendre(rst, pt_x, pt_y, pt_z, &
84 x_hat, y_hat, z_hat, &
85 resx, resy, resz, &
86 lx, el_ids, n_pt, tol, conv_pts) &
87 bind(c, name='opencl_find_rst_legendre')
88 use, intrinsic :: iso_c_binding
89 use num_types
90 implicit none
91 type(c_ptr), value :: rst
92 type(c_ptr), value :: pt_x, pt_y, pt_z
93 type(c_ptr), value :: x_hat, y_hat, z_hat
94 type(c_ptr), value :: resx, resy, resz
95 type(c_ptr), value :: el_ids, conv_pts
96 integer(c_int) :: lx, n_pt
97 real(c_rp) :: tol
98 end subroutine opencl_find_rst_legendre
99 end interface
100#endif
101
102contains
103
104 subroutine device_find_rst_legendre(rst_d, pt_x_d, pt_y_d, pt_z_d, &
105 x_hat_d, y_hat_d, z_hat_d, &
106 resx_d, resy_d, resz_d, &
107 lx, el_ids_d, n_pts, tol, conv_pts_d)
108 type(c_ptr) :: rst_d, pt_x_d, pt_y_d, pt_z_d
109 type(c_ptr) :: x_hat_d, y_hat_d, z_hat_d
110 type(c_ptr) :: resx_d, resy_d, resz_d
111 type(c_ptr) :: el_ids_d, conv_pts_d
112 integer :: lx, n_pts
113 real(kind=c_rp) :: tol
114
115#ifdef HAVE_HIP
116 call hip_find_rst_legendre(rst_d, pt_x_d, pt_y_d, pt_z_d, &
117 x_hat_d, y_hat_d, z_hat_d, &
118 resx_d, resy_d, resz_d, &
119 lx, el_ids_d, n_pts, tol, conv_pts_d)
120#elif HAVE_CUDA
121 call cuda_find_rst_legendre(rst_d, pt_x_d, pt_y_d, pt_z_d, &
122 x_hat_d, y_hat_d, z_hat_d, &
123 resx_d, resy_d, resz_d, &
124 lx, el_ids_d, n_pts, tol, conv_pts_d)
125#elif HAVE_OPENCL
126 call opencl_find_rst_legendre(rst_d, pt_x_d, pt_y_d, pt_z_d, &
127 x_hat_d, y_hat_d, z_hat_d, &
128 resx_d, resy_d, resz_d, &
129 lx, el_ids_d, n_pts, tol, conv_pts_d)
130#else
131 call neko_error('No device backend configured')
132#endif
133
134 end subroutine device_find_rst_legendre
135
void opencl_find_rst_legendre(void *rst, void *pt_x, void *pt_y, void *pt_z, void *x_hat, void *y_hat, void *z_hat, void *resx, void *resy, void *resz, int *lx, void *el_ids, int *n_pt, real *tol, void *conv_pts)
void cuda_find_rst_legendre(void *rst, void *pt_x, void *pt_y, void *pt_z, void *x_hat, void *y_hat, void *z_hat, void *resx, void *resy, void *resz, int *lx, void *el_ids, int *n_pt, real *tol, void *conv_pts)
subroutine, public device_find_rst_legendre(rst_d, pt_x_d, pt_y_d, pt_z_d, x_hat_d, y_hat_d, z_hat_d, resx_d, resy_d, resz_d, lx, el_ids_d, n_pts, tol, conv_pts_d)
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