xref: /libCEED/include/ceed/jit-source/hip/hip-shared-basis-tensor-templates.h (revision 02219a082eb38cf2d3edc97fbfe55fa395a4dc99)
1d275d636SJeremy L Thompson // Copyright (c) 2017-2025, 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 shared memory tensor product basis templates
10c0b5abf0SJeremy L Thompson #include <ceed/types.h>
119e201c85SYohann 
129e201c85SYohann //------------------------------------------------------------------------------
139e201c85SYohann // 1D
149e201c85SYohann //------------------------------------------------------------------------------
159e201c85SYohann 
169e201c85SYohann //------------------------------------------------------------------------------
179e201c85SYohann // 1D tensor contraction x
189e201c85SYohann //------------------------------------------------------------------------------
199e201c85SYohann template <int NUM_COMP, int P_1D, int Q_1D>
209e201c85SYohann inline __device__ void ContractX1d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
21d6c19ee8SJeremy L Thompson   __syncthreads();
229e201c85SYohann   data.slice[data.t_id_x] = *U;
239e201c85SYohann   __syncthreads();
249e201c85SYohann   *V = 0.0;
259e201c85SYohann   if (data.t_id_x < Q_1D) {
269e201c85SYohann     for (CeedInt i = 0; i < P_1D; i++) {
279e201c85SYohann       *V += B[i + data.t_id_x * P_1D] * data.slice[i];  // Contract x direction
289e201c85SYohann     }
299e201c85SYohann   }
309e201c85SYohann }
319e201c85SYohann 
329e201c85SYohann //------------------------------------------------------------------------------
339e201c85SYohann // 1D transpose tensor contraction x
349e201c85SYohann //------------------------------------------------------------------------------
359e201c85SYohann template <int NUM_COMP, int P_1D, int Q_1D>
369e201c85SYohann inline __device__ void ContractTransposeX1d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
37d6c19ee8SJeremy L Thompson   __syncthreads();
389e201c85SYohann   data.slice[data.t_id_x] = *U;
399e201c85SYohann   __syncthreads();
409e201c85SYohann   *V = 0.0;
419e201c85SYohann   if (data.t_id_x < P_1D) {
429e201c85SYohann     for (CeedInt i = 0; i < Q_1D; i++) {
439e201c85SYohann       *V += B[data.t_id_x + i * P_1D] * data.slice[i];  // Contract x direction
449e201c85SYohann     }
459e201c85SYohann   }
469e201c85SYohann }
479e201c85SYohann 
489e201c85SYohann //------------------------------------------------------------------------------
499e201c85SYohann // 1D interpolate to quadrature points
509e201c85SYohann //------------------------------------------------------------------------------
516b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
529e201c85SYohann inline __device__ void Interp1d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, CeedScalar *__restrict__ r_V) {
539e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
54db2becc9SJeremy L Thompson     ContractX1d<NUM_COMP, P_1D, Q_1D>(data, &r_U[comp], c_B, &r_V[comp]);
559e201c85SYohann   }
569e201c85SYohann }
579e201c85SYohann 
589e201c85SYohann //------------------------------------------------------------------------------
599e201c85SYohann // 1D interpolate transpose
609e201c85SYohann //------------------------------------------------------------------------------
616b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
622b730f8bSJeremy L Thompson inline __device__ void InterpTranspose1d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B,
632b730f8bSJeremy L Thompson                                          CeedScalar *__restrict__ r_V) {
649e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
65db2becc9SJeremy L Thompson     ContractTransposeX1d<NUM_COMP, P_1D, Q_1D>(data, &r_U[comp], c_B, &r_V[comp]);
669e201c85SYohann   }
679e201c85SYohann }
689e201c85SYohann 
699e201c85SYohann //------------------------------------------------------------------------------
709e201c85SYohann // 1D derivatives at quadrature points
719e201c85SYohann //------------------------------------------------------------------------------
726b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
732b730f8bSJeremy L Thompson inline __device__ void Grad1d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, const CeedScalar *c_G,
742b730f8bSJeremy L Thompson                               CeedScalar *__restrict__ r_V) {
759e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
76db2becc9SJeremy L Thompson     ContractX1d<NUM_COMP, P_1D, Q_1D>(data, &r_U[comp], c_G, &r_V[comp]);
779e201c85SYohann   }
789e201c85SYohann }
799e201c85SYohann 
809e201c85SYohann //------------------------------------------------------------------------------
819e201c85SYohann // 1D derivatives transpose
829e201c85SYohann //------------------------------------------------------------------------------
836b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
842b730f8bSJeremy L Thompson inline __device__ void GradTranspose1d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, const CeedScalar *c_G,
852b730f8bSJeremy L Thompson                                        CeedScalar *__restrict__ r_V) {
869e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
87db2becc9SJeremy L Thompson     ContractTransposeX1d<NUM_COMP, P_1D, Q_1D>(data, &r_U[comp], c_G, &r_V[comp]);
889e201c85SYohann   }
899e201c85SYohann }
909e201c85SYohann 
919e201c85SYohann //------------------------------------------------------------------------------
929e201c85SYohann // 1D quadrature weights
939e201c85SYohann //------------------------------------------------------------------------------
949b91271bSJeremy L Thompson template <int P_1D, int Q_1D>
959e201c85SYohann inline __device__ void Weight1d(SharedData_Hip &data, const CeedScalar *__restrict__ q_weight_1d, CeedScalar *w) {
969e201c85SYohann   *w = (data.t_id_x < Q_1D) ? q_weight_1d[data.t_id_x] : 0.0;
979e201c85SYohann }
989e201c85SYohann 
999e201c85SYohann //------------------------------------------------------------------------------
1009e201c85SYohann // 2D
1019e201c85SYohann //------------------------------------------------------------------------------
1029e201c85SYohann 
1039e201c85SYohann //------------------------------------------------------------------------------
1049e201c85SYohann // 2D tensor contraction x
1059e201c85SYohann //------------------------------------------------------------------------------
1066b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
1079e201c85SYohann inline __device__ void ContractX2d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
108d6c19ee8SJeremy L Thompson   __syncthreads();
1099e201c85SYohann   data.slice[data.t_id_x + data.t_id_y * T_1D] = *U;
1109e201c85SYohann   __syncthreads();
1119e201c85SYohann   *V = 0.0;
1129e201c85SYohann   if (data.t_id_x < Q_1D && data.t_id_y < P_1D) {
1139e201c85SYohann     for (CeedInt i = 0; i < P_1D; i++) {
1149e201c85SYohann       *V += B[i + data.t_id_x * P_1D] * data.slice[i + data.t_id_y * T_1D];  // Contract x direction
1159e201c85SYohann     }
1169e201c85SYohann   }
1179e201c85SYohann }
1189e201c85SYohann 
1199e201c85SYohann //------------------------------------------------------------------------------
1209e201c85SYohann // 2D tensor contract y
1219e201c85SYohann //------------------------------------------------------------------------------
1226b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
1239e201c85SYohann inline __device__ void ContractY2d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
124d6c19ee8SJeremy L Thompson   __syncthreads();
1259e201c85SYohann   data.slice[data.t_id_x + data.t_id_y * T_1D] = *U;
1269e201c85SYohann   __syncthreads();
1279e201c85SYohann   *V = 0.0;
1289e201c85SYohann   if (data.t_id_x < Q_1D && data.t_id_y < Q_1D) {
1299e201c85SYohann     for (CeedInt i = 0; i < P_1D; i++) {
1309e201c85SYohann       *V += B[i + data.t_id_y * P_1D] * data.slice[data.t_id_x + i * T_1D];  // Contract y direction
1319e201c85SYohann     }
1329e201c85SYohann   }
1339e201c85SYohann }
1349e201c85SYohann 
1359e201c85SYohann //------------------------------------------------------------------------------
1369e201c85SYohann // 2D transpose tensor contract y
1379e201c85SYohann //------------------------------------------------------------------------------
1386b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
1399e201c85SYohann inline __device__ void ContractTransposeY2d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
140d6c19ee8SJeremy L Thompson   __syncthreads();
1419e201c85SYohann   data.slice[data.t_id_x + data.t_id_y * T_1D] = *U;
1429e201c85SYohann   __syncthreads();
1439e201c85SYohann   *V = 0.0;
1449e201c85SYohann   if (data.t_id_x < Q_1D && data.t_id_y < P_1D) {
1459e201c85SYohann     for (CeedInt i = 0; i < Q_1D; i++) {
1469e201c85SYohann       *V += B[data.t_id_y + i * P_1D] * data.slice[data.t_id_x + i * T_1D];  // Contract y direction
1479e201c85SYohann     }
1489e201c85SYohann   }
1499e201c85SYohann }
1509e201c85SYohann 
1519e201c85SYohann //------------------------------------------------------------------------------
1529e201c85SYohann // 2D transpose tensor contract x
1539e201c85SYohann //------------------------------------------------------------------------------
1546b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
1559e201c85SYohann inline __device__ void ContractTransposeX2d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
156d6c19ee8SJeremy L Thompson   __syncthreads();
1579e201c85SYohann   data.slice[data.t_id_x + data.t_id_y * T_1D] = *U;
1589e201c85SYohann   __syncthreads();
1599e201c85SYohann   *V = 0.0;
1609e201c85SYohann   if (data.t_id_x < P_1D && data.t_id_y < P_1D) {
1619e201c85SYohann     for (CeedInt i = 0; i < Q_1D; i++) {
1629e201c85SYohann       *V += B[data.t_id_x + i * P_1D] * data.slice[i + data.t_id_y * T_1D];  // Contract x direction
1639e201c85SYohann     }
1649e201c85SYohann   }
1659e201c85SYohann }
1669e201c85SYohann 
1679e201c85SYohann //------------------------------------------------------------------------------
1689e201c85SYohann // 2D transpose tensor contract and add x
1699e201c85SYohann //------------------------------------------------------------------------------
1706b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
1719e201c85SYohann inline __device__ void ContractTransposeAddX2d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
172d6c19ee8SJeremy L Thompson   __syncthreads();
1739e201c85SYohann   data.slice[data.t_id_x + data.t_id_y * T_1D] = *U;
1749e201c85SYohann   __syncthreads();
1759e201c85SYohann   if (data.t_id_x < P_1D && data.t_id_y < P_1D) {
1769e201c85SYohann     for (CeedInt i = 0; i < Q_1D; i++) {
1779e201c85SYohann       *V += B[data.t_id_x + i * P_1D] * data.slice[i + data.t_id_y * T_1D];  // Contract x direction
1789e201c85SYohann     }
1799e201c85SYohann   }
1809e201c85SYohann }
1819e201c85SYohann 
1829e201c85SYohann //------------------------------------------------------------------------------
1839e201c85SYohann // 2D interpolate to quadrature points
1849e201c85SYohann //------------------------------------------------------------------------------
1856b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
1869e201c85SYohann inline __device__ void InterpTensor2d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, CeedScalar *__restrict__ r_V) {
1879e201c85SYohann   CeedScalar r_t[1];
1889e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
1896b92dc4bSJeremy L Thompson     ContractX2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp], c_B, r_t);
1906b92dc4bSJeremy L Thompson     ContractY2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t, c_B, &r_V[comp]);
1919e201c85SYohann   }
1929e201c85SYohann }
1939e201c85SYohann 
1949e201c85SYohann //------------------------------------------------------------------------------
1959e201c85SYohann // 2D interpolate transpose
1969e201c85SYohann //------------------------------------------------------------------------------
1976b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
1982b730f8bSJeremy L Thompson inline __device__ void InterpTransposeTensor2d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B,
1992b730f8bSJeremy L Thompson                                                CeedScalar *__restrict__ r_V) {
2009e201c85SYohann   CeedScalar r_t[1];
2019e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
2026b92dc4bSJeremy L Thompson     ContractTransposeY2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp], c_B, r_t);
2036b92dc4bSJeremy L Thompson     ContractTransposeX2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t, c_B, &r_V[comp]);
2049e201c85SYohann   }
2059e201c85SYohann }
2069e201c85SYohann 
2079e201c85SYohann //------------------------------------------------------------------------------
2089e201c85SYohann // 2D derivatives at quadrature points
2099e201c85SYohann //------------------------------------------------------------------------------
2106b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
2112b730f8bSJeremy L Thompson inline __device__ void GradTensor2d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, const CeedScalar *c_G,
2122b730f8bSJeremy L Thompson                                     CeedScalar *__restrict__ r_V) {
2139e201c85SYohann   CeedScalar r_t[1];
2149e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
2156b92dc4bSJeremy L Thompson     ContractX2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp], c_G, r_t);
2166b92dc4bSJeremy L Thompson     ContractY2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t, c_B, &r_V[comp + 0 * NUM_COMP]);
2176b92dc4bSJeremy L Thompson     ContractX2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp], c_B, r_t);
2186b92dc4bSJeremy L Thompson     ContractY2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t, c_G, &r_V[comp + 1 * NUM_COMP]);
2199e201c85SYohann   }
2209e201c85SYohann }
2219e201c85SYohann 
2229e201c85SYohann //------------------------------------------------------------------------------
2239e201c85SYohann // 2D derivatives transpose
2249e201c85SYohann //------------------------------------------------------------------------------
2256b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
2262b730f8bSJeremy L Thompson inline __device__ void GradTransposeTensor2d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, const CeedScalar *c_G,
2272b730f8bSJeremy L Thompson                                              CeedScalar *__restrict__ r_V) {
2289e201c85SYohann   CeedScalar r_t[1];
2299e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
2306b92dc4bSJeremy L Thompson     ContractTransposeY2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp + 0 * NUM_COMP], c_B, r_t);
2316b92dc4bSJeremy L Thompson     ContractTransposeX2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t, c_G, &r_V[comp]);
2326b92dc4bSJeremy L Thompson     ContractTransposeY2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp + 1 * NUM_COMP], c_G, r_t);
2336b92dc4bSJeremy L Thompson     ContractTransposeAddX2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t, c_B, &r_V[comp]);
2349e201c85SYohann   }
2359e201c85SYohann }
2369e201c85SYohann 
2379e201c85SYohann //------------------------------------------------------------------------------
238*02219a08SJeremy L Thompson // 2D derivatives at quadrature points, nodes and quadrature points collocated
239*02219a08SJeremy L Thompson //------------------------------------------------------------------------------
240*02219a08SJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
241*02219a08SJeremy L Thompson inline __device__ void GradTensorCollocatedNodes2d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
242*02219a08SJeremy L Thompson                                                    CeedScalar *__restrict__ r_V) {
243*02219a08SJeremy L Thompson   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
244*02219a08SJeremy L Thompson     ContractX2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp], c_G, &r_V[comp + 0 * NUM_COMP]);
245*02219a08SJeremy L Thompson     ContractY2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp], c_G, &r_V[comp + 1 * NUM_COMP]);
246*02219a08SJeremy L Thompson   }
247*02219a08SJeremy L Thompson }
248*02219a08SJeremy L Thompson 
249*02219a08SJeremy L Thompson //------------------------------------------------------------------------------
250*02219a08SJeremy L Thompson // 2D derivatives transpose, nodes and quadrature points collocated
251*02219a08SJeremy L Thompson //------------------------------------------------------------------------------
252*02219a08SJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
253*02219a08SJeremy L Thompson inline __device__ void GradTransposeTensorCollocatedNodes2d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
254*02219a08SJeremy L Thompson                                                             CeedScalar *__restrict__ r_V) {
255*02219a08SJeremy L Thompson   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
256*02219a08SJeremy L Thompson     ContractTransposeY2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp + 1 * NUM_COMP], c_G, &r_V[comp]);
257*02219a08SJeremy L Thompson     ContractTransposeAddX2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp + 0 * NUM_COMP], c_G, &r_V[comp]);
258*02219a08SJeremy L Thompson   }
259*02219a08SJeremy L Thompson }
260*02219a08SJeremy L Thompson 
261*02219a08SJeremy L Thompson //------------------------------------------------------------------------------
2629e201c85SYohann // 2D quadrature weights
2639e201c85SYohann //------------------------------------------------------------------------------
264ca595be6SJeremy L Thompson template <int P_1D, int Q_1D>
2659e201c85SYohann inline __device__ void WeightTensor2d(SharedData_Hip &data, const CeedScalar *__restrict__ q_weight_1d, CeedScalar *w) {
2662b730f8bSJeremy L Thompson   *w = (data.t_id_x < Q_1D && data.t_id_y < Q_1D) ? q_weight_1d[data.t_id_x] * q_weight_1d[data.t_id_y] : 0.0;
2679e201c85SYohann }
2689e201c85SYohann 
2699e201c85SYohann //------------------------------------------------------------------------------
2709e201c85SYohann // 3D
2719e201c85SYohann //------------------------------------------------------------------------------
2729e201c85SYohann 
2739e201c85SYohann //------------------------------------------------------------------------------
2749e201c85SYohann // 3D tensor contract x
2759e201c85SYohann //------------------------------------------------------------------------------
2766b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
2779e201c85SYohann inline __device__ void ContractX3d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
2789e201c85SYohann   CeedScalar r_B[P_1D];
2799e201c85SYohann   for (CeedInt i = 0; i < P_1D; i++) {
2809e201c85SYohann     r_B[i] = B[i + data.t_id_x * P_1D];
2819e201c85SYohann   }
2829e201c85SYohann 
2839e201c85SYohann   for (CeedInt k = 0; k < P_1D; k++) {
284d6c19ee8SJeremy L Thompson     __syncthreads();
2859e201c85SYohann     data.slice[data.t_id_x + data.t_id_y * T_1D] = U[k];
2869e201c85SYohann     __syncthreads();
2879e201c85SYohann     V[k] = 0.0;
2889e201c85SYohann     if (data.t_id_x < Q_1D && data.t_id_y < P_1D) {
2899e201c85SYohann       for (CeedInt i = 0; i < P_1D; i++) {
2909e201c85SYohann         V[k] += r_B[i] * data.slice[i + data.t_id_y * T_1D];  // Contract x direction
2919e201c85SYohann       }
2929e201c85SYohann     }
2939e201c85SYohann   }
2949e201c85SYohann }
2959e201c85SYohann 
2969e201c85SYohann //------------------------------------------------------------------------------
2979e201c85SYohann // 3D tensor contract y
2989e201c85SYohann //------------------------------------------------------------------------------
2996b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
3009e201c85SYohann inline __device__ void ContractY3d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
3019e201c85SYohann   CeedScalar r_B[P_1D];
3029e201c85SYohann   for (CeedInt i = 0; i < P_1D; i++) {
3039e201c85SYohann     r_B[i] = B[i + data.t_id_y * P_1D];
3049e201c85SYohann   }
3059e201c85SYohann 
3069e201c85SYohann   for (CeedInt k = 0; k < P_1D; k++) {
307d6c19ee8SJeremy L Thompson     __syncthreads();
3089e201c85SYohann     data.slice[data.t_id_x + data.t_id_y * T_1D] = U[k];
3099e201c85SYohann     __syncthreads();
3109e201c85SYohann     V[k] = 0.0;
3119e201c85SYohann     if (data.t_id_x < Q_1D && data.t_id_y < Q_1D) {
3129e201c85SYohann       for (CeedInt i = 0; i < P_1D; i++) {
3139e201c85SYohann         V[k] += r_B[i] * data.slice[data.t_id_x + i * T_1D];  // Contract y direction
3149e201c85SYohann       }
3159e201c85SYohann     }
3169e201c85SYohann   }
3179e201c85SYohann }
3189e201c85SYohann 
3199e201c85SYohann //------------------------------------------------------------------------------
3209e201c85SYohann // 3D tensor contract z
3219e201c85SYohann //------------------------------------------------------------------------------
3226b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
3239e201c85SYohann inline __device__ void ContractZ3d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
3249e201c85SYohann   for (CeedInt k = 0; k < Q_1D; k++) {
3259e201c85SYohann     V[k] = 0.0;
3269e201c85SYohann     if (data.t_id_x < Q_1D && data.t_id_y < Q_1D) {
3279e201c85SYohann       for (CeedInt i = 0; i < P_1D; i++) {
3289e201c85SYohann         V[k] += B[i + k * P_1D] * U[i];  // Contract z direction
3299e201c85SYohann       }
3309e201c85SYohann     }
3319e201c85SYohann   }
3329e201c85SYohann }
3339e201c85SYohann 
3349e201c85SYohann //------------------------------------------------------------------------------
3359e201c85SYohann // 3D transpose tensor contract z
3369e201c85SYohann //------------------------------------------------------------------------------
3376b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
3389e201c85SYohann inline __device__ void ContractTransposeZ3d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
3399e201c85SYohann   for (CeedInt k = 0; k < P_1D; k++) {
3409e201c85SYohann     V[k] = 0.0;
3419e201c85SYohann     if (data.t_id_x < Q_1D && data.t_id_y < Q_1D) {
3429e201c85SYohann       for (CeedInt i = 0; i < Q_1D; i++) {
3439e201c85SYohann         V[k] += B[k + i * P_1D] * U[i];  // Contract z direction
3449e201c85SYohann       }
3459e201c85SYohann     }
3469e201c85SYohann   }
3479e201c85SYohann }
3489e201c85SYohann 
3499e201c85SYohann //------------------------------------------------------------------------------
3509e201c85SYohann // 3D transpose tensor contract y
3519e201c85SYohann //------------------------------------------------------------------------------
3526b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
3539e201c85SYohann inline __device__ void ContractTransposeY3d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
3549e201c85SYohann   CeedScalar r_B[Q_1D];
3559e201c85SYohann   for (CeedInt i = 0; i < Q_1D; i++) {
3569e201c85SYohann     r_B[i] = B[data.t_id_y + i * P_1D];
3579e201c85SYohann   }
3589e201c85SYohann 
3599e201c85SYohann   for (CeedInt k = 0; k < P_1D; k++) {
360d6c19ee8SJeremy L Thompson     __syncthreads();
3619e201c85SYohann     data.slice[data.t_id_x + data.t_id_y * T_1D] = U[k];
3629e201c85SYohann     __syncthreads();
3639e201c85SYohann     V[k] = 0.0;
3649e201c85SYohann     if (data.t_id_x < Q_1D && data.t_id_y < P_1D) {
3659e201c85SYohann       for (CeedInt i = 0; i < Q_1D; i++) {
3669e201c85SYohann         V[k] += r_B[i] * data.slice[data.t_id_x + i * T_1D];  // Contract y direction
3679e201c85SYohann       }
3689e201c85SYohann     }
3699e201c85SYohann   }
3709e201c85SYohann }
3719e201c85SYohann 
3729e201c85SYohann //------------------------------------------------------------------------------
3739e201c85SYohann // 3D transpose tensor contract y
3749e201c85SYohann //------------------------------------------------------------------------------
3756b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
3769e201c85SYohann inline __device__ void ContractTransposeAddY3d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
3779e201c85SYohann   CeedScalar r_B[Q_1D];
3789e201c85SYohann   for (CeedInt i = 0; i < Q_1D; i++) {
3799e201c85SYohann     r_B[i] = B[data.t_id_y + i * P_1D];
3809e201c85SYohann   }
3819e201c85SYohann 
3829e201c85SYohann   for (CeedInt k = 0; k < P_1D; k++) {
383d6c19ee8SJeremy L Thompson     __syncthreads();
3849e201c85SYohann     data.slice[data.t_id_x + data.t_id_y * T_1D] = U[k];
3859e201c85SYohann     __syncthreads();
3869e201c85SYohann     if (data.t_id_x < Q_1D && data.t_id_y < P_1D) {
3879e201c85SYohann       for (CeedInt i = 0; i < Q_1D; i++) {
3889e201c85SYohann         V[k] += r_B[i] * data.slice[data.t_id_x + i * T_1D];  // Contract y direction
3899e201c85SYohann       }
3909e201c85SYohann     }
3919e201c85SYohann   }
3929e201c85SYohann }
3939e201c85SYohann 
3949e201c85SYohann //------------------------------------------------------------------------------
3959e201c85SYohann // 3D transpose tensor contract x
3969e201c85SYohann //------------------------------------------------------------------------------
3976b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
3989e201c85SYohann inline __device__ void ContractTransposeX3d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
3999e201c85SYohann   CeedScalar r_B[Q_1D];
4009e201c85SYohann   for (CeedInt i = 0; i < Q_1D; i++) {
4019e201c85SYohann     r_B[i] = B[data.t_id_x + i * P_1D];
4029e201c85SYohann   }
4039e201c85SYohann 
4049e201c85SYohann   for (CeedInt k = 0; k < P_1D; k++) {
405d6c19ee8SJeremy L Thompson     __syncthreads();
4069e201c85SYohann     data.slice[data.t_id_x + data.t_id_y * T_1D] = U[k];
4079e201c85SYohann     __syncthreads();
4089e201c85SYohann     V[k] = 0.0;
4099e201c85SYohann     if (data.t_id_x < P_1D && data.t_id_y < P_1D) {
4109e201c85SYohann       for (CeedInt i = 0; i < Q_1D; i++) {
4119e201c85SYohann         V[k] += r_B[i] * data.slice[i + data.t_id_y * T_1D];  // Contract x direction
4129e201c85SYohann       }
4139e201c85SYohann     }
4149e201c85SYohann   }
4159e201c85SYohann }
4169e201c85SYohann 
4179e201c85SYohann //------------------------------------------------------------------------------
4189e201c85SYohann // 3D transpose tensor contract add x
4199e201c85SYohann //------------------------------------------------------------------------------
4206b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
4219e201c85SYohann inline __device__ void ContractTransposeAddX3d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
4229e201c85SYohann   CeedScalar r_B[Q_1D];
4239e201c85SYohann   for (CeedInt i = 0; i < Q_1D; i++) {
4249e201c85SYohann     r_B[i] = B[data.t_id_x + i * P_1D];
4259e201c85SYohann   }
4269e201c85SYohann 
4279e201c85SYohann   for (CeedInt k = 0; k < P_1D; k++) {
428d6c19ee8SJeremy L Thompson     __syncthreads();
4299e201c85SYohann     data.slice[data.t_id_x + data.t_id_y * T_1D] = U[k];
4309e201c85SYohann     __syncthreads();
4319e201c85SYohann     if (data.t_id_x < P_1D && data.t_id_y < P_1D) {
4329e201c85SYohann       for (CeedInt i = 0; i < Q_1D; i++) {
4339e201c85SYohann         V[k] += r_B[i] * data.slice[i + data.t_id_y * T_1D];  // Contract x direction
4349e201c85SYohann       }
4359e201c85SYohann     }
4369e201c85SYohann   }
4379e201c85SYohann }
4389e201c85SYohann 
4399e201c85SYohann //------------------------------------------------------------------------------
4409e201c85SYohann // 3D interpolate to quadrature points
4419e201c85SYohann //------------------------------------------------------------------------------
4426b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
4439e201c85SYohann inline __device__ void InterpTensor3d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, CeedScalar *__restrict__ r_V) {
4449e201c85SYohann   CeedScalar r_t1[T_1D];
4459e201c85SYohann   CeedScalar r_t2[T_1D];
4469e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
4476b92dc4bSJeremy L Thompson     ContractX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * P_1D], c_B, r_t1);
4486b92dc4bSJeremy L Thompson     ContractY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_B, r_t2);
4496b92dc4bSJeremy L Thompson     ContractZ3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_B, &r_V[comp * Q_1D]);
4509e201c85SYohann   }
4519e201c85SYohann }
4529e201c85SYohann 
4539e201c85SYohann //------------------------------------------------------------------------------
4549e201c85SYohann // 3D interpolate transpose
4559e201c85SYohann //------------------------------------------------------------------------------
4566b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
4572b730f8bSJeremy L Thompson inline __device__ void InterpTransposeTensor3d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B,
4582b730f8bSJeremy L Thompson                                                CeedScalar *__restrict__ r_V) {
4599e201c85SYohann   CeedScalar r_t1[T_1D];
4609e201c85SYohann   CeedScalar r_t2[T_1D];
4619e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
4626b92dc4bSJeremy L Thompson     ContractTransposeZ3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * Q_1D], c_B, r_t1);
4636b92dc4bSJeremy L Thompson     ContractTransposeY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_B, r_t2);
4646b92dc4bSJeremy L Thompson     ContractTransposeX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_B, &r_V[comp * P_1D]);
4659e201c85SYohann   }
4669e201c85SYohann }
4679e201c85SYohann 
4689e201c85SYohann //------------------------------------------------------------------------------
4699e201c85SYohann // 3D derivatives at quadrature points
4709e201c85SYohann //------------------------------------------------------------------------------
4716b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
4722b730f8bSJeremy L Thompson inline __device__ void GradTensor3d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, const CeedScalar *c_G,
4732b730f8bSJeremy L Thompson                                     CeedScalar *__restrict__ r_V) {
4749e201c85SYohann   CeedScalar r_t1[T_1D];
4759e201c85SYohann   CeedScalar r_t2[T_1D];
4769e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
4776b92dc4bSJeremy L Thompson     ContractX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * P_1D], c_G, r_t1);
4786b92dc4bSJeremy L Thompson     ContractY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_B, r_t2);
4796b92dc4bSJeremy L Thompson     ContractZ3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_B, &r_V[comp * Q_1D + 0 * NUM_COMP * Q_1D]);
4806b92dc4bSJeremy L Thompson     ContractX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * P_1D], c_B, r_t1);
4816b92dc4bSJeremy L Thompson     ContractY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_G, r_t2);
4826b92dc4bSJeremy L Thompson     ContractZ3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_B, &r_V[comp * Q_1D + 1 * NUM_COMP * Q_1D]);
4836b92dc4bSJeremy L Thompson     ContractX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * P_1D], c_B, r_t1);
4846b92dc4bSJeremy L Thompson     ContractY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_B, r_t2);
4856b92dc4bSJeremy L Thompson     ContractZ3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_G, &r_V[comp * Q_1D + 2 * NUM_COMP * Q_1D]);
4869e201c85SYohann   }
4879e201c85SYohann }
4889e201c85SYohann 
4899e201c85SYohann //------------------------------------------------------------------------------
4909e201c85SYohann // 3D derivatives transpose
4919e201c85SYohann //------------------------------------------------------------------------------
4926b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
4932b730f8bSJeremy L Thompson inline __device__ void GradTransposeTensor3d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, const CeedScalar *c_G,
4942b730f8bSJeremy L Thompson                                              CeedScalar *__restrict__ r_V) {
4959e201c85SYohann   CeedScalar r_t1[T_1D];
4969e201c85SYohann   CeedScalar r_t2[T_1D];
4979e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
4986b92dc4bSJeremy L Thompson     ContractTransposeZ3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * Q_1D + 0 * NUM_COMP * Q_1D], c_B, r_t1);
4996b92dc4bSJeremy L Thompson     ContractTransposeY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_B, r_t2);
5006b92dc4bSJeremy L Thompson     ContractTransposeX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_G, &r_V[comp * P_1D]);
5016b92dc4bSJeremy L Thompson     ContractTransposeZ3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * Q_1D + 1 * NUM_COMP * Q_1D], c_B, r_t1);
5026b92dc4bSJeremy L Thompson     ContractTransposeY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_G, r_t2);
5036b92dc4bSJeremy L Thompson     ContractTransposeAddX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_B, &r_V[comp * P_1D]);
5046b92dc4bSJeremy L Thompson     ContractTransposeZ3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * Q_1D + 2 * NUM_COMP * Q_1D], c_G, r_t1);
5056b92dc4bSJeremy L Thompson     ContractTransposeY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_B, r_t2);
5066b92dc4bSJeremy L Thompson     ContractTransposeAddX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_B, &r_V[comp * P_1D]);
5079e201c85SYohann   }
5089e201c85SYohann }
5099e201c85SYohann 
5109e201c85SYohann //------------------------------------------------------------------------------
5119e201c85SYohann // 3D derivatives at quadrature points
5129e201c85SYohann //------------------------------------------------------------------------------
5136b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
5142b730f8bSJeremy L Thompson inline __device__ void GradTensorCollocated3d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, const CeedScalar *c_G,
5152b730f8bSJeremy L Thompson                                               CeedScalar *__restrict__ r_V) {
5169e201c85SYohann   CeedScalar r_t1[T_1D];
5179e201c85SYohann   CeedScalar r_t2[T_1D];
5189e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
5196b92dc4bSJeremy L Thompson     ContractX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * P_1D], c_B, r_t1);
5206b92dc4bSJeremy L Thompson     ContractY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_B, r_t2);
5216b92dc4bSJeremy L Thompson     ContractZ3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_B, r_t1);
5226b92dc4bSJeremy L Thompson     ContractX3d<NUM_COMP, Q_1D, Q_1D, T_1D>(data, r_t1, c_G, &r_V[comp * Q_1D + 0 * NUM_COMP * Q_1D]);
5236b92dc4bSJeremy L Thompson     ContractY3d<NUM_COMP, Q_1D, Q_1D, T_1D>(data, r_t1, c_G, &r_V[comp * Q_1D + 1 * NUM_COMP * Q_1D]);
5246b92dc4bSJeremy L Thompson     ContractZ3d<NUM_COMP, Q_1D, Q_1D, T_1D>(data, r_t1, c_G, &r_V[comp * Q_1D + 2 * NUM_COMP * Q_1D]);
5259e201c85SYohann   }
5269e201c85SYohann }
5279e201c85SYohann 
5289e201c85SYohann //------------------------------------------------------------------------------
5299e201c85SYohann // 3D derivatives transpose
5309e201c85SYohann //------------------------------------------------------------------------------
5316b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
5322b730f8bSJeremy L Thompson inline __device__ void GradTransposeTensorCollocated3d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B,
5332b730f8bSJeremy L Thompson                                                        const CeedScalar *c_G, CeedScalar *__restrict__ r_V) {
5349e201c85SYohann   CeedScalar r_t1[T_1D];
5359e201c85SYohann   CeedScalar r_t2[T_1D];
5369e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
5376b92dc4bSJeremy L Thompson     ContractTransposeZ3d<NUM_COMP, Q_1D, Q_1D, T_1D>(data, &r_U[comp * Q_1D + 2 * NUM_COMP * Q_1D], c_G, r_t2);
5386b92dc4bSJeremy L Thompson     ContractTransposeAddY3d<NUM_COMP, Q_1D, Q_1D, T_1D>(data, &r_U[comp * Q_1D + 1 * NUM_COMP * Q_1D], c_G, r_t2);
5396b92dc4bSJeremy L Thompson     ContractTransposeAddX3d<NUM_COMP, Q_1D, Q_1D, T_1D>(data, &r_U[comp * Q_1D + 0 * NUM_COMP * Q_1D], c_G, r_t2);
5406b92dc4bSJeremy L Thompson     ContractTransposeZ3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_B, r_t1);
5416b92dc4bSJeremy L Thompson     ContractTransposeY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_B, r_t2);
5426b92dc4bSJeremy L Thompson     ContractTransposeX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_B, &r_V[comp * P_1D]);
5439e201c85SYohann   }
5449e201c85SYohann }
5459e201c85SYohann 
5469e201c85SYohann //------------------------------------------------------------------------------
547*02219a08SJeremy L Thompson // 3D derivatives at quadrature points, nodes and quadrature points collocated
548*02219a08SJeremy L Thompson //------------------------------------------------------------------------------
549*02219a08SJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
550*02219a08SJeremy L Thompson inline __device__ void GradTensorCollocatedNodes3d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
551*02219a08SJeremy L Thompson                                                    CeedScalar *__restrict__ r_V) {
552*02219a08SJeremy L Thompson   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
553*02219a08SJeremy L Thompson     ContractX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * P_1D], c_G, &r_V[comp * Q_1D + 0 * NUM_COMP * Q_1D]);
554*02219a08SJeremy L Thompson     ContractY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * P_1D], c_G, &r_V[comp * Q_1D + 1 * NUM_COMP * Q_1D]);
555*02219a08SJeremy L Thompson     ContractZ3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * P_1D], c_G, &r_V[comp * Q_1D + 2 * NUM_COMP * Q_1D]);
556*02219a08SJeremy L Thompson   }
557*02219a08SJeremy L Thompson }
558*02219a08SJeremy L Thompson 
559*02219a08SJeremy L Thompson //------------------------------------------------------------------------------
560*02219a08SJeremy L Thompson // 3D derivatives transpose, nodes and quadrature points collocated
561*02219a08SJeremy L Thompson //------------------------------------------------------------------------------
562*02219a08SJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
563*02219a08SJeremy L Thompson inline __device__ void GradTransposeTensorCollocatedNodes3d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
564*02219a08SJeremy L Thompson                                                             CeedScalar *__restrict__ r_V) {
565*02219a08SJeremy L Thompson   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
566*02219a08SJeremy L Thompson     ContractTransposeZ3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * Q_1D + 2 * NUM_COMP * Q_1D], c_G, &r_V[comp * P_1D]);
567*02219a08SJeremy L Thompson     ContractTransposeAddY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * Q_1D + 1 * NUM_COMP * Q_1D], c_G, &r_V[comp * P_1D]);
568*02219a08SJeremy L Thompson     ContractTransposeAddX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * Q_1D + 0 * NUM_COMP * Q_1D], c_G, &r_V[comp * P_1D]);
569*02219a08SJeremy L Thompson   }
570*02219a08SJeremy L Thompson }
571*02219a08SJeremy L Thompson 
572*02219a08SJeremy L Thompson //------------------------------------------------------------------------------
5739e201c85SYohann // 3D quadrature weights
5749e201c85SYohann //------------------------------------------------------------------------------
5759b91271bSJeremy L Thompson template <int P_1D, int Q_1D>
5769e201c85SYohann inline __device__ void WeightTensor3d(SharedData_Hip &data, const CeedScalar *__restrict__ q_weight_1d, CeedScalar *w) {
5779e201c85SYohann   const bool       quad = (data.t_id_x < Q_1D && data.t_id_y < Q_1D);
5789e201c85SYohann   const CeedScalar pw   = quad ? q_weight_1d[data.t_id_x] * q_weight_1d[data.t_id_y] : 0.0;
5799e201c85SYohann   for (CeedInt q = 0; q < Q_1D; q++) {
5809e201c85SYohann     w[q] = quad ? pw * q_weight_1d[q] : 0.0;
5819e201c85SYohann   }
5829e201c85SYohann }
583