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