10b63de31SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 20b63de31SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 30b63de31SJeremy L Thompson // 40b63de31SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 50b63de31SJeremy L Thompson // 60b63de31SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 70b63de31SJeremy L Thompson 80b63de31SJeremy L Thompson /// @file 90b63de31SJeremy L Thompson /// Internal header for HIP offset element restriction kernels 100b63de31SJeremy L Thompson 110b63de31SJeremy L Thompson #include <ceed.h> 120b63de31SJeremy L Thompson 130b63de31SJeremy L Thompson //------------------------------------------------------------------------------ 140b63de31SJeremy L Thompson // E-vector -> L-vector, standard (with offsets) 150b63de31SJeremy L Thompson //------------------------------------------------------------------------------ 160b63de31SJeremy L Thompson #if !USE_DETERMINISTIC 170b63de31SJeremy L Thompson extern "C" __global__ void AtPointsTranspose(const CeedInt *__restrict__ indices, const CeedInt *__restrict__ points_per_elem, 180b63de31SJeremy L Thompson const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) { 190b63de31SJeremy L Thompson for (CeedInt node = blockIdx.x * blockDim.x + threadIdx.x; node < RSTR_NUM_ELEM * RSTR_ELEM_SIZE; node += blockDim.x * gridDim.x) { 200b63de31SJeremy L Thompson const CeedInt ind = indices[node]; 210b63de31SJeremy L Thompson const CeedInt loc_node = node % RSTR_ELEM_SIZE; 220b63de31SJeremy L Thompson const CeedInt elem = node / RSTR_ELEM_SIZE; 230b63de31SJeremy L Thompson 240b63de31SJeremy L Thompson if (loc_node >= points_per_elem[elem]) continue; 250b63de31SJeremy L Thompson for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) { 26*db2becc9SJeremy L Thompson atomicAdd(&v[ind + comp * RSTR_COMP_STRIDE], u[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE]); 270b63de31SJeremy L Thompson } 280b63de31SJeremy L Thompson } 290b63de31SJeremy L Thompson } 300b63de31SJeremy L Thompson #else 310b63de31SJeremy L Thompson extern "C" __global__ void AtPointsTranspose(const CeedInt *__restrict__ l_vec_indices, const CeedInt *__restrict__ t_indices, 320b63de31SJeremy L Thompson const CeedInt *__restrict__ points_per_elem, const CeedInt *__restrict__ t_offsets, 330b63de31SJeremy L Thompson const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) { 340b63de31SJeremy L Thompson CeedScalar value[RSTR_NUM_COMP]; 350b63de31SJeremy L Thompson 360b63de31SJeremy L Thompson for (CeedInt i = blockIdx.x * blockDim.x + threadIdx.x; i < RSTR_NUM_NODES; i += blockDim.x * gridDim.x) { 370b63de31SJeremy L Thompson const CeedInt ind = l_vec_indices[i]; 380b63de31SJeremy L Thompson const CeedInt range_1 = t_offsets[i]; 390b63de31SJeremy L Thompson const CeedInt range_N = t_offsets[i + 1]; 400b63de31SJeremy L Thompson 410b63de31SJeremy L Thompson for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) value[comp] = 0.0; 420b63de31SJeremy L Thompson 430b63de31SJeremy L Thompson for (CeedInt j = range_1; j < range_N; j++) { 440b63de31SJeremy L Thompson const CeedInt t_ind = t_indices[j]; 450b63de31SJeremy L Thompson const CeedInt loc_node = t_ind % RSTR_ELEM_SIZE; 460b63de31SJeremy L Thompson const CeedInt elem = t_ind / RSTR_ELEM_SIZE; 470b63de31SJeremy L Thompson 480b63de31SJeremy L Thompson if (loc_node >= points_per_elem[elem]) continue; 490b63de31SJeremy L Thompson for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) { 500b63de31SJeremy L Thompson value[comp] += u[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE]; 510b63de31SJeremy L Thompson } 520b63de31SJeremy L Thompson } 530b63de31SJeremy L Thompson 540b63de31SJeremy L Thompson for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) v[ind + comp * RSTR_COMP_STRIDE] += value[comp]; 550b63de31SJeremy L Thompson } 560b63de31SJeremy L Thompson } 570b63de31SJeremy L Thompson #endif 58