xref: /libCEED/rust/libceed-sys/c-src/include/ceed/jit-source/cuda/cuda-gen-templates.h (revision ca7355308a14805110a7d98dabf1f658d5ec16d5)
19e201c85SYohann // Copyright (c) 2017-2022, 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
1094b7b29bSJeremy L Thompson #ifndef CEED_CUDA_GEN_TEMPLATES_H
1194b7b29bSJeremy L Thompson #define CEED_CUDA_GEN_TEMPLATES_H
129e201c85SYohann 
139e201c85SYohann #include <ceed/types.h>
149e201c85SYohann 
159e201c85SYohann //------------------------------------------------------------------------------
169e201c85SYohann // Load matrices for basis actions
179e201c85SYohann //------------------------------------------------------------------------------
189e201c85SYohann template <int P, int Q>
199e201c85SYohann inline __device__ void loadMatrix(SharedData_Cuda &data, const CeedScalar *__restrict__ d_B, CeedScalar *B) {
202b730f8bSJeremy L Thompson   for (CeedInt i = data.t_id; i < P * Q; i += blockDim.x * blockDim.y * blockDim.z) B[i] = d_B[i];
219e201c85SYohann }
229e201c85SYohann 
239e201c85SYohann //------------------------------------------------------------------------------
249e201c85SYohann // 1D
259e201c85SYohann //------------------------------------------------------------------------------
269e201c85SYohann 
279e201c85SYohann //------------------------------------------------------------------------------
289e201c85SYohann // L-vector -> E-vector, offsets provided
299e201c85SYohann //------------------------------------------------------------------------------
30*ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
31*ca735530SJeremy L Thompson inline __device__ void readDofsOffset1d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
322b730f8bSJeremy L Thompson                                         const CeedScalar *__restrict__ d_u, CeedScalar *r_u) {
33*ca735530SJeremy L Thompson   if (data.t_id_x < P_1d) {
349e201c85SYohann     const CeedInt node = data.t_id_x;
35*ca735530SJeremy L Thompson     const CeedInt ind  = indices[node + elem * P_1d];
36*ca735530SJeremy L Thompson 
37*ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp];
389e201c85SYohann   }
399e201c85SYohann }
409e201c85SYohann 
419e201c85SYohann //------------------------------------------------------------------------------
429e201c85SYohann // L-vector -> E-vector, strided
439e201c85SYohann //------------------------------------------------------------------------------
44*ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
459e201c85SYohann inline __device__ void readDofsStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) {
46*ca735530SJeremy L Thompson   if (data.t_id_x < P_1d) {
479e201c85SYohann     const CeedInt node = data.t_id_x;
489e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
49*ca735530SJeremy L Thompson 
50*ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
519e201c85SYohann   }
529e201c85SYohann }
539e201c85SYohann 
549e201c85SYohann //------------------------------------------------------------------------------
559e201c85SYohann // E-vector -> L-vector, offsets provided
569e201c85SYohann //------------------------------------------------------------------------------
57*ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
58*ca735530SJeremy L Thompson inline __device__ void writeDofsOffset1d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
592b730f8bSJeremy L Thompson                                          const CeedScalar *r_v, CeedScalar *d_v) {
60*ca735530SJeremy L Thompson   if (data.t_id_x < P_1d) {
619e201c85SYohann     const CeedInt node = data.t_id_x;
62*ca735530SJeremy L Thompson     const CeedInt ind  = indices[node + elem * P_1d];
63*ca735530SJeremy L Thompson 
64*ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]);
659e201c85SYohann   }
669e201c85SYohann }
679e201c85SYohann 
689e201c85SYohann //------------------------------------------------------------------------------
699e201c85SYohann // E-vector -> L-vector, strided
709e201c85SYohann //------------------------------------------------------------------------------
71*ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
729e201c85SYohann inline __device__ void writeDofsStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *r_v, CeedScalar *d_v) {
73*ca735530SJeremy L Thompson   if (data.t_id_x < P_1d) {
749e201c85SYohann     const CeedInt node = data.t_id_x;
759e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
76*ca735530SJeremy L Thompson 
77*ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp];
789e201c85SYohann   }
799e201c85SYohann }
809e201c85SYohann 
819e201c85SYohann //------------------------------------------------------------------------------
829e201c85SYohann // 2D
839e201c85SYohann //------------------------------------------------------------------------------
849e201c85SYohann 
859e201c85SYohann //------------------------------------------------------------------------------
869e201c85SYohann // L-vector -> E-vector, offsets provided
879e201c85SYohann //------------------------------------------------------------------------------
88*ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
89*ca735530SJeremy L Thompson inline __device__ void readDofsOffset2d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
902b730f8bSJeremy L Thompson                                         const CeedScalar *__restrict__ d_u, CeedScalar *r_u) {
91*ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
92*ca735530SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
93*ca735530SJeremy L Thompson     const CeedInt ind  = indices[node + elem * P_1d * P_1d];
94*ca735530SJeremy L Thompson 
95*ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp];
969e201c85SYohann   }
979e201c85SYohann }
989e201c85SYohann 
999e201c85SYohann //------------------------------------------------------------------------------
1009e201c85SYohann // L-vector -> E-vector, strided
1019e201c85SYohann //------------------------------------------------------------------------------
102*ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
1039e201c85SYohann inline __device__ void readDofsStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) {
104*ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
105*ca735530SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
1069e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
107*ca735530SJeremy L Thompson 
108*ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
1099e201c85SYohann   }
1109e201c85SYohann }
1119e201c85SYohann 
1129e201c85SYohann //------------------------------------------------------------------------------
1139e201c85SYohann // E-vector -> L-vector, offsets provided
1149e201c85SYohann //------------------------------------------------------------------------------
115*ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
116*ca735530SJeremy L Thompson inline __device__ void writeDofsOffset2d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
1172b730f8bSJeremy L Thompson                                          const CeedScalar *r_v, CeedScalar *d_v) {
118*ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
119*ca735530SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
120*ca735530SJeremy L Thompson     const CeedInt ind  = indices[node + elem * P_1d * P_1d];
121*ca735530SJeremy L Thompson 
122*ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]);
1239e201c85SYohann   }
1249e201c85SYohann }
1259e201c85SYohann 
1269e201c85SYohann //------------------------------------------------------------------------------
1279e201c85SYohann // E-vector -> L-vector, strided
1289e201c85SYohann //------------------------------------------------------------------------------
129*ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
1309e201c85SYohann inline __device__ void writeDofsStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *r_v, CeedScalar *d_v) {
131*ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
132*ca735530SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
1339e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
134*ca735530SJeremy L Thompson 
135*ca735530SJeremy L Thompson     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 //------------------------------------------------------------------------------
1469e201c85SYohann // TODO: remove "Dofs" and "Quads" in the following function names?
1479e201c85SYohann //   - readDofsOffset3d -> readOffset3d ?
1489e201c85SYohann //   - readDofsStrided3d -> readStrided3d ?
1499e201c85SYohann //   - readSliceQuadsOffset3d -> readSliceOffset3d ?
1509e201c85SYohann //   - readSliceQuadsStrided3d -> readSliceStrided3d ?
1519e201c85SYohann //   - writeDofsOffset3d -> writeOffset3d ?
1529e201c85SYohann //   - writeDofsStrided3d -> writeStrided3d ?
153*ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
154*ca735530SJeremy L Thompson inline __device__ void readDofsOffset3d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
1552b730f8bSJeremy L Thompson                                         const CeedScalar *__restrict__ d_u, CeedScalar *r_u) {
156*ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
157*ca735530SJeremy L Thompson     for (CeedInt z = 0; z < P_1d; z++) {
158*ca735530SJeremy L Thompson       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
159*ca735530SJeremy L Thompson       const CeedInt ind  = indices[node + elem * P_1d * P_1d * P_1d];
160*ca735530SJeremy L Thompson 
161*ca735530SJeremy L Thompson       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 //------------------------------------------------------------------------------
168*ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
1699e201c85SYohann inline __device__ void readDofsStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) {
170*ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
171*ca735530SJeremy L Thompson     for (CeedInt z = 0; z < P_1d; z++) {
172*ca735530SJeremy L Thompson       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;
174*ca735530SJeremy L Thompson 
175*ca735530SJeremy L Thompson       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 //------------------------------------------------------------------------------
182*ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int Q_1d>
1832b730f8bSJeremy L Thompson inline __device__ void readSliceQuadsOffset3d(SharedData_Cuda &data, const CeedInt nquads, const CeedInt elem, const CeedInt q,
1842b730f8bSJeremy L Thompson                                               const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) {
185*ca735530SJeremy L Thompson   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
186*ca735530SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d;
187*ca735530SJeremy L Thompson     const CeedInt ind  = indices[node + elem * Q_1d * Q_1d * Q_1d];
188*ca735530SJeremy L Thompson 
189*ca735530SJeremy L Thompson     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 //------------------------------------------------------------------------------
196*ca735530SJeremy L Thompson template <int NUM_COMP, int Q_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
1972b730f8bSJeremy L Thompson inline __device__ void readSliceQuadsStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt q, const CeedScalar *__restrict__ d_u,
1982b730f8bSJeremy L Thompson                                                CeedScalar *r_u) {
199*ca735530SJeremy L Thompson   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
200*ca735530SJeremy L Thompson     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;
202*ca735530SJeremy L Thompson 
203*ca735530SJeremy L Thompson     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 //------------------------------------------------------------------------------
210*ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
211*ca735530SJeremy L Thompson inline __device__ void writeDofsOffset3d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
2122b730f8bSJeremy L Thompson                                          const CeedScalar *r_v, CeedScalar *d_v) {
213*ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
214*ca735530SJeremy L Thompson     for (CeedInt z = 0; z < P_1d; z++) {
215*ca735530SJeremy L Thompson       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
216*ca735530SJeremy L Thompson       const CeedInt ind  = indices[node + elem * P_1d * P_1d * P_1d];
217*ca735530SJeremy L Thompson 
218*ca735530SJeremy L Thompson       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 //------------------------------------------------------------------------------
225*ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
2269e201c85SYohann inline __device__ void writeDofsStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *r_v, CeedScalar *d_v) {
227*ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
228*ca735530SJeremy L Thompson     for (CeedInt z = 0; z < P_1d; z++) {
229*ca735530SJeremy L Thompson       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
2309e201c85SYohann       const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
231*ca735530SJeremy L Thompson 
232*ca735530SJeremy L Thompson       for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[z + comp * P_1d];
2339e201c85SYohann     }
2349e201c85SYohann }
2359e201c85SYohann 
2369e201c85SYohann //------------------------------------------------------------------------------
2379e201c85SYohann // 3D collocated derivatives computation
2389e201c85SYohann //------------------------------------------------------------------------------
239*ca735530SJeremy L Thompson template <int NUM_COMP, int Q_1d>
2402b730f8bSJeremy L Thompson inline __device__ void gradCollo3d(SharedData_Cuda &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
2412b730f8bSJeremy L Thompson                                    CeedScalar *__restrict__ r_V) {
242*ca735530SJeremy L Thompson   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
243*ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
244*ca735530SJeremy L Thompson       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[q + comp * Q_1d];
2459e201c85SYohann       __syncthreads();
2469e201c85SYohann       // X derivative
247*ca735530SJeremy L Thompson       r_V[comp + 0 * NUM_COMP] = 0.0;
248*ca735530SJeremy L Thompson       for (CeedInt i = 0; i < Q_1d; i++)
249*ca735530SJeremy 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];  // Contract x direction (X derivative)
2509e201c85SYohann       // Y derivative
251*ca735530SJeremy L Thompson       r_V[comp + 1 * NUM_COMP] = 0.0;
252*ca735530SJeremy L Thompson       for (CeedInt i = 0; i < Q_1d; i++)
253*ca735530SJeremy 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];  // Contract y direction (Y derivative)
2549e201c85SYohann       // Z derivative
255*ca735530SJeremy L Thompson       r_V[comp + 2 * NUM_COMP] = 0.0;
256*ca735530SJeremy L Thompson       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)
2579e201c85SYohann       __syncthreads();
2589e201c85SYohann     }
2599e201c85SYohann   }
2609e201c85SYohann }
2619e201c85SYohann 
2629e201c85SYohann //------------------------------------------------------------------------------
2639e201c85SYohann // 3D collocated derivatives transpose
2649e201c85SYohann //------------------------------------------------------------------------------
265*ca735530SJeremy L Thompson template <int NUM_COMP, int Q_1d>
2662b730f8bSJeremy L Thompson inline __device__ void gradColloTranspose3d(SharedData_Cuda &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
2672b730f8bSJeremy L Thompson                                             CeedScalar *__restrict__ r_V) {
268*ca735530SJeremy L Thompson   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
269*ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
2709e201c85SYohann       // X derivative
271*ca735530SJeremy L Thompson       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 0 * NUM_COMP];
2729e201c85SYohann       __syncthreads();
273*ca735530SJeremy L Thompson       for (CeedInt i = 0; i < Q_1d; i++)
274*ca735530SJeremy 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];  // Contract x direction (X derivative)
2759e201c85SYohann       __syncthreads();
2769e201c85SYohann       // Y derivative
277*ca735530SJeremy L Thompson       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 1 * NUM_COMP];
2789e201c85SYohann       __syncthreads();
279*ca735530SJeremy L Thompson       for (CeedInt i = 0; i < Q_1d; i++)
280*ca735530SJeremy 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];  // Contract y direction (Y derivative)
2819e201c85SYohann       __syncthreads();
2829e201c85SYohann       // Z derivative
283*ca735530SJeremy L Thompson       for (CeedInt i = 0; i < Q_1d; i++)
284*ca735530SJeremy L Thompson         r_V[i + comp * Q_1d] += c_G[i + q * Q_1d] * r_U[comp + 2 * NUM_COMP];  // PARTIAL contract z direction (Z derivative)
2859e201c85SYohann     }
2869e201c85SYohann   }
2879e201c85SYohann }
2889e201c85SYohann 
28994b7b29bSJeremy L Thompson #endif  // CEED_CUDA_GEN_TEMPLATES_H
290