xref: /libCEED/include/ceed/jit-source/hip/hip-shared-basis-tensor.h (revision db2becc9f302fe8eb3a32ace50ce3f3a5d42e6c4)
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 HIP shared memory tensor product basis
109e201c85SYohann 
119e201c85SYohann #include <ceed.h>
122b730f8bSJeremy L Thompson 
139e201c85SYohann #include "hip-shared-basis-read-write-templates.h"
149e201c85SYohann #include "hip-shared-basis-tensor-templates.h"
159e201c85SYohann 
169e201c85SYohann //------------------------------------------------------------------------------
179e201c85SYohann // Interp kernel by dim
189e201c85SYohann //------------------------------------------------------------------------------
192b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_INTERP_BLOCK_SIZE) __global__
202b730f8bSJeremy L Thompson     void Interp(const CeedInt num_elem, const CeedScalar *d_interp_1d, const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) {
219e201c85SYohann   extern __shared__ CeedScalar slice[];
22b2165e7aSSebastian Grimberg 
239e201c85SYohann   // load interp_1d into shared memory
249e201c85SYohann   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
259e201c85SYohann   loadMatrix<BASIS_P_1D * BASIS_Q_1D>(d_interp_1d, s_B);
269e201c85SYohann   __syncthreads();
279e201c85SYohann 
289e201c85SYohann   SharedData_Hip data;
299e201c85SYohann   data.t_id_x = threadIdx.x;
309e201c85SYohann   data.t_id_y = threadIdx.y;
319e201c85SYohann   data.t_id_z = threadIdx.z;
329e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
339e201c85SYohann   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
349e201c85SYohann 
359e201c85SYohann   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
369e201c85SYohann   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
379e201c85SYohann 
389e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
399e201c85SYohann     if (BASIS_DIM == 1) {
409e201c85SYohann       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, d_U, r_U);
419e201c85SYohann       Interp1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
429e201c85SYohann       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_V, d_V);
439e201c85SYohann     } else if (BASIS_DIM == 2) {
449e201c85SYohann       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);
459e201c85SYohann       InterpTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
469e201c85SYohann       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);
479e201c85SYohann     } else if (BASIS_DIM == 3) {
482b730f8bSJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
492b730f8bSJeremy L Thompson                                                        BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, d_U, r_U);
509e201c85SYohann       InterpTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
512b730f8bSJeremy L Thompson       WriteElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem,
522b730f8bSJeremy L Thompson                                                         BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, r_V, d_V);
539e201c85SYohann     }
549e201c85SYohann   }
559e201c85SYohann }
569e201c85SYohann 
572b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_INTERP_BLOCK_SIZE) __global__
582b730f8bSJeremy L Thompson     void InterpTranspose(const CeedInt num_elem, const CeedScalar *d_interp_1d, const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) {
599e201c85SYohann   extern __shared__ CeedScalar slice[];
60b2165e7aSSebastian Grimberg 
619e201c85SYohann   // load interp_1d into shared memory
629e201c85SYohann   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
639e201c85SYohann   loadMatrix<BASIS_P_1D * BASIS_Q_1D>(d_interp_1d, s_B);
649e201c85SYohann   __syncthreads();
659e201c85SYohann 
669e201c85SYohann   SharedData_Hip data;
679e201c85SYohann   data.t_id_x = threadIdx.x;
689e201c85SYohann   data.t_id_y = threadIdx.y;
699e201c85SYohann   data.t_id_z = threadIdx.z;
709e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
719e201c85SYohann   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
729e201c85SYohann 
739e201c85SYohann   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
749e201c85SYohann   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
759e201c85SYohann 
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);
799e201c85SYohann       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);
839e201c85SYohann       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);
889e201c85SYohann       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 
95*db2becc9SJeremy L Thompson extern "C" __launch_bounds__(BASIS_INTERP_BLOCK_SIZE) __global__
96*db2becc9SJeremy L Thompson     void InterpTransposeAdd(const CeedInt num_elem, const CeedScalar *d_interp_1d, const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) {
97*db2becc9SJeremy L Thompson   extern __shared__ CeedScalar slice[];
98*db2becc9SJeremy L Thompson 
99*db2becc9SJeremy L Thompson   // load interp_1d into shared memory
100*db2becc9SJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
101*db2becc9SJeremy L Thompson   loadMatrix<BASIS_P_1D * BASIS_Q_1D>(d_interp_1d, s_B);
102*db2becc9SJeremy L Thompson   __syncthreads();
103*db2becc9SJeremy L Thompson 
104*db2becc9SJeremy L Thompson   SharedData_Hip data;
105*db2becc9SJeremy L Thompson   data.t_id_x = threadIdx.x;
106*db2becc9SJeremy L Thompson   data.t_id_y = threadIdx.y;
107*db2becc9SJeremy L Thompson   data.t_id_z = threadIdx.z;
108*db2becc9SJeremy L Thompson   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
109*db2becc9SJeremy L Thompson   data.slice  = &slice[data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1)];
110*db2becc9SJeremy L Thompson 
111*db2becc9SJeremy L Thompson   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
112*db2becc9SJeremy L Thompson   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
113*db2becc9SJeremy L Thompson 
114*db2becc9SJeremy L Thompson   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
115*db2becc9SJeremy L Thompson     if (BASIS_DIM == 1) {
116*db2becc9SJeremy L Thompson       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U);
117*db2becc9SJeremy L Thompson       InterpTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
118*db2becc9SJeremy L Thompson       SumElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V);
119*db2becc9SJeremy L Thompson     } else if (BASIS_DIM == 2) {
120*db2becc9SJeremy 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);
121*db2becc9SJeremy L Thompson       InterpTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
122*db2becc9SJeremy 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);
123*db2becc9SJeremy L Thompson     } else if (BASIS_DIM == 3) {
124*db2becc9SJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem,
125*db2becc9SJeremy L Thompson                                                        BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
126*db2becc9SJeremy L Thompson       InterpTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
127*db2becc9SJeremy L Thompson       SumElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
128*db2becc9SJeremy L Thompson                                                       BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
129*db2becc9SJeremy L Thompson     }
130*db2becc9SJeremy L Thompson   }
131*db2becc9SJeremy L Thompson }
132*db2becc9SJeremy L Thompson 
1339e201c85SYohann //------------------------------------------------------------------------------
1349e201c85SYohann // Grad kernel by dim
1359e201c85SYohann //------------------------------------------------------------------------------
1362b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_GRAD_BLOCK_SIZE) __global__
1372b730f8bSJeremy L Thompson     void Grad(const CeedInt num_elem, const CeedScalar *d_interp_1d, const CeedScalar *d_grad_1d, const CeedScalar *__restrict__ d_U,
1389e201c85SYohann               CeedScalar *__restrict__ d_V) {
1399e201c85SYohann   extern __shared__ CeedScalar slice[];
140b2165e7aSSebastian Grimberg 
1419e201c85SYohann   // load interp_1d and grad_1d into shared memory
1429e201c85SYohann   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
1439e201c85SYohann   loadMatrix<BASIS_P_1D * BASIS_Q_1D>(d_interp_1d, s_B);
1449e201c85SYohann   __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)];
1459e201c85SYohann   loadMatrix<BASIS_Q_1D *(BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)>(d_grad_1d, s_G);
1469e201c85SYohann   __syncthreads();
1479e201c85SYohann 
1489e201c85SYohann   SharedData_Hip data;
1499e201c85SYohann   data.t_id_x = threadIdx.x;
1509e201c85SYohann   data.t_id_y = threadIdx.y;
1519e201c85SYohann   data.t_id_z = threadIdx.z;
1529e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
1539e201c85SYohann   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
1549e201c85SYohann 
1559e201c85SYohann   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
1569e201c85SYohann   CeedScalar r_V[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
1579e201c85SYohann 
1589e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
1599e201c85SYohann     if (BASIS_DIM == 1) {
1609e201c85SYohann       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, d_U, r_U);
1619e201c85SYohann       Grad1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
1629e201c85SYohann       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_V, d_V);
1639e201c85SYohann     } else if (BASIS_DIM == 2) {
1649e201c85SYohann       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);
1659e201c85SYohann       GradTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
1662b730f8bSJeremy 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,
1672b730f8bSJeremy L Thompson                                                                     d_V);
1689e201c85SYohann     } else if (BASIS_DIM == 3) {
1692b730f8bSJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
1702b730f8bSJeremy L Thompson                                                        BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, d_U, r_U);
1719e201c85SYohann       if (BASIS_HAS_COLLOCATED_GRAD) GradTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
1729e201c85SYohann       else GradTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
1732b730f8bSJeremy 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,
1742b730f8bSJeremy L Thompson                                                                     BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, r_V, d_V);
1759e201c85SYohann     }
1769e201c85SYohann   }
1779e201c85SYohann }
1789e201c85SYohann 
1792b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_GRAD_BLOCK_SIZE) __global__
1802b730f8bSJeremy L Thompson     void GradTranspose(const CeedInt num_elem, const CeedScalar *d_interp_1d, const CeedScalar *d_grad_1d, const CeedScalar *__restrict__ d_U,
1819e201c85SYohann                        CeedScalar *__restrict__ d_V) {
1829e201c85SYohann   extern __shared__ CeedScalar slice[];
183b2165e7aSSebastian Grimberg 
1849e201c85SYohann   // load interp_1d and grad_1d into shared memory
1859e201c85SYohann   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
1869e201c85SYohann   loadMatrix<BASIS_P_1D * BASIS_Q_1D>(d_interp_1d, s_B);
1879e201c85SYohann   __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)];
1889e201c85SYohann   loadMatrix<BASIS_Q_1D *(BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)>(d_grad_1d, s_G);
1899e201c85SYohann   __syncthreads();
1909e201c85SYohann 
1919e201c85SYohann   SharedData_Hip data;
1929e201c85SYohann   data.t_id_x = threadIdx.x;
1939e201c85SYohann   data.t_id_y = threadIdx.y;
1949e201c85SYohann   data.t_id_z = threadIdx.z;
1959e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
1969e201c85SYohann   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
1979e201c85SYohann 
1989e201c85SYohann   CeedScalar r_U[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
1999e201c85SYohann   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
2009e201c85SYohann 
2019e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
2029e201c85SYohann     if (BASIS_DIM == 1) {
2039e201c85SYohann       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U);
2049e201c85SYohann       GradTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
2059e201c85SYohann       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V);
2069e201c85SYohann     } else if (BASIS_DIM == 2) {
2072b730f8bSJeremy 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,
2082b730f8bSJeremy L Thompson                                                                    r_U);
2099e201c85SYohann       GradTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
2109e201c85SYohann       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);
2119e201c85SYohann     } else if (BASIS_DIM == 3) {
2122b730f8bSJeremy 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,
2132b730f8bSJeremy L Thompson                                                                    BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
2149e201c85SYohann       if (BASIS_HAS_COLLOCATED_GRAD) GradTransposeTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
2159e201c85SYohann       else GradTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
2162b730f8bSJeremy L Thompson       WriteElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
2172b730f8bSJeremy L Thompson                                                         BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
2189e201c85SYohann     }
2199e201c85SYohann   }
2209e201c85SYohann }
2219e201c85SYohann 
222*db2becc9SJeremy L Thompson extern "C" __launch_bounds__(BASIS_GRAD_BLOCK_SIZE) __global__
223*db2becc9SJeremy L Thompson     void GradTransposeAdd(const CeedInt num_elem, const CeedScalar *d_interp_1d, const CeedScalar *d_grad_1d, const CeedScalar *__restrict__ d_U,
224*db2becc9SJeremy L Thompson                           CeedScalar *__restrict__ d_V) {
225*db2becc9SJeremy L Thompson   extern __shared__ CeedScalar slice[];
226*db2becc9SJeremy L Thompson 
227*db2becc9SJeremy L Thompson   // load interp_1d and grad_1d into shared memory
228*db2becc9SJeremy L Thompson   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
229*db2becc9SJeremy L Thompson   loadMatrix<BASIS_P_1D * BASIS_Q_1D>(d_interp_1d, s_B);
230*db2becc9SJeremy L Thompson   __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)];
231*db2becc9SJeremy L Thompson   loadMatrix<BASIS_Q_1D *(BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)>(d_grad_1d, s_G);
232*db2becc9SJeremy L Thompson   __syncthreads();
233*db2becc9SJeremy L Thompson 
234*db2becc9SJeremy L Thompson   SharedData_Hip data;
235*db2becc9SJeremy L Thompson   data.t_id_x = threadIdx.x;
236*db2becc9SJeremy L Thompson   data.t_id_y = threadIdx.y;
237*db2becc9SJeremy L Thompson   data.t_id_z = threadIdx.z;
238*db2becc9SJeremy L Thompson   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
239*db2becc9SJeremy L Thompson   data.slice  = &slice[data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1)];
240*db2becc9SJeremy L Thompson 
241*db2becc9SJeremy L Thompson   CeedScalar r_U[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
242*db2becc9SJeremy L Thompson   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
243*db2becc9SJeremy L Thompson 
244*db2becc9SJeremy L Thompson   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
245*db2becc9SJeremy L Thompson     if (BASIS_DIM == 1) {
246*db2becc9SJeremy L Thompson       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U);
247*db2becc9SJeremy L Thompson       GradTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
248*db2becc9SJeremy L Thompson       SumElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V);
249*db2becc9SJeremy L Thompson     } else if (BASIS_DIM == 2) {
250*db2becc9SJeremy 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,
251*db2becc9SJeremy L Thompson                                                                    r_U);
252*db2becc9SJeremy L Thompson       GradTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
253*db2becc9SJeremy 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);
254*db2becc9SJeremy L Thompson     } else if (BASIS_DIM == 3) {
255*db2becc9SJeremy 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,
256*db2becc9SJeremy L Thompson                                                                    BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
257*db2becc9SJeremy 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);
258*db2becc9SJeremy L Thompson       else GradTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
259*db2becc9SJeremy L Thompson       SumElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
260*db2becc9SJeremy L Thompson                                                       BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
261*db2becc9SJeremy L Thompson     }
262*db2becc9SJeremy L Thompson   }
263*db2becc9SJeremy L Thompson }
264*db2becc9SJeremy L Thompson 
2659e201c85SYohann //------------------------------------------------------------------------------
2669e201c85SYohann // Weight kernels by dim
2679e201c85SYohann //------------------------------------------------------------------------------
2682b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_WEIGHT_BLOCK_SIZE) __global__
2692b730f8bSJeremy L Thompson     void Weight(const CeedInt num_elem, const CeedScalar *__restrict__ q_weight_1d, CeedScalar *__restrict__ d_W) {
2709e201c85SYohann   extern __shared__ CeedScalar slice[];
2719e201c85SYohann 
2729e201c85SYohann   SharedData_Hip 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 
2819e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
2829e201c85SYohann     if (BASIS_DIM == 1) {
2839e201c85SYohann       Weight1d<BASIS_Q_1D>(data, q_weight_1d, r_W);
2849e201c85SYohann       WriteElementStrided1d<1, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_W, d_W);
2859e201c85SYohann     } else if (BASIS_DIM == 2) {
2869e201c85SYohann       WeightTensor2d<BASIS_Q_1D>(data, q_weight_1d, r_W);
2879e201c85SYohann       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);
2889e201c85SYohann     } else if (BASIS_DIM == 3) {
2899e201c85SYohann       WeightTensor3d<BASIS_Q_1D>(data, q_weight_1d, r_W);
2902b730f8bSJeremy 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,
2912b730f8bSJeremy L Thompson                                            d_W);
2929e201c85SYohann     }
2939e201c85SYohann   }
2949e201c85SYohann }
295