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 10*94b7b29bSJeremy L Thompson #ifndef CEED_CUDA_SHARED_BASIS_READ_WRITE_TEMPLATES_H 11*94b7b29bSJeremy 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; 289e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 299e201c85SYohann r_u[comp] = d_u[ind + comp * strides_comp]; 309e201c85SYohann } 319e201c85SYohann } 329e201c85SYohann } 339e201c85SYohann 349e201c85SYohann //------------------------------------------------------------------------------ 359e201c85SYohann // Single element -> E-vector 369e201c85SYohann //------------------------------------------------------------------------------ 379e201c85SYohann template <int NUM_COMP, int P_1D> 382b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 392b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 409e201c85SYohann if (data.t_id_x < P_1D) { 419e201c85SYohann const CeedInt node = data.t_id_x; 429e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 439e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 449e201c85SYohann d_v[ind + comp * strides_comp] = r_v[comp]; 459e201c85SYohann } 469e201c85SYohann } 479e201c85SYohann } 489e201c85SYohann 499e201c85SYohann //------------------------------------------------------------------------------ 509e201c85SYohann // 2D 519e201c85SYohann //------------------------------------------------------------------------------ 529e201c85SYohann 539e201c85SYohann //------------------------------------------------------------------------------ 549e201c85SYohann // E-vector -> single element 559e201c85SYohann //------------------------------------------------------------------------------ 569e201c85SYohann template <int NUM_COMP, int P_1D> 572b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 582b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 599e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 609e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 619e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 629e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 639e201c85SYohann r_u[comp] = d_u[ind + comp * strides_comp]; 649e201c85SYohann } 659e201c85SYohann } 669e201c85SYohann } 679e201c85SYohann 689e201c85SYohann //------------------------------------------------------------------------------ 699e201c85SYohann // Single element -> E-vector 709e201c85SYohann //------------------------------------------------------------------------------ 719e201c85SYohann template <int NUM_COMP, int P_1D> 722b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 732b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 749e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 759e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 769e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 779e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 789e201c85SYohann d_v[ind + comp * strides_comp] = r_v[comp]; 799e201c85SYohann } 809e201c85SYohann } 819e201c85SYohann } 829e201c85SYohann 839e201c85SYohann //------------------------------------------------------------------------------ 849e201c85SYohann // 3D 859e201c85SYohann //------------------------------------------------------------------------------ 869e201c85SYohann 879e201c85SYohann //------------------------------------------------------------------------------ 889e201c85SYohann // E-vector -> single element 899e201c85SYohann //------------------------------------------------------------------------------ 909e201c85SYohann template <int NUM_COMP, int P_1D> 912b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 922b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 939e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 949e201c85SYohann for (CeedInt z = 0; z < P_1D; z++) { 959e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 969e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 979e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 989e201c85SYohann r_u[z + comp * P_1D] = d_u[ind + comp * strides_comp]; 999e201c85SYohann } 1009e201c85SYohann } 1019e201c85SYohann } 1029e201c85SYohann } 1039e201c85SYohann 1049e201c85SYohann //------------------------------------------------------------------------------ 1059e201c85SYohann // Single element -> E-vector 1069e201c85SYohann //------------------------------------------------------------------------------ 1079e201c85SYohann template <int NUM_COMP, int P_1D> 1082b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 1092b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 1109e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1119e201c85SYohann for (CeedInt z = 0; z < P_1D; z++) { 1129e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 1139e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 1149e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 1159e201c85SYohann d_v[ind + comp * strides_comp] = r_v[z + comp * P_1D]; 1169e201c85SYohann } 1179e201c85SYohann } 1189e201c85SYohann } 1199e201c85SYohann } 1209e201c85SYohann 1219e201c85SYohann //------------------------------------------------------------------------------ 1229e201c85SYohann 123*94b7b29bSJeremy L Thompson #endif // CEED_CUDA_SHARED_BASIS_READ_WRITE_TEMPLATES_H 124