xref: /libCEED/rust/libceed-sys/c-src/include/ceed/jit-source/hip/hip-gen-templates.h (revision c0b5abf0f23b15c4f0ada76f8abe9f8d2b6fa247)
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 HIP backend macro and type definitions for JiT source
109e201c85SYohann 
11*c0b5abf0SJeremy L Thompson #include <ceed/types.h>
129e201c85SYohann 
139e201c85SYohann //------------------------------------------------------------------------------
149e201c85SYohann // Load matrices for basis actions
159e201c85SYohann //------------------------------------------------------------------------------
169e201c85SYohann template <int P, int Q>
17672b0f2aSSebastian Grimberg inline __device__ void loadMatrix(SharedData_Hip &data, const CeedScalar *__restrict__ d_B, CeedScalar *B) {
182b730f8bSJeremy L Thompson   for (CeedInt i = data.t_id; i < P * Q; i += blockDim.x * blockDim.y * blockDim.z) B[i] = d_B[i];
199e201c85SYohann }
209e201c85SYohann 
219e201c85SYohann //------------------------------------------------------------------------------
229e201c85SYohann // 1D
239e201c85SYohann //------------------------------------------------------------------------------
249e201c85SYohann 
259e201c85SYohann //------------------------------------------------------------------------------
269e201c85SYohann // L-vector -> E-vector, offsets provided
279e201c85SYohann //------------------------------------------------------------------------------
28672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d>
29672b0f2aSSebastian Grimberg inline __device__ void readDofsOffset1d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
30672b0f2aSSebastian Grimberg                                         const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
31672b0f2aSSebastian Grimberg   if (data.t_id_x < P_1d) {
329e201c85SYohann     const CeedInt node = data.t_id_x;
33672b0f2aSSebastian Grimberg     const CeedInt ind  = indices[node + elem * P_1d];
34672b0f2aSSebastian Grimberg 
35672b0f2aSSebastian Grimberg     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp];
369e201c85SYohann   }
379e201c85SYohann }
389e201c85SYohann 
399e201c85SYohann //------------------------------------------------------------------------------
409e201c85SYohann // L-vector -> E-vector, strided
419e201c85SYohann //------------------------------------------------------------------------------
42672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
43672b0f2aSSebastian Grimberg inline __device__ void readDofsStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
44672b0f2aSSebastian Grimberg   if (data.t_id_x < P_1d) {
459e201c85SYohann     const CeedInt node = data.t_id_x;
469e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
47672b0f2aSSebastian Grimberg 
48672b0f2aSSebastian Grimberg     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
499e201c85SYohann   }
509e201c85SYohann }
519e201c85SYohann 
529e201c85SYohann //------------------------------------------------------------------------------
539e201c85SYohann // E-vector -> L-vector, offsets provided
549e201c85SYohann //------------------------------------------------------------------------------
55672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d>
56672b0f2aSSebastian Grimberg inline __device__ void writeDofsOffset1d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
57672b0f2aSSebastian Grimberg                                          const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) {
58672b0f2aSSebastian Grimberg   if (data.t_id_x < P_1d) {
599e201c85SYohann     const CeedInt node = data.t_id_x;
60672b0f2aSSebastian Grimberg     const CeedInt ind  = indices[node + elem * P_1d];
61672b0f2aSSebastian Grimberg 
62672b0f2aSSebastian Grimberg     for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]);
639e201c85SYohann   }
649e201c85SYohann }
659e201c85SYohann 
669e201c85SYohann //------------------------------------------------------------------------------
679e201c85SYohann // E-vector -> L-vector, strided
689e201c85SYohann //------------------------------------------------------------------------------
69672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
70672b0f2aSSebastian Grimberg inline __device__ void writeDofsStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ r_v,
71672b0f2aSSebastian Grimberg                                           CeedScalar *__restrict__ d_v) {
72672b0f2aSSebastian Grimberg   if (data.t_id_x < P_1d) {
739e201c85SYohann     const CeedInt node = data.t_id_x;
749e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
75672b0f2aSSebastian Grimberg 
76672b0f2aSSebastian Grimberg     for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp];
779e201c85SYohann   }
789e201c85SYohann }
799e201c85SYohann 
809e201c85SYohann //------------------------------------------------------------------------------
819e201c85SYohann // 2D
829e201c85SYohann //------------------------------------------------------------------------------
839e201c85SYohann 
849e201c85SYohann //------------------------------------------------------------------------------
859e201c85SYohann // L-vector -> E-vector, offsets provided
869e201c85SYohann //------------------------------------------------------------------------------
87672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d>
88672b0f2aSSebastian Grimberg inline __device__ void readDofsOffset2d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
89672b0f2aSSebastian Grimberg                                         const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
90672b0f2aSSebastian Grimberg   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
91672b0f2aSSebastian Grimberg     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
92672b0f2aSSebastian Grimberg     const CeedInt ind  = indices[node + elem * P_1d * P_1d];
93672b0f2aSSebastian Grimberg 
94672b0f2aSSebastian Grimberg     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp];
959e201c85SYohann   }
969e201c85SYohann }
979e201c85SYohann 
989e201c85SYohann //------------------------------------------------------------------------------
999e201c85SYohann // L-vector -> E-vector, strided
1009e201c85SYohann //------------------------------------------------------------------------------
101672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
102672b0f2aSSebastian Grimberg inline __device__ void readDofsStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
103672b0f2aSSebastian Grimberg   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
104672b0f2aSSebastian Grimberg     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
1059e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
106672b0f2aSSebastian Grimberg 
107672b0f2aSSebastian Grimberg     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
1089e201c85SYohann   }
1099e201c85SYohann }
1109e201c85SYohann 
1119e201c85SYohann //------------------------------------------------------------------------------
1129e201c85SYohann // E-vector -> L-vector, offsets provided
1139e201c85SYohann //------------------------------------------------------------------------------
114672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d>
115672b0f2aSSebastian Grimberg inline __device__ void writeDofsOffset2d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
116672b0f2aSSebastian Grimberg                                          const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) {
117672b0f2aSSebastian Grimberg   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
118672b0f2aSSebastian Grimberg     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
119672b0f2aSSebastian Grimberg     const CeedInt ind  = indices[node + elem * P_1d * P_1d];
120672b0f2aSSebastian Grimberg 
121672b0f2aSSebastian Grimberg     for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]);
1229e201c85SYohann   }
1239e201c85SYohann }
1249e201c85SYohann 
1259e201c85SYohann //------------------------------------------------------------------------------
1269e201c85SYohann // E-vector -> L-vector, strided
1279e201c85SYohann //------------------------------------------------------------------------------
128672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
129672b0f2aSSebastian Grimberg inline __device__ void writeDofsStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ r_v,
130672b0f2aSSebastian Grimberg                                           CeedScalar *__restrict__ d_v) {
131672b0f2aSSebastian Grimberg   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
132672b0f2aSSebastian Grimberg     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
1339e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
134672b0f2aSSebastian Grimberg 
135672b0f2aSSebastian Grimberg     for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp];
1369e201c85SYohann   }
1379e201c85SYohann }
1389e201c85SYohann 
1399e201c85SYohann //------------------------------------------------------------------------------
1409e201c85SYohann // 3D
1419e201c85SYohann //------------------------------------------------------------------------------
1429e201c85SYohann 
1439e201c85SYohann //------------------------------------------------------------------------------
1449e201c85SYohann // L-vector -> E-vector, offsets provided
1459e201c85SYohann //------------------------------------------------------------------------------
146672b0f2aSSebastian Grimberg // TODO: remove "Dofs" and "Quads" in the following function names?
147672b0f2aSSebastian Grimberg //   - readDofsOffset3d -> readOffset3d ?
148672b0f2aSSebastian Grimberg //   - readDofsStrided3d -> readStrided3d ?
149672b0f2aSSebastian Grimberg //   - readSliceQuadsOffset3d -> readSliceOffset3d ?
150672b0f2aSSebastian Grimberg //   - readSliceQuadsStrided3d -> readSliceStrided3d ?
151672b0f2aSSebastian Grimberg //   - writeDofsOffset3d -> writeOffset3d ?
152672b0f2aSSebastian Grimberg //   - writeDofsStrided3d -> writeStrided3d ?
153672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d>
154672b0f2aSSebastian Grimberg inline __device__ void readDofsOffset3d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
155672b0f2aSSebastian Grimberg                                         const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
156672b0f2aSSebastian Grimberg   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
157672b0f2aSSebastian Grimberg     for (CeedInt z = 0; z < P_1d; z++) {
158672b0f2aSSebastian Grimberg       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
159672b0f2aSSebastian Grimberg       const CeedInt ind  = indices[node + elem * P_1d * P_1d * P_1d];
160672b0f2aSSebastian Grimberg 
161672b0f2aSSebastian Grimberg       for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1d] = d_u[ind + COMP_STRIDE * comp];
1629e201c85SYohann     }
1639e201c85SYohann }
1649e201c85SYohann 
1659e201c85SYohann //------------------------------------------------------------------------------
1669e201c85SYohann // L-vector -> E-vector, strided
1679e201c85SYohann //------------------------------------------------------------------------------
168672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
169672b0f2aSSebastian Grimberg inline __device__ void readDofsStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
170672b0f2aSSebastian Grimberg   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
171672b0f2aSSebastian Grimberg     for (CeedInt z = 0; z < P_1d; z++) {
172672b0f2aSSebastian Grimberg       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
1739e201c85SYohann       const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
174672b0f2aSSebastian Grimberg 
175672b0f2aSSebastian Grimberg       for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1d] = d_u[ind + comp * STRIDES_COMP];
1769e201c85SYohann     }
1779e201c85SYohann }
1789e201c85SYohann 
1799e201c85SYohann //------------------------------------------------------------------------------
1809e201c85SYohann // E-vector -> Q-vector, offests provided
1819e201c85SYohann //------------------------------------------------------------------------------
182672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int Q_1d>
183672b0f2aSSebastian Grimberg inline __device__ void readSliceQuadsOffset3d(SharedData_Hip &data, const CeedInt nquads, const CeedInt elem, const CeedInt q,
184672b0f2aSSebastian Grimberg                                               const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
185672b0f2aSSebastian Grimberg   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
186672b0f2aSSebastian Grimberg     const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d;
187672b0f2aSSebastian Grimberg     const CeedInt ind  = indices[node + elem * Q_1d * Q_1d * Q_1d];
188672b0f2aSSebastian Grimberg 
189672b0f2aSSebastian Grimberg     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp];
1909e201c85SYohann   }
1919e201c85SYohann }
1929e201c85SYohann 
1939e201c85SYohann //------------------------------------------------------------------------------
1949e201c85SYohann // E-vector -> Q-vector, strided
1959e201c85SYohann //------------------------------------------------------------------------------
196672b0f2aSSebastian Grimberg template <int NUM_COMP, int Q_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
197672b0f2aSSebastian Grimberg inline __device__ void readSliceQuadsStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedInt q, const CeedScalar *__restrict__ d_u,
198672b0f2aSSebastian Grimberg                                                CeedScalar *__restrict__ r_u) {
199672b0f2aSSebastian Grimberg   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
200672b0f2aSSebastian Grimberg     const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d;
2019e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
202672b0f2aSSebastian Grimberg 
203672b0f2aSSebastian Grimberg     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
2049e201c85SYohann   }
2059e201c85SYohann }
2069e201c85SYohann 
2079e201c85SYohann //------------------------------------------------------------------------------
2089e201c85SYohann // E-vector -> L-vector, offsets provided
2099e201c85SYohann //------------------------------------------------------------------------------
210672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d>
211672b0f2aSSebastian Grimberg inline __device__ void writeDofsOffset3d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
212672b0f2aSSebastian Grimberg                                          const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) {
213672b0f2aSSebastian Grimberg   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
214672b0f2aSSebastian Grimberg     for (CeedInt z = 0; z < P_1d; z++) {
215672b0f2aSSebastian Grimberg       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
216672b0f2aSSebastian Grimberg       const CeedInt ind  = indices[node + elem * P_1d * P_1d * P_1d];
217672b0f2aSSebastian Grimberg 
218672b0f2aSSebastian Grimberg       for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[z + comp * P_1d]);
2199e201c85SYohann     }
2209e201c85SYohann }
2219e201c85SYohann 
2229e201c85SYohann //------------------------------------------------------------------------------
2239e201c85SYohann // E-vector -> L-vector, strided
2249e201c85SYohann //------------------------------------------------------------------------------
225672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
226672b0f2aSSebastian Grimberg inline __device__ void writeDofsStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ r_v,
227672b0f2aSSebastian Grimberg                                           CeedScalar *__restrict__ d_v) {
228672b0f2aSSebastian Grimberg   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
229672b0f2aSSebastian Grimberg     for (CeedInt z = 0; z < P_1d; z++) {
230672b0f2aSSebastian Grimberg       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
2319e201c85SYohann       const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
232672b0f2aSSebastian Grimberg 
233672b0f2aSSebastian Grimberg       for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[z + comp * P_1d];
2349e201c85SYohann     }
2359e201c85SYohann }
2369e201c85SYohann 
2379e201c85SYohann //------------------------------------------------------------------------------
2389e201c85SYohann // 3D collocated derivatives computation
2399e201c85SYohann //------------------------------------------------------------------------------
240672b0f2aSSebastian Grimberg template <int NUM_COMP, int Q_1d>
2412b730f8bSJeremy L Thompson inline __device__ void gradCollo3d(SharedData_Hip &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
2422b730f8bSJeremy L Thompson                                    CeedScalar *__restrict__ r_V) {
243672b0f2aSSebastian Grimberg   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
244672b0f2aSSebastian Grimberg     for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
245672b0f2aSSebastian Grimberg       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[q + comp * Q_1d];
2469e201c85SYohann       __syncthreads();
2479e201c85SYohann       // X derivative
248672b0f2aSSebastian Grimberg       r_V[comp + 0 * NUM_COMP] = 0.0;
249672b0f2aSSebastian Grimberg       for (CeedInt i = 0; i < Q_1d; i++)
250672b0f2aSSebastian Grimberg         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)
2519e201c85SYohann       // Y derivative
252672b0f2aSSebastian Grimberg       r_V[comp + 1 * NUM_COMP] = 0.0;
253672b0f2aSSebastian Grimberg       for (CeedInt i = 0; i < Q_1d; i++)
254672b0f2aSSebastian Grimberg         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)
2559e201c85SYohann       // Z derivative
256672b0f2aSSebastian Grimberg       r_V[comp + 2 * NUM_COMP] = 0.0;
257672b0f2aSSebastian Grimberg       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)
2589e201c85SYohann       __syncthreads();
2599e201c85SYohann     }
2609e201c85SYohann   }
2619e201c85SYohann }
2629e201c85SYohann 
2639e201c85SYohann //------------------------------------------------------------------------------
2649e201c85SYohann // 3D collocated derivatives transpose
2659e201c85SYohann //------------------------------------------------------------------------------
266672b0f2aSSebastian Grimberg template <int NUM_COMP, int Q_1d>
2672b730f8bSJeremy L Thompson inline __device__ void gradColloTranspose3d(SharedData_Hip &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
2682b730f8bSJeremy L Thompson                                             CeedScalar *__restrict__ r_V) {
269672b0f2aSSebastian Grimberg   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
270672b0f2aSSebastian Grimberg     for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
2719e201c85SYohann       // X derivative
272672b0f2aSSebastian Grimberg       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 0 * NUM_COMP];
2739e201c85SYohann       __syncthreads();
274672b0f2aSSebastian Grimberg       for (CeedInt i = 0; i < Q_1d; i++)
275672b0f2aSSebastian Grimberg         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)
2769e201c85SYohann       __syncthreads();
2779e201c85SYohann       // Y derivative
278672b0f2aSSebastian Grimberg       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 1 * NUM_COMP];
2799e201c85SYohann       __syncthreads();
280672b0f2aSSebastian Grimberg       for (CeedInt i = 0; i < Q_1d; i++)
281672b0f2aSSebastian Grimberg         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)
2829e201c85SYohann       __syncthreads();
2839e201c85SYohann       // Z derivative
284672b0f2aSSebastian Grimberg       for (CeedInt i = 0; i < Q_1d; i++)
285672b0f2aSSebastian Grimberg         r_V[i + comp * Q_1d] += c_G[i + q * Q_1d] * r_U[comp + 2 * NUM_COMP];  // PARTIAL contract z direction (Z derivative)
2869e201c85SYohann     }
2879e201c85SYohann   }
2889e201c85SYohann }
289