xref: /libCEED/include/ceed/jit-source/hip/hip-ref-basis-tensor.h (revision 3c9d155ad22979ea84178d1d41eea0c31217327b)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Internal header for HIP tensor product basis
10 #ifndef CEED_HIP_REF_BASIS_TENSOR_H
11 #define CEED_HIP_REF_BASIS_TENSOR_H
12 
13 #include <ceed.h>
14 
15 //------------------------------------------------------------------------------
16 // Tensor Basis Kernels
17 //------------------------------------------------------------------------------
18 
19 //------------------------------------------------------------------------------
20 // Interp
21 //------------------------------------------------------------------------------
22 extern "C" __global__ void Interp(const CeedInt num_elem, const CeedInt transpose, const CeedScalar *__restrict__ interp_1d,
23                                   const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) {
24   const CeedInt i = threadIdx.x;
25 
26   __shared__ CeedScalar s_mem[BASIS_Q_1D * BASIS_P_1D + 2 * BASIS_BUF_LEN];
27   CeedScalar           *s_interp_1d = s_mem;
28   CeedScalar           *s_buffer_1  = s_mem + BASIS_Q_1D * BASIS_P_1D;
29   CeedScalar           *s_buffer_2  = s_buffer_1 + BASIS_BUF_LEN;
30   for (CeedInt k = i; k < BASIS_Q_1D * BASIS_P_1D; k += blockDim.x) {
31     s_interp_1d[k] = interp_1d[k];
32   }
33 
34   const CeedInt P             = transpose ? BASIS_Q_1D : BASIS_P_1D;
35   const CeedInt Q             = transpose ? BASIS_P_1D : BASIS_Q_1D;
36   const CeedInt stride_0      = transpose ? 1 : BASIS_P_1D;
37   const CeedInt stride_1      = transpose ? BASIS_P_1D : 1;
38   const CeedInt u_stride      = transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES;
39   const CeedInt v_stride      = transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS;
40   const CeedInt u_comp_stride = num_elem * (transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES);
41   const CeedInt v_comp_stride = num_elem * (transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS);
42   const CeedInt u_size        = transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES;
43 
44   // Apply basis element by element
45   for (CeedInt elem = blockIdx.x; elem < num_elem; elem += gridDim.x) {
46     for (CeedInt comp = 0; comp < BASIS_NUM_COMP; comp++) {
47       const CeedScalar *cur_u = u + elem * u_stride + comp * u_comp_stride;
48       CeedScalar       *cur_v = v + elem * v_stride + comp * v_comp_stride;
49       for (CeedInt k = i; k < u_size; k += blockDim.x) {
50         s_buffer_1[k] = cur_u[k];
51       }
52       CeedInt pre  = u_size;
53       CeedInt post = 1;
54       for (CeedInt d = 0; d < BASIS_DIM; d++) {
55         __syncthreads();
56         // Update buffers used
57         pre /= P;
58         const CeedScalar *in  = d % 2 ? s_buffer_2 : s_buffer_1;
59         CeedScalar       *out = d == BASIS_DIM - 1 ? cur_v : (d % 2 ? s_buffer_1 : s_buffer_2);
60 
61         // Contract along middle index
62         const CeedInt writeLen = pre * post * Q;
63         for (CeedInt k = i; k < writeLen; k += blockDim.x) {
64           const CeedInt c = k % post;
65           const CeedInt j = (k / post) % Q;
66           const CeedInt a = k / (post * Q);
67 
68           CeedScalar vk = 0;
69           for (CeedInt b = 0; b < P; b++) vk += s_interp_1d[j * stride_0 + b * stride_1] * in[(a * P + b) * post + c];
70 
71           out[k] = vk;
72         }
73 
74         post *= Q;
75       }
76     }
77   }
78 }
79 
80 //------------------------------------------------------------------------------
81 // Grad
82 //------------------------------------------------------------------------------
83 extern "C" __global__ void Grad(const CeedInt num_elem, const CeedInt transpose, const CeedScalar *__restrict__ interp_1d,
84                                 const CeedScalar *__restrict__ grad_1d, const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) {
85   const CeedInt i = threadIdx.x;
86 
87   __shared__ CeedScalar s_mem[2 * (BASIS_Q_1D * BASIS_P_1D + BASIS_BUF_LEN)];
88   CeedScalar           *s_interp_1d = s_mem;
89   CeedScalar           *s_grad_1d   = s_interp_1d + BASIS_Q_1D * BASIS_P_1D;
90   CeedScalar           *s_buffer_1  = s_grad_1d + BASIS_Q_1D * BASIS_P_1D;
91   CeedScalar           *s_buffer_2  = s_buffer_1 + BASIS_BUF_LEN;
92   for (CeedInt k = i; k < BASIS_Q_1D * BASIS_P_1D; k += blockDim.x) {
93     s_interp_1d[k] = interp_1d[k];
94     s_grad_1d[k]   = grad_1d[k];
95   }
96 
97   const CeedInt P             = transpose ? BASIS_Q_1D : BASIS_P_1D;
98   const CeedInt Q             = transpose ? BASIS_P_1D : BASIS_Q_1D;
99   const CeedInt stride_0      = transpose ? 1 : BASIS_P_1D;
100   const CeedInt stride_1      = transpose ? BASIS_P_1D : 1;
101   const CeedInt u_stride      = transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES;
102   const CeedInt v_stride      = transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS;
103   const CeedInt u_comp_stride = num_elem * (transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES);
104   const CeedInt v_comp_stride = num_elem * (transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS);
105   const CeedInt u_dim_stride  = transpose ? num_elem * BASIS_NUM_QPTS * BASIS_NUM_COMP : 0;
106   const CeedInt v_dim_stride  = transpose ? 0 : num_elem * BASIS_NUM_QPTS * BASIS_NUM_COMP;
107 
108   // Apply basis element by element
109   for (CeedInt elem = blockIdx.x; elem < num_elem; elem += gridDim.x) {
110     for (CeedInt comp = 0; comp < BASIS_NUM_COMP; comp++) {
111       // dim*dim contractions for grad
112       for (CeedInt dim_1 = 0; dim_1 < BASIS_DIM; dim_1++) {
113         CeedInt           pre   = transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES;
114         CeedInt           post  = 1;
115         const CeedScalar *cur_u = u + elem * u_stride + dim_1 * u_dim_stride + comp * u_comp_stride;
116         CeedScalar       *cur_v = v + elem * v_stride + dim_1 * v_dim_stride + comp * v_comp_stride;
117         for (CeedInt dim_2 = 0; dim_2 < BASIS_DIM; dim_2++) {
118           __syncthreads();
119           // Update buffers used
120           pre /= P;
121           const CeedScalar *op  = dim_1 == dim_2 ? s_grad_1d : s_interp_1d;
122           const CeedScalar *in  = dim_2 == 0 ? cur_u : (dim_2 % 2 ? s_buffer_2 : s_buffer_1);
123           CeedScalar       *out = dim_2 == BASIS_DIM - 1 ? cur_v : (dim_2 % 2 ? s_buffer_1 : s_buffer_2);
124 
125           // Contract along middle index
126           const CeedInt writeLen = pre * post * Q;
127           for (CeedInt k = i; k < writeLen; k += blockDim.x) {
128             const CeedInt c   = k % post;
129             const CeedInt j   = (k / post) % Q;
130             const CeedInt a   = k / (post * Q);
131             CeedScalar    v_k = 0;
132             for (CeedInt b = 0; b < P; b++) v_k += op[j * stride_0 + b * stride_1] * in[(a * P + b) * post + c];
133 
134             if (transpose && dim_2 == BASIS_DIM - 1) out[k] += v_k;
135             else out[k] = v_k;
136           }
137 
138           post *= Q;
139         }
140       }
141     }
142   }
143 }
144 
145 //------------------------------------------------------------------------------
146 // 1D quadrature weights
147 //------------------------------------------------------------------------------
148 __device__ void Weight1d(const CeedInt num_elem, const CeedScalar *q_weight_1d, CeedScalar *w) {
149   const CeedInt i = threadIdx.x;
150   if (i < BASIS_Q_1D) {
151     const size_t elem = blockIdx.x;
152     if (elem < num_elem) w[elem * BASIS_Q_1D + i] = q_weight_1d[i];
153   }
154 }
155 
156 //------------------------------------------------------------------------------
157 // 2D quadrature weights
158 //------------------------------------------------------------------------------
159 __device__ void Weight2d(const CeedInt num_elem, const CeedScalar *q_weight_1d, CeedScalar *w) {
160   const CeedInt i = threadIdx.x;
161   const CeedInt j = threadIdx.y;
162   if (i < BASIS_Q_1D && j < BASIS_Q_1D) {
163     const size_t elem = blockIdx.x;
164     if (elem < num_elem) {
165       const size_t ind = (elem * BASIS_Q_1D + j) * BASIS_Q_1D + i;
166       w[ind]           = q_weight_1d[i] * q_weight_1d[j];
167     }
168   }
169 }
170 
171 //------------------------------------------------------------------------------
172 // 3D quadrature weights
173 //------------------------------------------------------------------------------
174 __device__ void Weight3d(const CeedInt num_elem, const CeedScalar *q_weight_1d, CeedScalar *w) {
175   const CeedInt i = threadIdx.x;
176   const CeedInt j = threadIdx.y;
177   if (i < BASIS_Q_1D && j < BASIS_Q_1D) {
178     const size_t elem = blockIdx.x;
179     if (elem < num_elem) {
180       for (CeedInt k = 0; k < BASIS_Q_1D; k++) {
181         const size_t ind = ((elem * BASIS_Q_1D + k) * BASIS_Q_1D + j) * BASIS_Q_1D + i;
182         w[ind]           = q_weight_1d[i] * q_weight_1d[j] * q_weight_1d[k];
183       }
184     }
185   }
186 }
187 
188 //------------------------------------------------------------------------------
189 // Quadrature weights
190 //------------------------------------------------------------------------------
191 extern "C" __global__ void Weight(const CeedInt num_elem, const CeedScalar *__restrict__ q_weight_1d, CeedScalar *__restrict__ v) {
192   if (BASIS_DIM == 1) Weight1d(num_elem, q_weight_1d, v);
193   else if (BASIS_DIM == 2) Weight2d(num_elem, q_weight_1d, v);
194   else if (BASIS_DIM == 3) Weight3d(num_elem, q_weight_1d, v);
195 }
196 
197 //------------------------------------------------------------------------------
198 
199 #endif  // CEED_HIP_REF_BASIS_TENSOR_H
200