xref: /libCEED/rust/libceed-sys/c-src/include/ceed/jit-source/cuda/cuda-gen-templates.h (revision 5aed82e4fa97acf4ba24a7f10a35f5303a6798e0)
1*5aed82e4SJeremy 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
1094b7b29bSJeremy L Thompson #ifndef CEED_CUDA_GEN_TEMPLATES_H
1194b7b29bSJeremy L Thompson #define CEED_CUDA_GEN_TEMPLATES_H
129e201c85SYohann 
131b7492f8SSebastian Grimberg #include <ceed.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 //------------------------------------------------------------------------------
30ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
31ca735530SJeremy L Thompson inline __device__ void readDofsOffset1d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
32672b0f2aSSebastian Grimberg                                         const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
33ca735530SJeremy L Thompson   if (data.t_id_x < P_1d) {
349e201c85SYohann     const CeedInt node = data.t_id_x;
35ca735530SJeremy L Thompson     const CeedInt ind  = indices[node + elem * P_1d];
36ca735530SJeremy L Thompson 
37ca735530SJeremy 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 //------------------------------------------------------------------------------
44ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
45672b0f2aSSebastian Grimberg inline __device__ void readDofsStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u,
46672b0f2aSSebastian Grimberg                                          CeedScalar *__restrict__ r_u) {
47ca735530SJeremy L Thompson   if (data.t_id_x < P_1d) {
489e201c85SYohann     const CeedInt node = data.t_id_x;
499e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
50ca735530SJeremy L Thompson 
51ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
529e201c85SYohann   }
539e201c85SYohann }
549e201c85SYohann 
559e201c85SYohann //------------------------------------------------------------------------------
569e201c85SYohann // E-vector -> L-vector, offsets provided
579e201c85SYohann //------------------------------------------------------------------------------
58ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
59ca735530SJeremy L Thompson inline __device__ void writeDofsOffset1d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
60672b0f2aSSebastian Grimberg                                          const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) {
61ca735530SJeremy L Thompson   if (data.t_id_x < P_1d) {
629e201c85SYohann     const CeedInt node = data.t_id_x;
63ca735530SJeremy L Thompson     const CeedInt ind  = indices[node + elem * P_1d];
64ca735530SJeremy L Thompson 
65ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]);
669e201c85SYohann   }
679e201c85SYohann }
689e201c85SYohann 
699e201c85SYohann //------------------------------------------------------------------------------
709e201c85SYohann // E-vector -> L-vector, strided
719e201c85SYohann //------------------------------------------------------------------------------
72ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
73672b0f2aSSebastian Grimberg inline __device__ void writeDofsStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ r_v,
74672b0f2aSSebastian Grimberg                                           CeedScalar *__restrict__ d_v) {
75ca735530SJeremy L Thompson   if (data.t_id_x < P_1d) {
769e201c85SYohann     const CeedInt node = data.t_id_x;
779e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
78ca735530SJeremy L Thompson 
79ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp];
809e201c85SYohann   }
819e201c85SYohann }
829e201c85SYohann 
839e201c85SYohann //------------------------------------------------------------------------------
849e201c85SYohann // 2D
859e201c85SYohann //------------------------------------------------------------------------------
869e201c85SYohann 
879e201c85SYohann //------------------------------------------------------------------------------
889e201c85SYohann // L-vector -> E-vector, offsets provided
899e201c85SYohann //------------------------------------------------------------------------------
90ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
91ca735530SJeremy L Thompson inline __device__ void readDofsOffset2d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
92672b0f2aSSebastian Grimberg                                         const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
93ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
94ca735530SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
95ca735530SJeremy L Thompson     const CeedInt ind  = indices[node + elem * P_1d * P_1d];
96ca735530SJeremy L Thompson 
97ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp];
989e201c85SYohann   }
999e201c85SYohann }
1009e201c85SYohann 
1019e201c85SYohann //------------------------------------------------------------------------------
1029e201c85SYohann // L-vector -> E-vector, strided
1039e201c85SYohann //------------------------------------------------------------------------------
104ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
105672b0f2aSSebastian Grimberg inline __device__ void readDofsStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u,
106672b0f2aSSebastian Grimberg                                          CeedScalar *__restrict__ r_u) {
107ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
108ca735530SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
1099e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
110ca735530SJeremy L Thompson 
111ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
1129e201c85SYohann   }
1139e201c85SYohann }
1149e201c85SYohann 
1159e201c85SYohann //------------------------------------------------------------------------------
1169e201c85SYohann // E-vector -> L-vector, offsets provided
1179e201c85SYohann //------------------------------------------------------------------------------
118ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
119ca735530SJeremy L Thompson inline __device__ void writeDofsOffset2d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
120672b0f2aSSebastian Grimberg                                          const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) {
121ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
122ca735530SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
123ca735530SJeremy L Thompson     const CeedInt ind  = indices[node + elem * P_1d * P_1d];
124ca735530SJeremy L Thompson 
125ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]);
1269e201c85SYohann   }
1279e201c85SYohann }
1289e201c85SYohann 
1299e201c85SYohann //------------------------------------------------------------------------------
1309e201c85SYohann // E-vector -> L-vector, strided
1319e201c85SYohann //------------------------------------------------------------------------------
132ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
133672b0f2aSSebastian Grimberg inline __device__ void writeDofsStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ r_v,
134672b0f2aSSebastian Grimberg                                           CeedScalar *__restrict__ d_v) {
135ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
136ca735530SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
1379e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
138ca735530SJeremy L Thompson 
139ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp];
1409e201c85SYohann   }
1419e201c85SYohann }
1429e201c85SYohann 
1439e201c85SYohann //------------------------------------------------------------------------------
1449e201c85SYohann // 3D
1459e201c85SYohann //------------------------------------------------------------------------------
1469e201c85SYohann 
1479e201c85SYohann //------------------------------------------------------------------------------
1489e201c85SYohann // L-vector -> E-vector, offsets provided
1499e201c85SYohann //------------------------------------------------------------------------------
1509e201c85SYohann // TODO: remove "Dofs" and "Quads" in the following function names?
1519e201c85SYohann //   - readDofsOffset3d -> readOffset3d ?
1529e201c85SYohann //   - readDofsStrided3d -> readStrided3d ?
1539e201c85SYohann //   - readSliceQuadsOffset3d -> readSliceOffset3d ?
1549e201c85SYohann //   - readSliceQuadsStrided3d -> readSliceStrided3d ?
1559e201c85SYohann //   - writeDofsOffset3d -> writeOffset3d ?
1569e201c85SYohann //   - writeDofsStrided3d -> writeStrided3d ?
157ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
158ca735530SJeremy L Thompson inline __device__ void readDofsOffset3d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
159672b0f2aSSebastian Grimberg                                         const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
160ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
161ca735530SJeremy L Thompson     for (CeedInt z = 0; z < P_1d; z++) {
162ca735530SJeremy L Thompson       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
163ca735530SJeremy L Thompson       const CeedInt ind  = indices[node + elem * P_1d * P_1d * P_1d];
164ca735530SJeremy L Thompson 
165ca735530SJeremy L Thompson       for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1d] = d_u[ind + COMP_STRIDE * comp];
1669e201c85SYohann     }
1679e201c85SYohann }
1689e201c85SYohann 
1699e201c85SYohann //------------------------------------------------------------------------------
1709e201c85SYohann // L-vector -> E-vector, strided
1719e201c85SYohann //------------------------------------------------------------------------------
172ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
173672b0f2aSSebastian Grimberg inline __device__ void readDofsStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u,
174672b0f2aSSebastian Grimberg                                          CeedScalar *__restrict__ r_u) {
175ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
176ca735530SJeremy L Thompson     for (CeedInt z = 0; z < P_1d; z++) {
177ca735530SJeremy L Thompson       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
1789e201c85SYohann       const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
179ca735530SJeremy L Thompson 
180ca735530SJeremy L Thompson       for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1d] = d_u[ind + comp * STRIDES_COMP];
1819e201c85SYohann     }
1829e201c85SYohann }
1839e201c85SYohann 
1849e201c85SYohann //------------------------------------------------------------------------------
1859e201c85SYohann // E-vector -> Q-vector, offests provided
1869e201c85SYohann //------------------------------------------------------------------------------
187ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int Q_1d>
1882b730f8bSJeremy L Thompson inline __device__ void readSliceQuadsOffset3d(SharedData_Cuda &data, const CeedInt nquads, const CeedInt elem, const CeedInt q,
189672b0f2aSSebastian Grimberg                                               const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
190ca735530SJeremy L Thompson   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
191ca735530SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d;
192ca735530SJeremy L Thompson     const CeedInt ind  = indices[node + elem * Q_1d * Q_1d * Q_1d];
193ca735530SJeremy L Thompson 
194ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp];
1959e201c85SYohann   }
1969e201c85SYohann }
1979e201c85SYohann 
1989e201c85SYohann //------------------------------------------------------------------------------
1999e201c85SYohann // E-vector -> Q-vector, strided
2009e201c85SYohann //------------------------------------------------------------------------------
201ca735530SJeremy L Thompson template <int NUM_COMP, int Q_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
2022b730f8bSJeremy L Thompson inline __device__ void readSliceQuadsStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt q, const CeedScalar *__restrict__ d_u,
203672b0f2aSSebastian Grimberg                                                CeedScalar *__restrict__ r_u) {
204ca735530SJeremy L Thompson   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
205ca735530SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d;
2069e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
207ca735530SJeremy L Thompson 
208ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
2099e201c85SYohann   }
2109e201c85SYohann }
2119e201c85SYohann 
2129e201c85SYohann //------------------------------------------------------------------------------
2139e201c85SYohann // E-vector -> L-vector, offsets provided
2149e201c85SYohann //------------------------------------------------------------------------------
215ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
216ca735530SJeremy L Thompson inline __device__ void writeDofsOffset3d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
217672b0f2aSSebastian Grimberg                                          const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) {
218ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
219ca735530SJeremy L Thompson     for (CeedInt z = 0; z < P_1d; z++) {
220ca735530SJeremy L Thompson       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
221ca735530SJeremy L Thompson       const CeedInt ind  = indices[node + elem * P_1d * P_1d * P_1d];
222ca735530SJeremy L Thompson 
223ca735530SJeremy L Thompson       for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[z + comp * P_1d]);
2249e201c85SYohann     }
2259e201c85SYohann }
2269e201c85SYohann 
2279e201c85SYohann //------------------------------------------------------------------------------
2289e201c85SYohann // E-vector -> L-vector, strided
2299e201c85SYohann //------------------------------------------------------------------------------
230ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
231672b0f2aSSebastian Grimberg inline __device__ void writeDofsStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ r_v,
232672b0f2aSSebastian Grimberg                                           CeedScalar *__restrict__ d_v) {
233ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
234ca735530SJeremy L Thompson     for (CeedInt z = 0; z < P_1d; z++) {
235ca735530SJeremy L Thompson       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
2369e201c85SYohann       const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
237ca735530SJeremy L Thompson 
238ca735530SJeremy L Thompson       for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[z + comp * P_1d];
2399e201c85SYohann     }
2409e201c85SYohann }
2419e201c85SYohann 
2429e201c85SYohann //------------------------------------------------------------------------------
2439e201c85SYohann // 3D collocated derivatives computation
2449e201c85SYohann //------------------------------------------------------------------------------
245ca735530SJeremy L Thompson template <int NUM_COMP, int Q_1d>
2462b730f8bSJeremy L Thompson inline __device__ void gradCollo3d(SharedData_Cuda &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
2472b730f8bSJeremy L Thompson                                    CeedScalar *__restrict__ r_V) {
248ca735530SJeremy L Thompson   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
249ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
250ca735530SJeremy L Thompson       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[q + comp * Q_1d];
2519e201c85SYohann       __syncthreads();
2529e201c85SYohann       // X derivative
253ca735530SJeremy L Thompson       r_V[comp + 0 * NUM_COMP] = 0.0;
254ca735530SJeremy L Thompson       for (CeedInt i = 0; i < Q_1d; i++)
255ca735530SJeremy 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)
2569e201c85SYohann       // Y derivative
257ca735530SJeremy L Thompson       r_V[comp + 1 * NUM_COMP] = 0.0;
258ca735530SJeremy L Thompson       for (CeedInt i = 0; i < Q_1d; i++)
259ca735530SJeremy 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)
2609e201c85SYohann       // Z derivative
261ca735530SJeremy L Thompson       r_V[comp + 2 * NUM_COMP] = 0.0;
262ca735530SJeremy 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)
2639e201c85SYohann       __syncthreads();
2649e201c85SYohann     }
2659e201c85SYohann   }
2669e201c85SYohann }
2679e201c85SYohann 
2689e201c85SYohann //------------------------------------------------------------------------------
2699e201c85SYohann // 3D collocated derivatives transpose
2709e201c85SYohann //------------------------------------------------------------------------------
271ca735530SJeremy L Thompson template <int NUM_COMP, int Q_1d>
2722b730f8bSJeremy L Thompson inline __device__ void gradColloTranspose3d(SharedData_Cuda &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
2732b730f8bSJeremy L Thompson                                             CeedScalar *__restrict__ r_V) {
274ca735530SJeremy L Thompson   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
275ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
2769e201c85SYohann       // X derivative
277ca735530SJeremy L Thompson       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 0 * NUM_COMP];
2789e201c85SYohann       __syncthreads();
279ca735530SJeremy L Thompson       for (CeedInt i = 0; i < Q_1d; i++)
280ca735530SJeremy 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)
2819e201c85SYohann       __syncthreads();
2829e201c85SYohann       // Y derivative
283ca735530SJeremy L Thompson       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 1 * NUM_COMP];
2849e201c85SYohann       __syncthreads();
285ca735530SJeremy L Thompson       for (CeedInt i = 0; i < Q_1d; i++)
286ca735530SJeremy 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)
2879e201c85SYohann       __syncthreads();
2889e201c85SYohann       // Z derivative
289ca735530SJeremy L Thompson       for (CeedInt i = 0; i < Q_1d; i++)
290ca735530SJeremy 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)
2919e201c85SYohann     }
2929e201c85SYohann   }
2939e201c85SYohann }
2949e201c85SYohann 
29594b7b29bSJeremy L Thompson #endif  // CEED_CUDA_GEN_TEMPLATES_H
296