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 10c0b5abf0SJeremy L Thompson #include <ceed/types.h> 119e201c85SYohann 129e201c85SYohann //------------------------------------------------------------------------------ 139e201c85SYohann // 1D 149e201c85SYohann //------------------------------------------------------------------------------ 159e201c85SYohann 169e201c85SYohann //------------------------------------------------------------------------------ 179e201c85SYohann // E-vector -> single element 189e201c85SYohann //------------------------------------------------------------------------------ 199e201c85SYohann template <int NUM_COMP, int P_1D> 202b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 212b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 229e201c85SYohann if (data.t_id_x < P_1D) { 239e201c85SYohann const CeedInt node = data.t_id_x; 249e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 25ca735530SJeremy L Thompson 269e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 279e201c85SYohann r_u[comp] = d_u[ind + comp * strides_comp]; 289e201c85SYohann } 299e201c85SYohann } 309e201c85SYohann } 319e201c85SYohann 329e201c85SYohann //------------------------------------------------------------------------------ 339e201c85SYohann // Single element -> E-vector 349e201c85SYohann //------------------------------------------------------------------------------ 359e201c85SYohann template <int NUM_COMP, int P_1D> 362b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 372b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 389e201c85SYohann if (data.t_id_x < P_1D) { 399e201c85SYohann const CeedInt node = data.t_id_x; 409e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 41ca735530SJeremy L Thompson 429e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 439e201c85SYohann d_v[ind + comp * strides_comp] = r_v[comp]; 449e201c85SYohann } 459e201c85SYohann } 469e201c85SYohann } 479e201c85SYohann 48db2becc9SJeremy L Thompson template <int NUM_COMP, int P_1D> 49db2becc9SJeremy L Thompson inline __device__ void SumElementStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 50db2becc9SJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 51db2becc9SJeremy L Thompson if (data.t_id_x < P_1D) { 52db2becc9SJeremy L Thompson const CeedInt node = data.t_id_x; 53db2becc9SJeremy L Thompson const CeedInt ind = node * strides_node + elem * strides_elem; 54db2becc9SJeremy L Thompson 55db2becc9SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 56db2becc9SJeremy L Thompson d_v[ind + comp * strides_comp] += r_v[comp]; 57db2becc9SJeremy L Thompson } 58db2becc9SJeremy L Thompson } 59db2becc9SJeremy L Thompson } 60db2becc9SJeremy L Thompson 619e201c85SYohann //------------------------------------------------------------------------------ 629e201c85SYohann // 2D 639e201c85SYohann //------------------------------------------------------------------------------ 649e201c85SYohann 659e201c85SYohann //------------------------------------------------------------------------------ 669e201c85SYohann // E-vector -> single element 679e201c85SYohann //------------------------------------------------------------------------------ 689e201c85SYohann template <int NUM_COMP, int P_1D> 692b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 702b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 719e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 729e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 739e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 74ca735530SJeremy L Thompson 759e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 769e201c85SYohann r_u[comp] = d_u[ind + comp * strides_comp]; 779e201c85SYohann } 789e201c85SYohann } 799e201c85SYohann } 809e201c85SYohann 819e201c85SYohann //------------------------------------------------------------------------------ 829e201c85SYohann // Single element -> E-vector 839e201c85SYohann //------------------------------------------------------------------------------ 849e201c85SYohann template <int NUM_COMP, int P_1D> 852b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 862b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 879e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 889e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 899e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 90ca735530SJeremy L Thompson 919e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 929e201c85SYohann d_v[ind + comp * strides_comp] = r_v[comp]; 939e201c85SYohann } 949e201c85SYohann } 959e201c85SYohann } 969e201c85SYohann 97db2becc9SJeremy L Thompson template <int NUM_COMP, int P_1D> 98db2becc9SJeremy L Thompson inline __device__ void SumElementStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 99db2becc9SJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 100db2becc9SJeremy L Thompson if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 101db2becc9SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 102db2becc9SJeremy L Thompson const CeedInt ind = node * strides_node + elem * strides_elem; 103db2becc9SJeremy L Thompson 104db2becc9SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 105db2becc9SJeremy L Thompson d_v[ind + comp * strides_comp] += r_v[comp]; 106db2becc9SJeremy L Thompson } 107db2becc9SJeremy L Thompson } 108db2becc9SJeremy L Thompson } 109db2becc9SJeremy L Thompson 1109e201c85SYohann //------------------------------------------------------------------------------ 1119e201c85SYohann // 3D 1129e201c85SYohann //------------------------------------------------------------------------------ 1139e201c85SYohann 1149e201c85SYohann //------------------------------------------------------------------------------ 1159e201c85SYohann // E-vector -> single element 1169e201c85SYohann //------------------------------------------------------------------------------ 1179e201c85SYohann template <int NUM_COMP, int P_1D> 1182b730f8bSJeremy L Thompson inline __device__ void ReadElementStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 1192b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 1209e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1219e201c85SYohann for (CeedInt z = 0; z < P_1D; z++) { 1229e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 1239e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 124ca735530SJeremy L Thompson 1259e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 1269e201c85SYohann r_u[z + comp * P_1D] = d_u[ind + comp * strides_comp]; 1279e201c85SYohann } 1289e201c85SYohann } 1299e201c85SYohann } 1309e201c85SYohann } 1319e201c85SYohann 1329e201c85SYohann //------------------------------------------------------------------------------ 1339e201c85SYohann // Single element -> E-vector 1349e201c85SYohann //------------------------------------------------------------------------------ 1359e201c85SYohann template <int NUM_COMP, int P_1D> 1362b730f8bSJeremy L Thompson inline __device__ void WriteElementStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 1372b730f8bSJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 1389e201c85SYohann if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1399e201c85SYohann for (CeedInt z = 0; z < P_1D; z++) { 1409e201c85SYohann const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 1419e201c85SYohann const CeedInt ind = node * strides_node + elem * strides_elem; 142ca735530SJeremy L Thompson 1439e201c85SYohann for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 1449e201c85SYohann d_v[ind + comp * strides_comp] = r_v[z + comp * P_1D]; 1459e201c85SYohann } 1469e201c85SYohann } 1479e201c85SYohann } 1489e201c85SYohann } 149db2becc9SJeremy L Thompson 150db2becc9SJeremy L Thompson template <int NUM_COMP, int P_1D> 151db2becc9SJeremy L Thompson inline __device__ void SumElementStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt strides_node, const CeedInt strides_comp, 152db2becc9SJeremy L Thompson const CeedInt strides_elem, const CeedScalar *r_v, CeedScalar *d_v) { 153db2becc9SJeremy L Thompson if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 154db2becc9SJeremy L Thompson for (CeedInt z = 0; z < P_1D; z++) { 155db2becc9SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 156db2becc9SJeremy L Thompson const CeedInt ind = node * strides_node + elem * strides_elem; 157db2becc9SJeremy L Thompson 158db2becc9SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 159db2becc9SJeremy L Thompson d_v[ind + comp * strides_comp] += r_v[z + comp * P_1D]; 160db2becc9SJeremy L Thompson } 161db2becc9SJeremy L Thompson } 162db2becc9SJeremy L Thompson } 163db2becc9SJeremy L Thompson } 164*b6a2eb79SJeremy L Thompson 165*b6a2eb79SJeremy L Thompson //------------------------------------------------------------------------------ 166*b6a2eb79SJeremy L Thompson // AtPoints 167*b6a2eb79SJeremy L Thompson //------------------------------------------------------------------------------ 168*b6a2eb79SJeremy L Thompson 169*b6a2eb79SJeremy L Thompson //------------------------------------------------------------------------------ 170*b6a2eb79SJeremy L Thompson // E-vector -> single point 171*b6a2eb79SJeremy L Thompson //------------------------------------------------------------------------------ 172*b6a2eb79SJeremy L Thompson template <int NUM_COMP, int NUM_PTS> 173*b6a2eb79SJeremy L Thompson inline __device__ void ReadPoint(SharedData_Cuda &data, const CeedInt elem, const CeedInt p, const CeedInt points_in_elem, 174*b6a2eb79SJeremy L Thompson const CeedInt strides_point, const CeedInt strides_comp, const CeedInt strides_elem, 175*b6a2eb79SJeremy L Thompson const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 176*b6a2eb79SJeremy L Thompson const CeedInt ind = (p % NUM_PTS) * strides_point + elem * strides_elem; 177*b6a2eb79SJeremy L Thompson 178*b6a2eb79SJeremy L Thompson if (p < points_in_elem) { 179*b6a2eb79SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 180*b6a2eb79SJeremy L Thompson r_u[comp] = d_u[ind + comp * strides_comp]; 181*b6a2eb79SJeremy L Thompson } 182*b6a2eb79SJeremy L Thompson } else { 183*b6a2eb79SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 184*b6a2eb79SJeremy L Thompson r_u[comp] = 0.0; 185*b6a2eb79SJeremy L Thompson } 186*b6a2eb79SJeremy L Thompson } 187*b6a2eb79SJeremy L Thompson } 188*b6a2eb79SJeremy L Thompson 189*b6a2eb79SJeremy L Thompson //------------------------------------------------------------------------------ 190*b6a2eb79SJeremy L Thompson // Single point -> E-vector 191*b6a2eb79SJeremy L Thompson //------------------------------------------------------------------------------ 192*b6a2eb79SJeremy L Thompson template <int NUM_COMP, int NUM_PTS> 193*b6a2eb79SJeremy L Thompson inline __device__ void WritePoint(SharedData_Cuda &data, const CeedInt elem, const CeedInt p, const CeedInt points_in_elem, 194*b6a2eb79SJeremy L Thompson const CeedInt strides_point, const CeedInt strides_comp, const CeedInt strides_elem, const CeedScalar *r_v, 195*b6a2eb79SJeremy L Thompson CeedScalar *d_v) { 196*b6a2eb79SJeremy L Thompson if (p < points_in_elem) { 197*b6a2eb79SJeremy L Thompson const CeedInt ind = (p % NUM_PTS) * strides_point + elem * strides_elem; 198*b6a2eb79SJeremy L Thompson 199*b6a2eb79SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 200*b6a2eb79SJeremy L Thompson d_v[ind + comp * strides_comp] = r_v[comp]; 201*b6a2eb79SJeremy L Thompson } 202*b6a2eb79SJeremy L Thompson } 203*b6a2eb79SJeremy L Thompson } 204