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 CUDA shared memory basis read/write templates 1094b7b29bSJeremy L Thompson #ifndef CEED_CUDA_SHARED_BASIS_READ_WRITE_TEMPLATES_H 1194b7b29bSJeremy L Thompson #define CEED_CUDA_SHARED_BASIS_READ_WRITE_TEMPLATES_H 129e201c85SYohann 139e201c85SYohann #include <ceed.h> 149e201c85SYohann 159e201c85SYohann //------------------------------------------------------------------------------ 169e201c85SYohann // 1D 179e201c85SYohann //------------------------------------------------------------------------------ 189e201c85SYohann 199e201c85SYohann //------------------------------------------------------------------------------ 209e201c85SYohann // E-vector -> single element 219e201c85SYohann //------------------------------------------------------------------------------ 229e201c85SYohann template <int NUM_COMP, int P_1D> 232b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 242b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 259e201c85SYohann if (data.t_id_x < P_1D) { 269e201c85SYohann const CeedInt node = data.t_id_x; 279e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 28*ca735530SJeremy L Thompson 299e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 309e201c85SYohann r_u[comp] = d_u[ind + comp * strides_comp]; 319e201c85SYohann } 329e201c85SYohann } 339e201c85SYohann } 349e201c85SYohann 359e201c85SYohann //------------------------------------------------------------------------------ 369e201c85SYohann // Single element -> E-vector 379e201c85SYohann //------------------------------------------------------------------------------ 389e201c85SYohann template <int NUM_COMP, int P_1D> 392b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 402b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 419e201c85SYohann if (data.t_id_x < P_1D) { 429e201c85SYohann const CeedInt node = data.t_id_x; 439e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 44*ca735530SJeremy L Thompson 459e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 469e201c85SYohann d_v[ind + comp * strides_comp] = r_v[comp]; 479e201c85SYohann } 489e201c85SYohann } 499e201c85SYohann } 509e201c85SYohann 519e201c85SYohann //------------------------------------------------------------------------------ 529e201c85SYohann // 2D 539e201c85SYohann //------------------------------------------------------------------------------ 549e201c85SYohann 559e201c85SYohann //------------------------------------------------------------------------------ 569e201c85SYohann // E-vector -> single element 579e201c85SYohann //------------------------------------------------------------------------------ 589e201c85SYohann template <int NUM_COMP, int P_1D> 592b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 602b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 619e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 629e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 639e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 64*ca735530SJeremy L Thompson 659e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 669e201c85SYohann r_u[comp] = d_u[ind + comp * strides_comp]; 679e201c85SYohann } 689e201c85SYohann } 699e201c85SYohann } 709e201c85SYohann 719e201c85SYohann //------------------------------------------------------------------------------ 729e201c85SYohann // Single element -> E-vector 739e201c85SYohann //------------------------------------------------------------------------------ 749e201c85SYohann template <int NUM_COMP, int P_1D> 752b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 762b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 779e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 789e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 799e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 80*ca735530SJeremy L Thompson 819e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 829e201c85SYohann d_v[ind + comp * strides_comp] = r_v[comp]; 839e201c85SYohann } 849e201c85SYohann } 859e201c85SYohann } 869e201c85SYohann 879e201c85SYohann //------------------------------------------------------------------------------ 889e201c85SYohann // 3D 899e201c85SYohann //------------------------------------------------------------------------------ 909e201c85SYohann 919e201c85SYohann //------------------------------------------------------------------------------ 929e201c85SYohann // E-vector -> single element 939e201c85SYohann //------------------------------------------------------------------------------ 949e201c85SYohann template <int NUM_COMP, int P_1D> 952b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 962b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 979e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 989e201c85SYohann for (CeedInt z = 0; z < P_1D; z++) { 999e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 1009e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 101*ca735530SJeremy L Thompson 1029e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 1039e201c85SYohann r_u[z + comp * P_1D] = d_u[ind + comp * strides_comp]; 1049e201c85SYohann } 1059e201c85SYohann } 1069e201c85SYohann } 1079e201c85SYohann } 1089e201c85SYohann 1099e201c85SYohann //------------------------------------------------------------------------------ 1109e201c85SYohann // Single element -> E-vector 1119e201c85SYohann //------------------------------------------------------------------------------ 1129e201c85SYohann template <int NUM_COMP, int P_1D> 1132b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 1142b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 1159e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1169e201c85SYohann for (CeedInt z = 0; z < P_1D; z++) { 1179e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 1189e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 119*ca735530SJeremy L Thompson 1209e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 1219e201c85SYohann d_v[ind + comp * strides_comp] = r_v[z + comp * P_1D]; 1229e201c85SYohann } 1239e201c85SYohann } 1249e201c85SYohann } 1259e201c85SYohann } 1269e201c85SYohann 1279e201c85SYohann //------------------------------------------------------------------------------ 1289e201c85SYohann 12994b7b29bSJeremy L Thompson #endif // CEED_CUDA_SHARED_BASIS_READ_WRITE_TEMPLATES_H 130