xref: /libCEED/include/ceed/jit-source/hip/hip-gen-templates.h (revision c2cc34eeef4dafbcd1d00ef770c45974c5ea3da2)
1 // Copyright (c) 2017-2024, 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 HIP backend macro and type definitions for JiT source
10 
11 #include <ceed.h>
12 
13 //------------------------------------------------------------------------------
14 // Load matrices for basis actions
15 //------------------------------------------------------------------------------
16 template <int P, int Q>
17 inline __device__ void loadMatrix(SharedData_Hip &data, const CeedScalar *__restrict__ d_B, CeedScalar *B) {
18   for (CeedInt i = data.t_id; i < P * Q; i += blockDim.x * blockDim.y * blockDim.z) B[i] = d_B[i];
19 }
20 
21 //------------------------------------------------------------------------------
22 // 1D
23 //------------------------------------------------------------------------------
24 
25 //------------------------------------------------------------------------------
26 // L-vector -> E-vector, offsets provided
27 //------------------------------------------------------------------------------
28 template <int NUM_COMP, int COMP_STRIDE, int P_1d>
29 inline __device__ void readDofsOffset1d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
30                                         const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
31   if (data.t_id_x < P_1d) {
32     const CeedInt node = data.t_id_x;
33     const CeedInt ind  = indices[node + elem * P_1d];
34 
35     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp];
36   }
37 }
38 
39 //------------------------------------------------------------------------------
40 // L-vector -> E-vector, strided
41 //------------------------------------------------------------------------------
42 template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
43 inline __device__ void readDofsStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
44   if (data.t_id_x < P_1d) {
45     const CeedInt node = data.t_id_x;
46     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
47 
48     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
49   }
50 }
51 
52 //------------------------------------------------------------------------------
53 // E-vector -> L-vector, offsets provided
54 //------------------------------------------------------------------------------
55 template <int NUM_COMP, int COMP_STRIDE, int P_1d>
56 inline __device__ void writeDofsOffset1d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
57                                          const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) {
58   if (data.t_id_x < P_1d) {
59     const CeedInt node = data.t_id_x;
60     const CeedInt ind  = indices[node + elem * P_1d];
61 
62     for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]);
63   }
64 }
65 
66 //------------------------------------------------------------------------------
67 // E-vector -> L-vector, strided
68 //------------------------------------------------------------------------------
69 template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
70 inline __device__ void writeDofsStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ r_v,
71                                           CeedScalar *__restrict__ d_v) {
72   if (data.t_id_x < P_1d) {
73     const CeedInt node = data.t_id_x;
74     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
75 
76     for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp];
77   }
78 }
79 
80 //------------------------------------------------------------------------------
81 // 2D
82 //------------------------------------------------------------------------------
83 
84 //------------------------------------------------------------------------------
85 // L-vector -> E-vector, offsets provided
86 //------------------------------------------------------------------------------
87 template <int NUM_COMP, int COMP_STRIDE, int P_1d>
88 inline __device__ void readDofsOffset2d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
89                                         const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
90   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
91     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
92     const CeedInt ind  = indices[node + elem * P_1d * P_1d];
93 
94     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp];
95   }
96 }
97 
98 //------------------------------------------------------------------------------
99 // L-vector -> E-vector, strided
100 //------------------------------------------------------------------------------
101 template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
102 inline __device__ void readDofsStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
103   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
104     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
105     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
106 
107     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
108   }
109 }
110 
111 //------------------------------------------------------------------------------
112 // E-vector -> L-vector, offsets provided
113 //------------------------------------------------------------------------------
114 template <int NUM_COMP, int COMP_STRIDE, int P_1d>
115 inline __device__ void writeDofsOffset2d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
116                                          const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) {
117   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
118     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
119     const CeedInt ind  = indices[node + elem * P_1d * P_1d];
120 
121     for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]);
122   }
123 }
124 
125 //------------------------------------------------------------------------------
126 // E-vector -> L-vector, strided
127 //------------------------------------------------------------------------------
128 template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
129 inline __device__ void writeDofsStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ r_v,
130                                           CeedScalar *__restrict__ d_v) {
131   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
132     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
133     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
134 
135     for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp];
136   }
137 }
138 
139 //------------------------------------------------------------------------------
140 // 3D
141 //------------------------------------------------------------------------------
142 
143 //------------------------------------------------------------------------------
144 // L-vector -> E-vector, offsets provided
145 //------------------------------------------------------------------------------
146 // TODO: remove "Dofs" and "Quads" in the following function names?
147 //   - readDofsOffset3d -> readOffset3d ?
148 //   - readDofsStrided3d -> readStrided3d ?
149 //   - readSliceQuadsOffset3d -> readSliceOffset3d ?
150 //   - readSliceQuadsStrided3d -> readSliceStrided3d ?
151 //   - writeDofsOffset3d -> writeOffset3d ?
152 //   - writeDofsStrided3d -> writeStrided3d ?
153 template <int NUM_COMP, int COMP_STRIDE, int P_1d>
154 inline __device__ void readDofsOffset3d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
155                                         const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
156   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
157     for (CeedInt z = 0; z < P_1d; z++) {
158       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
159       const CeedInt ind  = indices[node + elem * P_1d * P_1d * P_1d];
160 
161       for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1d] = d_u[ind + COMP_STRIDE * comp];
162     }
163 }
164 
165 //------------------------------------------------------------------------------
166 // L-vector -> E-vector, strided
167 //------------------------------------------------------------------------------
168 template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
169 inline __device__ void readDofsStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
170   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
171     for (CeedInt z = 0; z < P_1d; z++) {
172       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
173       const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
174 
175       for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1d] = d_u[ind + comp * STRIDES_COMP];
176     }
177 }
178 
179 //------------------------------------------------------------------------------
180 // E-vector -> Q-vector, offests provided
181 //------------------------------------------------------------------------------
182 template <int NUM_COMP, int COMP_STRIDE, int Q_1d>
183 inline __device__ void readSliceQuadsOffset3d(SharedData_Hip &data, const CeedInt nquads, const CeedInt elem, const CeedInt q,
184                                               const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
185   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
186     const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d;
187     const CeedInt ind  = indices[node + elem * Q_1d * Q_1d * Q_1d];
188 
189     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp];
190   }
191 }
192 
193 //------------------------------------------------------------------------------
194 // E-vector -> Q-vector, strided
195 //------------------------------------------------------------------------------
196 template <int NUM_COMP, int Q_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
197 inline __device__ void readSliceQuadsStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedInt q, const CeedScalar *__restrict__ d_u,
198                                                CeedScalar *__restrict__ r_u) {
199   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
200     const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d;
201     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
202 
203     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
204   }
205 }
206 
207 //------------------------------------------------------------------------------
208 // E-vector -> L-vector, offsets provided
209 //------------------------------------------------------------------------------
210 template <int NUM_COMP, int COMP_STRIDE, int P_1d>
211 inline __device__ void writeDofsOffset3d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
212                                          const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) {
213   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
214     for (CeedInt z = 0; z < P_1d; z++) {
215       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
216       const CeedInt ind  = indices[node + elem * P_1d * P_1d * P_1d];
217 
218       for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[z + comp * P_1d]);
219     }
220 }
221 
222 //------------------------------------------------------------------------------
223 // E-vector -> L-vector, strided
224 //------------------------------------------------------------------------------
225 template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
226 inline __device__ void writeDofsStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ r_v,
227                                           CeedScalar *__restrict__ d_v) {
228   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
229     for (CeedInt z = 0; z < P_1d; z++) {
230       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
231       const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
232 
233       for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[z + comp * P_1d];
234     }
235 }
236 
237 //------------------------------------------------------------------------------
238 // 3D collocated derivatives computation
239 //------------------------------------------------------------------------------
240 template <int NUM_COMP, int Q_1d>
241 inline __device__ void gradCollo3d(SharedData_Hip &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
242                                    CeedScalar *__restrict__ r_V) {
243   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
244     for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
245       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[q + comp * Q_1d];
246       __syncthreads();
247       // X derivative
248       r_V[comp + 0 * NUM_COMP] = 0.0;
249       for (CeedInt i = 0; i < Q_1d; i++)
250         r_V[comp + 0 * NUM_COMP] += c_G[i + data.t_id_x * Q_1d] * data.slice[i + data.t_id_y * T_1D];  // Contract x direction (X derivative)
251       // Y derivative
252       r_V[comp + 1 * NUM_COMP] = 0.0;
253       for (CeedInt i = 0; i < Q_1d; i++)
254         r_V[comp + 1 * NUM_COMP] += c_G[i + data.t_id_y * Q_1d] * data.slice[data.t_id_x + i * T_1D];  // Contract y direction (Y derivative)
255       // Z derivative
256       r_V[comp + 2 * NUM_COMP] = 0.0;
257       for (CeedInt i = 0; i < Q_1d; i++) r_V[comp + 2 * NUM_COMP] += c_G[i + q * Q_1d] * r_U[i + comp * Q_1d];  // Contract z direction (Z derivative)
258       __syncthreads();
259     }
260   }
261 }
262 
263 //------------------------------------------------------------------------------
264 // 3D collocated derivatives transpose
265 //------------------------------------------------------------------------------
266 template <int NUM_COMP, int Q_1d>
267 inline __device__ void gradColloTranspose3d(SharedData_Hip &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
268                                             CeedScalar *__restrict__ r_V) {
269   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
270     for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
271       // X derivative
272       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 0 * NUM_COMP];
273       __syncthreads();
274       for (CeedInt i = 0; i < Q_1d; i++)
275         r_V[q + comp * Q_1d] += c_G[data.t_id_x + i * Q_1d] * data.slice[i + data.t_id_y * T_1D];  // Contract x direction (X derivative)
276       __syncthreads();
277       // Y derivative
278       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 1 * NUM_COMP];
279       __syncthreads();
280       for (CeedInt i = 0; i < Q_1d; i++)
281         r_V[q + comp * Q_1d] += c_G[data.t_id_y + i * Q_1d] * data.slice[data.t_id_x + i * T_1D];  // Contract y direction (Y derivative)
282       __syncthreads();
283       // Z derivative
284       for (CeedInt i = 0; i < Q_1d; i++)
285         r_V[i + comp * Q_1d] += c_G[i + q * Q_1d] * r_U[comp + 2 * NUM_COMP];  // PARTIAL contract z direction (Z derivative)
286     }
287   }
288 }
289