1*5aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2cf8cbdd6SSebastian Grimberg // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3cf8cbdd6SSebastian Grimberg // 4cf8cbdd6SSebastian Grimberg // SPDX-License-Identifier: BSD-2-Clause 5cf8cbdd6SSebastian Grimberg // 6cf8cbdd6SSebastian Grimberg // This file is part of CEED: http://github.com/ceed 7cf8cbdd6SSebastian Grimberg 8cf8cbdd6SSebastian Grimberg /// @file 9cf8cbdd6SSebastian Grimberg /// Internal header for CUDA curl-oriented element restriction kernels 10cf8cbdd6SSebastian Grimberg #ifndef CEED_CUDA_REF_RESTRICTION_CURL_ORIENTED_H 11cf8cbdd6SSebastian Grimberg #define CEED_CUDA_REF_RESTRICTION_CURL_ORIENTED_H 12cf8cbdd6SSebastian Grimberg 13cf8cbdd6SSebastian Grimberg #include <ceed.h> 14cf8cbdd6SSebastian Grimberg 15cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------ 16cf8cbdd6SSebastian Grimberg // L-vector -> E-vector, curl-oriented 17cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------ 18cf8cbdd6SSebastian Grimberg extern "C" __global__ void CurlOrientedNoTranspose(const CeedInt *__restrict__ indices, const CeedInt8 *__restrict__ curl_orients, 19cf8cbdd6SSebastian Grimberg const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) { 20cf8cbdd6SSebastian Grimberg for (CeedInt node = blockIdx.x * blockDim.x + threadIdx.x; node < RSTR_NUM_ELEM * RSTR_ELEM_SIZE; node += blockDim.x * gridDim.x) { 21cf8cbdd6SSebastian Grimberg const CeedInt loc_node = node % RSTR_ELEM_SIZE; 22cf8cbdd6SSebastian Grimberg const CeedInt elem = node / RSTR_ELEM_SIZE; 23cf8cbdd6SSebastian Grimberg const CeedInt ind_dl = loc_node > 0 ? indices[node - 1] : 0; 24cf8cbdd6SSebastian Grimberg const CeedInt ind_d = indices[node]; 25cf8cbdd6SSebastian Grimberg const CeedInt ind_du = loc_node < (RSTR_ELEM_SIZE - 1) ? indices[node + 1] : 0; 26cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_dl = curl_orients[3 * node + 0]; 27cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_d = curl_orients[3 * node + 1]; 28cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_du = curl_orients[3 * node + 2]; 29cf8cbdd6SSebastian Grimberg 30cf8cbdd6SSebastian Grimberg for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) { 31cf8cbdd6SSebastian Grimberg CeedScalar value = 0.0; 32cf8cbdd6SSebastian Grimberg value += loc_node > 0 ? u[ind_dl + comp * RSTR_COMP_STRIDE] * curl_orient_dl : 0.0; 33cf8cbdd6SSebastian Grimberg value += u[ind_d + comp * RSTR_COMP_STRIDE] * curl_orient_d; 34cf8cbdd6SSebastian Grimberg value += loc_node < (RSTR_ELEM_SIZE - 1) ? u[ind_du + comp * RSTR_COMP_STRIDE] * curl_orient_du : 0.0; 35cf8cbdd6SSebastian Grimberg v[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] = value; 36cf8cbdd6SSebastian Grimberg } 37cf8cbdd6SSebastian Grimberg } 38cf8cbdd6SSebastian Grimberg } 39cf8cbdd6SSebastian Grimberg 40cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------ 41cf8cbdd6SSebastian Grimberg // L-vector -> E-vector, unsigned curl-oriented 42cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------ 43cf8cbdd6SSebastian Grimberg extern "C" __global__ void CurlOrientedUnsignedNoTranspose(const CeedInt *__restrict__ indices, const CeedInt8 *__restrict__ curl_orients, 44cf8cbdd6SSebastian Grimberg const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) { 45cf8cbdd6SSebastian Grimberg for (CeedInt node = blockIdx.x * blockDim.x + threadIdx.x; node < RSTR_NUM_ELEM * RSTR_ELEM_SIZE; node += blockDim.x * gridDim.x) { 46cf8cbdd6SSebastian Grimberg const CeedInt loc_node = node % RSTR_ELEM_SIZE; 47cf8cbdd6SSebastian Grimberg const CeedInt elem = node / RSTR_ELEM_SIZE; 48cf8cbdd6SSebastian Grimberg const CeedInt ind_dl = loc_node > 0 ? indices[node - 1] : 0; 49cf8cbdd6SSebastian Grimberg const CeedInt ind_d = indices[node]; 50cf8cbdd6SSebastian Grimberg const CeedInt ind_du = loc_node < (RSTR_ELEM_SIZE - 1) ? indices[node + 1] : 0; 51cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_dl = abs(curl_orients[3 * node + 0]); 52cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_d = abs(curl_orients[3 * node + 1]); 53cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_du = abs(curl_orients[3 * node + 2]); 54cf8cbdd6SSebastian Grimberg 55cf8cbdd6SSebastian Grimberg for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) { 56cf8cbdd6SSebastian Grimberg CeedScalar value = 0.0; 57cf8cbdd6SSebastian Grimberg value += loc_node > 0 ? u[ind_dl + comp * RSTR_COMP_STRIDE] * curl_orient_dl : 0.0; 58cf8cbdd6SSebastian Grimberg value += u[ind_d + comp * RSTR_COMP_STRIDE] * curl_orient_d; 59cf8cbdd6SSebastian Grimberg value += loc_node < (RSTR_ELEM_SIZE - 1) ? u[ind_du + comp * RSTR_COMP_STRIDE] * curl_orient_du : 0.0; 60cf8cbdd6SSebastian Grimberg v[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] = value; 61cf8cbdd6SSebastian Grimberg } 62cf8cbdd6SSebastian Grimberg } 63cf8cbdd6SSebastian Grimberg } 64cf8cbdd6SSebastian Grimberg 65cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------ 66cf8cbdd6SSebastian Grimberg // E-vector -> L-vector, curl-oriented 67cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------ 68cf8cbdd6SSebastian Grimberg #if !USE_DETERMINISTIC 69cf8cbdd6SSebastian Grimberg extern "C" __global__ void CurlOrientedTranspose(const CeedInt *__restrict__ indices, const CeedInt8 *__restrict__ curl_orients, 70cf8cbdd6SSebastian Grimberg const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) { 71cf8cbdd6SSebastian Grimberg for (CeedInt node = blockIdx.x * blockDim.x + threadIdx.x; node < RSTR_NUM_ELEM * RSTR_ELEM_SIZE; node += blockDim.x * gridDim.x) { 72cf8cbdd6SSebastian Grimberg const CeedInt ind = indices[node]; 73cf8cbdd6SSebastian Grimberg const CeedInt loc_node = node % RSTR_ELEM_SIZE; 74cf8cbdd6SSebastian Grimberg const CeedInt elem = node / RSTR_ELEM_SIZE; 75cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_du = loc_node > 0 ? curl_orients[3 * node - 1] : 0.0; 76cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_d = curl_orients[3 * node + 1]; 77cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_dl = loc_node < (RSTR_ELEM_SIZE - 1) ? curl_orients[3 * node + 3] : 0.0; 78cf8cbdd6SSebastian Grimberg 79cf8cbdd6SSebastian Grimberg for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) { 80cf8cbdd6SSebastian Grimberg CeedScalar value = 0.0; 81cf8cbdd6SSebastian Grimberg value += loc_node > 0 ? u[loc_node - 1 + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_du : 0.0; 82cf8cbdd6SSebastian Grimberg value += u[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_d; 83cf8cbdd6SSebastian Grimberg value += 84cf8cbdd6SSebastian Grimberg loc_node < (RSTR_ELEM_SIZE - 1) ? u[loc_node + 1 + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_dl : 0.0; 85cf8cbdd6SSebastian Grimberg atomicAdd(v + ind + comp * RSTR_COMP_STRIDE, value); 86cf8cbdd6SSebastian Grimberg } 87cf8cbdd6SSebastian Grimberg } 88cf8cbdd6SSebastian Grimberg } 89cf8cbdd6SSebastian Grimberg #else 90cf8cbdd6SSebastian Grimberg extern "C" __global__ void CurlOrientedTranspose(const CeedInt *__restrict__ l_vec_indices, const CeedInt *__restrict__ t_indices, 91cf8cbdd6SSebastian Grimberg const CeedInt *__restrict__ t_offsets, const CeedInt8 *__restrict__ curl_orients, 92cf8cbdd6SSebastian Grimberg const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) { 93cf8cbdd6SSebastian Grimberg CeedScalar value[RSTR_NUM_COMP]; 94cf8cbdd6SSebastian Grimberg 95cf8cbdd6SSebastian Grimberg for (CeedInt i = blockIdx.x * blockDim.x + threadIdx.x; i < RSTR_NUM_NODES; i += blockDim.x * gridDim.x) { 96cf8cbdd6SSebastian Grimberg const CeedInt ind = l_vec_indices[i]; 97cf8cbdd6SSebastian Grimberg const CeedInt range_1 = t_offsets[i]; 98cf8cbdd6SSebastian Grimberg const CeedInt range_N = t_offsets[i + 1]; 99cf8cbdd6SSebastian Grimberg 100cf8cbdd6SSebastian Grimberg for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) value[comp] = 0.0; 101cf8cbdd6SSebastian Grimberg 102cf8cbdd6SSebastian Grimberg for (CeedInt j = range_1; j < range_N; j++) { 103cf8cbdd6SSebastian Grimberg const CeedInt t_ind = t_indices[j]; 104cf8cbdd6SSebastian Grimberg const CeedInt loc_node = t_ind % RSTR_ELEM_SIZE; 105cf8cbdd6SSebastian Grimberg const CeedInt elem = t_ind / RSTR_ELEM_SIZE; 106cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_du = loc_node > 0 ? curl_orients[3 * t_ind - 1] : 0.0; 107cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_d = curl_orients[3 * t_ind + 1]; 108cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_dl = loc_node < (RSTR_ELEM_SIZE - 1) ? curl_orients[3 * t_ind + 3] : 0.0; 109cf8cbdd6SSebastian Grimberg 110cf8cbdd6SSebastian Grimberg for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) { 111cf8cbdd6SSebastian Grimberg value[comp] += loc_node > 0 ? u[loc_node - 1 + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_du : 0.0; 112cf8cbdd6SSebastian Grimberg value[comp] += u[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_d; 113cf8cbdd6SSebastian Grimberg value[comp] += 114cf8cbdd6SSebastian Grimberg loc_node < (RSTR_ELEM_SIZE - 1) ? u[loc_node + 1 + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_dl : 0.0; 115cf8cbdd6SSebastian Grimberg } 116cf8cbdd6SSebastian Grimberg } 117cf8cbdd6SSebastian Grimberg 118cf8cbdd6SSebastian Grimberg for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) v[ind + comp * RSTR_COMP_STRIDE] += value[comp]; 119cf8cbdd6SSebastian Grimberg } 120cf8cbdd6SSebastian Grimberg } 121cf8cbdd6SSebastian Grimberg #endif 122cf8cbdd6SSebastian Grimberg 123cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------ 124cf8cbdd6SSebastian Grimberg // E-vector -> L-vector, unsigned curl-oriented 125cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------ 126cf8cbdd6SSebastian Grimberg #if !USE_DETERMINISTIC 127cf8cbdd6SSebastian Grimberg extern "C" __global__ void CurlOrientedUnsignedTranspose(const CeedInt *__restrict__ indices, const CeedInt8 *__restrict__ curl_orients, 128cf8cbdd6SSebastian Grimberg const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) { 129cf8cbdd6SSebastian Grimberg for (CeedInt node = blockIdx.x * blockDim.x + threadIdx.x; node < RSTR_NUM_ELEM * RSTR_ELEM_SIZE; node += blockDim.x * gridDim.x) { 130cf8cbdd6SSebastian Grimberg const CeedInt loc_node = node % RSTR_ELEM_SIZE; 131cf8cbdd6SSebastian Grimberg const CeedInt elem = node / RSTR_ELEM_SIZE; 132cf8cbdd6SSebastian Grimberg const CeedInt ind = indices[node]; 133cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_du = loc_node > 0 ? abs(curl_orients[3 * node - 1]) : 0.0; 134cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_d = abs(curl_orients[3 * node + 1]); 135cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_dl = loc_node < (RSTR_ELEM_SIZE - 1) ? abs(curl_orients[3 * node + 3]) : 0.0; 136cf8cbdd6SSebastian Grimberg 137cf8cbdd6SSebastian Grimberg for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) { 138cf8cbdd6SSebastian Grimberg CeedScalar value = 0.0; 139cf8cbdd6SSebastian Grimberg value += loc_node > 0 ? u[loc_node - 1 + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_du : 0.0; 140cf8cbdd6SSebastian Grimberg value += u[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_d; 141cf8cbdd6SSebastian Grimberg value += 142cf8cbdd6SSebastian Grimberg loc_node < (RSTR_ELEM_SIZE - 1) ? u[loc_node + 1 + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_dl : 0.0; 143cf8cbdd6SSebastian Grimberg atomicAdd(v + ind + comp * RSTR_COMP_STRIDE, value); 144cf8cbdd6SSebastian Grimberg } 145cf8cbdd6SSebastian Grimberg } 146cf8cbdd6SSebastian Grimberg } 147cf8cbdd6SSebastian Grimberg #else 148cf8cbdd6SSebastian Grimberg extern "C" __global__ void CurlOrientedUnsignedTranspose(const CeedInt *__restrict__ l_vec_indices, const CeedInt *__restrict__ t_indices, 149cf8cbdd6SSebastian Grimberg const CeedInt *__restrict__ t_offsets, const CeedInt8 *__restrict__ curl_orients, 150cf8cbdd6SSebastian Grimberg const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) { 151cf8cbdd6SSebastian Grimberg CeedScalar value[RSTR_NUM_COMP]; 152cf8cbdd6SSebastian Grimberg 153cf8cbdd6SSebastian Grimberg for (CeedInt i = blockIdx.x * blockDim.x + threadIdx.x; i < RSTR_NUM_NODES; i += blockDim.x * gridDim.x) { 154cf8cbdd6SSebastian Grimberg const CeedInt ind = l_vec_indices[i]; 155cf8cbdd6SSebastian Grimberg const CeedInt range_1 = t_offsets[i]; 156cf8cbdd6SSebastian Grimberg const CeedInt range_N = t_offsets[i + 1]; 157cf8cbdd6SSebastian Grimberg 158cf8cbdd6SSebastian Grimberg for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) value[comp] = 0.0; 159cf8cbdd6SSebastian Grimberg 160cf8cbdd6SSebastian Grimberg for (CeedInt j = range_1; j < range_N; j++) { 161cf8cbdd6SSebastian Grimberg const CeedInt t_ind = t_indices[j]; 162cf8cbdd6SSebastian Grimberg const CeedInt loc_node = t_ind % RSTR_ELEM_SIZE; 163cf8cbdd6SSebastian Grimberg const CeedInt elem = t_ind / RSTR_ELEM_SIZE; 164cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_du = loc_node > 0 ? abs(curl_orients[3 * t_ind - 1]) : 0.0; 165cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_d = abs(curl_orients[3 * t_ind + 1]); 166cf8cbdd6SSebastian Grimberg const CeedInt8 curl_orient_dl = loc_node < (RSTR_ELEM_SIZE - 1) ? abs(curl_orients[3 * t_ind + 3]) : 0.0; 167cf8cbdd6SSebastian Grimberg 168cf8cbdd6SSebastian Grimberg for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) { 169cf8cbdd6SSebastian Grimberg value[comp] += loc_node > 0 ? u[loc_node - 1 + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_du : 0.0; 170cf8cbdd6SSebastian Grimberg value[comp] += u[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_d; 171cf8cbdd6SSebastian Grimberg value[comp] += 172cf8cbdd6SSebastian Grimberg loc_node < (RSTR_ELEM_SIZE - 1) ? u[loc_node + 1 + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_dl : 0.0; 173cf8cbdd6SSebastian Grimberg } 174cf8cbdd6SSebastian Grimberg } 175cf8cbdd6SSebastian Grimberg 176cf8cbdd6SSebastian Grimberg for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) v[ind + comp * RSTR_COMP_STRIDE] += value[comp]; 177cf8cbdd6SSebastian Grimberg } 178cf8cbdd6SSebastian Grimberg } 179cf8cbdd6SSebastian Grimberg #endif 180cf8cbdd6SSebastian Grimberg 181cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------ 182cf8cbdd6SSebastian Grimberg 183cf8cbdd6SSebastian Grimberg #endif // CEED_CUDA_REF_RESTRICTION_CURL_ORIENTED_H 184