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