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 CUDA shared memory basis read/write templates 109e201c85SYohann 11*c0b5abf0SJeremy L Thompson #include <ceed/types.h> 129e201c85SYohann 139e201c85SYohann //------------------------------------------------------------------------------ 149e201c85SYohann // 1D 159e201c85SYohann //------------------------------------------------------------------------------ 169e201c85SYohann 179e201c85SYohann //------------------------------------------------------------------------------ 189e201c85SYohann // E-vector -> single element 199e201c85SYohann //------------------------------------------------------------------------------ 209e201c85SYohann template <int NUM_COMP, int P_1D> 212b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 222b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 239e201c85SYohann if (data.t_id_x < P_1D) { 249e201c85SYohann const CeedInt node = data.t_id_x; 259e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 26ca735530SJeremy L Thompson 279e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 289e201c85SYohann r_u[comp] = d_u[ind + comp * strides_comp]; 299e201c85SYohann } 309e201c85SYohann } 319e201c85SYohann } 329e201c85SYohann 339e201c85SYohann //------------------------------------------------------------------------------ 349e201c85SYohann // Single element -> E-vector 359e201c85SYohann //------------------------------------------------------------------------------ 369e201c85SYohann template <int NUM_COMP, int P_1D> 372b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 382b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 399e201c85SYohann if (data.t_id_x < P_1D) { 409e201c85SYohann const CeedInt node = data.t_id_x; 419e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 42ca735530SJeremy L Thompson 439e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 449e201c85SYohann d_v[ind + comp * strides_comp] = r_v[comp]; 459e201c85SYohann } 469e201c85SYohann } 479e201c85SYohann } 489e201c85SYohann 49db2becc9SJeremy L Thompson template <int NUM_COMP, int P_1D> 50db2becc9SJeremy L Thompson inline __device__ void SumElementStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 51db2becc9SJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 52db2becc9SJeremy L Thompson if (data.t_id_x < P_1D) { 53db2becc9SJeremy L Thompson const CeedInt node = data.t_id_x; 54db2becc9SJeremy L Thompson const CeedInt ind = node * strides_node + elem * strides_elem; 55db2becc9SJeremy L Thompson 56db2becc9SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 57db2becc9SJeremy L Thompson d_v[ind + comp * strides_comp] += r_v[comp]; 58db2becc9SJeremy L Thompson } 59db2becc9SJeremy L Thompson } 60db2becc9SJeremy L Thompson } 61db2becc9SJeremy L Thompson 629e201c85SYohann //------------------------------------------------------------------------------ 639e201c85SYohann // 2D 649e201c85SYohann //------------------------------------------------------------------------------ 659e201c85SYohann 669e201c85SYohann //------------------------------------------------------------------------------ 679e201c85SYohann // E-vector -> single element 689e201c85SYohann //------------------------------------------------------------------------------ 699e201c85SYohann template <int NUM_COMP, int P_1D> 702b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 712b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 729e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 739e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 749e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 75ca735530SJeremy L Thompson 769e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 779e201c85SYohann r_u[comp] = d_u[ind + comp * strides_comp]; 789e201c85SYohann } 799e201c85SYohann } 809e201c85SYohann } 819e201c85SYohann 829e201c85SYohann //------------------------------------------------------------------------------ 839e201c85SYohann // Single element -> E-vector 849e201c85SYohann //------------------------------------------------------------------------------ 859e201c85SYohann template <int NUM_COMP, int P_1D> 862b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 872b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 889e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 899e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 909e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 91ca735530SJeremy L Thompson 929e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 939e201c85SYohann d_v[ind + comp * strides_comp] = r_v[comp]; 949e201c85SYohann } 959e201c85SYohann } 969e201c85SYohann } 979e201c85SYohann 98db2becc9SJeremy L Thompson template <int NUM_COMP, int P_1D> 99db2becc9SJeremy L Thompson inline __device__ void SumElementStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 100db2becc9SJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 101db2becc9SJeremy L Thompson if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 102db2becc9SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 103db2becc9SJeremy L Thompson const CeedInt ind = node * strides_node + elem * strides_elem; 104db2becc9SJeremy L Thompson 105db2becc9SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 106db2becc9SJeremy L Thompson d_v[ind + comp * strides_comp] += r_v[comp]; 107db2becc9SJeremy L Thompson } 108db2becc9SJeremy L Thompson } 109db2becc9SJeremy L Thompson } 110db2becc9SJeremy L Thompson 1119e201c85SYohann //------------------------------------------------------------------------------ 1129e201c85SYohann // 3D 1139e201c85SYohann //------------------------------------------------------------------------------ 1149e201c85SYohann 1159e201c85SYohann //------------------------------------------------------------------------------ 1169e201c85SYohann // E-vector -> single element 1179e201c85SYohann //------------------------------------------------------------------------------ 1189e201c85SYohann template <int NUM_COMP, int P_1D> 1192b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 1202b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 1219e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1229e201c85SYohann for (CeedInt z = 0; z < P_1D; z++) { 1239e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 1249e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 125ca735530SJeremy L Thompson 1269e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 1279e201c85SYohann r_u[z + comp * P_1D] = d_u[ind + comp * strides_comp]; 1289e201c85SYohann } 1299e201c85SYohann } 1309e201c85SYohann } 1319e201c85SYohann } 1329e201c85SYohann 1339e201c85SYohann //------------------------------------------------------------------------------ 1349e201c85SYohann // Single element -> E-vector 1359e201c85SYohann //------------------------------------------------------------------------------ 1369e201c85SYohann template <int NUM_COMP, int P_1D> 1372b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 1382b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 1399e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1409e201c85SYohann for (CeedInt z = 0; z < P_1D; z++) { 1419e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 1429e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 143ca735530SJeremy L Thompson 1449e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 1459e201c85SYohann d_v[ind + comp * strides_comp] = r_v[z + comp * P_1D]; 1469e201c85SYohann } 1479e201c85SYohann } 1489e201c85SYohann } 1499e201c85SYohann } 150db2becc9SJeremy L Thompson 151db2becc9SJeremy L Thompson template <int NUM_COMP, int P_1D> 152db2becc9SJeremy L Thompson inline __device__ void SumElementStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 153db2becc9SJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 154db2becc9SJeremy L Thompson if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 155db2becc9SJeremy L Thompson for (CeedInt z = 0; z < P_1D; z++) { 156db2becc9SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 157db2becc9SJeremy L Thompson const CeedInt ind = node * strides_node + elem * strides_elem; 158db2becc9SJeremy L Thompson 159db2becc9SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 160db2becc9SJeremy L Thompson d_v[ind + comp * strides_comp] += r_v[z + comp * P_1D]; 161db2becc9SJeremy L Thompson } 162db2becc9SJeremy L Thompson } 163db2becc9SJeremy L Thompson } 164db2becc9SJeremy L Thompson } 165