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 basis read/write templates 10c0b5abf0SJeremy L Thompson #include <ceed/types.h> 119e201c85SYohann 129e201c85SYohann //------------------------------------------------------------------------------ 139e201c85SYohann // Helper function: load matrices for basis actions 149e201c85SYohann //------------------------------------------------------------------------------ 159e201c85SYohann template <int SIZE> 169e201c85SYohann inline __device__ void loadMatrix(const CeedScalar *d_B, CeedScalar *B) { 179e201c85SYohann CeedInt tid = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 18672b0f2aSSebastian Grimberg 192b730f8bSJeremy L Thompson for (CeedInt i = tid; i < SIZE; i += blockDim.x * blockDim.y * blockDim.z) B[i] = d_B[i]; 209e201c85SYohann } 219e201c85SYohann 229e201c85SYohann //------------------------------------------------------------------------------ 239e201c85SYohann // 1D 249e201c85SYohann //------------------------------------------------------------------------------ 259e201c85SYohann 269e201c85SYohann //------------------------------------------------------------------------------ 279e201c85SYohann // E-vector -> single element 289e201c85SYohann //------------------------------------------------------------------------------ 299e201c85SYohann template <int NUM_COMP, int P_1D> 302b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 312b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 329e201c85SYohann if (data.t_id_x < P_1D) { 339e201c85SYohann const CeedInt node = data.t_id_x; 349e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 35672b0f2aSSebastian Grimberg 369e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 379e201c85SYohann r_u[comp] = d_u[ind + comp * strides_comp]; 389e201c85SYohann } 399e201c85SYohann } 409e201c85SYohann } 419e201c85SYohann 429e201c85SYohann //------------------------------------------------------------------------------ 439e201c85SYohann // Single element -> E-vector 449e201c85SYohann //------------------------------------------------------------------------------ 459e201c85SYohann template <int NUM_COMP, int P_1D> 462b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 472b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 489e201c85SYohann if (data.t_id_x < P_1D) { 499e201c85SYohann const CeedInt node = data.t_id_x; 509e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 51672b0f2aSSebastian Grimberg 529e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 539e201c85SYohann d_v[ind + comp * strides_comp] = r_v[comp]; 549e201c85SYohann } 559e201c85SYohann } 569e201c85SYohann } 579e201c85SYohann 58db2becc9SJeremy L Thompson template <int NUM_COMP, int P_1D> 59db2becc9SJeremy L Thompson inline __device__ void SumElementStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 60db2becc9SJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 61db2becc9SJeremy L Thompson if (data.t_id_x < P_1D) { 62db2becc9SJeremy L Thompson const CeedInt node = data.t_id_x; 63db2becc9SJeremy L Thompson const CeedInt ind = node * strides_node + elem * strides_elem; 64db2becc9SJeremy L Thompson 65db2becc9SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 66db2becc9SJeremy L Thompson d_v[ind + comp * strides_comp] += r_v[comp]; 67db2becc9SJeremy L Thompson } 68db2becc9SJeremy L Thompson } 69db2becc9SJeremy L Thompson } 70db2becc9SJeremy L Thompson 719e201c85SYohann //------------------------------------------------------------------------------ 729e201c85SYohann // 2D 739e201c85SYohann //------------------------------------------------------------------------------ 749e201c85SYohann 759e201c85SYohann //------------------------------------------------------------------------------ 769e201c85SYohann // E-vector -> single element 779e201c85SYohann //------------------------------------------------------------------------------ 789e201c85SYohann template <int NUM_COMP, int P_1D> 792b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 802b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 819e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 829e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 839e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 84672b0f2aSSebastian Grimberg 859e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 869e201c85SYohann r_u[comp] = d_u[ind + comp * strides_comp]; 879e201c85SYohann } 889e201c85SYohann } 899e201c85SYohann } 909e201c85SYohann 919e201c85SYohann //------------------------------------------------------------------------------ 929e201c85SYohann // Single element -> E-vector 939e201c85SYohann //------------------------------------------------------------------------------ 949e201c85SYohann template <int NUM_COMP, int P_1D> 952b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 962b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 979e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 989e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 999e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 100672b0f2aSSebastian Grimberg 1019e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 1029e201c85SYohann d_v[ind + comp * strides_comp] = r_v[comp]; 1039e201c85SYohann } 1049e201c85SYohann } 1059e201c85SYohann } 1069e201c85SYohann 107db2becc9SJeremy L Thompson template <int NUM_COMP, int P_1D> 108db2becc9SJeremy L Thompson inline __device__ void SumElementStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 109db2becc9SJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 110db2becc9SJeremy L Thompson if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 111db2becc9SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 112db2becc9SJeremy L Thompson const CeedInt ind = node * strides_node + elem * strides_elem; 113db2becc9SJeremy L Thompson 114db2becc9SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 115db2becc9SJeremy L Thompson d_v[ind + comp * strides_comp] += r_v[comp]; 116db2becc9SJeremy L Thompson } 117db2becc9SJeremy L Thompson } 118db2becc9SJeremy L Thompson } 119db2becc9SJeremy L Thompson 1209e201c85SYohann //------------------------------------------------------------------------------ 1219e201c85SYohann // 3D 1229e201c85SYohann //------------------------------------------------------------------------------ 1239e201c85SYohann 1249e201c85SYohann //------------------------------------------------------------------------------ 1259e201c85SYohann // E-vector -> single element 1269e201c85SYohann //------------------------------------------------------------------------------ 1279e201c85SYohann template <int NUM_COMP, int P_1D> 1282b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 1292b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 1309e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1319e201c85SYohann for (CeedInt z = 0; z < P_1D; z++) { 1329e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 1339e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 134672b0f2aSSebastian Grimberg 1359e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 1369e201c85SYohann r_u[z + comp * P_1D] = d_u[ind + comp * strides_comp]; 1379e201c85SYohann } 1389e201c85SYohann } 1399e201c85SYohann } 1409e201c85SYohann } 1419e201c85SYohann 1429e201c85SYohann //------------------------------------------------------------------------------ 1439e201c85SYohann // Single element -> E-vector 1449e201c85SYohann //------------------------------------------------------------------------------ 1459e201c85SYohann template <int NUM_COMP, int P_1D> 1462b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 1472b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 1489e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1499e201c85SYohann for (CeedInt z = 0; z < P_1D; z++) { 1509e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 1519e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 152672b0f2aSSebastian Grimberg 1539e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 1549e201c85SYohann d_v[ind + comp * strides_comp] = r_v[z + comp * P_1D]; 1559e201c85SYohann } 1569e201c85SYohann } 1579e201c85SYohann } 1589e201c85SYohann } 159db2becc9SJeremy L Thompson 160db2becc9SJeremy L Thompson template <int NUM_COMP, int P_1D> 161db2becc9SJeremy L Thompson inline __device__ void SumElementStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 162db2becc9SJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 163db2becc9SJeremy L Thompson if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 164db2becc9SJeremy L Thompson for (CeedInt z = 0; z < P_1D; z++) { 165db2becc9SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 166db2becc9SJeremy L Thompson const CeedInt ind = node * strides_node + elem * strides_elem; 167db2becc9SJeremy L Thompson 168db2becc9SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 169db2becc9SJeremy L Thompson d_v[ind + comp * strides_comp] += r_v[z + comp * P_1D]; 170db2becc9SJeremy L Thompson } 171db2becc9SJeremy L Thompson } 172db2becc9SJeremy L Thompson } 173db2becc9SJeremy L Thompson } 174*b6a2eb79SJeremy L Thompson 175*b6a2eb79SJeremy L Thompson //------------------------------------------------------------------------------ 176*b6a2eb79SJeremy L Thompson // AtPoints 177*b6a2eb79SJeremy L Thompson //------------------------------------------------------------------------------ 178*b6a2eb79SJeremy L Thompson 179*b6a2eb79SJeremy L Thompson //------------------------------------------------------------------------------ 180*b6a2eb79SJeremy L Thompson // E-vector -> single point 181*b6a2eb79SJeremy L Thompson //------------------------------------------------------------------------------ 182*b6a2eb79SJeremy L Thompson template <int NUM_COMP, int NUM_PTS> 183*b6a2eb79SJeremy L Thompson inline __device__ void ReadPoint(SharedData_Hip &data, const CeedInt elem, const CeedInt p, const CeedInt points_in_elem, const CeedInt strides_point, 184*b6a2eb79SJeremy L Thompson const CeedInt strides_comp, const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 185*b6a2eb79SJeremy L Thompson const CeedInt ind = (p % NUM_PTS) * strides_point + elem * strides_elem; 186*b6a2eb79SJeremy L Thompson 187*b6a2eb79SJeremy L Thompson if (p < points_in_elem) { 188*b6a2eb79SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 189*b6a2eb79SJeremy L Thompson r_u[comp] = d_u[ind + comp * strides_comp]; 190*b6a2eb79SJeremy L Thompson } 191*b6a2eb79SJeremy L Thompson } else { 192*b6a2eb79SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 193*b6a2eb79SJeremy L Thompson r_u[comp] = 0.0; 194*b6a2eb79SJeremy L Thompson } 195*b6a2eb79SJeremy L Thompson } 196*b6a2eb79SJeremy L Thompson } 197*b6a2eb79SJeremy L Thompson 198*b6a2eb79SJeremy L Thompson //------------------------------------------------------------------------------ 199*b6a2eb79SJeremy L Thompson // Single point -> E-vector 200*b6a2eb79SJeremy L Thompson //------------------------------------------------------------------------------ 201*b6a2eb79SJeremy L Thompson template <int NUM_COMP, int NUM_PTS> 202*b6a2eb79SJeremy L Thompson inline __device__ void WritePoint(SharedData_Hip &data, const CeedInt elem, const CeedInt p, const CeedInt points_in_elem, 203*b6a2eb79SJeremy L Thompson const CeedInt strides_point, const CeedInt strides_comp, const CeedInt strides_elem, const CeedScalar *r_v, 204*b6a2eb79SJeremy L Thompson CeedScalar *d_v) { 205*b6a2eb79SJeremy L Thompson if (p < points_in_elem) { 206*b6a2eb79SJeremy L Thompson const CeedInt ind = (p % NUM_PTS) * strides_point + elem * strides_elem; 207*b6a2eb79SJeremy L Thompson 208*b6a2eb79SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 209*b6a2eb79SJeremy L Thompson d_v[ind + comp * strides_comp] = r_v[comp]; 210*b6a2eb79SJeremy L Thompson } 211*b6a2eb79SJeremy L Thompson } 212*b6a2eb79SJeremy L Thompson } 213