xref: /libCEED/include/ceed/jit-source/cuda/cuda-gen-templates.h (revision 99421279ebf997f2fb157e1e18f23b9dc112bbac)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
29e201c85SYohann // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
39e201c85SYohann //
49e201c85SYohann // SPDX-License-Identifier: BSD-2-Clause
59e201c85SYohann //
69e201c85SYohann // This file is part of CEED:  http://github.com/ceed
79e201c85SYohann 
89e201c85SYohann /// @file
99e201c85SYohann /// Internal header for CUDA backend macro and type definitions for JiT source
10c0b5abf0SJeremy L Thompson #include <ceed/types.h>
119e201c85SYohann 
129e201c85SYohann //------------------------------------------------------------------------------
139e201c85SYohann // Load matrices for basis actions
149e201c85SYohann //------------------------------------------------------------------------------
159e201c85SYohann template <int P, int Q>
16f815fac9SJeremy L Thompson inline __device__ void LoadMatrix(SharedData_Cuda &data, const CeedScalar *__restrict__ d_B, CeedScalar *B) {
172b730f8bSJeremy L Thompson   for (CeedInt i = data.t_id; i < P * Q; i += blockDim.x * blockDim.y * blockDim.z) B[i] = d_B[i];
189e201c85SYohann }
199e201c85SYohann 
209e201c85SYohann //------------------------------------------------------------------------------
218b97b69aSJeremy L Thompson // AtPoints
228b97b69aSJeremy L Thompson //------------------------------------------------------------------------------
238b97b69aSJeremy L Thompson 
248b97b69aSJeremy L Thompson //------------------------------------------------------------------------------
258b97b69aSJeremy L Thompson // L-vector -> single point
268b97b69aSJeremy L Thompson //------------------------------------------------------------------------------
278b97b69aSJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int NUM_PTS>
288b97b69aSJeremy L Thompson inline __device__ void ReadPoint(SharedData_Cuda &data, const CeedInt elem, const CeedInt p, const CeedInt points_in_elem,
298b97b69aSJeremy L Thompson                                  const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) {
308b97b69aSJeremy L Thompson   const CeedInt ind = indices[p + elem * NUM_PTS];
318b97b69aSJeremy L Thompson 
328b97b69aSJeremy L Thompson   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
338b97b69aSJeremy L Thompson     r_u[comp] = d_u[ind + comp * COMP_STRIDE];
348b97b69aSJeremy L Thompson   }
358b97b69aSJeremy L Thompson }
368b97b69aSJeremy L Thompson 
378b97b69aSJeremy L Thompson //------------------------------------------------------------------------------
388b97b69aSJeremy L Thompson // Single point -> L-vector
398b97b69aSJeremy L Thompson //------------------------------------------------------------------------------
408b97b69aSJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int NUM_PTS>
418b97b69aSJeremy L Thompson inline __device__ void WritePoint(SharedData_Cuda &data, const CeedInt elem, const CeedInt p, const CeedInt points_in_elem,
428b97b69aSJeremy L Thompson                                   const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ r_u, CeedScalar *d_u) {
438b97b69aSJeremy L Thompson   if (p < points_in_elem) {
448b97b69aSJeremy L Thompson     const CeedInt ind = indices[p + elem * NUM_PTS];
458b97b69aSJeremy L Thompson 
468b97b69aSJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
478b97b69aSJeremy L Thompson       d_u[ind + comp * COMP_STRIDE] += r_u[comp];
488b97b69aSJeremy L Thompson     }
498b97b69aSJeremy L Thompson   }
508b97b69aSJeremy L Thompson }
518b97b69aSJeremy L Thompson 
528b97b69aSJeremy L Thompson //------------------------------------------------------------------------------
539e201c85SYohann // 1D
549e201c85SYohann //------------------------------------------------------------------------------
559e201c85SYohann 
569e201c85SYohann //------------------------------------------------------------------------------
579e201c85SYohann // L-vector -> E-vector, offsets provided
589e201c85SYohann //------------------------------------------------------------------------------
59*99421279SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1D>
60f815fac9SJeremy L Thompson inline __device__ void ReadLVecStandard1d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
61672b0f2aSSebastian Grimberg                                           const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
62*99421279SJeremy L Thompson   if (data.t_id_x < P_1D) {
639e201c85SYohann     const CeedInt node = data.t_id_x;
64*99421279SJeremy L Thompson     const CeedInt ind  = indices[node + elem * P_1D];
65ca735530SJeremy L Thompson 
66ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp];
679e201c85SYohann   }
689e201c85SYohann }
699e201c85SYohann 
709e201c85SYohann //------------------------------------------------------------------------------
719e201c85SYohann // L-vector -> E-vector, strided
729e201c85SYohann //------------------------------------------------------------------------------
73*99421279SJeremy L Thompson template <int NUM_COMP, int P_1D, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
74f815fac9SJeremy L Thompson inline __device__ void ReadLVecStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u,
75672b0f2aSSebastian Grimberg                                          CeedScalar *__restrict__ r_u) {
76*99421279SJeremy L Thompson   if (data.t_id_x < P_1D) {
779e201c85SYohann     const CeedInt node = data.t_id_x;
789e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
79ca735530SJeremy L Thompson 
80ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
819e201c85SYohann   }
829e201c85SYohann }
839e201c85SYohann 
849e201c85SYohann //------------------------------------------------------------------------------
859e201c85SYohann // E-vector -> L-vector, offsets provided
869e201c85SYohann //------------------------------------------------------------------------------
87*99421279SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1D>
88f815fac9SJeremy L Thompson inline __device__ void WriteLVecStandard1d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
89672b0f2aSSebastian Grimberg                                            const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) {
90*99421279SJeremy L Thompson   if (data.t_id_x < P_1D) {
919e201c85SYohann     const CeedInt node = data.t_id_x;
92*99421279SJeremy L Thompson     const CeedInt ind  = indices[node + elem * P_1D];
93ca735530SJeremy L Thompson 
94ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]);
959e201c85SYohann   }
969e201c85SYohann }
979e201c85SYohann 
989e201c85SYohann //------------------------------------------------------------------------------
999e201c85SYohann // E-vector -> L-vector, strided
1009e201c85SYohann //------------------------------------------------------------------------------
101*99421279SJeremy L Thompson template <int NUM_COMP, int P_1D, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
102f815fac9SJeremy L Thompson inline __device__ void WriteLVecStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ r_v,
103672b0f2aSSebastian Grimberg                                           CeedScalar *__restrict__ d_v) {
104*99421279SJeremy L Thompson   if (data.t_id_x < P_1D) {
1059e201c85SYohann     const CeedInt node = data.t_id_x;
1069e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
107ca735530SJeremy L Thompson 
108ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp];
1099e201c85SYohann   }
1109e201c85SYohann }
1119e201c85SYohann 
1129e201c85SYohann //------------------------------------------------------------------------------
1139e201c85SYohann // 2D
1149e201c85SYohann //------------------------------------------------------------------------------
1159e201c85SYohann 
1169e201c85SYohann //------------------------------------------------------------------------------
1179e201c85SYohann // L-vector -> E-vector, offsets provided
1189e201c85SYohann //------------------------------------------------------------------------------
119*99421279SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1D>
120f815fac9SJeremy L Thompson inline __device__ void ReadLVecStandard2d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
121672b0f2aSSebastian Grimberg                                           const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
122*99421279SJeremy L Thompson   if (data.t_id_x < P_1D && data.t_id_y < P_1D) {
123*99421279SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * P_1D;
124*99421279SJeremy L Thompson     const CeedInt ind  = indices[node + elem * P_1D * P_1D];
125ca735530SJeremy L Thompson 
126ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp];
1279e201c85SYohann   }
1289e201c85SYohann }
1299e201c85SYohann 
1309e201c85SYohann //------------------------------------------------------------------------------
1319e201c85SYohann // L-vector -> E-vector, strided
1329e201c85SYohann //------------------------------------------------------------------------------
133*99421279SJeremy L Thompson template <int NUM_COMP, int P_1D, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
134f815fac9SJeremy L Thompson inline __device__ void ReadLVecStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u,
135672b0f2aSSebastian Grimberg                                          CeedScalar *__restrict__ r_u) {
136*99421279SJeremy L Thompson   if (data.t_id_x < P_1D && data.t_id_y < P_1D) {
137*99421279SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * P_1D;
1389e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
139ca735530SJeremy L Thompson 
140ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
1419e201c85SYohann   }
1429e201c85SYohann }
1439e201c85SYohann 
1449e201c85SYohann //------------------------------------------------------------------------------
1459e201c85SYohann // E-vector -> L-vector, offsets provided
1469e201c85SYohann //------------------------------------------------------------------------------
147*99421279SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1D>
148f815fac9SJeremy L Thompson inline __device__ void WriteLVecStandard2d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
149672b0f2aSSebastian Grimberg                                            const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) {
150*99421279SJeremy L Thompson   if (data.t_id_x < P_1D && data.t_id_y < P_1D) {
151*99421279SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * P_1D;
152*99421279SJeremy L Thompson     const CeedInt ind  = indices[node + elem * P_1D * P_1D];
153ca735530SJeremy L Thompson 
154ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]);
1559e201c85SYohann   }
1569e201c85SYohann }
1579e201c85SYohann 
1589e201c85SYohann //------------------------------------------------------------------------------
1599e201c85SYohann // E-vector -> L-vector, strided
1609e201c85SYohann //------------------------------------------------------------------------------
161*99421279SJeremy L Thompson template <int NUM_COMP, int P_1D, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
162f815fac9SJeremy L Thompson inline __device__ void WriteLVecStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ r_v,
163672b0f2aSSebastian Grimberg                                           CeedScalar *__restrict__ d_v) {
164*99421279SJeremy L Thompson   if (data.t_id_x < P_1D && data.t_id_y < P_1D) {
165*99421279SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * P_1D;
1669e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
167ca735530SJeremy L Thompson 
168ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp];
1699e201c85SYohann   }
1709e201c85SYohann }
1719e201c85SYohann 
1729e201c85SYohann //------------------------------------------------------------------------------
1739e201c85SYohann // 3D
1749e201c85SYohann //------------------------------------------------------------------------------
1759e201c85SYohann 
1769e201c85SYohann //------------------------------------------------------------------------------
1779e201c85SYohann // L-vector -> E-vector, offsets provided
1789e201c85SYohann //------------------------------------------------------------------------------
179*99421279SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1D>
180f815fac9SJeremy L Thompson inline __device__ void ReadLVecStandard3d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
181672b0f2aSSebastian Grimberg                                           const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
182*99421279SJeremy L Thompson   if (data.t_id_x < P_1D && data.t_id_y < P_1D) {
183*99421279SJeremy L Thompson     for (CeedInt z = 0; z < P_1D; z++) {
184*99421279SJeremy L Thompson       const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D;
185*99421279SJeremy L Thompson       const CeedInt ind  = indices[node + elem * P_1D * P_1D * P_1D];
186ca735530SJeremy L Thompson 
187*99421279SJeremy L Thompson       for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1D] = d_u[ind + COMP_STRIDE * comp];
1889e201c85SYohann     }
1899e201c85SYohann   }
190688b5473SJeremy L Thompson }
1919e201c85SYohann 
1929e201c85SYohann //------------------------------------------------------------------------------
1939e201c85SYohann // L-vector -> E-vector, strided
1949e201c85SYohann //------------------------------------------------------------------------------
195*99421279SJeremy L Thompson template <int NUM_COMP, int P_1D, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
196f815fac9SJeremy L Thompson inline __device__ void ReadLVecStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u,
197672b0f2aSSebastian Grimberg                                          CeedScalar *__restrict__ r_u) {
198*99421279SJeremy L Thompson   if (data.t_id_x < P_1D && data.t_id_y < P_1D) {
199*99421279SJeremy L Thompson     for (CeedInt z = 0; z < P_1D; z++) {
200*99421279SJeremy L Thompson       const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D;
2019e201c85SYohann       const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
202ca735530SJeremy L Thompson 
203*99421279SJeremy L Thompson       for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1D] = d_u[ind + comp * STRIDES_COMP];
2049e201c85SYohann     }
2059e201c85SYohann   }
206688b5473SJeremy L Thompson }
2079e201c85SYohann 
2089e201c85SYohann //------------------------------------------------------------------------------
2099e201c85SYohann // E-vector -> Q-vector, offests provided
2109e201c85SYohann //------------------------------------------------------------------------------
211*99421279SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int Q_1D>
212f815fac9SJeremy L Thompson inline __device__ void ReadEVecSliceStandard3d(SharedData_Cuda &data, const CeedInt nquads, const CeedInt elem, const CeedInt q,
213f815fac9SJeremy L Thompson                                                const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ d_u,
214f815fac9SJeremy L Thompson                                                CeedScalar *__restrict__ r_u) {
215*99421279SJeremy L Thompson   if (data.t_id_x < Q_1D && data.t_id_y < Q_1D) {
216*99421279SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * Q_1D + q * Q_1D * Q_1D;
217*99421279SJeremy L Thompson     const CeedInt ind  = indices[node + elem * Q_1D * Q_1D * Q_1D];
218ca735530SJeremy L Thompson 
219ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp];
2209e201c85SYohann   }
2219e201c85SYohann }
2229e201c85SYohann 
2239e201c85SYohann //------------------------------------------------------------------------------
2249e201c85SYohann // E-vector -> Q-vector, strided
2259e201c85SYohann //------------------------------------------------------------------------------
226*99421279SJeremy L Thompson template <int NUM_COMP, int Q_1D, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
227f815fac9SJeremy L Thompson inline __device__ void ReadEVecSliceStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt q, const CeedScalar *__restrict__ d_u,
228672b0f2aSSebastian Grimberg                                               CeedScalar *__restrict__ r_u) {
229*99421279SJeremy L Thompson   if (data.t_id_x < Q_1D && data.t_id_y < Q_1D) {
230*99421279SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * Q_1D + q * Q_1D * Q_1D;
2319e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
232ca735530SJeremy L Thompson 
233ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
2349e201c85SYohann   }
2359e201c85SYohann }
2369e201c85SYohann 
2379e201c85SYohann //------------------------------------------------------------------------------
2389e201c85SYohann // E-vector -> L-vector, offsets provided
2399e201c85SYohann //------------------------------------------------------------------------------
240*99421279SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1D>
241f815fac9SJeremy L Thompson inline __device__ void WriteLVecStandard3d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
242672b0f2aSSebastian Grimberg                                            const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) {
243*99421279SJeremy L Thompson   if (data.t_id_x < P_1D && data.t_id_y < P_1D) {
244*99421279SJeremy L Thompson     for (CeedInt z = 0; z < P_1D; z++) {
245*99421279SJeremy L Thompson       const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D;
246*99421279SJeremy L Thompson       const CeedInt ind  = indices[node + elem * P_1D * P_1D * P_1D];
247ca735530SJeremy L Thompson 
248*99421279SJeremy L Thompson       for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[z + comp * P_1D]);
2499e201c85SYohann     }
2509e201c85SYohann   }
251688b5473SJeremy L Thompson }
2529e201c85SYohann 
2539e201c85SYohann //------------------------------------------------------------------------------
2549e201c85SYohann // E-vector -> L-vector, strided
2559e201c85SYohann //------------------------------------------------------------------------------
256*99421279SJeremy L Thompson template <int NUM_COMP, int P_1D, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
257f815fac9SJeremy L Thompson inline __device__ void WriteLVecStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ r_v,
258672b0f2aSSebastian Grimberg                                           CeedScalar *__restrict__ d_v) {
259*99421279SJeremy L Thompson   if (data.t_id_x < P_1D && data.t_id_y < P_1D) {
260*99421279SJeremy L Thompson     for (CeedInt z = 0; z < P_1D; z++) {
261*99421279SJeremy L Thompson       const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D;
2629e201c85SYohann       const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
263ca735530SJeremy L Thompson 
264*99421279SJeremy L Thompson       for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[z + comp * P_1D];
2659e201c85SYohann     }
2669e201c85SYohann   }
267688b5473SJeremy L Thompson }
2689e201c85SYohann 
2699e201c85SYohann //------------------------------------------------------------------------------
2709e201c85SYohann // 3D collocated derivatives computation
2719e201c85SYohann //------------------------------------------------------------------------------
272*99421279SJeremy L Thompson template <int NUM_COMP, int Q_1D, int T_1D>
273f815fac9SJeremy L Thompson inline __device__ void GradColloSlice3d(SharedData_Cuda &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
2742b730f8bSJeremy L Thompson                                         CeedScalar *__restrict__ r_V) {
275*99421279SJeremy L Thompson   if (data.t_id_x < Q_1D && data.t_id_y < Q_1D) {
276ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
277*99421279SJeremy L Thompson       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[q + comp * Q_1D];
2789e201c85SYohann       __syncthreads();
2799e201c85SYohann       // X derivative
280ca735530SJeremy L Thompson       r_V[comp + 0 * NUM_COMP] = 0.0;
281*99421279SJeremy L Thompson       for (CeedInt i = 0; i < Q_1D; i++) {
282*99421279SJeremy L Thompson         r_V[comp + 0 * NUM_COMP] += c_G[i + data.t_id_x * Q_1D] * data.slice[i + data.t_id_y * T_1D];
283688b5473SJeremy L Thompson       }
2849e201c85SYohann       // Y derivative
285ca735530SJeremy L Thompson       r_V[comp + 1 * NUM_COMP] = 0.0;
286*99421279SJeremy L Thompson       for (CeedInt i = 0; i < Q_1D; i++) {
287*99421279SJeremy L Thompson         r_V[comp + 1 * NUM_COMP] += c_G[i + data.t_id_y * Q_1D] * data.slice[data.t_id_x + i * T_1D];
288688b5473SJeremy L Thompson       }
2899e201c85SYohann       // Z derivative
290ca735530SJeremy L Thompson       r_V[comp + 2 * NUM_COMP] = 0.0;
291*99421279SJeremy L Thompson       for (CeedInt i = 0; i < Q_1D; i++) {
292*99421279SJeremy L Thompson         r_V[comp + 2 * NUM_COMP] += c_G[i + q * Q_1D] * r_U[i + comp * Q_1D];
293688b5473SJeremy L Thompson       }
2949e201c85SYohann       __syncthreads();
2959e201c85SYohann     }
2969e201c85SYohann   }
2979e201c85SYohann }
2989e201c85SYohann 
2999e201c85SYohann //------------------------------------------------------------------------------
3009e201c85SYohann // 3D collocated derivatives transpose
3019e201c85SYohann //------------------------------------------------------------------------------
302*99421279SJeremy L Thompson template <int NUM_COMP, int Q_1D, int T_1D>
303f815fac9SJeremy L Thompson inline __device__ void GradColloSliceTranspose3d(SharedData_Cuda &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
3042b730f8bSJeremy L Thompson                                                  CeedScalar *__restrict__ r_V) {
305*99421279SJeremy L Thompson   if (data.t_id_x < Q_1D && data.t_id_y < Q_1D) {
306ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
307ca735530SJeremy L Thompson       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 0 * NUM_COMP];
3089e201c85SYohann       __syncthreads();
309688b5473SJeremy L Thompson       // X derivative
310*99421279SJeremy L Thompson       for (CeedInt i = 0; i < Q_1D; i++) {
311*99421279SJeremy L Thompson         r_V[q + comp * Q_1D] += c_G[data.t_id_x + i * Q_1D] * data.slice[i + data.t_id_y * T_1D];
312688b5473SJeremy L Thompson       }
3139e201c85SYohann       __syncthreads();
3149e201c85SYohann       // Y derivative
315ca735530SJeremy L Thompson       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 1 * NUM_COMP];
3169e201c85SYohann       __syncthreads();
317*99421279SJeremy L Thompson       for (CeedInt i = 0; i < Q_1D; i++) {
318*99421279SJeremy L Thompson         r_V[q + comp * Q_1D] += c_G[data.t_id_y + i * Q_1D] * data.slice[data.t_id_x + i * T_1D];
319688b5473SJeremy L Thompson       }
3209e201c85SYohann       __syncthreads();
3219e201c85SYohann       // Z derivative
322*99421279SJeremy L Thompson       for (CeedInt i = 0; i < Q_1D; i++) {
323*99421279SJeremy L Thompson         r_V[i + comp * Q_1D] += c_G[i + q * Q_1D] * r_U[comp + 2 * NUM_COMP];
324688b5473SJeremy L Thompson       }
3259e201c85SYohann     }
3269e201c85SYohann   }
3279e201c85SYohann }
328