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> 142b730f8bSJeremy 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 //------------------------------------------------------------------------------ 212b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_INTERP_BLOCK_SIZE) __global__ 222b730f8bSJeremy 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[]; 24*b2165e7aSSebastian Grimberg 259e201c85SYohann // load interp_1d into shared memory 269e201c85SYohann __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D]; 279e201c85SYohann loadMatrix<BASIS_P_1D * BASIS_Q_1D>(d_interp_1d, s_B); 289e201c85SYohann __syncthreads(); 299e201c85SYohann 309e201c85SYohann SharedData_Hip data; 319e201c85SYohann data.t_id_x = threadIdx.x; 329e201c85SYohann data.t_id_y = threadIdx.y; 339e201c85SYohann data.t_id_z = threadIdx.z; 349e201c85SYohann data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 359e201c85SYohann data.slice = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1); 369e201c85SYohann 379e201c85SYohann CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)]; 389e201c85SYohann CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)]; 399e201c85SYohann 409e201c85SYohann for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 419e201c85SYohann if (BASIS_DIM == 1) { 429e201c85SYohann ReadElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, d_U, r_U); 439e201c85SYohann Interp1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V); 449e201c85SYohann WriteElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_V, d_V); 459e201c85SYohann } else if (BASIS_DIM == 2) { 469e201c85SYohann 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); 479e201c85SYohann InterpTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V); 489e201c85SYohann 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); 499e201c85SYohann } else if (BASIS_DIM == 3) { 502b730f8bSJeremy L Thompson ReadElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem, 512b730f8bSJeremy L Thompson BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, d_U, r_U); 529e201c85SYohann InterpTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V); 532b730f8bSJeremy L Thompson WriteElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem, 542b730f8bSJeremy L Thompson BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, r_V, d_V); 559e201c85SYohann } 569e201c85SYohann } 579e201c85SYohann } 589e201c85SYohann 592b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_INTERP_BLOCK_SIZE) __global__ 602b730f8bSJeremy L Thompson void InterpTranspose(const CeedInt num_elem, const CeedScalar *d_interp_1d, const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) { 619e201c85SYohann extern __shared__ CeedScalar slice[]; 62*b2165e7aSSebastian Grimberg 639e201c85SYohann // load interp_1d into shared memory 649e201c85SYohann __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D]; 659e201c85SYohann loadMatrix<BASIS_P_1D * BASIS_Q_1D>(d_interp_1d, s_B); 669e201c85SYohann __syncthreads(); 679e201c85SYohann 689e201c85SYohann SharedData_Hip data; 699e201c85SYohann data.t_id_x = threadIdx.x; 709e201c85SYohann data.t_id_y = threadIdx.y; 719e201c85SYohann data.t_id_z = threadIdx.z; 729e201c85SYohann data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 739e201c85SYohann data.slice = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1); 749e201c85SYohann 759e201c85SYohann CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)]; 769e201c85SYohann CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)]; 779e201c85SYohann 789e201c85SYohann for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 799e201c85SYohann if (BASIS_DIM == 1) { 809e201c85SYohann ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U); 819e201c85SYohann InterpTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V); 829e201c85SYohann WriteElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V); 839e201c85SYohann } else if (BASIS_DIM == 2) { 849e201c85SYohann 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); 859e201c85SYohann InterpTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V); 869e201c85SYohann 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); 879e201c85SYohann } else if (BASIS_DIM == 3) { 882b730f8bSJeremy L Thompson ReadElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem, 892b730f8bSJeremy L Thompson BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U); 909e201c85SYohann InterpTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V); 912b730f8bSJeremy L Thompson WriteElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem, 922b730f8bSJeremy L Thompson BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V); 939e201c85SYohann } 949e201c85SYohann } 959e201c85SYohann } 969e201c85SYohann 979e201c85SYohann //------------------------------------------------------------------------------ 989e201c85SYohann // Grad kernel by dim 999e201c85SYohann //------------------------------------------------------------------------------ 1002b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_GRAD_BLOCK_SIZE) __global__ 1012b730f8bSJeremy L Thompson void Grad(const CeedInt num_elem, const CeedScalar *d_interp_1d, const CeedScalar *d_grad_1d, const CeedScalar *__restrict__ d_U, 1029e201c85SYohann CeedScalar *__restrict__ d_V) { 1039e201c85SYohann extern __shared__ CeedScalar slice[]; 104*b2165e7aSSebastian Grimberg 1059e201c85SYohann // load interp_1d and grad_1d into shared memory 1069e201c85SYohann __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D]; 1079e201c85SYohann loadMatrix<BASIS_P_1D * BASIS_Q_1D>(d_interp_1d, s_B); 1089e201c85SYohann __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)]; 1099e201c85SYohann loadMatrix<BASIS_Q_1D *(BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)>(d_grad_1d, s_G); 1109e201c85SYohann __syncthreads(); 1119e201c85SYohann 1129e201c85SYohann SharedData_Hip data; 1139e201c85SYohann data.t_id_x = threadIdx.x; 1149e201c85SYohann data.t_id_y = threadIdx.y; 1159e201c85SYohann data.t_id_z = threadIdx.z; 1169e201c85SYohann data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 1179e201c85SYohann data.slice = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1); 1189e201c85SYohann 1199e201c85SYohann CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)]; 1209e201c85SYohann CeedScalar r_V[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)]; 1219e201c85SYohann 1229e201c85SYohann for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 1239e201c85SYohann if (BASIS_DIM == 1) { 1249e201c85SYohann ReadElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, d_U, r_U); 1259e201c85SYohann Grad1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 1269e201c85SYohann WriteElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_V, d_V); 1279e201c85SYohann } else if (BASIS_DIM == 2) { 1289e201c85SYohann 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); 1299e201c85SYohann GradTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 1302b730f8bSJeremy 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, 1312b730f8bSJeremy L Thompson d_V); 1329e201c85SYohann } else if (BASIS_DIM == 3) { 1332b730f8bSJeremy L Thompson ReadElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem, 1342b730f8bSJeremy L Thompson BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, d_U, r_U); 1359e201c85SYohann if (BASIS_HAS_COLLOCATED_GRAD) GradTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 1369e201c85SYohann else GradTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 1372b730f8bSJeremy 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, 1382b730f8bSJeremy L Thompson BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, r_V, d_V); 1399e201c85SYohann } 1409e201c85SYohann } 1419e201c85SYohann } 1429e201c85SYohann 1432b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_GRAD_BLOCK_SIZE) __global__ 1442b730f8bSJeremy L Thompson void GradTranspose(const CeedInt num_elem, const CeedScalar *d_interp_1d, const CeedScalar *d_grad_1d, const CeedScalar *__restrict__ d_U, 1459e201c85SYohann CeedScalar *__restrict__ d_V) { 1469e201c85SYohann extern __shared__ CeedScalar slice[]; 147*b2165e7aSSebastian Grimberg 1489e201c85SYohann // load interp_1d and grad_1d into shared memory 1499e201c85SYohann __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D]; 1509e201c85SYohann loadMatrix<BASIS_P_1D * BASIS_Q_1D>(d_interp_1d, s_B); 1519e201c85SYohann __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)]; 1529e201c85SYohann loadMatrix<BASIS_Q_1D *(BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)>(d_grad_1d, s_G); 1539e201c85SYohann __syncthreads(); 1549e201c85SYohann 1559e201c85SYohann SharedData_Hip data; 1569e201c85SYohann data.t_id_x = threadIdx.x; 1579e201c85SYohann data.t_id_y = threadIdx.y; 1589e201c85SYohann data.t_id_z = threadIdx.z; 1599e201c85SYohann data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 1609e201c85SYohann data.slice = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1); 1619e201c85SYohann 1629e201c85SYohann CeedScalar r_U[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)]; 1639e201c85SYohann CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)]; 1649e201c85SYohann 1659e201c85SYohann for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 1669e201c85SYohann if (BASIS_DIM == 1) { 1679e201c85SYohann ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U); 1689e201c85SYohann GradTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 1699e201c85SYohann WriteElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V); 1709e201c85SYohann } else if (BASIS_DIM == 2) { 1712b730f8bSJeremy 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, 1722b730f8bSJeremy L Thompson r_U); 1739e201c85SYohann GradTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 1749e201c85SYohann 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); 1759e201c85SYohann } else if (BASIS_DIM == 3) { 1762b730f8bSJeremy 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, 1772b730f8bSJeremy L Thompson BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U); 1789e201c85SYohann if (BASIS_HAS_COLLOCATED_GRAD) GradTransposeTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 1799e201c85SYohann else GradTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 1802b730f8bSJeremy L Thompson WriteElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem, 1812b730f8bSJeremy L Thompson BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V); 1829e201c85SYohann } 1839e201c85SYohann } 1849e201c85SYohann } 1859e201c85SYohann 1869e201c85SYohann //------------------------------------------------------------------------------ 1879e201c85SYohann // Weight kernels by dim 1889e201c85SYohann //------------------------------------------------------------------------------ 1892b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_WEIGHT_BLOCK_SIZE) __global__ 1902b730f8bSJeremy L Thompson void Weight(const CeedInt num_elem, const CeedScalar *__restrict__ q_weight_1d, CeedScalar *__restrict__ d_W) { 1919e201c85SYohann extern __shared__ CeedScalar slice[]; 1929e201c85SYohann 1939e201c85SYohann SharedData_Hip data; 1949e201c85SYohann data.t_id_x = threadIdx.x; 1959e201c85SYohann data.t_id_y = threadIdx.y; 1969e201c85SYohann data.t_id_z = threadIdx.z; 1979e201c85SYohann data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 1989e201c85SYohann data.slice = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1); 1999e201c85SYohann 2009e201c85SYohann CeedScalar r_W[BASIS_DIM > 2 ? BASIS_Q_1D : 1]; 2019e201c85SYohann 2029e201c85SYohann for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 2039e201c85SYohann if (BASIS_DIM == 1) { 2049e201c85SYohann Weight1d<BASIS_Q_1D>(data, q_weight_1d, r_W); 2059e201c85SYohann WriteElementStrided1d<1, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_W, d_W); 2069e201c85SYohann } else if (BASIS_DIM == 2) { 2079e201c85SYohann WeightTensor2d<BASIS_Q_1D>(data, q_weight_1d, r_W); 2089e201c85SYohann 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); 2099e201c85SYohann } else if (BASIS_DIM == 3) { 2109e201c85SYohann WeightTensor3d<BASIS_Q_1D>(data, q_weight_1d, r_W); 2112b730f8bSJeremy 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, 2122b730f8bSJeremy L Thompson d_W); 2139e201c85SYohann } 2149e201c85SYohann } 2159e201c85SYohann } 2169e201c85SYohann 2179e201c85SYohann //------------------------------------------------------------------------------ 2189e201c85SYohann 2199e201c85SYohann #endif 220