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 10c0b5abf0SJeremy L Thompson #include <ceed/types.h> 112b730f8bSJeremy L Thompson 129e201c85SYohann #include "hip-shared-basis-read-write-templates.h" 139e201c85SYohann #include "hip-shared-basis-tensor-templates.h" 149e201c85SYohann 159e201c85SYohann //------------------------------------------------------------------------------ 169e201c85SYohann // Interp kernel by dim 179e201c85SYohann //------------------------------------------------------------------------------ 182b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_INTERP_BLOCK_SIZE) __global__ 19*aa4002adSJeremy L Thompson void Interp(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) { 209e201c85SYohann extern __shared__ CeedScalar slice[]; 21b2165e7aSSebastian Grimberg 229e201c85SYohann SharedData_Hip data; 239e201c85SYohann data.t_id_x = threadIdx.x; 249e201c85SYohann data.t_id_y = threadIdx.y; 259e201c85SYohann data.t_id_z = threadIdx.z; 269e201c85SYohann data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 279e201c85SYohann data.slice = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1); 289e201c85SYohann 299e201c85SYohann CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)]; 309e201c85SYohann CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)]; 319e201c85SYohann 32*aa4002adSJeremy L Thompson // load interp_1d into shared memory 33*aa4002adSJeremy L Thompson __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D]; 34*aa4002adSJeremy L Thompson LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B); 35*aa4002adSJeremy L Thompson __syncthreads(); 36*aa4002adSJeremy L Thompson 37*aa4002adSJeremy L Thompson // Apply basis element by element 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__ 58*aa4002adSJeremy L Thompson void InterpTranspose(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) { 599e201c85SYohann extern __shared__ CeedScalar slice[]; 60b2165e7aSSebastian Grimberg 619e201c85SYohann SharedData_Hip data; 629e201c85SYohann data.t_id_x = threadIdx.x; 639e201c85SYohann data.t_id_y = threadIdx.y; 649e201c85SYohann data.t_id_z = threadIdx.z; 659e201c85SYohann data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 669e201c85SYohann data.slice = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1); 679e201c85SYohann 689e201c85SYohann CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)]; 699e201c85SYohann CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)]; 709e201c85SYohann 71*aa4002adSJeremy L Thompson // load interp_1d into shared memory 72*aa4002adSJeremy L Thompson __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D]; 73*aa4002adSJeremy L Thompson LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B); 74*aa4002adSJeremy L Thompson __syncthreads(); 75*aa4002adSJeremy L Thompson 76*aa4002adSJeremy L Thompson // Apply basis element by element 779e201c85SYohann for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 789e201c85SYohann if (BASIS_DIM == 1) { 799e201c85SYohann ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U); 809e201c85SYohann InterpTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V); 819e201c85SYohann WriteElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V); 829e201c85SYohann } else if (BASIS_DIM == 2) { 839e201c85SYohann 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); 849e201c85SYohann InterpTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V); 859e201c85SYohann 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); 869e201c85SYohann } else if (BASIS_DIM == 3) { 872b730f8bSJeremy L Thompson ReadElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem, 882b730f8bSJeremy L Thompson BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U); 899e201c85SYohann InterpTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V); 902b730f8bSJeremy L Thompson WriteElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem, 912b730f8bSJeremy L Thompson BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V); 929e201c85SYohann } 939e201c85SYohann } 949e201c85SYohann } 959e201c85SYohann 96db2becc9SJeremy L Thompson extern "C" __launch_bounds__(BASIS_INTERP_BLOCK_SIZE) __global__ 97*aa4002adSJeremy L Thompson void InterpTransposeAdd(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) { 98db2becc9SJeremy L Thompson extern __shared__ CeedScalar slice[]; 99db2becc9SJeremy L Thompson 100db2becc9SJeremy L Thompson SharedData_Hip data; 101db2becc9SJeremy L Thompson data.t_id_x = threadIdx.x; 102db2becc9SJeremy L Thompson data.t_id_y = threadIdx.y; 103db2becc9SJeremy L Thompson data.t_id_z = threadIdx.z; 104db2becc9SJeremy L Thompson data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 105db2becc9SJeremy L Thompson data.slice = &slice[data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1)]; 106db2becc9SJeremy L Thompson 107db2becc9SJeremy L Thompson CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)]; 108db2becc9SJeremy L Thompson CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)]; 109db2becc9SJeremy L Thompson 110*aa4002adSJeremy L Thompson // load interp_1d into shared memory 111*aa4002adSJeremy L Thompson __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D]; 112*aa4002adSJeremy L Thompson LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B); 113*aa4002adSJeremy L Thompson __syncthreads(); 114*aa4002adSJeremy L Thompson 115*aa4002adSJeremy L Thompson // Apply basis element by element 116db2becc9SJeremy L Thompson for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 117db2becc9SJeremy L Thompson if (BASIS_DIM == 1) { 118db2becc9SJeremy L Thompson ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U); 119db2becc9SJeremy L Thompson InterpTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V); 120db2becc9SJeremy L Thompson SumElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V); 121db2becc9SJeremy L Thompson } else if (BASIS_DIM == 2) { 122db2becc9SJeremy 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); 123db2becc9SJeremy L Thompson InterpTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V); 124db2becc9SJeremy 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); 125db2becc9SJeremy L Thompson } else if (BASIS_DIM == 3) { 126db2becc9SJeremy L Thompson ReadElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem, 127db2becc9SJeremy L Thompson BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U); 128db2becc9SJeremy L Thompson InterpTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, r_V); 129db2becc9SJeremy L Thompson SumElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem, 130db2becc9SJeremy L Thompson BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V); 131db2becc9SJeremy L Thompson } 132db2becc9SJeremy L Thompson } 133db2becc9SJeremy L Thompson } 134db2becc9SJeremy L Thompson 1359e201c85SYohann //------------------------------------------------------------------------------ 1369e201c85SYohann // Grad kernel by dim 1379e201c85SYohann //------------------------------------------------------------------------------ 138*aa4002adSJeremy L Thompson extern "C" __launch_bounds__(BASIS_GRAD_BLOCK_SIZE) __global__ void Grad(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *c_G, 139*aa4002adSJeremy L Thompson const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) { 1409e201c85SYohann extern __shared__ CeedScalar slice[]; 141b2165e7aSSebastian Grimberg 1429e201c85SYohann SharedData_Hip data; 1439e201c85SYohann data.t_id_x = threadIdx.x; 1449e201c85SYohann data.t_id_y = threadIdx.y; 1459e201c85SYohann data.t_id_z = threadIdx.z; 1469e201c85SYohann data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 1479e201c85SYohann data.slice = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1); 1489e201c85SYohann 1499e201c85SYohann CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)]; 1509e201c85SYohann CeedScalar r_V[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)]; 1519e201c85SYohann 152*aa4002adSJeremy L Thompson // load interp_1d and grad_1d into shared memory 153*aa4002adSJeremy L Thompson __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D]; 154*aa4002adSJeremy L Thompson LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B); 155*aa4002adSJeremy L Thompson __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)]; 156*aa4002adSJeremy L Thompson LoadMatrix<BASIS_Q_1D, BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D>(data, c_G, s_G); 157*aa4002adSJeremy L Thompson __syncthreads(); 158*aa4002adSJeremy L Thompson 159*aa4002adSJeremy L Thompson // Apply basis element by element 1609e201c85SYohann for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 1619e201c85SYohann if (BASIS_DIM == 1) { 1629e201c85SYohann ReadElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, d_U, r_U); 1639e201c85SYohann Grad1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 1649e201c85SYohann WriteElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_V, d_V); 1659e201c85SYohann } else if (BASIS_DIM == 2) { 1669e201c85SYohann 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); 1679e201c85SYohann GradTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 1682b730f8bSJeremy 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, 1692b730f8bSJeremy L Thompson d_V); 1709e201c85SYohann } else if (BASIS_DIM == 3) { 1712b730f8bSJeremy L Thompson ReadElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem, 1722b730f8bSJeremy L Thompson BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, d_U, r_U); 1739e201c85SYohann if (BASIS_HAS_COLLOCATED_GRAD) GradTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 1749e201c85SYohann else GradTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 1752b730f8bSJeremy 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, 1762b730f8bSJeremy L Thompson BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, r_V, d_V); 1779e201c85SYohann } 1789e201c85SYohann } 1799e201c85SYohann } 1809e201c85SYohann 1812b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_GRAD_BLOCK_SIZE) __global__ 182*aa4002adSJeremy L Thompson void GradTranspose(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *c_G, const CeedScalar *__restrict__ d_U, 1839e201c85SYohann CeedScalar *__restrict__ d_V) { 1849e201c85SYohann extern __shared__ CeedScalar slice[]; 185b2165e7aSSebastian Grimberg 1869e201c85SYohann SharedData_Hip data; 1879e201c85SYohann data.t_id_x = threadIdx.x; 1889e201c85SYohann data.t_id_y = threadIdx.y; 1899e201c85SYohann data.t_id_z = threadIdx.z; 1909e201c85SYohann data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 1919e201c85SYohann data.slice = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1); 1929e201c85SYohann 1939e201c85SYohann CeedScalar r_U[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)]; 1949e201c85SYohann CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)]; 1959e201c85SYohann 196*aa4002adSJeremy L Thompson // load interp_1d and grad_1d into shared memory 197*aa4002adSJeremy L Thompson __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D]; 198*aa4002adSJeremy L Thompson LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B); 199*aa4002adSJeremy L Thompson __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)]; 200*aa4002adSJeremy L Thompson LoadMatrix<BASIS_Q_1D, BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D>(data, c_G, s_G); 201*aa4002adSJeremy L Thompson __syncthreads(); 202*aa4002adSJeremy L Thompson 203*aa4002adSJeremy L Thompson // Apply basis element by element 2049e201c85SYohann for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 2059e201c85SYohann if (BASIS_DIM == 1) { 2069e201c85SYohann ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U); 2079e201c85SYohann GradTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 2089e201c85SYohann WriteElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V); 2099e201c85SYohann } else if (BASIS_DIM == 2) { 2102b730f8bSJeremy 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, 2112b730f8bSJeremy L Thompson r_U); 2129e201c85SYohann GradTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 2139e201c85SYohann 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); 2149e201c85SYohann } else if (BASIS_DIM == 3) { 2152b730f8bSJeremy 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, 2162b730f8bSJeremy L Thompson BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U); 2179e201c85SYohann if (BASIS_HAS_COLLOCATED_GRAD) GradTransposeTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 2189e201c85SYohann else GradTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 2192b730f8bSJeremy L Thompson WriteElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem, 2202b730f8bSJeremy L Thompson BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V); 2219e201c85SYohann } 2229e201c85SYohann } 2239e201c85SYohann } 2249e201c85SYohann 225db2becc9SJeremy L Thompson extern "C" __launch_bounds__(BASIS_GRAD_BLOCK_SIZE) __global__ 226*aa4002adSJeremy L Thompson void GradTransposeAdd(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *c_G, const CeedScalar *__restrict__ d_U, 227db2becc9SJeremy L Thompson CeedScalar *__restrict__ d_V) { 228db2becc9SJeremy L Thompson extern __shared__ CeedScalar slice[]; 229db2becc9SJeremy L Thompson 230db2becc9SJeremy L Thompson SharedData_Hip data; 231db2becc9SJeremy L Thompson data.t_id_x = threadIdx.x; 232db2becc9SJeremy L Thompson data.t_id_y = threadIdx.y; 233db2becc9SJeremy L Thompson data.t_id_z = threadIdx.z; 234db2becc9SJeremy L Thompson data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 235db2becc9SJeremy L Thompson data.slice = &slice[data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1)]; 236db2becc9SJeremy L Thompson 237db2becc9SJeremy L Thompson CeedScalar r_U[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)]; 238db2becc9SJeremy L Thompson CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)]; 239db2becc9SJeremy L Thompson 240*aa4002adSJeremy L Thompson // load interp_1d and grad_1d into shared memory 241*aa4002adSJeremy L Thompson __shared__ CeedScalar s_B[BASIS_P_1D * BASIS_Q_1D]; 242*aa4002adSJeremy L Thompson LoadMatrix<BASIS_P_1D, BASIS_Q_1D>(data, c_B, s_B); 243*aa4002adSJeremy L Thompson __shared__ CeedScalar s_G[BASIS_Q_1D * (BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D)]; 244*aa4002adSJeremy L Thompson LoadMatrix<BASIS_Q_1D, BASIS_HAS_COLLOCATED_GRAD ? BASIS_Q_1D : BASIS_P_1D>(data, c_G, s_G); 245*aa4002adSJeremy L Thompson __syncthreads(); 246*aa4002adSJeremy L Thompson 247*aa4002adSJeremy L Thompson // Apply basis element by element 248db2becc9SJeremy L Thompson for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 249db2becc9SJeremy L Thompson if (BASIS_DIM == 1) { 250db2becc9SJeremy L Thompson ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U); 251db2becc9SJeremy L Thompson GradTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 252db2becc9SJeremy L Thompson SumElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V); 253db2becc9SJeremy L Thompson } else if (BASIS_DIM == 2) { 254db2becc9SJeremy 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, 255db2becc9SJeremy L Thompson r_U); 256db2becc9SJeremy L Thompson GradTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 257db2becc9SJeremy 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); 258db2becc9SJeremy L Thompson } else if (BASIS_DIM == 3) { 259db2becc9SJeremy 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, 260db2becc9SJeremy L Thompson BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U); 261db2becc9SJeremy 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); 262db2becc9SJeremy L Thompson else GradTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, s_B, s_G, r_V); 263db2becc9SJeremy L Thompson SumElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem, 264db2becc9SJeremy L Thompson BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V); 265db2becc9SJeremy L Thompson } 266db2becc9SJeremy L Thompson } 267db2becc9SJeremy L Thompson } 268db2becc9SJeremy L Thompson 2699e201c85SYohann //------------------------------------------------------------------------------ 2709e201c85SYohann // Weight kernels by dim 2719e201c85SYohann //------------------------------------------------------------------------------ 2722b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BASIS_WEIGHT_BLOCK_SIZE) __global__ 2732b730f8bSJeremy L Thompson void Weight(const CeedInt num_elem, const CeedScalar *__restrict__ q_weight_1d, CeedScalar *__restrict__ d_W) { 2749e201c85SYohann extern __shared__ CeedScalar slice[]; 2759e201c85SYohann 2769e201c85SYohann SharedData_Hip data; 2779e201c85SYohann data.t_id_x = threadIdx.x; 2789e201c85SYohann data.t_id_y = threadIdx.y; 2799e201c85SYohann data.t_id_z = threadIdx.z; 2809e201c85SYohann data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 2819e201c85SYohann data.slice = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1); 2829e201c85SYohann 2839e201c85SYohann CeedScalar r_W[BASIS_DIM > 2 ? BASIS_Q_1D : 1]; 2849e201c85SYohann 2859e201c85SYohann for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 2869e201c85SYohann if (BASIS_DIM == 1) { 2879e201c85SYohann Weight1d<BASIS_Q_1D>(data, q_weight_1d, r_W); 2889e201c85SYohann WriteElementStrided1d<1, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_W, d_W); 2899e201c85SYohann } else if (BASIS_DIM == 2) { 2909e201c85SYohann WeightTensor2d<BASIS_Q_1D>(data, q_weight_1d, r_W); 2919e201c85SYohann 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); 2929e201c85SYohann } else if (BASIS_DIM == 3) { 2939e201c85SYohann WeightTensor3d<BASIS_Q_1D>(data, q_weight_1d, r_W); 2942b730f8bSJeremy 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, 2952b730f8bSJeremy L Thompson d_W); 2969e201c85SYohann } 2979e201c85SYohann } 2989e201c85SYohann } 299