xref: /libCEED/rust/libceed-sys/c-src/include/ceed/jit-source/cuda/cuda-gen-templates.h (revision f815fac990b20019e227e6950cc74a96439f9eba)
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>
16*f815fac9SJeremy 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 //------------------------------------------------------------------------------
219e201c85SYohann // 1D
229e201c85SYohann //------------------------------------------------------------------------------
239e201c85SYohann 
249e201c85SYohann //------------------------------------------------------------------------------
259e201c85SYohann // L-vector -> E-vector, offsets provided
269e201c85SYohann //------------------------------------------------------------------------------
27ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
28*f815fac9SJeremy L Thompson inline __device__ void ReadLVecStandard1d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
29672b0f2aSSebastian Grimberg                                           const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
30ca735530SJeremy L Thompson   if (data.t_id_x < P_1d) {
319e201c85SYohann     const CeedInt node = data.t_id_x;
32ca735530SJeremy L Thompson     const CeedInt ind  = indices[node + elem * P_1d];
33ca735530SJeremy L Thompson 
34ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp];
359e201c85SYohann   }
369e201c85SYohann }
379e201c85SYohann 
389e201c85SYohann //------------------------------------------------------------------------------
399e201c85SYohann // L-vector -> E-vector, strided
409e201c85SYohann //------------------------------------------------------------------------------
41ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
42*f815fac9SJeremy L Thompson inline __device__ void ReadLVecStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u,
43672b0f2aSSebastian Grimberg                                          CeedScalar *__restrict__ r_u) {
44ca735530SJeremy L Thompson   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;
47ca735530SJeremy L Thompson 
48ca735530SJeremy L Thompson     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 //------------------------------------------------------------------------------
55ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
56*f815fac9SJeremy L Thompson inline __device__ void WriteLVecStandard1d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
57672b0f2aSSebastian Grimberg                                            const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) {
58ca735530SJeremy L Thompson   if (data.t_id_x < P_1d) {
599e201c85SYohann     const CeedInt node = data.t_id_x;
60ca735530SJeremy L Thompson     const CeedInt ind  = indices[node + elem * P_1d];
61ca735530SJeremy L Thompson 
62ca735530SJeremy L Thompson     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 //------------------------------------------------------------------------------
69ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
70*f815fac9SJeremy L Thompson inline __device__ void WriteLVecStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ r_v,
71672b0f2aSSebastian Grimberg                                           CeedScalar *__restrict__ d_v) {
72ca735530SJeremy L Thompson   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;
75ca735530SJeremy L Thompson 
76ca735530SJeremy L Thompson     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 //------------------------------------------------------------------------------
87ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
88*f815fac9SJeremy L Thompson inline __device__ void ReadLVecStandard2d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
89672b0f2aSSebastian Grimberg                                           const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
90ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
91ca735530SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
92ca735530SJeremy L Thompson     const CeedInt ind  = indices[node + elem * P_1d * P_1d];
93ca735530SJeremy L Thompson 
94ca735530SJeremy L Thompson     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 //------------------------------------------------------------------------------
101ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
102*f815fac9SJeremy L Thompson inline __device__ void ReadLVecStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u,
103672b0f2aSSebastian Grimberg                                          CeedScalar *__restrict__ r_u) {
104ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
105ca735530SJeremy 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;
107ca735530SJeremy L Thompson 
108ca735530SJeremy 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 //------------------------------------------------------------------------------
115ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
116*f815fac9SJeremy L Thompson inline __device__ void WriteLVecStandard2d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
117672b0f2aSSebastian Grimberg                                            const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) {
118ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
119ca735530SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
120ca735530SJeremy L Thompson     const CeedInt ind  = indices[node + elem * P_1d * P_1d];
121ca735530SJeremy L Thompson 
122ca735530SJeremy 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 //------------------------------------------------------------------------------
129ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
130*f815fac9SJeremy L Thompson inline __device__ void WriteLVecStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ r_v,
131672b0f2aSSebastian Grimberg                                           CeedScalar *__restrict__ d_v) {
132ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d) {
133ca735530SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * P_1d;
1349e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
135ca735530SJeremy L Thompson 
136ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp];
1379e201c85SYohann   }
1389e201c85SYohann }
1399e201c85SYohann 
1409e201c85SYohann //------------------------------------------------------------------------------
1419e201c85SYohann // 3D
1429e201c85SYohann //------------------------------------------------------------------------------
1439e201c85SYohann 
1449e201c85SYohann //------------------------------------------------------------------------------
1459e201c85SYohann // L-vector -> E-vector, offsets provided
1469e201c85SYohann //------------------------------------------------------------------------------
147ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
148*f815fac9SJeremy L Thompson inline __device__ void ReadLVecStandard3d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
149672b0f2aSSebastian Grimberg                                           const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) {
150ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
151ca735530SJeremy L Thompson     for (CeedInt z = 0; z < P_1d; z++) {
152ca735530SJeremy L Thompson       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
153ca735530SJeremy L Thompson       const CeedInt ind  = indices[node + elem * P_1d * P_1d * P_1d];
154ca735530SJeremy L Thompson 
155ca735530SJeremy L Thompson       for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1d] = d_u[ind + COMP_STRIDE * comp];
1569e201c85SYohann     }
1579e201c85SYohann }
1589e201c85SYohann 
1599e201c85SYohann //------------------------------------------------------------------------------
1609e201c85SYohann // L-vector -> E-vector, strided
1619e201c85SYohann //------------------------------------------------------------------------------
162ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
163*f815fac9SJeremy L Thompson inline __device__ void ReadLVecStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u,
164672b0f2aSSebastian Grimberg                                          CeedScalar *__restrict__ r_u) {
165ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
166ca735530SJeremy L Thompson     for (CeedInt z = 0; z < P_1d; z++) {
167ca735530SJeremy L Thompson       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
1689e201c85SYohann       const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
169ca735530SJeremy L Thompson 
170ca735530SJeremy L Thompson       for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1d] = d_u[ind + comp * STRIDES_COMP];
1719e201c85SYohann     }
1729e201c85SYohann }
1739e201c85SYohann 
1749e201c85SYohann //------------------------------------------------------------------------------
1759e201c85SYohann // E-vector -> Q-vector, offests provided
1769e201c85SYohann //------------------------------------------------------------------------------
177ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int Q_1d>
178*f815fac9SJeremy L Thompson inline __device__ void ReadEVecSliceStandard3d(SharedData_Cuda &data, const CeedInt nquads, const CeedInt elem, const CeedInt q,
179*f815fac9SJeremy L Thompson                                                const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ d_u,
180*f815fac9SJeremy L Thompson                                                CeedScalar *__restrict__ r_u) {
181ca735530SJeremy L Thompson   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
182ca735530SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d;
183ca735530SJeremy L Thompson     const CeedInt ind  = indices[node + elem * Q_1d * Q_1d * Q_1d];
184ca735530SJeremy L Thompson 
185ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp];
1869e201c85SYohann   }
1879e201c85SYohann }
1889e201c85SYohann 
1899e201c85SYohann //------------------------------------------------------------------------------
1909e201c85SYohann // E-vector -> Q-vector, strided
1919e201c85SYohann //------------------------------------------------------------------------------
192ca735530SJeremy L Thompson template <int NUM_COMP, int Q_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
193*f815fac9SJeremy L Thompson inline __device__ void ReadEVecSliceStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt q, const CeedScalar *__restrict__ d_u,
194672b0f2aSSebastian Grimberg                                               CeedScalar *__restrict__ r_u) {
195ca735530SJeremy L Thompson   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
196ca735530SJeremy L Thompson     const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d;
1979e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
198ca735530SJeremy L Thompson 
199ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
2009e201c85SYohann   }
2019e201c85SYohann }
2029e201c85SYohann 
2039e201c85SYohann //------------------------------------------------------------------------------
2049e201c85SYohann // E-vector -> L-vector, offsets provided
2059e201c85SYohann //------------------------------------------------------------------------------
206ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d>
207*f815fac9SJeremy L Thompson inline __device__ void WriteLVecStandard3d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices,
208672b0f2aSSebastian Grimberg                                            const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) {
209ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
210ca735530SJeremy L Thompson     for (CeedInt z = 0; z < P_1d; z++) {
211ca735530SJeremy L Thompson       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
212ca735530SJeremy L Thompson       const CeedInt ind  = indices[node + elem * P_1d * P_1d * P_1d];
213ca735530SJeremy L Thompson 
214ca735530SJeremy L Thompson       for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[z + comp * P_1d]);
2159e201c85SYohann     }
2169e201c85SYohann }
2179e201c85SYohann 
2189e201c85SYohann //------------------------------------------------------------------------------
2199e201c85SYohann // E-vector -> L-vector, strided
2209e201c85SYohann //------------------------------------------------------------------------------
221ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
222*f815fac9SJeremy L Thompson inline __device__ void WriteLVecStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ r_v,
223672b0f2aSSebastian Grimberg                                           CeedScalar *__restrict__ d_v) {
224ca735530SJeremy L Thompson   if (data.t_id_x < P_1d && data.t_id_y < P_1d)
225ca735530SJeremy L Thompson     for (CeedInt z = 0; z < P_1d; z++) {
226ca735530SJeremy L Thompson       const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d;
2279e201c85SYohann       const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
228ca735530SJeremy L Thompson 
229ca735530SJeremy L Thompson       for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[z + comp * P_1d];
2309e201c85SYohann     }
2319e201c85SYohann }
2329e201c85SYohann 
2339e201c85SYohann //------------------------------------------------------------------------------
2349e201c85SYohann // 3D collocated derivatives computation
2359e201c85SYohann //------------------------------------------------------------------------------
236ca735530SJeremy L Thompson template <int NUM_COMP, int Q_1d>
237*f815fac9SJeremy L Thompson inline __device__ void GradColloSlice3d(SharedData_Cuda &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
2382b730f8bSJeremy L Thompson                                         CeedScalar *__restrict__ r_V) {
239ca735530SJeremy L Thompson   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
240ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
241ca735530SJeremy L Thompson       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[q + comp * Q_1d];
2429e201c85SYohann       __syncthreads();
2439e201c85SYohann       // X derivative
244ca735530SJeremy L Thompson       r_V[comp + 0 * NUM_COMP] = 0.0;
245ca735530SJeremy L Thompson       for (CeedInt i = 0; i < Q_1d; i++)
246ca735530SJeremy 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)
2479e201c85SYohann       // Y derivative
248ca735530SJeremy L Thompson       r_V[comp + 1 * NUM_COMP] = 0.0;
249ca735530SJeremy L Thompson       for (CeedInt i = 0; i < Q_1d; i++)
250ca735530SJeremy 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)
2519e201c85SYohann       // Z derivative
252ca735530SJeremy L Thompson       r_V[comp + 2 * NUM_COMP] = 0.0;
253ca735530SJeremy 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)
2549e201c85SYohann       __syncthreads();
2559e201c85SYohann     }
2569e201c85SYohann   }
2579e201c85SYohann }
2589e201c85SYohann 
2599e201c85SYohann //------------------------------------------------------------------------------
2609e201c85SYohann // 3D collocated derivatives transpose
2619e201c85SYohann //------------------------------------------------------------------------------
262ca735530SJeremy L Thompson template <int NUM_COMP, int Q_1d>
263*f815fac9SJeremy L Thompson inline __device__ void GradColloSliceTranspose3d(SharedData_Cuda &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
2642b730f8bSJeremy L Thompson                                                  CeedScalar *__restrict__ r_V) {
265ca735530SJeremy L Thompson   if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) {
266ca735530SJeremy L Thompson     for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
2679e201c85SYohann       // X derivative
268ca735530SJeremy L Thompson       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 0 * NUM_COMP];
2699e201c85SYohann       __syncthreads();
270ca735530SJeremy L Thompson       for (CeedInt i = 0; i < Q_1d; i++)
271ca735530SJeremy 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)
2729e201c85SYohann       __syncthreads();
2739e201c85SYohann       // Y derivative
274ca735530SJeremy L Thompson       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 1 * NUM_COMP];
2759e201c85SYohann       __syncthreads();
276ca735530SJeremy L Thompson       for (CeedInt i = 0; i < Q_1d; i++)
277ca735530SJeremy 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)
2789e201c85SYohann       __syncthreads();
2799e201c85SYohann       // Z derivative
280ca735530SJeremy L Thompson       for (CeedInt i = 0; i < Q_1d; i++)
281ca735530SJeremy 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)
2829e201c85SYohann     }
2839e201c85SYohann   }
2849e201c85SYohann }
285