xref: /libCEED/include/ceed/jit-source/cuda/cuda-shared-basis-tensor.h (revision d4cc18453651bd0f94c1a2e078b2646a92dafdcc)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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 //------------------------------------------------------------------------------
1621292910SJeremy L Thompson // Interp kernels by dim
179e201c85SYohann //------------------------------------------------------------------------------
Interp(const CeedInt num_elem,const CeedScalar * c_B,const CeedScalar * __restrict__ d_U,CeedScalar * __restrict__ d_V)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;
2699421279SJeremy L Thompson   data.slice  = slice + data.t_id_z * BASIS_T_1D * (BASIS_DIM > 1 ? BASIS_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 
31aa4002adSJeremy L Thompson   // load interp_1d into shared memory
32aa4002adSJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
33aa4002adSJeremy L Thompson   LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B);
34aa4002adSJeremy L Thompson   __syncthreads();
35aa4002adSJeremy L Thompson 
36aa4002adSJeremy 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);
4099421279SJeremy L Thompson       Interp1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_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);
4499421279SJeremy L Thompson       InterpTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_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);
4999421279SJeremy L Thompson       InterpTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_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 
InterpCollocated(const CeedInt num_elem,const CeedScalar * c_B,const CeedScalar * __restrict__ d_U,CeedScalar * __restrict__ d_V)5621292910SJeremy L Thompson extern "C" __global__ void InterpCollocated(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *__restrict__ d_U,
5721292910SJeremy L Thompson                                             CeedScalar *__restrict__ d_V) {
5821292910SJeremy L Thompson   extern __shared__ CeedScalar slice[];
5921292910SJeremy L Thompson 
6021292910SJeremy L Thompson   SharedData_Cuda data;
6121292910SJeremy L Thompson   data.t_id_x = threadIdx.x;
6221292910SJeremy L Thompson   data.t_id_y = threadIdx.y;
6321292910SJeremy L Thompson   data.t_id_z = threadIdx.z;
6421292910SJeremy L Thompson   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
6521292910SJeremy L Thompson   data.slice  = slice + data.t_id_z * BASIS_T_1D * (BASIS_DIM > 1 ? BASIS_T_1D : 1);
6621292910SJeremy L Thompson 
6721292910SJeremy L Thompson   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
6821292910SJeremy L Thompson 
6921292910SJeremy L Thompson   // Apply basis element by element
7021292910SJeremy L Thompson   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
7121292910SJeremy L Thompson     if (BASIS_DIM == 1) {
7221292910SJeremy L Thompson       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, d_U, r_U);
7321292910SJeremy L Thompson       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_U, d_V);
7421292910SJeremy L Thompson     } else if (BASIS_DIM == 2) {
7521292910SJeremy 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);
7621292910SJeremy L Thompson       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_U, d_V);
7721292910SJeremy L Thompson     } else if (BASIS_DIM == 3) {
7821292910SJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
7921292910SJeremy L Thompson                                                        BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, d_U, r_U);
8021292910SJeremy L Thompson       WriteElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem,
8121292910SJeremy L Thompson                                                         BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, r_U, d_V);
8221292910SJeremy L Thompson     }
8321292910SJeremy L Thompson   }
8421292910SJeremy L Thompson }
8521292910SJeremy L Thompson 
InterpTranspose(const CeedInt num_elem,const CeedScalar * c_B,const CeedScalar * __restrict__ d_U,CeedScalar * __restrict__ d_V)862b730f8bSJeremy L Thompson extern "C" __global__ void InterpTranspose(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *__restrict__ d_U,
879e201c85SYohann                                            CeedScalar *__restrict__ d_V) {
889e201c85SYohann   extern __shared__ CeedScalar slice[];
899e201c85SYohann 
909e201c85SYohann   SharedData_Cuda data;
919e201c85SYohann   data.t_id_x = threadIdx.x;
929e201c85SYohann   data.t_id_y = threadIdx.y;
939e201c85SYohann   data.t_id_z = threadIdx.z;
949e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
9599421279SJeremy L Thompson   data.slice  = slice + data.t_id_z * BASIS_T_1D * (BASIS_DIM > 1 ? BASIS_T_1D : 1);
969e201c85SYohann 
979e201c85SYohann   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
989e201c85SYohann   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
999e201c85SYohann 
100aa4002adSJeremy L Thompson   // load interp_1d into shared memory
101aa4002adSJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
102aa4002adSJeremy L Thompson   LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B);
103aa4002adSJeremy L Thompson   __syncthreads();
104aa4002adSJeremy L Thompson 
105aa4002adSJeremy L Thompson   // Apply basis element by element
1069e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
1079e201c85SYohann     if (BASIS_DIM == 1) {
1089e201c85SYohann       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U);
10999421279SJeremy L Thompson       InterpTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, r_V);
1109e201c85SYohann       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V);
1119e201c85SYohann     } else if (BASIS_DIM == 2) {
1129e201c85SYohann       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);
11399421279SJeremy L Thompson       InterpTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, r_V);
1149e201c85SYohann       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);
1159e201c85SYohann     } else if (BASIS_DIM == 3) {
1162b730f8bSJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem,
1172b730f8bSJeremy L Thompson                                                        BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
11899421279SJeremy L Thompson       InterpTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, r_V);
1192b730f8bSJeremy L Thompson       WriteElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
1202b730f8bSJeremy L Thompson                                                         BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
1219e201c85SYohann     }
1229e201c85SYohann   }
1239e201c85SYohann }
1249e201c85SYohann 
InterpCollocatedTranspose(const CeedInt num_elem,const CeedScalar * c_B,const CeedScalar * __restrict__ d_U,CeedScalar * __restrict__ d_V)12521292910SJeremy L Thompson extern "C" __global__ void InterpCollocatedTranspose(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *__restrict__ d_U,
12621292910SJeremy L Thompson                                                      CeedScalar *__restrict__ d_V) {
12721292910SJeremy L Thompson   extern __shared__ CeedScalar slice[];
12821292910SJeremy L Thompson 
12921292910SJeremy L Thompson   SharedData_Cuda data;
13021292910SJeremy L Thompson   data.t_id_x = threadIdx.x;
13121292910SJeremy L Thompson   data.t_id_y = threadIdx.y;
13221292910SJeremy L Thompson   data.t_id_z = threadIdx.z;
13321292910SJeremy L Thompson   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
13421292910SJeremy L Thompson   data.slice  = slice + data.t_id_z * BASIS_T_1D * (BASIS_DIM > 1 ? BASIS_T_1D : 1);
13521292910SJeremy L Thompson 
13621292910SJeremy L Thompson   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
13721292910SJeremy L Thompson 
13821292910SJeremy L Thompson   // Apply basis element by element
13921292910SJeremy L Thompson   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
14021292910SJeremy L Thompson     if (BASIS_DIM == 1) {
14121292910SJeremy L Thompson       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U);
14221292910SJeremy L Thompson       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_U, d_V);
14321292910SJeremy L Thompson     } else if (BASIS_DIM == 2) {
14421292910SJeremy 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);
14521292910SJeremy L Thompson       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_U, d_V);
14621292910SJeremy L Thompson     } else if (BASIS_DIM == 3) {
14721292910SJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem,
14821292910SJeremy L Thompson                                                        BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
14921292910SJeremy L Thompson       WriteElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
15021292910SJeremy L Thompson                                                         BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_U, d_V);
15121292910SJeremy L Thompson     }
15221292910SJeremy L Thompson   }
15321292910SJeremy L Thompson }
15421292910SJeremy L Thompson 
InterpTransposeAdd(const CeedInt num_elem,const CeedScalar * c_B,const CeedScalar * __restrict__ d_U,CeedScalar * __restrict__ d_V)155db2becc9SJeremy L Thompson extern "C" __global__ void InterpTransposeAdd(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *__restrict__ d_U,
156db2becc9SJeremy L Thompson                                               CeedScalar *__restrict__ d_V) {
157db2becc9SJeremy L Thompson   extern __shared__ CeedScalar slice[];
158db2becc9SJeremy L Thompson 
159db2becc9SJeremy L Thompson   SharedData_Cuda data;
160db2becc9SJeremy L Thompson   data.t_id_x = threadIdx.x;
161db2becc9SJeremy L Thompson   data.t_id_y = threadIdx.y;
162db2becc9SJeremy L Thompson   data.t_id_z = threadIdx.z;
163db2becc9SJeremy L Thompson   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
16499421279SJeremy L Thompson   data.slice  = slice + data.t_id_z * BASIS_T_1D * (BASIS_DIM > 1 ? BASIS_T_1D : 1);
165db2becc9SJeremy L Thompson 
166db2becc9SJeremy L Thompson   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
167db2becc9SJeremy L Thompson   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
168db2becc9SJeremy L Thompson 
169aa4002adSJeremy L Thompson   // load interp_1d into shared memory
170aa4002adSJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
171aa4002adSJeremy L Thompson   LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B);
172aa4002adSJeremy L Thompson   __syncthreads();
173aa4002adSJeremy L Thompson 
174aa4002adSJeremy L Thompson   // Apply basis element by element
175db2becc9SJeremy L Thompson   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
176db2becc9SJeremy L Thompson     if (BASIS_DIM == 1) {
177db2becc9SJeremy L Thompson       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U);
17899421279SJeremy L Thompson       InterpTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, r_V);
179db2becc9SJeremy L Thompson       SumElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V);
180db2becc9SJeremy L Thompson     } else if (BASIS_DIM == 2) {
181db2becc9SJeremy 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);
18299421279SJeremy L Thompson       InterpTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, r_V);
183db2becc9SJeremy 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);
184db2becc9SJeremy L Thompson     } else if (BASIS_DIM == 3) {
185db2becc9SJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem,
186db2becc9SJeremy L Thompson                                                        BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
18799421279SJeremy L Thompson       InterpTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, r_V);
188db2becc9SJeremy L Thompson       SumElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
189db2becc9SJeremy L Thompson                                                       BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
190db2becc9SJeremy L Thompson     }
191db2becc9SJeremy L Thompson   }
192db2becc9SJeremy L Thompson }
193db2becc9SJeremy L Thompson 
InterpCollocatedTransposeAdd(const CeedInt num_elem,const CeedScalar * c_B,const CeedScalar * __restrict__ d_U,CeedScalar * __restrict__ d_V)19421292910SJeremy L Thompson extern "C" __global__ void InterpCollocatedTransposeAdd(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *__restrict__ d_U,
19521292910SJeremy L Thompson                                                         CeedScalar *__restrict__ d_V) {
19621292910SJeremy L Thompson   extern __shared__ CeedScalar slice[];
19721292910SJeremy L Thompson 
19821292910SJeremy L Thompson   SharedData_Cuda data;
19921292910SJeremy L Thompson   data.t_id_x = threadIdx.x;
20021292910SJeremy L Thompson   data.t_id_y = threadIdx.y;
20121292910SJeremy L Thompson   data.t_id_z = threadIdx.z;
20221292910SJeremy L Thompson   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
20321292910SJeremy L Thompson   data.slice  = slice + data.t_id_z * BASIS_T_1D * (BASIS_DIM > 1 ? BASIS_T_1D : 1);
20421292910SJeremy L Thompson 
20521292910SJeremy L Thompson   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
20621292910SJeremy L Thompson 
20721292910SJeremy L Thompson   // Apply basis element by element
20821292910SJeremy L Thompson   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
20921292910SJeremy L Thompson     if (BASIS_DIM == 1) {
21021292910SJeremy L Thompson       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U);
21121292910SJeremy L Thompson       SumElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_U, d_V);
21221292910SJeremy L Thompson     } else if (BASIS_DIM == 2) {
21321292910SJeremy 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);
21421292910SJeremy 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_U, d_V);
21521292910SJeremy L Thompson     } else if (BASIS_DIM == 3) {
21621292910SJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem,
21721292910SJeremy L Thompson                                                        BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
21821292910SJeremy L Thompson       SumElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
21921292910SJeremy L Thompson                                                       BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_U, d_V);
22021292910SJeremy L Thompson     }
22121292910SJeremy L Thompson   }
22221292910SJeremy L Thompson }
22321292910SJeremy L Thompson 
2249e201c85SYohann //------------------------------------------------------------------------------
2259e201c85SYohann // Grad kernel by dim
2269e201c85SYohann //------------------------------------------------------------------------------
Grad(const CeedInt num_elem,const CeedScalar * c_B,const CeedScalar * c_G,const CeedScalar * __restrict__ d_U,CeedScalar * __restrict__ d_V)2272b730f8bSJeremy L Thompson extern "C" __global__ void Grad(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *c_G, const CeedScalar *__restrict__ d_U,
2289e201c85SYohann                                 CeedScalar *__restrict__ d_V) {
2299e201c85SYohann   extern __shared__ CeedScalar slice[];
2309e201c85SYohann 
2319e201c85SYohann   SharedData_Cuda data;
2329e201c85SYohann   data.t_id_x = threadIdx.x;
2339e201c85SYohann   data.t_id_y = threadIdx.y;
2349e201c85SYohann   data.t_id_z = threadIdx.z;
2359e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
23699421279SJeremy L Thompson   data.slice  = slice + data.t_id_z * BASIS_T_1D * (BASIS_DIM > 1 ? BASIS_T_1D : 1);
2379e201c85SYohann 
2389e201c85SYohann   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
2399e201c85SYohann   CeedScalar r_V[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
2409e201c85SYohann 
241aa4002adSJeremy L Thompson   // load interp_1d and grad_1d into shared memory
242aa4002adSJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
243aa4002adSJeremy L Thompson   LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B);
244aa4002adSJeremy L Thompson   __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)];
245aa4002adSJeremy L Thompson   LoadMatrix<BASIS_Q_1D, BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D>(data, c_G, s_G);
246aa4002adSJeremy L Thompson   __syncthreads();
247aa4002adSJeremy L Thompson 
248aa4002adSJeremy L Thompson   // Apply basis element by element
2499e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
2509e201c85SYohann     if (BASIS_DIM == 1) {
2519e201c85SYohann       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, d_U, r_U);
25299421279SJeremy L Thompson       Grad1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, s_G, r_V);
2539e201c85SYohann       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_V, d_V);
2549e201c85SYohann     } else if (BASIS_DIM == 2) {
2559e201c85SYohann       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);
25699421279SJeremy L Thompson       GradTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, s_G, r_V);
2572b730f8bSJeremy 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,
2582b730f8bSJeremy L Thompson                                                                     d_V);
2599e201c85SYohann     } else if (BASIS_DIM == 3) {
2602b730f8bSJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
2612b730f8bSJeremy L Thompson                                                        BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, d_U, r_U);
26299421279SJeremy L Thompson       if (BASIS_HAS_COLLOCATED_GRAD) GradTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, s_G, r_V);
26399421279SJeremy L Thompson       else GradTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, s_G, r_V);
2642b730f8bSJeremy 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,
2652b730f8bSJeremy L Thompson                                                                     BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, r_V, d_V);
2669e201c85SYohann     }
2679e201c85SYohann   }
2689e201c85SYohann }
2699e201c85SYohann 
GradCollocated(const CeedInt num_elem,const CeedScalar * c_B,const CeedScalar * c_G,const CeedScalar * __restrict__ d_U,CeedScalar * __restrict__ d_V)27021292910SJeremy L Thompson extern "C" __global__ void GradCollocated(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *c_G, const CeedScalar *__restrict__ d_U,
27121292910SJeremy L Thompson                                           CeedScalar *__restrict__ d_V) {
27221292910SJeremy L Thompson   extern __shared__ CeedScalar slice[];
27321292910SJeremy L Thompson 
27421292910SJeremy L Thompson   SharedData_Cuda data;
27521292910SJeremy L Thompson   data.t_id_x = threadIdx.x;
27621292910SJeremy L Thompson   data.t_id_y = threadIdx.y;
27721292910SJeremy L Thompson   data.t_id_z = threadIdx.z;
27821292910SJeremy L Thompson   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
27921292910SJeremy L Thompson   data.slice  = slice + data.t_id_z * BASIS_T_1D * (BASIS_DIM > 1 ? BASIS_T_1D : 1);
28021292910SJeremy L Thompson 
28121292910SJeremy L Thompson   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
28221292910SJeremy L Thompson   CeedScalar r_V[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
28321292910SJeremy L Thompson 
28421292910SJeremy L Thompson   // load grad_1d into shared memory
28521292910SJeremy L Thompson   __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)];
28621292910SJeremy L Thompson   LoadMatrix<BASIS_Q_1D, BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D>(data, c_G, s_G);
28721292910SJeremy L Thompson   __syncthreads();
28821292910SJeremy L Thompson 
28921292910SJeremy L Thompson   // Apply basis element by element
29021292910SJeremy L Thompson   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
29121292910SJeremy L Thompson     if (BASIS_DIM == 1) {
29221292910SJeremy L Thompson       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, d_U, r_U);
29321292910SJeremy L Thompson       Grad1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, NULL, s_G, r_V);
29421292910SJeremy L Thompson       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_V, d_V);
29521292910SJeremy L Thompson     } else if (BASIS_DIM == 2) {
29621292910SJeremy 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);
2970ccda8ebSJeremy L Thompson       GradTensorCollocatedNodes2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, NULL, s_G, r_V);
29821292910SJeremy 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,
29921292910SJeremy L Thompson                                                                     d_V);
30021292910SJeremy L Thompson     } else if (BASIS_DIM == 3) {
30121292910SJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
30221292910SJeremy L Thompson                                                        BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, d_U, r_U);
3030ccda8ebSJeremy L Thompson       GradTensorCollocatedNodes3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, NULL, s_G, r_V);
30421292910SJeremy 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,
30521292910SJeremy L Thompson                                                                     BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, r_V, d_V);
30621292910SJeremy L Thompson     }
30721292910SJeremy L Thompson   }
30821292910SJeremy L Thompson }
30921292910SJeremy L Thompson 
GradTranspose(const CeedInt num_elem,const CeedScalar * c_B,const CeedScalar * c_G,const CeedScalar * __restrict__ d_U,CeedScalar * __restrict__ d_V)3102b730f8bSJeremy L Thompson extern "C" __global__ void GradTranspose(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *c_G, const CeedScalar *__restrict__ d_U,
3119e201c85SYohann                                          CeedScalar *__restrict__ d_V) {
3129e201c85SYohann   extern __shared__ CeedScalar slice[];
3139e201c85SYohann 
3149e201c85SYohann   SharedData_Cuda data;
3159e201c85SYohann   data.t_id_x = threadIdx.x;
3169e201c85SYohann   data.t_id_y = threadIdx.y;
3179e201c85SYohann   data.t_id_z = threadIdx.z;
3189e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
31999421279SJeremy L Thompson   data.slice  = slice + data.t_id_z * BASIS_T_1D * (BASIS_DIM > 1 ? BASIS_T_1D : 1);
3209e201c85SYohann 
3219e201c85SYohann   CeedScalar r_U[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
3229e201c85SYohann   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
3239e201c85SYohann 
324aa4002adSJeremy L Thompson   // load interp_1d and grad_1d into shared memory
325aa4002adSJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
326aa4002adSJeremy L Thompson   LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B);
327aa4002adSJeremy L Thompson   __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)];
328aa4002adSJeremy L Thompson   LoadMatrix<BASIS_Q_1D, BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D>(data, c_G, s_G);
329aa4002adSJeremy L Thompson   __syncthreads();
330aa4002adSJeremy L Thompson 
331aa4002adSJeremy L Thompson   // Apply basis element by element
3329e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
3339e201c85SYohann     if (BASIS_DIM == 1) {
3349e201c85SYohann       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U);
33599421279SJeremy L Thompson       GradTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, s_G, r_V);
3369e201c85SYohann       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V);
3379e201c85SYohann     } else if (BASIS_DIM == 2) {
3382b730f8bSJeremy 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,
3392b730f8bSJeremy L Thompson                                                                    r_U);
34099421279SJeremy L Thompson       GradTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, s_G, r_V);
3419e201c85SYohann       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);
3429e201c85SYohann     } else if (BASIS_DIM == 3) {
3432b730f8bSJeremy 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,
3442b730f8bSJeremy L Thompson                                                                    BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
34599421279SJeremy L Thompson       if (BASIS_HAS_COLLOCATED_GRAD) GradTransposeTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, s_G, r_V);
34699421279SJeremy L Thompson       else GradTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, s_G, r_V);
3472b730f8bSJeremy L Thompson       WriteElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
3482b730f8bSJeremy L Thompson                                                         BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
3499e201c85SYohann     }
3509e201c85SYohann   }
3519e201c85SYohann }
3529e201c85SYohann 
GradCollocatedTranspose(const CeedInt num_elem,const CeedScalar * c_B,const CeedScalar * c_G,const CeedScalar * __restrict__ d_U,CeedScalar * __restrict__ d_V)35321292910SJeremy L Thompson extern "C" __global__ void GradCollocatedTranspose(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *c_G,
35421292910SJeremy L Thompson                                                    const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) {
35521292910SJeremy L Thompson   extern __shared__ CeedScalar slice[];
35621292910SJeremy L Thompson 
35721292910SJeremy L Thompson   SharedData_Cuda data;
35821292910SJeremy L Thompson   data.t_id_x = threadIdx.x;
35921292910SJeremy L Thompson   data.t_id_y = threadIdx.y;
36021292910SJeremy L Thompson   data.t_id_z = threadIdx.z;
36121292910SJeremy L Thompson   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
36221292910SJeremy L Thompson   data.slice  = slice + data.t_id_z * BASIS_T_1D * (BASIS_DIM > 1 ? BASIS_T_1D : 1);
36321292910SJeremy L Thompson 
36421292910SJeremy L Thompson   CeedScalar r_U[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
36521292910SJeremy L Thompson   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
36621292910SJeremy L Thompson 
36721292910SJeremy L Thompson   // load grad_1d into shared memory
36821292910SJeremy L Thompson   __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)];
36921292910SJeremy L Thompson   LoadMatrix<BASIS_Q_1D, BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D>(data, c_G, s_G);
37021292910SJeremy L Thompson   __syncthreads();
37121292910SJeremy L Thompson 
37221292910SJeremy L Thompson   // Apply basis element by element
37321292910SJeremy L Thompson   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
37421292910SJeremy L Thompson     if (BASIS_DIM == 1) {
37521292910SJeremy L Thompson       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U);
37621292910SJeremy L Thompson       GradTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, NULL, s_G, r_V);
37721292910SJeremy L Thompson       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V);
37821292910SJeremy L Thompson     } else if (BASIS_DIM == 2) {
37921292910SJeremy 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,
38021292910SJeremy L Thompson                                                                    r_U);
3810ccda8ebSJeremy L Thompson       GradTransposeTensorCollocatedNodes2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, NULL, s_G, r_V);
38221292910SJeremy L Thompson       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);
38321292910SJeremy L Thompson     } else if (BASIS_DIM == 3) {
38421292910SJeremy 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,
38521292910SJeremy L Thompson                                                                    BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
3860ccda8ebSJeremy L Thompson       GradTransposeTensorCollocatedNodes3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, NULL, s_G, r_V);
38721292910SJeremy L Thompson       WriteElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
38821292910SJeremy L Thompson                                                         BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
38921292910SJeremy L Thompson     }
39021292910SJeremy L Thompson   }
39121292910SJeremy L Thompson }
39221292910SJeremy L Thompson 
GradTransposeAdd(const CeedInt num_elem,const CeedScalar * c_B,const CeedScalar * c_G,const CeedScalar * __restrict__ d_U,CeedScalar * __restrict__ d_V)393db2becc9SJeremy L Thompson extern "C" __global__ void GradTransposeAdd(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *c_G, const CeedScalar *__restrict__ d_U,
394db2becc9SJeremy L Thompson                                             CeedScalar *__restrict__ d_V) {
395db2becc9SJeremy L Thompson   extern __shared__ CeedScalar slice[];
396db2becc9SJeremy L Thompson 
397db2becc9SJeremy L Thompson   SharedData_Cuda data;
398db2becc9SJeremy L Thompson   data.t_id_x = threadIdx.x;
399db2becc9SJeremy L Thompson   data.t_id_y = threadIdx.y;
400db2becc9SJeremy L Thompson   data.t_id_z = threadIdx.z;
401db2becc9SJeremy L Thompson   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
40299421279SJeremy L Thompson   data.slice  = slice + data.t_id_z * BASIS_T_1D * (BASIS_DIM > 1 ? BASIS_T_1D : 1);
403db2becc9SJeremy L Thompson 
404db2becc9SJeremy L Thompson   CeedScalar r_U[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
405db2becc9SJeremy L Thompson   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
406db2becc9SJeremy L Thompson 
407aa4002adSJeremy L Thompson   // load interp_1d and grad_1d into shared memory
408aa4002adSJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
409aa4002adSJeremy L Thompson   LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B);
410aa4002adSJeremy L Thompson   __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)];
411aa4002adSJeremy L Thompson   LoadMatrix<BASIS_Q_1D, BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D>(data, c_G, s_G);
412aa4002adSJeremy L Thompson   __syncthreads();
413aa4002adSJeremy L Thompson 
414aa4002adSJeremy L Thompson   // Apply basis element by element
415db2becc9SJeremy L Thompson   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
416db2becc9SJeremy L Thompson     if (BASIS_DIM == 1) {
417db2becc9SJeremy L Thompson       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U);
41899421279SJeremy L Thompson       GradTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, s_G, r_V);
419db2becc9SJeremy L Thompson       SumElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V);
420db2becc9SJeremy L Thompson     } else if (BASIS_DIM == 2) {
421db2becc9SJeremy 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,
422db2becc9SJeremy L Thompson                                                                    r_U);
42399421279SJeremy L Thompson       GradTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, s_G, r_V);
424db2becc9SJeremy 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);
425db2becc9SJeremy L Thompson     } else if (BASIS_DIM == 3) {
426db2becc9SJeremy 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,
427db2becc9SJeremy L Thompson                                                                    BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
42899421279SJeremy L Thompson       if (BASIS_HAS_COLLOCATED_GRAD) GradTransposeTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, s_G, r_V);
42999421279SJeremy L Thompson       else GradTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, s_B, s_G, r_V);
430db2becc9SJeremy L Thompson       SumElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
431db2becc9SJeremy L Thompson                                                       BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
432db2becc9SJeremy L Thompson     }
433db2becc9SJeremy L Thompson   }
434db2becc9SJeremy L Thompson }
435db2becc9SJeremy L Thompson 
GradCollocatedTransposeAdd(const CeedInt num_elem,const CeedScalar * c_B,const CeedScalar * c_G,const CeedScalar * __restrict__ d_U,CeedScalar * __restrict__ d_V)43621292910SJeremy L Thompson extern "C" __global__ void GradCollocatedTransposeAdd(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *c_G,
43721292910SJeremy L Thompson                                                       const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) {
43821292910SJeremy L Thompson   extern __shared__ CeedScalar slice[];
43921292910SJeremy L Thompson 
44021292910SJeremy L Thompson   SharedData_Cuda data;
44121292910SJeremy L Thompson   data.t_id_x = threadIdx.x;
44221292910SJeremy L Thompson   data.t_id_y = threadIdx.y;
44321292910SJeremy L Thompson   data.t_id_z = threadIdx.z;
44421292910SJeremy L Thompson   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
44521292910SJeremy L Thompson   data.slice  = slice + data.t_id_z * BASIS_T_1D * (BASIS_DIM > 1 ? BASIS_T_1D : 1);
44621292910SJeremy L Thompson 
44721292910SJeremy L Thompson   CeedScalar r_U[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
44821292910SJeremy L Thompson   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
44921292910SJeremy L Thompson 
45021292910SJeremy L Thompson   // load grad_1d into shared memory
45121292910SJeremy L Thompson   __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)];
45221292910SJeremy L Thompson   LoadMatrix<BASIS_Q_1D, BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D>(data, c_G, s_G);
45321292910SJeremy L Thompson   __syncthreads();
45421292910SJeremy L Thompson 
45521292910SJeremy L Thompson   // Apply basis element by element
45621292910SJeremy L Thompson   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
45721292910SJeremy L Thompson     if (BASIS_DIM == 1) {
45821292910SJeremy L Thompson       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U);
45921292910SJeremy L Thompson       GradTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, NULL, s_G, r_V);
46021292910SJeremy L Thompson       SumElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V);
46121292910SJeremy L Thompson     } else if (BASIS_DIM == 2) {
46221292910SJeremy 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,
46321292910SJeremy L Thompson                                                                    r_U);
4640ccda8ebSJeremy L Thompson       GradTransposeTensorCollocatedNodes2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, NULL, s_G, r_V);
46521292910SJeremy 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);
46621292910SJeremy L Thompson     } else if (BASIS_DIM == 3) {
46721292910SJeremy 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,
46821292910SJeremy L Thompson                                                                    BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
4690ccda8ebSJeremy L Thompson       GradTransposeTensorCollocatedNodes3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D, BASIS_T_1D>(data, r_U, NULL, s_G, r_V);
47021292910SJeremy L Thompson       SumElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
47121292910SJeremy L Thompson                                                       BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
47221292910SJeremy L Thompson     }
47321292910SJeremy L Thompson   }
47421292910SJeremy L Thompson }
47521292910SJeremy L Thompson 
4769e201c85SYohann //------------------------------------------------------------------------------
4779e201c85SYohann // Weight kernels by dim
4789e201c85SYohann //------------------------------------------------------------------------------
Weight(const CeedInt num_elem,const CeedScalar * __restrict__ q_weight_1d,CeedScalar * __restrict__ d_W)4792b730f8bSJeremy L Thompson extern "C" __global__ void Weight(const CeedInt num_elem, const CeedScalar *__restrict__ q_weight_1d, CeedScalar *__restrict__ d_W) {
4809e201c85SYohann   extern __shared__ CeedScalar slice[];
4819e201c85SYohann 
4829e201c85SYohann   SharedData_Cuda data;
4839e201c85SYohann   data.t_id_x = threadIdx.x;
4849e201c85SYohann   data.t_id_y = threadIdx.y;
4859e201c85SYohann   data.t_id_z = threadIdx.z;
4869e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
48799421279SJeremy L Thompson   data.slice  = slice + data.t_id_z * BASIS_T_1D * (BASIS_DIM > 1 ? BASIS_T_1D : 1);
4889e201c85SYohann 
4899e201c85SYohann   CeedScalar r_W[BASIS_DIM > 2 ? BASIS_Q_1D : 1];
4909e201c85SYohann 
491aa4002adSJeremy L Thompson   // Apply basis element by element
4929e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
4939e201c85SYohann     if (BASIS_DIM == 1) {
494343e3094SJeremy L Thompson       Weight1d<BASIS_P_1D, BASIS_Q_1D>(data, q_weight_1d, r_W);
4959e201c85SYohann       WriteElementStrided1d<1, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_W, d_W);
4969e201c85SYohann     } else if (BASIS_DIM == 2) {
497343e3094SJeremy L Thompson       WeightTensor2d<BASIS_P_1D, BASIS_Q_1D>(data, q_weight_1d, r_W);
4989e201c85SYohann       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);
4999e201c85SYohann     } else if (BASIS_DIM == 3) {
500343e3094SJeremy L Thompson       WeightTensor3d<BASIS_P_1D, BASIS_Q_1D>(data, q_weight_1d, r_W);
5012b730f8bSJeremy 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,
5022b730f8bSJeremy L Thompson                                            d_W);
5039e201c85SYohann     }
5049e201c85SYohann   }
5059e201c85SYohann }
506