xref: /libCEED/include/ceed/jit-source/cuda/cuda-shared-basis-tensor-at-points.h (revision aa4002ad1a9f7f3442c9f6afb79353f990ebef22)
19e1d4b82SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
29e1d4b82SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
39e1d4b82SJeremy L Thompson //
49e1d4b82SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
59e1d4b82SJeremy L Thompson //
69e1d4b82SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
79e1d4b82SJeremy L Thompson 
89e1d4b82SJeremy L Thompson /// @file
99e1d4b82SJeremy L Thompson /// Internal header for CUDA tensor product basis with AtPoints evaluation
109e1d4b82SJeremy L Thompson #include <ceed/types.h>
119e1d4b82SJeremy L Thompson 
129e1d4b82SJeremy L Thompson #include "cuda-shared-basis-read-write-templates.h"
139e1d4b82SJeremy L Thompson #include "cuda-shared-basis-tensor-at-points-templates.h"
149e1d4b82SJeremy L Thompson #include "cuda-shared-basis-tensor-templates.h"
159e1d4b82SJeremy L Thompson 
169e1d4b82SJeremy L Thompson //------------------------------------------------------------------------------
179e1d4b82SJeremy L Thompson // Tensor Basis Kernels AtPoints
189e1d4b82SJeremy L Thompson //------------------------------------------------------------------------------
199e1d4b82SJeremy L Thompson 
209e1d4b82SJeremy L Thompson //------------------------------------------------------------------------------
219e1d4b82SJeremy L Thompson // Interp
229e1d4b82SJeremy L Thompson //------------------------------------------------------------------------------
239e1d4b82SJeremy L Thompson extern "C" __global__ void InterpAtPoints(const CeedInt num_elem, const CeedScalar *__restrict__ c_B, const CeedInt *__restrict__ points_per_elem,
249e1d4b82SJeremy L Thompson                                           const CeedScalar *__restrict__ d_X, const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) {
259e1d4b82SJeremy L Thompson   extern __shared__ CeedScalar slice[];
269e1d4b82SJeremy L Thompson 
279e1d4b82SJeremy L Thompson   SharedData_Cuda data;
289e1d4b82SJeremy L Thompson   data.t_id_x = threadIdx.x;
299e1d4b82SJeremy L Thompson   data.t_id_y = threadIdx.y;
309e1d4b82SJeremy L Thompson   data.t_id_z = threadIdx.z;
319e1d4b82SJeremy L Thompson   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
329e1d4b82SJeremy L Thompson   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
339e1d4b82SJeremy L Thompson 
34b6a2eb79SJeremy L Thompson   CeedScalar r_X[BASIS_DIM];
359e1d4b82SJeremy L Thompson   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
369e1d4b82SJeremy L Thompson   CeedScalar r_C[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
379e1d4b82SJeremy L Thompson   CeedScalar r_V[BASIS_NUM_COMP];
389e1d4b82SJeremy L Thompson 
39*aa4002adSJeremy L Thompson   // load interp_1d into shared memory
40*aa4002adSJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
41*aa4002adSJeremy L Thompson   LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B);
42*aa4002adSJeremy L Thompson   __syncthreads();
43*aa4002adSJeremy L Thompson 
449e1d4b82SJeremy L Thompson   // Apply basis element by element
459e1d4b82SJeremy L Thompson   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
469e1d4b82SJeremy L Thompson     // Map to coefficients
479e1d4b82SJeremy L Thompson     if (BASIS_DIM == 1) {
489e1d4b82SJeremy L Thompson       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, d_U, r_U);
49*aa4002adSJeremy L Thompson       Interp1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_C);
509e1d4b82SJeremy L Thompson     } else if (BASIS_DIM == 2) {
519e1d4b82SJeremy L Thompson       ReadElementStrided2d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * num_elem, BASIS_P_1D * BASIS_P_1D, d_U, r_U);
52*aa4002adSJeremy L Thompson       InterpTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_C);
539e1d4b82SJeremy L Thompson     } else if (BASIS_DIM == 3) {
549e1d4b82SJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
559e1d4b82SJeremy L Thompson                                                        BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, d_U, r_U);
56*aa4002adSJeremy L Thompson       InterpTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_C);
579e1d4b82SJeremy L Thompson     }
589e1d4b82SJeremy L Thompson 
599e1d4b82SJeremy L Thompson     // Map to points
60b6a2eb79SJeremy L Thompson     const CeedInt point_loop_bound = (blockDim.x * blockDim.y) * ceil(1.0 * BASIS_NUM_PTS / (blockDim.x * blockDim.y));
61b6a2eb79SJeremy L Thompson 
62b6a2eb79SJeremy L Thompson     for (CeedInt i = threadIdx.x + threadIdx.y * blockDim.x; i < point_loop_bound; i += blockDim.x * blockDim.y) {
63b6a2eb79SJeremy L Thompson       const CeedInt p = i % BASIS_NUM_PTS;
64b6a2eb79SJeremy L Thompson 
65b6a2eb79SJeremy L Thompson       ReadPoint<BASIS_DIM, BASIS_NUM_PTS>(data, elem, p, BASIS_NUM_PTS, 1, num_elem * BASIS_NUM_PTS, BASIS_NUM_PTS, d_X, r_X);
66b6a2eb79SJeremy L Thompson       if (BASIS_DIM == 1) {
67b6a2eb79SJeremy L Thompson         InterpAtPoints1d<BASIS_NUM_COMP, BASIS_NUM_PTS, BASIS_Q_1D>(data, i, r_C, r_X, r_V);
68b6a2eb79SJeremy L Thompson       } else if (BASIS_DIM == 2) {
69b6a2eb79SJeremy L Thompson         InterpAtPoints2d<BASIS_NUM_COMP, BASIS_NUM_PTS, BASIS_Q_1D>(data, i, r_C, r_X, r_V);
70b6a2eb79SJeremy L Thompson       } else if (BASIS_DIM == 3) {
71b6a2eb79SJeremy L Thompson         InterpAtPoints3d<BASIS_NUM_COMP, BASIS_NUM_PTS, BASIS_Q_1D>(data, i, r_C, r_X, r_V);
72b6a2eb79SJeremy L Thompson       }
73b6a2eb79SJeremy L Thompson       WritePoint<BASIS_NUM_COMP, BASIS_NUM_PTS>(data, elem, p, BASIS_NUM_PTS, 1, num_elem * BASIS_NUM_PTS, BASIS_NUM_PTS, r_V, d_V);
74b6a2eb79SJeremy L Thompson     }
759e1d4b82SJeremy L Thompson   }
769e1d4b82SJeremy L Thompson }
779e1d4b82SJeremy L Thompson 
789e1d4b82SJeremy L Thompson extern "C" __global__ void InterpTransposeAtPoints(const CeedInt num_elem, const CeedScalar *__restrict__ c_B,
799e1d4b82SJeremy L Thompson                                                    const CeedInt *__restrict__ points_per_elem, const CeedScalar *__restrict__ d_X,
809e1d4b82SJeremy L Thompson                                                    const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) {
819e1d4b82SJeremy L Thompson   extern __shared__ CeedScalar slice[];
829e1d4b82SJeremy L Thompson 
839e1d4b82SJeremy L Thompson   SharedData_Cuda data;
849e1d4b82SJeremy L Thompson   data.t_id_x = threadIdx.x;
859e1d4b82SJeremy L Thompson   data.t_id_y = threadIdx.y;
869e1d4b82SJeremy L Thompson   data.t_id_z = threadIdx.z;
879e1d4b82SJeremy L Thompson   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
889e1d4b82SJeremy L Thompson   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
899e1d4b82SJeremy L Thompson 
90b6a2eb79SJeremy L Thompson   CeedScalar r_X[BASIS_DIM];
919e1d4b82SJeremy L Thompson   CeedScalar r_U[BASIS_NUM_COMP];
929e1d4b82SJeremy L Thompson   CeedScalar r_C[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
939e1d4b82SJeremy L Thompson   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
949e1d4b82SJeremy L Thompson 
95*aa4002adSJeremy L Thompson   // load interp_1d into shared memory
96*aa4002adSJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
97*aa4002adSJeremy L Thompson   LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B);
98*aa4002adSJeremy L Thompson   __syncthreads();
99*aa4002adSJeremy L Thompson 
1009e1d4b82SJeremy L Thompson   // Apply basis element by element
1019e1d4b82SJeremy L Thompson   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
102b6a2eb79SJeremy L Thompson     // Clear register
103b6a2eb79SJeremy L Thompson     for (CeedInt i = 0; i < BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1); i++) r_C[i] = 0.0;
104b6a2eb79SJeremy L Thompson 
1059e1d4b82SJeremy L Thompson     // Map from points
106b6a2eb79SJeremy L Thompson     const CeedInt point_loop_bound = (blockDim.x * blockDim.y) * ceil(1.0 * BASIS_NUM_PTS / (blockDim.x * blockDim.y));
107b6a2eb79SJeremy L Thompson 
108b6a2eb79SJeremy L Thompson     for (CeedInt i = threadIdx.x + threadIdx.y * blockDim.x; i < point_loop_bound; i += blockDim.x * blockDim.y) {
109b6a2eb79SJeremy L Thompson       const CeedInt p = i % BASIS_NUM_PTS;
110b6a2eb79SJeremy L Thompson 
111b6a2eb79SJeremy L Thompson       ReadPoint<BASIS_DIM, BASIS_NUM_PTS>(data, elem, p, BASIS_NUM_PTS, 1, num_elem * BASIS_NUM_PTS, BASIS_NUM_PTS, d_X, r_X);
112b6a2eb79SJeremy L Thompson       ReadPoint<BASIS_NUM_COMP, BASIS_NUM_PTS>(data, elem, i, points_per_elem[elem], 1, num_elem * BASIS_NUM_PTS, BASIS_NUM_PTS, d_U, r_U);
113b6a2eb79SJeremy L Thompson       if (BASIS_DIM == 1) {
114b6a2eb79SJeremy L Thompson         InterpTransposeAtPoints1d<BASIS_NUM_COMP, BASIS_NUM_PTS, BASIS_Q_1D>(data, i, r_U, r_X, r_C);
115b6a2eb79SJeremy L Thompson       } else if (BASIS_DIM == 2) {
116b6a2eb79SJeremy L Thompson         InterpTransposeAtPoints2d<BASIS_NUM_COMP, BASIS_NUM_PTS, BASIS_Q_1D>(data, i, r_U, r_X, r_C);
117b6a2eb79SJeremy L Thompson       } else if (BASIS_DIM == 3) {
118b6a2eb79SJeremy L Thompson         InterpTransposeAtPoints3d<BASIS_NUM_COMP, BASIS_NUM_PTS, BASIS_Q_1D>(data, i, r_U, r_X, r_C);
119b6a2eb79SJeremy L Thompson       }
120b6a2eb79SJeremy L Thompson     }
121b6a2eb79SJeremy L Thompson     __syncthreads();
1229e1d4b82SJeremy L Thompson 
1239e1d4b82SJeremy L Thompson     // Map from coefficients
1249e1d4b82SJeremy L Thompson     if (BASIS_DIM == 1) {
125*aa4002adSJeremy L Thompson       InterpTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_C, s_B, r_V);
1269e1d4b82SJeremy L Thompson       SumElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V);
1279e1d4b82SJeremy L Thompson     } else if (BASIS_DIM == 2) {
128*aa4002adSJeremy L Thompson       InterpTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_C, s_B, r_V);
1299e1d4b82SJeremy L Thompson       SumElementStrided2d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * num_elem, BASIS_P_1D * BASIS_P_1D, r_V, d_V);
1309e1d4b82SJeremy L Thompson     } else if (BASIS_DIM == 3) {
131*aa4002adSJeremy L Thompson       InterpTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_C, s_B, r_V);
1329e1d4b82SJeremy L Thompson       SumElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
1339e1d4b82SJeremy L Thompson                                                       BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
1349e1d4b82SJeremy L Thompson     }
1359e1d4b82SJeremy L Thompson   }
1369e1d4b82SJeremy L Thompson }
1379e1d4b82SJeremy L Thompson 
1389e1d4b82SJeremy L Thompson //------------------------------------------------------------------------------
1399e1d4b82SJeremy L Thompson // Grad
1409e1d4b82SJeremy L Thompson //------------------------------------------------------------------------------
1419e1d4b82SJeremy L Thompson extern "C" __global__ void GradAtPoints(const CeedInt num_elem, const CeedScalar *__restrict__ c_B, const CeedInt *__restrict__ points_per_elem,
1429e1d4b82SJeremy L Thompson                                         const CeedScalar *__restrict__ d_X, const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) {
1439e1d4b82SJeremy L Thompson   extern __shared__ CeedScalar slice[];
1449e1d4b82SJeremy L Thompson 
1459e1d4b82SJeremy L Thompson   SharedData_Cuda data;
1469e1d4b82SJeremy L Thompson   data.t_id_x = threadIdx.x;
1479e1d4b82SJeremy L Thompson   data.t_id_y = threadIdx.y;
1489e1d4b82SJeremy L Thompson   data.t_id_z = threadIdx.z;
1499e1d4b82SJeremy L Thompson   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
1509e1d4b82SJeremy L Thompson   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
1519e1d4b82SJeremy L Thompson 
152b6a2eb79SJeremy L Thompson   CeedScalar r_X[BASIS_DIM];
1539e1d4b82SJeremy L Thompson   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
1549e1d4b82SJeremy L Thompson   CeedScalar r_C[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
1559e1d4b82SJeremy L Thompson   CeedScalar r_V[BASIS_NUM_COMP * BASIS_DIM];
1569e1d4b82SJeremy L Thompson 
157*aa4002adSJeremy L Thompson   // load interp_1d into shared memory
158*aa4002adSJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
159*aa4002adSJeremy L Thompson   LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B);
160*aa4002adSJeremy L Thompson   __syncthreads();
161*aa4002adSJeremy L Thompson 
1629e1d4b82SJeremy L Thompson   // Apply basis element by element
1639e1d4b82SJeremy L Thompson   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
1649e1d4b82SJeremy L Thompson     // Map to coefficients
1659e1d4b82SJeremy L Thompson     if (BASIS_DIM == 1) {
1669e1d4b82SJeremy L Thompson       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, d_U, r_U);
167*aa4002adSJeremy L Thompson       Interp1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_C);
1689e1d4b82SJeremy L Thompson     } else if (BASIS_DIM == 2) {
1699e1d4b82SJeremy L Thompson       ReadElementStrided2d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * num_elem, BASIS_P_1D * BASIS_P_1D, d_U, r_U);
170*aa4002adSJeremy L Thompson       InterpTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_C);
1719e1d4b82SJeremy L Thompson     } else if (BASIS_DIM == 3) {
1729e1d4b82SJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
1739e1d4b82SJeremy L Thompson                                                        BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, d_U, r_U);
174*aa4002adSJeremy L Thompson       InterpTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_C);
1759e1d4b82SJeremy L Thompson     }
1769e1d4b82SJeremy L Thompson 
1779e1d4b82SJeremy L Thompson     // Map to points
178b6a2eb79SJeremy L Thompson     const CeedInt point_loop_bound = (blockDim.x * blockDim.y) * ceil(1.0 * BASIS_NUM_PTS / (blockDim.x * blockDim.y));
179b6a2eb79SJeremy L Thompson 
180b6a2eb79SJeremy L Thompson     for (CeedInt i = threadIdx.x + threadIdx.y * blockDim.x; i < point_loop_bound; i += blockDim.x * blockDim.y) {
181b6a2eb79SJeremy L Thompson       const CeedInt p = i % BASIS_NUM_PTS;
182b6a2eb79SJeremy L Thompson 
183b6a2eb79SJeremy L Thompson       ReadPoint<BASIS_DIM, BASIS_NUM_PTS>(data, elem, p, BASIS_NUM_PTS, 1, num_elem * BASIS_NUM_PTS, BASIS_NUM_PTS, d_X, r_X);
184b6a2eb79SJeremy L Thompson       if (BASIS_DIM == 1) {
185b6a2eb79SJeremy L Thompson         GradAtPoints1d<BASIS_NUM_COMP, BASIS_NUM_PTS, BASIS_Q_1D>(data, i, r_C, r_X, r_V);
186b6a2eb79SJeremy L Thompson       } else if (BASIS_DIM == 2) {
187b6a2eb79SJeremy L Thompson         GradAtPoints2d<BASIS_NUM_COMP, BASIS_NUM_PTS, BASIS_Q_1D>(data, i, r_C, r_X, r_V);
188b6a2eb79SJeremy L Thompson       } else if (BASIS_DIM == 3) {
189b6a2eb79SJeremy L Thompson         GradAtPoints3d<BASIS_NUM_COMP, BASIS_NUM_PTS, BASIS_Q_1D>(data, i, r_C, r_X, r_V);
190b6a2eb79SJeremy L Thompson       }
191b6a2eb79SJeremy L Thompson       WritePoint<BASIS_NUM_COMP * BASIS_DIM, BASIS_NUM_PTS>(data, elem, p, BASIS_NUM_PTS, 1, num_elem * BASIS_NUM_PTS, BASIS_NUM_PTS, r_V, d_V);
192b6a2eb79SJeremy L Thompson     }
1939e1d4b82SJeremy L Thompson   }
1949e1d4b82SJeremy L Thompson }
1959e1d4b82SJeremy L Thompson 
1969e1d4b82SJeremy L Thompson extern "C" __global__ void GradTransposeAtPoints(const CeedInt num_elem, const CeedScalar *__restrict__ c_B,
1979e1d4b82SJeremy L Thompson                                                  const CeedInt *__restrict__ points_per_elem, const CeedScalar *__restrict__ d_X,
1989e1d4b82SJeremy L Thompson                                                  const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) {
1999e1d4b82SJeremy L Thompson   extern __shared__ CeedScalar slice[];
2009e1d4b82SJeremy L Thompson 
2019e1d4b82SJeremy L Thompson   SharedData_Cuda data;
2029e1d4b82SJeremy L Thompson   data.t_id_x = threadIdx.x;
2039e1d4b82SJeremy L Thompson   data.t_id_y = threadIdx.y;
2049e1d4b82SJeremy L Thompson   data.t_id_z = threadIdx.z;
2059e1d4b82SJeremy L Thompson   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
2069e1d4b82SJeremy L Thompson   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
2079e1d4b82SJeremy L Thompson 
208b6a2eb79SJeremy L Thompson   CeedScalar r_X[BASIS_DIM];
2099e1d4b82SJeremy L Thompson   CeedScalar r_U[BASIS_NUM_COMP * BASIS_DIM];
2109e1d4b82SJeremy L Thompson   CeedScalar r_C[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
2119e1d4b82SJeremy L Thompson   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
2129e1d4b82SJeremy L Thompson 
213*aa4002adSJeremy L Thompson   // load interp_1d into shared memory
214*aa4002adSJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
215*aa4002adSJeremy L Thompson   LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B);
216*aa4002adSJeremy L Thompson   __syncthreads();
217*aa4002adSJeremy L Thompson 
2189e1d4b82SJeremy L Thompson   // Apply basis element by element
2199e1d4b82SJeremy L Thompson   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
220b6a2eb79SJeremy L Thompson     // Clear register
221b6a2eb79SJeremy L Thompson     for (CeedInt i = 0; i < BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1); i++) r_C[i] = 0.0;
222b6a2eb79SJeremy L Thompson 
2239e1d4b82SJeremy L Thompson     // Map from points
224b6a2eb79SJeremy L Thompson     const CeedInt point_loop_bound = (blockDim.x * blockDim.y) * ceil(1.0 * BASIS_NUM_PTS / (blockDim.x * blockDim.y));
225b6a2eb79SJeremy L Thompson 
226b6a2eb79SJeremy L Thompson     for (CeedInt i = threadIdx.x + threadIdx.y * blockDim.x; i < point_loop_bound; i += blockDim.x * blockDim.y) {
227b6a2eb79SJeremy L Thompson       const CeedInt p = i % BASIS_NUM_PTS;
228b6a2eb79SJeremy L Thompson 
229b6a2eb79SJeremy L Thompson       ReadPoint<BASIS_DIM, BASIS_NUM_PTS>(data, elem, p, BASIS_NUM_PTS, 1, num_elem * BASIS_NUM_PTS, BASIS_NUM_PTS, d_X, r_X);
230b6a2eb79SJeremy L Thompson       ReadPoint<BASIS_NUM_COMP * BASIS_DIM, BASIS_NUM_PTS>(data, elem, i, points_per_elem[elem], 1, num_elem * BASIS_NUM_PTS, BASIS_NUM_PTS, d_U,
231b6a2eb79SJeremy L Thompson                                                            r_U);
232b6a2eb79SJeremy L Thompson       if (BASIS_DIM == 1) {
233b6a2eb79SJeremy L Thompson         GradTransposeAtPoints1d<BASIS_NUM_COMP, BASIS_NUM_PTS, BASIS_Q_1D>(data, i, r_U, r_X, r_C);
234b6a2eb79SJeremy L Thompson       } else if (BASIS_DIM == 2) {
235b6a2eb79SJeremy L Thompson         GradTransposeAtPoints2d<BASIS_NUM_COMP, BASIS_NUM_PTS, BASIS_Q_1D>(data, i, r_U, r_X, r_C);
236b6a2eb79SJeremy L Thompson       } else if (BASIS_DIM == 3) {
237b6a2eb79SJeremy L Thompson         GradTransposeAtPoints3d<BASIS_NUM_COMP, BASIS_NUM_PTS, BASIS_Q_1D>(data, i, r_U, r_X, r_C);
238b6a2eb79SJeremy L Thompson       }
239b6a2eb79SJeremy L Thompson     }
240b6a2eb79SJeremy L Thompson     __syncthreads();
2419e1d4b82SJeremy L Thompson 
2429e1d4b82SJeremy L Thompson     // Map from coefficients
2439e1d4b82SJeremy L Thompson     if (BASIS_DIM == 1) {
244*aa4002adSJeremy L Thompson       InterpTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_C, s_B, r_V);
2459e1d4b82SJeremy L Thompson       SumElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V);
2469e1d4b82SJeremy L Thompson     } else if (BASIS_DIM == 2) {
247*aa4002adSJeremy L Thompson       InterpTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_C, s_B, r_V);
2489e1d4b82SJeremy L Thompson       SumElementStrided2d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * num_elem, BASIS_P_1D * BASIS_P_1D, r_V, d_V);
2499e1d4b82SJeremy L Thompson     } else if (BASIS_DIM == 3) {
250*aa4002adSJeremy L Thompson       InterpTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_C, s_B, r_V);
2519e1d4b82SJeremy L Thompson       SumElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
2529e1d4b82SJeremy L Thompson                                                       BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
2539e1d4b82SJeremy L Thompson     }
2549e1d4b82SJeremy L Thompson   }
2559e1d4b82SJeremy L Thompson }
256