xref: /libCEED/include/ceed/jit-source/hip/hip-shared-basis-tensor.h (revision 2b730f8b5a9c809740a0b3b302db43a719c636b1)
19e201c85SYohann // Copyright (c) 2017-2022, 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 #ifndef _ceed_hip_shared_basis_tensor_h
119e201c85SYohann #define _ceed_hip_shared_basis_tensor_h
129e201c85SYohann 
139e201c85SYohann #include <ceed.h>
14*2b730f8bSJeremy L Thompson 
159e201c85SYohann #include "hip-shared-basis-read-write-templates.h"
169e201c85SYohann #include "hip-shared-basis-tensor-templates.h"
179e201c85SYohann 
189e201c85SYohann //------------------------------------------------------------------------------
199e201c85SYohann // Interp kernel by dim
209e201c85SYohann //------------------------------------------------------------------------------
21*2b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_INTERP_BLOCK_SIZE) __global__
22*2b730f8bSJeremy L Thompson     void Interp(const CeedInt num_elem, const CeedScalar *d_interp_1d, const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) {
239e201c85SYohann   extern __shared__ CeedScalar slice[];
249e201c85SYohann   // load interp_1d into shared memory
259e201c85SYohann   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
269e201c85SYohann   loadMatrix<BASIS_P_1D * BASIS_Q_1D>(d_interp_1d, s_B);
279e201c85SYohann   __syncthreads();
289e201c85SYohann 
299e201c85SYohann   SharedData_Hip data;
309e201c85SYohann   data.t_id_x = threadIdx.x;
319e201c85SYohann   data.t_id_y = threadIdx.y;
329e201c85SYohann   data.t_id_z = threadIdx.z;
339e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
349e201c85SYohann   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
359e201c85SYohann 
369e201c85SYohann   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
379e201c85SYohann   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
389e201c85SYohann 
399e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
409e201c85SYohann     if (BASIS_DIM == 1) {
419e201c85SYohann       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, d_U, r_U);
429e201c85SYohann       Interp1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
439e201c85SYohann       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_V, d_V);
449e201c85SYohann     } else if (BASIS_DIM == 2) {
459e201c85SYohann       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);
469e201c85SYohann       InterpTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
479e201c85SYohann       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);
489e201c85SYohann     } else if (BASIS_DIM == 3) {
49*2b730f8bSJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
50*2b730f8bSJeremy L Thompson                                                        BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, d_U, r_U);
519e201c85SYohann       InterpTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V);
52*2b730f8bSJeremy L Thompson       WriteElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem,
53*2b730f8bSJeremy L Thompson                                                         BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, r_V, d_V);
549e201c85SYohann     }
559e201c85SYohann   }
569e201c85SYohann }
579e201c85SYohann 
58*2b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_INTERP_BLOCK_SIZE) __global__
59*2b730f8bSJeremy L Thompson     void InterpTranspose(const CeedInt num_elem, const CeedScalar *d_interp_1d, const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) {
609e201c85SYohann   extern __shared__ CeedScalar slice[];
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) {
86*2b730f8bSJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem,
87*2b730f8bSJeremy 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);
89*2b730f8bSJeremy L Thompson       WriteElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
90*2b730f8bSJeremy L Thompson                                                         BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
919e201c85SYohann     }
929e201c85SYohann   }
939e201c85SYohann }
949e201c85SYohann 
959e201c85SYohann //------------------------------------------------------------------------------
969e201c85SYohann // Grad kernel by dim
979e201c85SYohann //------------------------------------------------------------------------------
98*2b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_GRAD_BLOCK_SIZE) __global__
99*2b730f8bSJeremy L Thompson     void Grad(const CeedInt num_elem, const CeedScalar *d_interp_1d, const CeedScalar *d_grad_1d, const CeedScalar *__restrict__ d_U,
1009e201c85SYohann               CeedScalar *__restrict__ d_V) {
1019e201c85SYohann   extern __shared__ CeedScalar slice[];
1029e201c85SYohann   // load interp_1d and grad_1d into shared memory
1039e201c85SYohann   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
1049e201c85SYohann   loadMatrix<BASIS_P_1D * BASIS_Q_1D>(d_interp_1d, s_B);
1059e201c85SYohann   __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)];
1069e201c85SYohann   loadMatrix<BASIS_Q_1D *(BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)>(d_grad_1d, s_G);
1079e201c85SYohann   __syncthreads();
1089e201c85SYohann 
1099e201c85SYohann   SharedData_Hip data;
1109e201c85SYohann   data.t_id_x = threadIdx.x;
1119e201c85SYohann   data.t_id_y = threadIdx.y;
1129e201c85SYohann   data.t_id_z = threadIdx.z;
1139e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
1149e201c85SYohann   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
1159e201c85SYohann 
1169e201c85SYohann   CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
1179e201c85SYohann   CeedScalar r_V[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
1189e201c85SYohann 
1199e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
1209e201c85SYohann     if (BASIS_DIM == 1) {
1219e201c85SYohann       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, d_U, r_U);
1229e201c85SYohann       Grad1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
1239e201c85SYohann       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_V, d_V);
1249e201c85SYohann     } else if (BASIS_DIM == 2) {
1259e201c85SYohann       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);
1269e201c85SYohann       GradTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
127*2b730f8bSJeremy 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,
128*2b730f8bSJeremy L Thompson                                                                     d_V);
1299e201c85SYohann     } else if (BASIS_DIM == 3) {
130*2b730f8bSJeremy L Thompson       ReadElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
131*2b730f8bSJeremy L Thompson                                                        BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, d_U, r_U);
1329e201c85SYohann       if (BASIS_HAS_COLLOCATED_GRAD) GradTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
1339e201c85SYohann       else GradTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
134*2b730f8bSJeremy 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,
135*2b730f8bSJeremy L Thompson                                                                     BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, r_V, d_V);
1369e201c85SYohann     }
1379e201c85SYohann   }
1389e201c85SYohann }
1399e201c85SYohann 
140*2b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_GRAD_BLOCK_SIZE) __global__
141*2b730f8bSJeremy L Thompson     void GradTranspose(const CeedInt num_elem, const CeedScalar *d_interp_1d, const CeedScalar *d_grad_1d, const CeedScalar *__restrict__ d_U,
1429e201c85SYohann                        CeedScalar *__restrict__ d_V) {
1439e201c85SYohann   extern __shared__ CeedScalar slice[];
1449e201c85SYohann   // load interp_1d and grad_1d into shared memory
1459e201c85SYohann   __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D];
1469e201c85SYohann   loadMatrix<BASIS_P_1D * BASIS_Q_1D>(d_interp_1d, s_B);
1479e201c85SYohann   __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)];
1489e201c85SYohann   loadMatrix<BASIS_Q_1D *(BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)>(d_grad_1d, s_G);
1499e201c85SYohann   __syncthreads();
1509e201c85SYohann 
1519e201c85SYohann   SharedData_Hip data;
1529e201c85SYohann   data.t_id_x = threadIdx.x;
1539e201c85SYohann   data.t_id_y = threadIdx.y;
1549e201c85SYohann   data.t_id_z = threadIdx.z;
1559e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
1569e201c85SYohann   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
1579e201c85SYohann 
1589e201c85SYohann   CeedScalar r_U[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)];
1599e201c85SYohann   CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)];
1609e201c85SYohann 
1619e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
1629e201c85SYohann     if (BASIS_DIM == 1) {
1639e201c85SYohann       ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U);
1649e201c85SYohann       GradTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
1659e201c85SYohann       WriteElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V);
1669e201c85SYohann     } else if (BASIS_DIM == 2) {
167*2b730f8bSJeremy 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,
168*2b730f8bSJeremy L Thompson                                                                    r_U);
1699e201c85SYohann       GradTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
1709e201c85SYohann       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);
1719e201c85SYohann     } else if (BASIS_DIM == 3) {
172*2b730f8bSJeremy 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,
173*2b730f8bSJeremy L Thompson                                                                    BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U);
1749e201c85SYohann       if (BASIS_HAS_COLLOCATED_GRAD) GradTransposeTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
1759e201c85SYohann       else GradTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V);
176*2b730f8bSJeremy L Thompson       WriteElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem,
177*2b730f8bSJeremy L Thompson                                                         BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V);
1789e201c85SYohann     }
1799e201c85SYohann   }
1809e201c85SYohann }
1819e201c85SYohann 
1829e201c85SYohann //------------------------------------------------------------------------------
1839e201c85SYohann // Weight kernels by dim
1849e201c85SYohann //------------------------------------------------------------------------------
185*2b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_WEIGHT_BLOCK_SIZE) __global__
186*2b730f8bSJeremy L Thompson     void Weight(const CeedInt num_elem, const CeedScalar *__restrict__ q_weight_1d, CeedScalar *__restrict__ d_W) {
1879e201c85SYohann   extern __shared__ CeedScalar slice[];
1889e201c85SYohann 
1899e201c85SYohann   SharedData_Hip data;
1909e201c85SYohann   data.t_id_x = threadIdx.x;
1919e201c85SYohann   data.t_id_y = threadIdx.y;
1929e201c85SYohann   data.t_id_z = threadIdx.z;
1939e201c85SYohann   data.t_id   = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x;
1949e201c85SYohann   data.slice  = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1);
1959e201c85SYohann 
1969e201c85SYohann   CeedScalar r_W[BASIS_DIM > 2 ? BASIS_Q_1D : 1];
1979e201c85SYohann 
1989e201c85SYohann   for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) {
1999e201c85SYohann     if (BASIS_DIM == 1) {
2009e201c85SYohann       Weight1d<BASIS_Q_1D>(data, q_weight_1d, r_W);
2019e201c85SYohann       WriteElementStrided1d<1, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_W, d_W);
2029e201c85SYohann     } else if (BASIS_DIM == 2) {
2039e201c85SYohann       WeightTensor2d<BASIS_Q_1D>(data, q_weight_1d, r_W);
2049e201c85SYohann       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);
2059e201c85SYohann     } else if (BASIS_DIM == 3) {
2069e201c85SYohann       WeightTensor3d<BASIS_Q_1D>(data, q_weight_1d, r_W);
207*2b730f8bSJeremy 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,
208*2b730f8bSJeremy L Thompson                                            d_W);
2099e201c85SYohann     }
2109e201c85SYohann   }
2119e201c85SYohann }
2129e201c85SYohann 
2139e201c85SYohann //------------------------------------------------------------------------------
2149e201c85SYohann 
2159e201c85SYohann #endif
216