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