1 // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 /// @file 9 /// Internal header for HIP offset element restriction kernels 10 #include <ceed/types.h> 11 12 //------------------------------------------------------------------------------ 13 // E-vector -> L-vector, standard (with offsets) 14 //------------------------------------------------------------------------------ 15 #if !USE_DETERMINISTIC 16 extern "C" __global__ void AtPointsTranspose(const CeedInt *__restrict__ indices, const CeedInt *__restrict__ points_per_elem, 17 const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) { 18 for (CeedInt node = blockIdx.x * blockDim.x + threadIdx.x; node < RSTR_NUM_ELEM * RSTR_ELEM_SIZE; node += blockDim.x * gridDim.x) { 19 const CeedInt ind = indices[node]; 20 const CeedInt loc_node = node % RSTR_ELEM_SIZE; 21 const CeedInt elem = node / RSTR_ELEM_SIZE; 22 23 if (loc_node >= points_per_elem[elem]) continue; 24 for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) { 25 atomicAdd(&v[ind + comp * RSTR_COMP_STRIDE], u[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE]); 26 } 27 } 28 } 29 #else 30 extern "C" __global__ void AtPointsTranspose(const CeedInt *__restrict__ l_vec_indices, const CeedInt *__restrict__ t_indices, 31 const CeedInt *__restrict__ points_per_elem, const CeedInt *__restrict__ t_offsets, 32 const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) { 33 CeedScalar value[RSTR_NUM_COMP]; 34 35 for (CeedInt i = blockIdx.x * blockDim.x + threadIdx.x; i < RSTR_NUM_NODES; i += blockDim.x * gridDim.x) { 36 const CeedInt ind = l_vec_indices[i]; 37 const CeedInt range_1 = t_offsets[i]; 38 const CeedInt range_N = t_offsets[i + 1]; 39 40 for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) value[comp] = 0.0; 41 42 for (CeedInt j = range_1; j < range_N; j++) { 43 const CeedInt t_ind = t_indices[j]; 44 const CeedInt loc_node = t_ind % RSTR_ELEM_SIZE; 45 const CeedInt elem = t_ind / RSTR_ELEM_SIZE; 46 47 if (loc_node >= points_per_elem[elem]) continue; 48 for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) { 49 value[comp] += u[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE]; 50 } 51 } 52 53 for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) v[ind + comp * RSTR_COMP_STRIDE] += value[comp]; 54 } 55 } 56 #endif 57