xref: /libCEED/include/ceed/jit-source/cuda/cuda-shared-basis-tensor.h (revision aa4002ad1a9f7f3442c9f6afb79353f990ebef22)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
29e201c85SYohann // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
39e201c85SYohann //
49e201c85SYohann // SPDX-License-Identifier: BSD-2-Clause
59e201c85SYohann //
69e201c85SYohann // This file is part of CEED:  http://github.com/ceed
79e201c85SYohann 
89e201c85SYohann /// @file
99e201c85SYohann /// Internal header for CUDA shared memory tensor product basis
10c0b5abf0SJeremy L Thompson #include <ceed/types.h>
112b730f8bSJeremy L Thompson 
129e201c85SYohann #include "cuda-shared-basis-read-write-templates.h"
139e201c85SYohann #include "cuda-shared-basis-tensor-templates.h"
149e201c85SYohann 
159e201c85SYohann //------------------------------------------------------------------------------
169e201c85SYohann // Interp kernel by dim
179e201c85SYohann //------------------------------------------------------------------------------
182b730f8bSJeremy L Thompson extern "C" __global__ void Interp(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) {
199e201c85SYohann   extern __shared__ CeedScalar slice[];
209e201c85SYohann 
219e201c85SYohann   SharedData_Cuda data;
229e201c85SYohann   data.t_id_x = threadIdx.x;
239e201c85SYohann   data.t_id_y = threadIdx.y;
249e201c85SYohann   data.t_id_z = threadIdx.z;
259e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
269e201c85SYohann   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
279e201c85SYohann 
289e201c85SYohann   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
299e201c85SYohann   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
309e201c85SYohann 
31*aa4002adSJeremy L Thompson   // load interp_1d into shared memory
32*aa4002adSJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
33*aa4002adSJeremy L Thompson   LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B);
34*aa4002adSJeremy L Thompson   __syncthreads();
35*aa4002adSJeremy L Thompson 
36*aa4002adSJeremy L Thompson   // Apply basis element by element
379e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
389e201c85SYohann     if (BASIS_DIM == 1) {
399e201c85SYohann       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, d_U, r_U);
40*aa4002adSJeremy L Thompson       Interp1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
419e201c85SYohann       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_V, d_V);
429e201c85SYohann     } else if (BASIS_DIM == 2) {
439e201c85SYohann       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);
44*aa4002adSJeremy L Thompson       InterpTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
459e201c85SYohann       WriteElementStrided2d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * num_elem, BASIS_Q_1D * BASIS_Q_1D, r_V, d_V);
469e201c85SYohann     } else if (BASIS_DIM == 3) {
472b730f8bSJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
482b730f8bSJeremy L Thompson                                                        BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, d_U, r_U);
49*aa4002adSJeremy L Thompson       InterpTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
502b730f8bSJeremy L Thompson       WriteElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem,
512b730f8bSJeremy L Thompson                                                         BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, r_V, d_V);
529e201c85SYohann     }
539e201c85SYohann   }
549e201c85SYohann }
559e201c85SYohann 
562b730f8bSJeremy L Thompson extern "C" __global__ void InterpTranspose(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *__restrict__ d_U,
579e201c85SYohann                                            CeedScalar *__restrict__ d_V) {
589e201c85SYohann   extern __shared__ CeedScalar slice[];
599e201c85SYohann 
609e201c85SYohann   SharedData_Cuda data;
619e201c85SYohann   data.t_id_x = threadIdx.x;
629e201c85SYohann   data.t_id_y = threadIdx.y;
639e201c85SYohann   data.t_id_z = threadIdx.z;
649e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
659e201c85SYohann   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
669e201c85SYohann 
679e201c85SYohann   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
689e201c85SYohann   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
699e201c85SYohann 
70*aa4002adSJeremy L Thompson   // load interp_1d into shared memory
71*aa4002adSJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
72*aa4002adSJeremy L Thompson   LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B);
73*aa4002adSJeremy L Thompson   __syncthreads();
74*aa4002adSJeremy L Thompson 
75*aa4002adSJeremy L Thompson   // Apply basis element by element
769e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
779e201c85SYohann     if (BASIS_DIM == 1) {
789e201c85SYohann       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U);
79*aa4002adSJeremy L Thompson       InterpTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
809e201c85SYohann       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V);
819e201c85SYohann     } else if (BASIS_DIM == 2) {
829e201c85SYohann       ReadElementStrided2d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * num_elem, BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
83*aa4002adSJeremy L Thompson       InterpTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
849e201c85SYohann       WriteElementStrided2d<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);
859e201c85SYohann     } else if (BASIS_DIM == 3) {
862b730f8bSJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem,
872b730f8bSJeremy L Thompson                                                        BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
88*aa4002adSJeremy L Thompson       InterpTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
892b730f8bSJeremy L Thompson       WriteElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
902b730f8bSJeremy L Thompson                                                         BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
919e201c85SYohann     }
929e201c85SYohann   }
939e201c85SYohann }
949e201c85SYohann 
95db2becc9SJeremy L Thompson extern "C" __global__ void InterpTransposeAdd(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *__restrict__ d_U,
96db2becc9SJeremy L Thompson                                               CeedScalar *__restrict__ d_V) {
97db2becc9SJeremy L Thompson   extern __shared__ CeedScalar slice[];
98db2becc9SJeremy L Thompson 
99db2becc9SJeremy L Thompson   SharedData_Cuda data;
100db2becc9SJeremy L Thompson   data.t_id_x = threadIdx.x;
101db2becc9SJeremy L Thompson   data.t_id_y = threadIdx.y;
102db2becc9SJeremy L Thompson   data.t_id_z = threadIdx.z;
103db2becc9SJeremy L Thompson   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
104db2becc9SJeremy L Thompson   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
105db2becc9SJeremy L Thompson 
106db2becc9SJeremy L Thompson   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
107db2becc9SJeremy L Thompson   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
108db2becc9SJeremy L Thompson 
109*aa4002adSJeremy L Thompson   // load interp_1d into shared memory
110*aa4002adSJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
111*aa4002adSJeremy L Thompson   LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B);
112*aa4002adSJeremy L Thompson   __syncthreads();
113*aa4002adSJeremy L Thompson 
114*aa4002adSJeremy L Thompson   // Apply basis element by element
115db2becc9SJeremy L Thompson   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
116db2becc9SJeremy L Thompson     if (BASIS_DIM == 1) {
117db2becc9SJeremy L Thompson       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U);
118*aa4002adSJeremy L Thompson       InterpTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
119db2becc9SJeremy L Thompson       SumElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V);
120db2becc9SJeremy L Thompson     } else if (BASIS_DIM == 2) {
121db2becc9SJeremy L Thompson       ReadElementStrided2d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * num_elem, BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
122*aa4002adSJeremy L Thompson       InterpTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
123db2becc9SJeremy 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);
124db2becc9SJeremy L Thompson     } else if (BASIS_DIM == 3) {
125db2becc9SJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem,
126db2becc9SJeremy L Thompson                                                        BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
127*aa4002adSJeremy L Thompson       InterpTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
128db2becc9SJeremy L Thompson       SumElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
129db2becc9SJeremy L Thompson                                                       BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
130db2becc9SJeremy L Thompson     }
131db2becc9SJeremy L Thompson   }
132db2becc9SJeremy L Thompson }
133db2becc9SJeremy L Thompson 
1349e201c85SYohann //------------------------------------------------------------------------------
1359e201c85SYohann // Grad kernel by dim
1369e201c85SYohann //------------------------------------------------------------------------------
1372b730f8bSJeremy L Thompson extern "C" __global__ void Grad(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *c_G, const CeedScalar *__restrict__ d_U,
1389e201c85SYohann                                 CeedScalar *__restrict__ d_V) {
1399e201c85SYohann   extern __shared__ CeedScalar slice[];
1409e201c85SYohann 
1419e201c85SYohann   SharedData_Cuda data;
1429e201c85SYohann   data.t_id_x = threadIdx.x;
1439e201c85SYohann   data.t_id_y = threadIdx.y;
1449e201c85SYohann   data.t_id_z = threadIdx.z;
1459e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
1469e201c85SYohann   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
1479e201c85SYohann 
1489e201c85SYohann   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
1499e201c85SYohann   CeedScalar r_V[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
1509e201c85SYohann 
151*aa4002adSJeremy L Thompson   // load interp_1d and grad_1d into shared memory
152*aa4002adSJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
153*aa4002adSJeremy L Thompson   LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B);
154*aa4002adSJeremy L Thompson   __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)];
155*aa4002adSJeremy L Thompson   LoadMatrix<BASIS_Q_1D, BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D>(data, c_G, s_G);
156*aa4002adSJeremy L Thompson   __syncthreads();
157*aa4002adSJeremy L Thompson 
158*aa4002adSJeremy L Thompson   // Apply basis element by element
1599e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
1609e201c85SYohann     if (BASIS_DIM == 1) {
1619e201c85SYohann       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, d_U, r_U);
162*aa4002adSJeremy L Thompson       Grad1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
1639e201c85SYohann       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_V, d_V);
1649e201c85SYohann     } else if (BASIS_DIM == 2) {
1659e201c85SYohann       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);
166*aa4002adSJeremy L Thompson       GradTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
1672b730f8bSJeremy L Thompson       WriteElementStrided2d<BASIS_NUM_COMP * BASIS_DIM, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * num_elem, BASIS_Q_1D * BASIS_Q_1D, r_V,
1682b730f8bSJeremy L Thompson                                                                     d_V);
1699e201c85SYohann     } else if (BASIS_DIM == 3) {
1702b730f8bSJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
1712b730f8bSJeremy L Thompson                                                        BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, d_U, r_U);
172*aa4002adSJeremy L Thompson       if (BASIS_HAS_COLLOCATED_GRAD) GradTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
173*aa4002adSJeremy L Thompson       else GradTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
1742b730f8bSJeremy L Thompson       WriteElementStrided3d<BASIS_NUM_COMP * BASIS_DIM, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem,
1752b730f8bSJeremy L Thompson                                                                     BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, r_V, d_V);
1769e201c85SYohann     }
1779e201c85SYohann   }
1789e201c85SYohann }
1799e201c85SYohann 
1802b730f8bSJeremy L Thompson extern "C" __global__ void GradTranspose(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *c_G, const CeedScalar *__restrict__ d_U,
1819e201c85SYohann                                          CeedScalar *__restrict__ d_V) {
1829e201c85SYohann   extern __shared__ CeedScalar slice[];
1839e201c85SYohann 
1849e201c85SYohann   SharedData_Cuda data;
1859e201c85SYohann   data.t_id_x = threadIdx.x;
1869e201c85SYohann   data.t_id_y = threadIdx.y;
1879e201c85SYohann   data.t_id_z = threadIdx.z;
1889e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
1899e201c85SYohann   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
1909e201c85SYohann 
1919e201c85SYohann   CeedScalar r_U[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
1929e201c85SYohann   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
1939e201c85SYohann 
194*aa4002adSJeremy L Thompson   // load interp_1d and grad_1d into shared memory
195*aa4002adSJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
196*aa4002adSJeremy L Thompson   LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B);
197*aa4002adSJeremy L Thompson   __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)];
198*aa4002adSJeremy L Thompson   LoadMatrix<BASIS_Q_1D, BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D>(data, c_G, s_G);
199*aa4002adSJeremy L Thompson   __syncthreads();
200*aa4002adSJeremy L Thompson 
201*aa4002adSJeremy L Thompson   // Apply basis element by element
2029e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
2039e201c85SYohann     if (BASIS_DIM == 1) {
2049e201c85SYohann       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U);
205*aa4002adSJeremy L Thompson       GradTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
2069e201c85SYohann       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V);
2079e201c85SYohann     } else if (BASIS_DIM == 2) {
2082b730f8bSJeremy L Thompson       ReadElementStrided2d<BASIS_NUM_COMP * BASIS_DIM, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * num_elem, BASIS_Q_1D * BASIS_Q_1D, d_U,
2092b730f8bSJeremy L Thompson                                                                    r_U);
210*aa4002adSJeremy L Thompson       GradTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
2119e201c85SYohann       WriteElementStrided2d<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);
2129e201c85SYohann     } else if (BASIS_DIM == 3) {
2132b730f8bSJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP * BASIS_DIM, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem,
2142b730f8bSJeremy L Thompson                                                                    BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
215*aa4002adSJeremy L Thompson       if (BASIS_HAS_COLLOCATED_GRAD) GradTransposeTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
216*aa4002adSJeremy L Thompson       else GradTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
2172b730f8bSJeremy L Thompson       WriteElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
2182b730f8bSJeremy L Thompson                                                         BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
2199e201c85SYohann     }
2209e201c85SYohann   }
2219e201c85SYohann }
2229e201c85SYohann 
223db2becc9SJeremy L Thompson extern "C" __global__ void GradTransposeAdd(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *c_G, const CeedScalar *__restrict__ d_U,
224db2becc9SJeremy L Thompson                                             CeedScalar *__restrict__ d_V) {
225db2becc9SJeremy L Thompson   extern __shared__ CeedScalar slice[];
226db2becc9SJeremy L Thompson 
227db2becc9SJeremy L Thompson   SharedData_Cuda data;
228db2becc9SJeremy L Thompson   data.t_id_x = threadIdx.x;
229db2becc9SJeremy L Thompson   data.t_id_y = threadIdx.y;
230db2becc9SJeremy L Thompson   data.t_id_z = threadIdx.z;
231db2becc9SJeremy L Thompson   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
232db2becc9SJeremy L Thompson   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
233db2becc9SJeremy L Thompson 
234db2becc9SJeremy L Thompson   CeedScalar r_U[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
235db2becc9SJeremy L Thompson   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
236db2becc9SJeremy L Thompson 
237*aa4002adSJeremy L Thompson   // load interp_1d and grad_1d into shared memory
238*aa4002adSJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
239*aa4002adSJeremy L Thompson   LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B);
240*aa4002adSJeremy L Thompson   __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)];
241*aa4002adSJeremy L Thompson   LoadMatrix<BASIS_Q_1D, BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D>(data, c_G, s_G);
242*aa4002adSJeremy L Thompson   __syncthreads();
243*aa4002adSJeremy L Thompson 
244*aa4002adSJeremy L Thompson   // Apply basis element by element
245db2becc9SJeremy L Thompson   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
246db2becc9SJeremy L Thompson     if (BASIS_DIM == 1) {
247db2becc9SJeremy L Thompson       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U);
248*aa4002adSJeremy L Thompson       GradTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
249db2becc9SJeremy L Thompson       SumElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V);
250db2becc9SJeremy L Thompson     } else if (BASIS_DIM == 2) {
251db2becc9SJeremy L Thompson       ReadElementStrided2d<BASIS_NUM_COMP * BASIS_DIM, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * num_elem, BASIS_Q_1D * BASIS_Q_1D, d_U,
252db2becc9SJeremy L Thompson                                                                    r_U);
253*aa4002adSJeremy L Thompson       GradTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
254db2becc9SJeremy 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);
255db2becc9SJeremy L Thompson     } else if (BASIS_DIM == 3) {
256db2becc9SJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP * BASIS_DIM, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem,
257db2becc9SJeremy L Thompson                                                                    BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
258*aa4002adSJeremy L Thompson       if (BASIS_HAS_COLLOCATED_GRAD) GradTransposeTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
259*aa4002adSJeremy L Thompson       else GradTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
260db2becc9SJeremy L Thompson       SumElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
261db2becc9SJeremy L Thompson                                                       BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
262db2becc9SJeremy L Thompson     }
263db2becc9SJeremy L Thompson   }
264db2becc9SJeremy L Thompson }
265db2becc9SJeremy L Thompson 
2669e201c85SYohann //------------------------------------------------------------------------------
2679e201c85SYohann // Weight kernels by dim
2689e201c85SYohann //------------------------------------------------------------------------------
2692b730f8bSJeremy L Thompson extern "C" __global__ void Weight(const CeedInt num_elem, const CeedScalar *__restrict__ q_weight_1d, CeedScalar *__restrict__ d_W) {
2709e201c85SYohann   extern __shared__ CeedScalar slice[];
2719e201c85SYohann 
2729e201c85SYohann   SharedData_Cuda data;
2739e201c85SYohann   data.t_id_x = threadIdx.x;
2749e201c85SYohann   data.t_id_y = threadIdx.y;
2759e201c85SYohann   data.t_id_z = threadIdx.z;
2769e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
2779e201c85SYohann   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
2789e201c85SYohann 
2799e201c85SYohann   CeedScalar r_W[BASIS_DIM > 2 ? BASIS_Q_1D : 1];
2809e201c85SYohann 
281*aa4002adSJeremy L Thompson   // Apply basis element by element
2829e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
2839e201c85SYohann     if (BASIS_DIM == 1) {
2849e201c85SYohann       Weight1d<BASIS_Q_1D>(data, q_weight_1d, r_W);
2859e201c85SYohann       WriteElementStrided1d<1, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_W, d_W);
2869e201c85SYohann     } else if (BASIS_DIM == 2) {
2879e201c85SYohann       WeightTensor2d<BASIS_Q_1D>(data, q_weight_1d, r_W);
2889e201c85SYohann       WriteElementStrided2d<1, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * num_elem, BASIS_Q_1D * BASIS_Q_1D, r_W, d_W);
2899e201c85SYohann     } else if (BASIS_DIM == 3) {
2909e201c85SYohann       WeightTensor3d<BASIS_Q_1D>(data, q_weight_1d, r_W);
2912b730f8bSJeremy L Thompson       WriteElementStrided3d<1, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, r_W,
2922b730f8bSJeremy L Thompson                                            d_W);
2939e201c85SYohann     }
2949e201c85SYohann   }
2959e201c85SYohann }
296