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 109e201c85SYohann 11*c0b5abf0SJeremy L Thompson #include <ceed/types.h> 129e201c85SYohann 139e201c85SYohann //------------------------------------------------------------------------------ 149e201c85SYohann // Helper function: load matrices for basis actions 159e201c85SYohann //------------------------------------------------------------------------------ 169e201c85SYohann template <int SIZE> 179e201c85SYohann inline __device__ void loadMatrix(const CeedScalar *d_B, CeedScalar *B) { 189e201c85SYohann CeedInt tid = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 19672b0f2aSSebastian Grimberg 202b730f8bSJeremy L Thompson for (CeedInt i = tid; i < SIZE; i += blockDim.x * blockDim.y * blockDim.z) B[i] = d_B[i]; 219e201c85SYohann } 229e201c85SYohann 239e201c85SYohann //------------------------------------------------------------------------------ 249e201c85SYohann // 1D 259e201c85SYohann //------------------------------------------------------------------------------ 269e201c85SYohann 279e201c85SYohann //------------------------------------------------------------------------------ 289e201c85SYohann // E-vector -> single element 299e201c85SYohann //------------------------------------------------------------------------------ 309e201c85SYohann template <int NUM_COMP, int P_1D> 312b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 322b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 339e201c85SYohann if (data.t_id_x < P_1D) { 349e201c85SYohann const CeedInt node = data.t_id_x; 359e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 36672b0f2aSSebastian Grimberg 379e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 389e201c85SYohann r_u[comp] = d_u[ind + comp * strides_comp]; 399e201c85SYohann } 409e201c85SYohann } 419e201c85SYohann } 429e201c85SYohann 439e201c85SYohann //------------------------------------------------------------------------------ 449e201c85SYohann // Single element -> E-vector 459e201c85SYohann //------------------------------------------------------------------------------ 469e201c85SYohann template <int NUM_COMP, int P_1D> 472b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 482b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 499e201c85SYohann if (data.t_id_x < P_1D) { 509e201c85SYohann const CeedInt node = data.t_id_x; 519e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 52672b0f2aSSebastian Grimberg 539e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 549e201c85SYohann d_v[ind + comp * strides_comp] = r_v[comp]; 559e201c85SYohann } 569e201c85SYohann } 579e201c85SYohann } 589e201c85SYohann 59db2becc9SJeremy L Thompson template <int NUM_COMP, int P_1D> 60db2becc9SJeremy L Thompson inline __device__ void SumElementStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 61db2becc9SJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 62db2becc9SJeremy L Thompson if (data.t_id_x < P_1D) { 63db2becc9SJeremy L Thompson const CeedInt node = data.t_id_x; 64db2becc9SJeremy L Thompson const CeedInt ind = node * strides_node + elem * strides_elem; 65db2becc9SJeremy L Thompson 66db2becc9SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 67db2becc9SJeremy L Thompson d_v[ind + comp * strides_comp] += r_v[comp]; 68db2becc9SJeremy L Thompson } 69db2becc9SJeremy L Thompson } 70db2becc9SJeremy L Thompson } 71db2becc9SJeremy L Thompson 729e201c85SYohann //------------------------------------------------------------------------------ 739e201c85SYohann // 2D 749e201c85SYohann //------------------------------------------------------------------------------ 759e201c85SYohann 769e201c85SYohann //------------------------------------------------------------------------------ 779e201c85SYohann // E-vector -> single element 789e201c85SYohann //------------------------------------------------------------------------------ 799e201c85SYohann template <int NUM_COMP, int P_1D> 802b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 812b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 829e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 839e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 849e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 85672b0f2aSSebastian Grimberg 869e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 879e201c85SYohann r_u[comp] = d_u[ind + comp * strides_comp]; 889e201c85SYohann } 899e201c85SYohann } 909e201c85SYohann } 919e201c85SYohann 929e201c85SYohann //------------------------------------------------------------------------------ 939e201c85SYohann // Single element -> E-vector 949e201c85SYohann //------------------------------------------------------------------------------ 959e201c85SYohann template <int NUM_COMP, int P_1D> 962b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 972b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 989e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 999e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 1009e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 101672b0f2aSSebastian Grimberg 1029e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 1039e201c85SYohann d_v[ind + comp * strides_comp] = r_v[comp]; 1049e201c85SYohann } 1059e201c85SYohann } 1069e201c85SYohann } 1079e201c85SYohann 108db2becc9SJeremy L Thompson template <int NUM_COMP, int P_1D> 109db2becc9SJeremy L Thompson inline __device__ void SumElementStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 110db2becc9SJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 111db2becc9SJeremy L Thompson if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 112db2becc9SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 113db2becc9SJeremy L Thompson const CeedInt ind = node * strides_node + elem * strides_elem; 114db2becc9SJeremy L Thompson 115db2becc9SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 116db2becc9SJeremy L Thompson d_v[ind + comp * strides_comp] += r_v[comp]; 117db2becc9SJeremy L Thompson } 118db2becc9SJeremy L Thompson } 119db2becc9SJeremy L Thompson } 120db2becc9SJeremy L Thompson 1219e201c85SYohann //------------------------------------------------------------------------------ 1229e201c85SYohann // 3D 1239e201c85SYohann //------------------------------------------------------------------------------ 1249e201c85SYohann 1259e201c85SYohann //------------------------------------------------------------------------------ 1269e201c85SYohann // E-vector -> single element 1279e201c85SYohann //------------------------------------------------------------------------------ 1289e201c85SYohann template <int NUM_COMP, int P_1D> 1292b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 1302b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 1319e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1329e201c85SYohann for (CeedInt z = 0; z < P_1D; z++) { 1339e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 1349e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 135672b0f2aSSebastian Grimberg 1369e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 1379e201c85SYohann r_u[z + comp * P_1D] = d_u[ind + comp * strides_comp]; 1389e201c85SYohann } 1399e201c85SYohann } 1409e201c85SYohann } 1419e201c85SYohann } 1429e201c85SYohann 1439e201c85SYohann //------------------------------------------------------------------------------ 1449e201c85SYohann // Single element -> E-vector 1459e201c85SYohann //------------------------------------------------------------------------------ 1469e201c85SYohann template <int NUM_COMP, int P_1D> 1472b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 1482b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 1499e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1509e201c85SYohann for (CeedInt z = 0; z < P_1D; z++) { 1519e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 1529e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 153672b0f2aSSebastian Grimberg 1549e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 1559e201c85SYohann d_v[ind + comp * strides_comp] = r_v[z + comp * P_1D]; 1569e201c85SYohann } 1579e201c85SYohann } 1589e201c85SYohann } 1599e201c85SYohann } 160db2becc9SJeremy L Thompson 161db2becc9SJeremy L Thompson template <int NUM_COMP, int P_1D> 162db2becc9SJeremy L Thompson inline __device__ void SumElementStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 163db2becc9SJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 164db2becc9SJeremy L Thompson if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 165db2becc9SJeremy L Thompson for (CeedInt z = 0; z < P_1D; z++) { 166db2becc9SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 167db2becc9SJeremy L Thompson const CeedInt ind = node * strides_node + elem * strides_elem; 168db2becc9SJeremy L Thompson 169db2becc9SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 170db2becc9SJeremy L Thompson d_v[ind + comp * strides_comp] += r_v[z + comp * P_1D]; 171db2becc9SJeremy L Thompson } 172db2becc9SJeremy L Thompson } 173db2becc9SJeremy L Thompson } 174db2becc9SJeremy L Thompson } 175