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 basis read/write templates 109e201c85SYohann #ifndef _ceed_hip_shared_basis_read_write_templates_h 119e201c85SYohann #define _ceed_hip_shared_basis_read_write_templates_h 129e201c85SYohann 139e201c85SYohann #include <ceed.h> 149e201c85SYohann 159e201c85SYohann //------------------------------------------------------------------------------ 169e201c85SYohann // Helper function: load matrices for basis actions 179e201c85SYohann //------------------------------------------------------------------------------ 189e201c85SYohann template <int SIZE> 199e201c85SYohann inline __device__ void loadMatrix(const CeedScalar *d_B, CeedScalar *B) { 209e201c85SYohann CeedInt tid = threadIdx.x + threadIdx.y * blockDim.x + threadIdx.z * blockDim.y * blockDim.x; 21*2b730f8bSJeremy L Thompson for (CeedInt i = tid; i < SIZE; i += blockDim.x * blockDim.y * blockDim.z) B[i] = d_B[i]; 229e201c85SYohann } 239e201c85SYohann 249e201c85SYohann //------------------------------------------------------------------------------ 259e201c85SYohann // 1D 269e201c85SYohann //------------------------------------------------------------------------------ 279e201c85SYohann 289e201c85SYohann //------------------------------------------------------------------------------ 299e201c85SYohann // E-vector -> single element 309e201c85SYohann //------------------------------------------------------------------------------ 319e201c85SYohann template <int NUM_COMP, int P_1D> 32*2b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 33*2b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 349e201c85SYohann if (data.t_id_x < P_1D) { 359e201c85SYohann const CeedInt node = data.t_id_x; 369e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 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> 47*2b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 48*2b730f8bSJeremy 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; 529e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 539e201c85SYohann d_v[ind + comp * strides_comp] = r_v[comp]; 549e201c85SYohann } 559e201c85SYohann } 569e201c85SYohann } 579e201c85SYohann 589e201c85SYohann //------------------------------------------------------------------------------ 599e201c85SYohann // 2D 609e201c85SYohann //------------------------------------------------------------------------------ 619e201c85SYohann 629e201c85SYohann //------------------------------------------------------------------------------ 639e201c85SYohann // E-vector -> single element 649e201c85SYohann //------------------------------------------------------------------------------ 659e201c85SYohann template <int NUM_COMP, int P_1D> 66*2b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 67*2b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 689e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 699e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 709e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 719e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 729e201c85SYohann r_u[comp] = d_u[ind + comp * strides_comp]; 739e201c85SYohann } 749e201c85SYohann } 759e201c85SYohann } 769e201c85SYohann 779e201c85SYohann //------------------------------------------------------------------------------ 789e201c85SYohann // Single element -> E-vector 799e201c85SYohann //------------------------------------------------------------------------------ 809e201c85SYohann template <int NUM_COMP, int P_1D> 81*2b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 82*2b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 839e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 849e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 859e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 869e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 879e201c85SYohann d_v[ind + comp * strides_comp] = r_v[comp]; 889e201c85SYohann } 899e201c85SYohann } 909e201c85SYohann } 919e201c85SYohann 929e201c85SYohann //------------------------------------------------------------------------------ 939e201c85SYohann // 3D 949e201c85SYohann //------------------------------------------------------------------------------ 959e201c85SYohann 969e201c85SYohann //------------------------------------------------------------------------------ 979e201c85SYohann // E-vector -> single element 989e201c85SYohann //------------------------------------------------------------------------------ 999e201c85SYohann template <int NUM_COMP, int P_1D> 100*2b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 101*2b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 1029e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1039e201c85SYohann for (CeedInt z = 0; z < P_1D; z++) { 1049e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 1059e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 1069e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 1079e201c85SYohann r_u[z + comp * P_1D] = d_u[ind + comp * strides_comp]; 1089e201c85SYohann } 1099e201c85SYohann } 1109e201c85SYohann } 1119e201c85SYohann } 1129e201c85SYohann 1139e201c85SYohann //------------------------------------------------------------------------------ 1149e201c85SYohann // Single element -> E-vector 1159e201c85SYohann //------------------------------------------------------------------------------ 1169e201c85SYohann template <int NUM_COMP, int P_1D> 117*2b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 118*2b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 1199e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1209e201c85SYohann for (CeedInt z = 0; z < P_1D; z++) { 1219e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 1229e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 1239e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 1249e201c85SYohann d_v[ind + comp * strides_comp] = r_v[z + comp * P_1D]; 1259e201c85SYohann } 1269e201c85SYohann } 1279e201c85SYohann } 1289e201c85SYohann } 1299e201c85SYohann 1309e201c85SYohann //------------------------------------------------------------------------------ 1319e201c85SYohann 1329e201c85SYohann #endif 133