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 CUDA shared memory tensor product basis 109e201c85SYohann 11*c0b5abf0SJeremy L Thompson #include <ceed/types.h> 122b730f8bSJeremy L Thompson 139e201c85SYohann #include "cuda-shared-basis-read-write-templates.h" 149e201c85SYohann #include "cuda-shared-basis-tensor-templates.h" 159e201c85SYohann 169e201c85SYohann //------------------------------------------------------------------------------ 179e201c85SYohann // Interp kernel by dim 189e201c85SYohann //------------------------------------------------------------------------------ 192b730f8bSJeremy L Thompson extern "C" __global__ void Interp(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) { 209e201c85SYohann extern __shared__ CeedScalar slice[]; 219e201c85SYohann 229e201c85SYohann SharedData_Cuda 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 329e201c85SYohann for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 339e201c85SYohann if (BASIS_DIM == 1) { 349e201c85SYohann ReadElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, d_U, r_U); 359e201c85SYohann Interp1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, r_V); 369e201c85SYohann WriteElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_V, d_V); 379e201c85SYohann } else if (BASIS_DIM == 2) { 389e201c85SYohann 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); 399e201c85SYohann InterpTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, r_V); 409e201c85SYohann 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); 419e201c85SYohann } else if (BASIS_DIM == 3) { 422b730f8bSJeremy L Thompson ReadElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem, 432b730f8bSJeremy L Thompson BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, d_U, r_U); 449e201c85SYohann InterpTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, r_V); 452b730f8bSJeremy L Thompson WriteElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem, 462b730f8bSJeremy L Thompson BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, r_V, d_V); 479e201c85SYohann } 489e201c85SYohann } 499e201c85SYohann } 509e201c85SYohann 512b730f8bSJeremy L Thompson extern "C" __global__ void InterpTranspose(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *__restrict__ d_U, 529e201c85SYohann CeedScalar *__restrict__ d_V) { 539e201c85SYohann extern __shared__ CeedScalar slice[]; 549e201c85SYohann 559e201c85SYohann SharedData_Cuda data; 569e201c85SYohann data.t_id_x = threadIdx.x; 579e201c85SYohann data.t_id_y = threadIdx.y; 589e201c85SYohann data.t_id_z = threadIdx.z; 599e201c85SYohann data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 609e201c85SYohann data.slice = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1); 619e201c85SYohann 629e201c85SYohann CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)]; 639e201c85SYohann CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)]; 649e201c85SYohann 659e201c85SYohann for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 669e201c85SYohann if (BASIS_DIM == 1) { 679e201c85SYohann ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U); 689e201c85SYohann InterpTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, r_V); 699e201c85SYohann WriteElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V); 709e201c85SYohann } else if (BASIS_DIM == 2) { 719e201c85SYohann 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); 729e201c85SYohann InterpTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, r_V); 739e201c85SYohann 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); 749e201c85SYohann } else if (BASIS_DIM == 3) { 752b730f8bSJeremy L Thompson ReadElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem, 762b730f8bSJeremy L Thompson BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U); 779e201c85SYohann InterpTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, r_V); 782b730f8bSJeremy L Thompson WriteElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem, 792b730f8bSJeremy L Thompson BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V); 809e201c85SYohann } 819e201c85SYohann } 829e201c85SYohann } 839e201c85SYohann 84db2becc9SJeremy L Thompson extern "C" __global__ void InterpTransposeAdd(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *__restrict__ d_U, 85db2becc9SJeremy L Thompson CeedScalar *__restrict__ d_V) { 86db2becc9SJeremy L Thompson extern __shared__ CeedScalar slice[]; 87db2becc9SJeremy L Thompson 88db2becc9SJeremy L Thompson SharedData_Cuda data; 89db2becc9SJeremy L Thompson data.t_id_x = threadIdx.x; 90db2becc9SJeremy L Thompson data.t_id_y = threadIdx.y; 91db2becc9SJeremy L Thompson data.t_id_z = threadIdx.z; 92db2becc9SJeremy L Thompson data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 93db2becc9SJeremy L Thompson data.slice = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1); 94db2becc9SJeremy L Thompson 95db2becc9SJeremy L Thompson CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)]; 96db2becc9SJeremy L Thompson CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)]; 97db2becc9SJeremy L Thompson 98db2becc9SJeremy L Thompson for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 99db2becc9SJeremy L Thompson if (BASIS_DIM == 1) { 100db2becc9SJeremy L Thompson ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U); 101db2becc9SJeremy L Thompson InterpTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, r_V); 102db2becc9SJeremy L Thompson SumElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V); 103db2becc9SJeremy L Thompson } else if (BASIS_DIM == 2) { 104db2becc9SJeremy 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); 105db2becc9SJeremy L Thompson InterpTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, r_V); 106db2becc9SJeremy 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); 107db2becc9SJeremy L Thompson } else if (BASIS_DIM == 3) { 108db2becc9SJeremy L Thompson ReadElementStrided3d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D * num_elem, 109db2becc9SJeremy L Thompson BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U); 110db2becc9SJeremy L Thompson InterpTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, r_V); 111db2becc9SJeremy L Thompson SumElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem, 112db2becc9SJeremy L Thompson BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V); 113db2becc9SJeremy L Thompson } 114db2becc9SJeremy L Thompson } 115db2becc9SJeremy L Thompson } 116db2becc9SJeremy L Thompson 1179e201c85SYohann //------------------------------------------------------------------------------ 1189e201c85SYohann // Grad kernel by dim 1199e201c85SYohann //------------------------------------------------------------------------------ 1202b730f8bSJeremy L Thompson extern "C" __global__ void Grad(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *c_G, const CeedScalar *__restrict__ d_U, 1219e201c85SYohann CeedScalar *__restrict__ d_V) { 1229e201c85SYohann extern __shared__ CeedScalar slice[]; 1239e201c85SYohann 1249e201c85SYohann SharedData_Cuda data; 1259e201c85SYohann data.t_id_x = threadIdx.x; 1269e201c85SYohann data.t_id_y = threadIdx.y; 1279e201c85SYohann data.t_id_z = threadIdx.z; 1289e201c85SYohann data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 1299e201c85SYohann data.slice = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1); 1309e201c85SYohann 1319e201c85SYohann CeedScalar r_U[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)]; 1329e201c85SYohann CeedScalar r_V[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)]; 1339e201c85SYohann 1349e201c85SYohann for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 1359e201c85SYohann if (BASIS_DIM == 1) { 1369e201c85SYohann ReadElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, d_U, r_U); 1379e201c85SYohann Grad1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, c_G, r_V); 1389e201c85SYohann WriteElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_V, d_V); 1399e201c85SYohann } else if (BASIS_DIM == 2) { 1409e201c85SYohann 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); 1419e201c85SYohann GradTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, c_G, r_V); 1422b730f8bSJeremy 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, 1432b730f8bSJeremy L Thompson d_V); 1449e201c85SYohann } else if (BASIS_DIM == 3) { 1452b730f8bSJeremy L Thompson ReadElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem, 1462b730f8bSJeremy L Thompson BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, d_U, r_U); 1479e201c85SYohann if (BASIS_HAS_COLLOCATED_GRAD) GradTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, c_G, r_V); 1489e201c85SYohann else GradTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, c_G, r_V); 1492b730f8bSJeremy 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, 1502b730f8bSJeremy L Thompson BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, r_V, d_V); 1519e201c85SYohann } 1529e201c85SYohann } 1539e201c85SYohann } 1549e201c85SYohann 1552b730f8bSJeremy L Thompson extern "C" __global__ void GradTranspose(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *c_G, const CeedScalar *__restrict__ d_U, 1569e201c85SYohann CeedScalar *__restrict__ d_V) { 1579e201c85SYohann extern __shared__ CeedScalar slice[]; 1589e201c85SYohann 1599e201c85SYohann SharedData_Cuda data; 1609e201c85SYohann data.t_id_x = threadIdx.x; 1619e201c85SYohann data.t_id_y = threadIdx.y; 1629e201c85SYohann data.t_id_z = threadIdx.z; 1639e201c85SYohann data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 1649e201c85SYohann data.slice = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1); 1659e201c85SYohann 1669e201c85SYohann CeedScalar r_U[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)]; 1679e201c85SYohann CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)]; 1689e201c85SYohann 1699e201c85SYohann for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 1709e201c85SYohann if (BASIS_DIM == 1) { 1719e201c85SYohann ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U); 1729e201c85SYohann GradTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, c_G, r_V); 1739e201c85SYohann WriteElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V); 1749e201c85SYohann } else if (BASIS_DIM == 2) { 1752b730f8bSJeremy 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, 1762b730f8bSJeremy L Thompson r_U); 1779e201c85SYohann GradTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, c_G, r_V); 1789e201c85SYohann 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); 1799e201c85SYohann } else if (BASIS_DIM == 3) { 1802b730f8bSJeremy 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, 1812b730f8bSJeremy L Thompson BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U); 1829e201c85SYohann if (BASIS_HAS_COLLOCATED_GRAD) GradTransposeTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, c_G, r_V); 1839e201c85SYohann else GradTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, c_G, r_V); 1842b730f8bSJeremy L Thompson WriteElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem, 1852b730f8bSJeremy L Thompson BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V); 1869e201c85SYohann } 1879e201c85SYohann } 1889e201c85SYohann } 1899e201c85SYohann 190db2becc9SJeremy L Thompson extern "C" __global__ void GradTransposeAdd(const CeedInt num_elem, const CeedScalar *c_B, const CeedScalar *c_G, const CeedScalar *__restrict__ d_U, 191db2becc9SJeremy L Thompson CeedScalar *__restrict__ d_V) { 192db2becc9SJeremy L Thompson extern __shared__ CeedScalar slice[]; 193db2becc9SJeremy L Thompson 194db2becc9SJeremy L Thompson SharedData_Cuda data; 195db2becc9SJeremy L Thompson data.t_id_x = threadIdx.x; 196db2becc9SJeremy L Thompson data.t_id_y = threadIdx.y; 197db2becc9SJeremy L Thompson data.t_id_z = threadIdx.z; 198db2becc9SJeremy L Thompson data.t_id = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 199db2becc9SJeremy L Thompson data.slice = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1); 200db2becc9SJeremy L Thompson 201db2becc9SJeremy L Thompson CeedScalar r_U[BASIS_NUM_COMP * BASIS_DIM * (BASIS_DIM > 2 ? BASIS_Q_1D : 1)]; 202db2becc9SJeremy L Thompson CeedScalar r_V[BASIS_NUM_COMP * (BASIS_DIM > 2 ? BASIS_P_1D : 1)]; 203db2becc9SJeremy L Thompson 204db2becc9SJeremy L Thompson for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 205db2becc9SJeremy L Thompson if (BASIS_DIM == 1) { 206db2becc9SJeremy L Thompson ReadElementStrided1d<BASIS_NUM_COMP, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, d_U, r_U); 207db2becc9SJeremy L Thompson GradTranspose1d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, c_G, r_V); 208db2becc9SJeremy L Thompson SumElementStrided1d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * num_elem, BASIS_P_1D, r_V, d_V); 209db2becc9SJeremy L Thompson } else if (BASIS_DIM == 2) { 210db2becc9SJeremy 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, 211db2becc9SJeremy L Thompson r_U); 212db2becc9SJeremy L Thompson GradTransposeTensor2d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, c_G, r_V); 213db2becc9SJeremy 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); 214db2becc9SJeremy L Thompson } else if (BASIS_DIM == 3) { 215db2becc9SJeremy 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, 216db2becc9SJeremy L Thompson BASIS_Q_1D * BASIS_Q_1D * BASIS_Q_1D, d_U, r_U); 217db2becc9SJeremy L Thompson if (BASIS_HAS_COLLOCATED_GRAD) GradTransposeTensorCollocated3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, c_G, r_V); 218db2becc9SJeremy L Thompson else GradTransposeTensor3d<BASIS_NUM_COMP, BASIS_P_1D, BASIS_Q_1D>(data, r_U, c_B, c_G, r_V); 219db2becc9SJeremy L Thompson SumElementStrided3d<BASIS_NUM_COMP, BASIS_P_1D>(data, elem, 1, BASIS_P_1D * BASIS_P_1D * BASIS_P_1D * num_elem, 220db2becc9SJeremy L Thompson BASIS_P_1D * BASIS_P_1D * BASIS_P_1D, r_V, d_V); 221db2becc9SJeremy L Thompson } 222db2becc9SJeremy L Thompson } 223db2becc9SJeremy L Thompson } 224db2becc9SJeremy L Thompson 2259e201c85SYohann //------------------------------------------------------------------------------ 2269e201c85SYohann // Weight kernels by dim 2279e201c85SYohann //------------------------------------------------------------------------------ 2282b730f8bSJeremy L Thompson extern "C" __global__ void Weight(const CeedInt num_elem, const CeedScalar *__restrict__ q_weight_1d, CeedScalar *__restrict__ d_W) { 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; 2369e201c85SYohann data.slice = slice + data.t_id_z * T_1D * (BASIS_DIM > 1 ? T_1D : 1); 2379e201c85SYohann 2389e201c85SYohann CeedScalar r_W[BASIS_DIM > 2 ? BASIS_Q_1D : 1]; 2399e201c85SYohann 2409e201c85SYohann for (CeedInt elem = blockIdx.x * blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x * blockDim.z) { 2419e201c85SYohann if (BASIS_DIM == 1) { 2429e201c85SYohann Weight1d<BASIS_Q_1D>(data, q_weight_1d, r_W); 2439e201c85SYohann WriteElementStrided1d<1, BASIS_Q_1D>(data, elem, 1, BASIS_Q_1D * num_elem, BASIS_Q_1D, r_W, d_W); 2449e201c85SYohann } else if (BASIS_DIM == 2) { 2459e201c85SYohann WeightTensor2d<BASIS_Q_1D>(data, q_weight_1d, r_W); 2469e201c85SYohann 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); 2479e201c85SYohann } else if (BASIS_DIM == 3) { 2489e201c85SYohann WeightTensor3d<BASIS_Q_1D>(data, q_weight_1d, r_W); 2492b730f8bSJeremy 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, 2502b730f8bSJeremy L Thompson d_W); 2519e201c85SYohann } 2529e201c85SYohann } 2539e201c85SYohann } 254