xref: /libCEED/include/ceed/jit-source/hip/hip-ref-restriction-curl-oriented.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 HIP curl-oriented element restriction kernels
10cf8cbdd6SSebastian Grimberg 
11*c0b5abf0SJeremy L Thompson #include <ceed/types.h>
12cf8cbdd6SSebastian Grimberg 
13cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------
14cf8cbdd6SSebastian Grimberg // L-vector -> E-vector, curl-oriented
15cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------
16cf8cbdd6SSebastian Grimberg extern "C" __global__ void CurlOrientedNoTranspose(const CeedInt *__restrict__ indices, const CeedInt8 *__restrict__ curl_orients,
17cf8cbdd6SSebastian Grimberg                                                    const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) {
18cf8cbdd6SSebastian Grimberg   for (CeedInt node = blockIdx.x * blockDim.x + threadIdx.x; node < RSTR_NUM_ELEM * RSTR_ELEM_SIZE; node += blockDim.x * gridDim.x) {
19cf8cbdd6SSebastian Grimberg     const CeedInt  loc_node       = node % RSTR_ELEM_SIZE;
20cf8cbdd6SSebastian Grimberg     const CeedInt  elem           = node / RSTR_ELEM_SIZE;
21cf8cbdd6SSebastian Grimberg     const CeedInt  ind_dl         = loc_node > 0 ? indices[node - 1] : 0;
22cf8cbdd6SSebastian Grimberg     const CeedInt  ind_d          = indices[node];
23cf8cbdd6SSebastian Grimberg     const CeedInt  ind_du         = loc_node < (RSTR_ELEM_SIZE - 1) ? indices[node + 1] : 0;
24cf8cbdd6SSebastian Grimberg     const CeedInt8 curl_orient_dl = curl_orients[3 * node + 0];
25cf8cbdd6SSebastian Grimberg     const CeedInt8 curl_orient_d  = curl_orients[3 * node + 1];
26cf8cbdd6SSebastian Grimberg     const CeedInt8 curl_orient_du = curl_orients[3 * node + 2];
27cf8cbdd6SSebastian Grimberg 
28cf8cbdd6SSebastian Grimberg     for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) {
29cf8cbdd6SSebastian Grimberg       CeedScalar value = 0.0;
30cf8cbdd6SSebastian Grimberg       value += loc_node > 0 ? u[ind_dl + comp * RSTR_COMP_STRIDE] * curl_orient_dl : 0.0;
31cf8cbdd6SSebastian Grimberg       value += u[ind_d + comp * RSTR_COMP_STRIDE] * curl_orient_d;
32cf8cbdd6SSebastian Grimberg       value += loc_node < (RSTR_ELEM_SIZE - 1) ? u[ind_du + comp * RSTR_COMP_STRIDE] * curl_orient_du : 0.0;
33cf8cbdd6SSebastian Grimberg       v[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] = value;
34cf8cbdd6SSebastian Grimberg     }
35cf8cbdd6SSebastian Grimberg   }
36cf8cbdd6SSebastian Grimberg }
37cf8cbdd6SSebastian Grimberg 
38cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------
39cf8cbdd6SSebastian Grimberg // L-vector -> E-vector, unsigned curl-oriented
40cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------
41cf8cbdd6SSebastian Grimberg extern "C" __global__ void CurlOrientedUnsignedNoTranspose(const CeedInt *__restrict__ indices, const CeedInt8 *__restrict__ curl_orients,
42cf8cbdd6SSebastian Grimberg                                                            const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) {
43cf8cbdd6SSebastian Grimberg   for (CeedInt node = blockIdx.x * blockDim.x + threadIdx.x; node < RSTR_NUM_ELEM * RSTR_ELEM_SIZE; node += blockDim.x * gridDim.x) {
44cf8cbdd6SSebastian Grimberg     const CeedInt  loc_node       = node % RSTR_ELEM_SIZE;
45cf8cbdd6SSebastian Grimberg     const CeedInt  elem           = node / RSTR_ELEM_SIZE;
46cf8cbdd6SSebastian Grimberg     const CeedInt  ind_dl         = loc_node > 0 ? indices[node - 1] : 0;
47cf8cbdd6SSebastian Grimberg     const CeedInt  ind_d          = indices[node];
48cf8cbdd6SSebastian Grimberg     const CeedInt  ind_du         = loc_node < (RSTR_ELEM_SIZE - 1) ? indices[node + 1] : 0;
49cf8cbdd6SSebastian Grimberg     const CeedInt8 curl_orient_dl = abs(curl_orients[3 * node + 0]);
50cf8cbdd6SSebastian Grimberg     const CeedInt8 curl_orient_d  = abs(curl_orients[3 * node + 1]);
51cf8cbdd6SSebastian Grimberg     const CeedInt8 curl_orient_du = abs(curl_orients[3 * node + 2]);
52cf8cbdd6SSebastian Grimberg 
53cf8cbdd6SSebastian Grimberg     for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) {
54cf8cbdd6SSebastian Grimberg       CeedScalar value = 0.0;
55cf8cbdd6SSebastian Grimberg       value += loc_node > 0 ? u[ind_dl + comp * RSTR_COMP_STRIDE] * curl_orient_dl : 0.0;
56cf8cbdd6SSebastian Grimberg       value += u[ind_d + comp * RSTR_COMP_STRIDE] * curl_orient_d;
57cf8cbdd6SSebastian Grimberg       value += loc_node < (RSTR_ELEM_SIZE - 1) ? u[ind_du + comp * RSTR_COMP_STRIDE] * curl_orient_du : 0.0;
58cf8cbdd6SSebastian Grimberg       v[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] = value;
59cf8cbdd6SSebastian Grimberg     }
60cf8cbdd6SSebastian Grimberg   }
61cf8cbdd6SSebastian Grimberg }
62cf8cbdd6SSebastian Grimberg 
63cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------
64cf8cbdd6SSebastian Grimberg // E-vector -> L-vector, curl-oriented
65cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------
66cf8cbdd6SSebastian Grimberg #if !USE_DETERMINISTIC
67cf8cbdd6SSebastian Grimberg extern "C" __global__ void CurlOrientedTranspose(const CeedInt *__restrict__ indices, const CeedInt8 *__restrict__ curl_orients,
68cf8cbdd6SSebastian Grimberg                                                  const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) {
69cf8cbdd6SSebastian Grimberg   for (CeedInt node = blockIdx.x * blockDim.x + threadIdx.x; node < RSTR_NUM_ELEM * RSTR_ELEM_SIZE; node += blockDim.x * gridDim.x) {
70cf8cbdd6SSebastian Grimberg     const CeedInt  ind            = indices[node];
71cf8cbdd6SSebastian Grimberg     const CeedInt  loc_node       = node % RSTR_ELEM_SIZE;
72cf8cbdd6SSebastian Grimberg     const CeedInt  elem           = node / RSTR_ELEM_SIZE;
73cf8cbdd6SSebastian Grimberg     const CeedInt8 curl_orient_du = loc_node > 0 ? curl_orients[3 * node - 1] : 0.0;
74cf8cbdd6SSebastian Grimberg     const CeedInt8 curl_orient_d  = curl_orients[3 * node + 1];
75cf8cbdd6SSebastian Grimberg     const CeedInt8 curl_orient_dl = loc_node < (RSTR_ELEM_SIZE - 1) ? curl_orients[3 * node + 3] : 0.0;
76cf8cbdd6SSebastian Grimberg 
77cf8cbdd6SSebastian Grimberg     for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) {
78cf8cbdd6SSebastian Grimberg       CeedScalar value = 0.0;
79cf8cbdd6SSebastian Grimberg       value += loc_node > 0 ? u[loc_node - 1 + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_du : 0.0;
80cf8cbdd6SSebastian Grimberg       value += u[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_d;
81cf8cbdd6SSebastian Grimberg       value +=
82cf8cbdd6SSebastian Grimberg           loc_node < (RSTR_ELEM_SIZE - 1) ? u[loc_node + 1 + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_dl : 0.0;
83db2becc9SJeremy L Thompson       atomicAdd(&v[ind + comp * RSTR_COMP_STRIDE], value);
84cf8cbdd6SSebastian Grimberg     }
85cf8cbdd6SSebastian Grimberg   }
86cf8cbdd6SSebastian Grimberg }
87cf8cbdd6SSebastian Grimberg #else
88cf8cbdd6SSebastian Grimberg extern "C" __global__ void CurlOrientedTranspose(const CeedInt *__restrict__ l_vec_indices, const CeedInt *__restrict__ t_indices,
89cf8cbdd6SSebastian Grimberg                                                  const CeedInt *__restrict__ t_offsets, const CeedInt8 *__restrict__ curl_orients,
90cf8cbdd6SSebastian Grimberg                                                  const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) {
91cf8cbdd6SSebastian Grimberg   CeedScalar value[RSTR_NUM_COMP];
92cf8cbdd6SSebastian Grimberg 
93cf8cbdd6SSebastian Grimberg   for (CeedInt i = blockIdx.x * blockDim.x + threadIdx.x; i < RSTR_NUM_NODES; i += blockDim.x * gridDim.x) {
94cf8cbdd6SSebastian Grimberg     const CeedInt ind     = l_vec_indices[i];
95cf8cbdd6SSebastian Grimberg     const CeedInt range_1 = t_offsets[i];
96cf8cbdd6SSebastian Grimberg     const CeedInt range_N = t_offsets[i + 1];
97cf8cbdd6SSebastian Grimberg 
98cf8cbdd6SSebastian Grimberg     for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) value[comp] = 0.0;
99cf8cbdd6SSebastian Grimberg 
100cf8cbdd6SSebastian Grimberg     for (CeedInt j = range_1; j < range_N; j++) {
101cf8cbdd6SSebastian Grimberg       const CeedInt  t_ind          = t_indices[j];
102cf8cbdd6SSebastian Grimberg       const CeedInt  loc_node       = t_ind % RSTR_ELEM_SIZE;
103cf8cbdd6SSebastian Grimberg       const CeedInt  elem           = t_ind / RSTR_ELEM_SIZE;
104cf8cbdd6SSebastian Grimberg       const CeedInt8 curl_orient_du = loc_node > 0 ? curl_orients[3 * t_ind - 1] : 0.0;
105cf8cbdd6SSebastian Grimberg       const CeedInt8 curl_orient_d  = curl_orients[3 * t_ind + 1];
106cf8cbdd6SSebastian Grimberg       const CeedInt8 curl_orient_dl = loc_node < (RSTR_ELEM_SIZE - 1) ? curl_orients[3 * t_ind + 3] : 0.0;
107cf8cbdd6SSebastian Grimberg 
108cf8cbdd6SSebastian Grimberg       for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) {
109cf8cbdd6SSebastian Grimberg         value[comp] += loc_node > 0 ? u[loc_node - 1 + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_du : 0.0;
110cf8cbdd6SSebastian Grimberg         value[comp] += u[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_d;
111cf8cbdd6SSebastian Grimberg         value[comp] +=
112cf8cbdd6SSebastian Grimberg             loc_node < (RSTR_ELEM_SIZE - 1) ? u[loc_node + 1 + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_dl : 0.0;
113cf8cbdd6SSebastian Grimberg       }
114cf8cbdd6SSebastian Grimberg     }
115cf8cbdd6SSebastian Grimberg 
116cf8cbdd6SSebastian Grimberg     for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) v[ind + comp * RSTR_COMP_STRIDE] += value[comp];
117cf8cbdd6SSebastian Grimberg   }
118cf8cbdd6SSebastian Grimberg }
119cf8cbdd6SSebastian Grimberg #endif
120cf8cbdd6SSebastian Grimberg 
121cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------
122cf8cbdd6SSebastian Grimberg // E-vector -> L-vector, unsigned curl-oriented
123cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------
124cf8cbdd6SSebastian Grimberg #if !USE_DETERMINISTIC
125cf8cbdd6SSebastian Grimberg extern "C" __global__ void CurlOrientedUnsignedTranspose(const CeedInt *__restrict__ indices, const CeedInt8 *__restrict__ curl_orients,
126cf8cbdd6SSebastian Grimberg                                                          const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) {
127cf8cbdd6SSebastian Grimberg   for (CeedInt node = blockIdx.x * blockDim.x + threadIdx.x; node < RSTR_NUM_ELEM * RSTR_ELEM_SIZE; node += blockDim.x * gridDim.x) {
128cf8cbdd6SSebastian Grimberg     const CeedInt  loc_node       = node % RSTR_ELEM_SIZE;
129cf8cbdd6SSebastian Grimberg     const CeedInt  elem           = node / RSTR_ELEM_SIZE;
130cf8cbdd6SSebastian Grimberg     const CeedInt  ind            = indices[node];
131cf8cbdd6SSebastian Grimberg     const CeedInt8 curl_orient_du = loc_node > 0 ? abs(curl_orients[3 * node - 1]) : 0.0;
132cf8cbdd6SSebastian Grimberg     const CeedInt8 curl_orient_d  = abs(curl_orients[3 * node + 1]);
133cf8cbdd6SSebastian Grimberg     const CeedInt8 curl_orient_dl = loc_node < (RSTR_ELEM_SIZE - 1) ? abs(curl_orients[3 * node + 3]) : 0.0;
134cf8cbdd6SSebastian Grimberg 
135cf8cbdd6SSebastian Grimberg     for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) {
136cf8cbdd6SSebastian Grimberg       CeedScalar value = 0.0;
137cf8cbdd6SSebastian Grimberg       value += loc_node > 0 ? u[loc_node - 1 + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_du : 0.0;
138cf8cbdd6SSebastian Grimberg       value += u[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_d;
139cf8cbdd6SSebastian Grimberg       value +=
140cf8cbdd6SSebastian Grimberg           loc_node < (RSTR_ELEM_SIZE - 1) ? u[loc_node + 1 + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_dl : 0.0;
141db2becc9SJeremy L Thompson       atomicAdd(&v[ind + comp * RSTR_COMP_STRIDE], value);
142cf8cbdd6SSebastian Grimberg     }
143cf8cbdd6SSebastian Grimberg   }
144cf8cbdd6SSebastian Grimberg }
145cf8cbdd6SSebastian Grimberg #else
146cf8cbdd6SSebastian Grimberg extern "C" __global__ void CurlOrientedUnsignedTranspose(const CeedInt *__restrict__ l_vec_indices, const CeedInt *__restrict__ t_indices,
147cf8cbdd6SSebastian Grimberg                                                          const CeedInt *__restrict__ t_offsets, const CeedInt8 *__restrict__ curl_orients,
148cf8cbdd6SSebastian Grimberg                                                          const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) {
149cf8cbdd6SSebastian Grimberg   CeedScalar value[RSTR_NUM_COMP];
150cf8cbdd6SSebastian Grimberg 
151cf8cbdd6SSebastian Grimberg   for (CeedInt i = blockIdx.x * blockDim.x + threadIdx.x; i < RSTR_NUM_NODES; i += blockDim.x * gridDim.x) {
152cf8cbdd6SSebastian Grimberg     const CeedInt ind     = l_vec_indices[i];
153cf8cbdd6SSebastian Grimberg     const CeedInt range_1 = t_offsets[i];
154cf8cbdd6SSebastian Grimberg     const CeedInt range_N = t_offsets[i + 1];
155cf8cbdd6SSebastian Grimberg 
156cf8cbdd6SSebastian Grimberg     for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) value[comp] = 0.0;
157cf8cbdd6SSebastian Grimberg 
158cf8cbdd6SSebastian Grimberg     for (CeedInt j = range_1; j < range_N; j++) {
159cf8cbdd6SSebastian Grimberg       const CeedInt  t_ind          = t_indices[j];
160cf8cbdd6SSebastian Grimberg       const CeedInt  loc_node       = t_ind % RSTR_ELEM_SIZE;
161cf8cbdd6SSebastian Grimberg       const CeedInt  elem           = t_ind / RSTR_ELEM_SIZE;
162cf8cbdd6SSebastian Grimberg       const CeedInt8 curl_orient_du = loc_node > 0 ? abs(curl_orients[3 * t_ind - 1]) : 0.0;
163cf8cbdd6SSebastian Grimberg       const CeedInt8 curl_orient_d  = abs(curl_orients[3 * t_ind + 1]);
164cf8cbdd6SSebastian Grimberg       const CeedInt8 curl_orient_dl = loc_node < (RSTR_ELEM_SIZE - 1) ? abs(curl_orients[3 * t_ind + 3]) : 0.0;
165cf8cbdd6SSebastian Grimberg 
166cf8cbdd6SSebastian Grimberg       for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) {
167cf8cbdd6SSebastian Grimberg         value[comp] += loc_node > 0 ? u[loc_node - 1 + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_du : 0.0;
168cf8cbdd6SSebastian Grimberg         value[comp] += u[loc_node + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_d;
169cf8cbdd6SSebastian Grimberg         value[comp] +=
170cf8cbdd6SSebastian Grimberg             loc_node < (RSTR_ELEM_SIZE - 1) ? u[loc_node + 1 + comp * RSTR_ELEM_SIZE * RSTR_NUM_ELEM + elem * RSTR_ELEM_SIZE] * curl_orient_dl : 0.0;
171cf8cbdd6SSebastian Grimberg       }
172cf8cbdd6SSebastian Grimberg     }
173cf8cbdd6SSebastian Grimberg 
174cf8cbdd6SSebastian Grimberg     for (CeedInt comp = 0; comp < RSTR_NUM_COMP; comp++) v[ind + comp * RSTR_COMP_STRIDE] += value[comp];
175cf8cbdd6SSebastian Grimberg   }
176cf8cbdd6SSebastian Grimberg }
177cf8cbdd6SSebastian Grimberg #endif
178