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 CUDA strided element restriction kernels 10 #include <ceed/types.h> 11 12 //------------------------------------------------------------------------------ 13 // L-vector -> E-vector, strided 14 //------------------------------------------------------------------------------ StridedNoTranspose(const CeedScalar * __restrict__ u,CeedScalar * __restrict__ v)15extern "C" __global__ void StridedNoTranspose(const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) { 16 for (CeedInt node = blockIdx.x * blockDim.x + threadIdx.x; node < RSTR_NUM_ELEM * RSTR_ELEM_SIZE; node += blockDim.x * gridDim.x) { 17 const CeedInt loc_node = node % RSTR_ELEM_SIZE; 18 const CeedInt elem = node / RSTR_ELEM_SIZE; 19 20 for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) { 21 v[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] = 22 u[loc_node * RSTR_STRIDE_NODES + comp * RSTR_STRIDE_COMP + elem * RSTR_STRIDE_ELEM]; 23 } 24 } 25 } 26 27 //------------------------------------------------------------------------------ 28 // E-vector -> L-vector, strided 29 //------------------------------------------------------------------------------ StridedTranspose(const CeedScalar * __restrict__ u,CeedScalar * __restrict__ v)30extern "C" __global__ void StridedTranspose(const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) { 31 for (CeedInt node = blockIdx.x * blockDim.x + threadIdx.x; node < RSTR_NUM_ELEM * RSTR_ELEM_SIZE; node += blockDim.x * gridDim.x) { 32 const CeedInt loc_node = node % RSTR_ELEM_SIZE; 33 const CeedInt elem = node / RSTR_ELEM_SIZE; 34 35 for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) { 36 v[loc_node * RSTR_STRIDE_NODES + comp * RSTR_STRIDE_COMP + elem * RSTR_STRIDE_ELEM] += 37 u[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE]; 38 } 39 } 40 } 41