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