xref: /libCEED/include/ceed/jit-source/hip/hip-ref-restriction-offset.h (revision d4cc18453651bd0f94c1a2e078b2646a92dafdcc)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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 HIP offset element restriction kernels
10c0b5abf0SJeremy L Thompson #include <ceed/types.h>
11cf8cbdd6SSebastian Grimberg 
12cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------
13cf8cbdd6SSebastian Grimberg // L-vector -> E-vector, standard (with offsets)
14cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------
OffsetNoTranspose(const CeedInt * __restrict__ indices,const CeedScalar * __restrict__ u,CeedScalar * __restrict__ v)15cf8cbdd6SSebastian Grimberg extern "C" __global__ void OffsetNoTranspose(const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) {
16cf8cbdd6SSebastian Grimberg   for (CeedInt node = blockIdx.x * blockDim.x + threadIdx.x; node < RSTR_NUM_ELEM * RSTR_ELEM_SIZE; node += blockDim.x * gridDim.x) {
17cf8cbdd6SSebastian Grimberg     const CeedInt ind      = indices[node];
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] = u[ind + comp * RSTR_COMP_STRIDE];
23cf8cbdd6SSebastian Grimberg     }
24cf8cbdd6SSebastian Grimberg   }
25cf8cbdd6SSebastian Grimberg }
26cf8cbdd6SSebastian Grimberg 
27cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------
28cf8cbdd6SSebastian Grimberg // E-vector -> L-vector, standard (with offsets)
29cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------
30cf8cbdd6SSebastian Grimberg #if !USE_DETERMINISTIC
OffsetTranspose(const CeedInt * __restrict__ indices,const CeedScalar * __restrict__ u,CeedScalar * __restrict__ v)31cf8cbdd6SSebastian Grimberg extern "C" __global__ void OffsetTranspose(const CeedInt *__restrict__ indices, 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 ind      = indices[node];
34cf8cbdd6SSebastian Grimberg     const CeedInt loc_node = node % RSTR_ELEM_SIZE;
35cf8cbdd6SSebastian Grimberg     const CeedInt elem     = node / RSTR_ELEM_SIZE;
36cf8cbdd6SSebastian Grimberg 
37cf8cbdd6SSebastian Grimberg     for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) {
38db2becc9SJeremy L Thompson       atomicAdd(&v[ind + comp * RSTR_COMP_STRIDE], u[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE]);
39cf8cbdd6SSebastian Grimberg     }
40cf8cbdd6SSebastian Grimberg   }
41cf8cbdd6SSebastian Grimberg }
42cf8cbdd6SSebastian Grimberg #else
OffsetTranspose(const CeedInt * __restrict__ l_vec_indices,const CeedInt * __restrict__ t_indices,const CeedInt * __restrict__ t_offsets,const CeedScalar * __restrict__ u,CeedScalar * __restrict__ v)43cf8cbdd6SSebastian Grimberg extern "C" __global__ void OffsetTranspose(const CeedInt *__restrict__ l_vec_indices, const CeedInt *__restrict__ t_indices,
44cf8cbdd6SSebastian Grimberg                                            const CeedInt *__restrict__ t_offsets, const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) {
45cf8cbdd6SSebastian Grimberg   CeedScalar value[RSTR_NUM_COMP];
46cf8cbdd6SSebastian Grimberg 
47cf8cbdd6SSebastian Grimberg   for (CeedInt i = blockIdx.x * blockDim.x + threadIdx.x; i < RSTR_NUM_NODES; i += blockDim.x * gridDim.x) {
48cf8cbdd6SSebastian Grimberg     const CeedInt ind     = l_vec_indices[i];
49cf8cbdd6SSebastian Grimberg     const CeedInt range_1 = t_offsets[i];
50cf8cbdd6SSebastian Grimberg     const CeedInt range_N = t_offsets[i + 1];
51cf8cbdd6SSebastian Grimberg 
52cf8cbdd6SSebastian Grimberg     for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) value[comp] = 0.0;
53cf8cbdd6SSebastian Grimberg 
54cf8cbdd6SSebastian Grimberg     for (CeedInt j = range_1; j < range_N; j++) {
55cf8cbdd6SSebastian Grimberg       const CeedInt t_ind    = t_indices[j];
56cf8cbdd6SSebastian Grimberg       const CeedInt loc_node = t_ind % RSTR_ELEM_SIZE;
57cf8cbdd6SSebastian Grimberg       const CeedInt elem     = t_ind / RSTR_ELEM_SIZE;
58cf8cbdd6SSebastian Grimberg 
59cf8cbdd6SSebastian Grimberg       for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) {
60cf8cbdd6SSebastian Grimberg         value[comp] += u[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE];
61cf8cbdd6SSebastian Grimberg       }
62cf8cbdd6SSebastian Grimberg     }
63cf8cbdd6SSebastian Grimberg 
64cf8cbdd6SSebastian Grimberg     for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) v[ind + comp * RSTR_COMP_STRIDE] += value[comp];
65cf8cbdd6SSebastian Grimberg   }
66cf8cbdd6SSebastian Grimberg }
67cf8cbdd6SSebastian Grimberg #endif
68