xref: /libCEED/include/ceed/jit-source/hip/hip-shared-basis-tensor-templates.h (revision 6b92dc4b38fb4a2ad1e306f024e5961058ad12b3)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
29e201c85SYohann // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
39e201c85SYohann //
49e201c85SYohann // SPDX-License-Identifier: BSD-2-Clause
59e201c85SYohann //
69e201c85SYohann // This file is part of CEED:  http://github.com/ceed
79e201c85SYohann 
89e201c85SYohann /// @file
99e201c85SYohann /// Internal header for HIP 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) {
219e201c85SYohann   data.slice[data.t_id_x] = *U;
229e201c85SYohann   __syncthreads();
239e201c85SYohann   *V = 0.0;
249e201c85SYohann   if (data.t_id_x < Q_1D) {
259e201c85SYohann     for (CeedInt i = 0; i < P_1D; i++) {
269e201c85SYohann       *V += B[i + data.t_id_x * P_1D] * data.slice[i];  // Contract x direction
279e201c85SYohann     }
289e201c85SYohann   }
299e201c85SYohann   __syncthreads();
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) {
379e201c85SYohann   data.slice[data.t_id_x] = *U;
389e201c85SYohann   __syncthreads();
399e201c85SYohann   *V = 0.0;
409e201c85SYohann   if (data.t_id_x < P_1D) {
419e201c85SYohann     for (CeedInt i = 0; i < Q_1D; i++) {
429e201c85SYohann       *V += B[data.t_id_x + i * P_1D] * data.slice[i];  // Contract x direction
439e201c85SYohann     }
449e201c85SYohann   }
459e201c85SYohann   __syncthreads();
469e201c85SYohann }
479e201c85SYohann 
489e201c85SYohann //------------------------------------------------------------------------------
499e201c85SYohann // 1D interpolate to quadrature points
509e201c85SYohann //------------------------------------------------------------------------------
51*6b92dc4bSJeremy 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 //------------------------------------------------------------------------------
61*6b92dc4bSJeremy 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 //------------------------------------------------------------------------------
72*6b92dc4bSJeremy 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 //------------------------------------------------------------------------------
83*6b92dc4bSJeremy 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 //------------------------------------------------------------------------------
949e201c85SYohann template <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 //------------------------------------------------------------------------------
106*6b92dc4bSJeremy 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) {
1089e201c85SYohann   data.slice[data.t_id_x + data.t_id_y * T_1D] = *U;
1099e201c85SYohann   __syncthreads();
1109e201c85SYohann   *V = 0.0;
1119e201c85SYohann   if (data.t_id_x < Q_1D && data.t_id_y < P_1D) {
1129e201c85SYohann     for (CeedInt i = 0; i < P_1D; i++) {
1139e201c85SYohann       *V += B[i + data.t_id_x * P_1D] * data.slice[i + data.t_id_y * T_1D];  // Contract x direction
1149e201c85SYohann     }
1159e201c85SYohann   }
1169e201c85SYohann   __syncthreads();
1179e201c85SYohann }
1189e201c85SYohann 
1199e201c85SYohann //------------------------------------------------------------------------------
1209e201c85SYohann // 2D tensor contract y
1219e201c85SYohann //------------------------------------------------------------------------------
122*6b92dc4bSJeremy 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) {
1249e201c85SYohann   data.slice[data.t_id_x + data.t_id_y * T_1D] = *U;
1259e201c85SYohann   __syncthreads();
1269e201c85SYohann   *V = 0.0;
1279e201c85SYohann   if (data.t_id_x < Q_1D && data.t_id_y < Q_1D) {
1289e201c85SYohann     for (CeedInt i = 0; i < P_1D; i++) {
1299e201c85SYohann       *V += B[i + data.t_id_y * P_1D] * data.slice[data.t_id_x + i * T_1D];  // Contract y direction
1309e201c85SYohann     }
1319e201c85SYohann   }
1329e201c85SYohann   __syncthreads();
1339e201c85SYohann }
1349e201c85SYohann 
1359e201c85SYohann //------------------------------------------------------------------------------
1369e201c85SYohann // 2D transpose tensor contract y
1379e201c85SYohann //------------------------------------------------------------------------------
138*6b92dc4bSJeremy 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) {
1409e201c85SYohann   data.slice[data.t_id_x + data.t_id_y * T_1D] = *U;
1419e201c85SYohann   __syncthreads();
1429e201c85SYohann   *V = 0.0;
1439e201c85SYohann   if (data.t_id_x < Q_1D && data.t_id_y < P_1D) {
1449e201c85SYohann     for (CeedInt i = 0; i < Q_1D; i++) {
1459e201c85SYohann       *V += B[data.t_id_y + i * P_1D] * data.slice[data.t_id_x + i * T_1D];  // Contract y direction
1469e201c85SYohann     }
1479e201c85SYohann   }
1489e201c85SYohann   __syncthreads();
1499e201c85SYohann }
1509e201c85SYohann 
1519e201c85SYohann //------------------------------------------------------------------------------
1529e201c85SYohann // 2D transpose tensor contract x
1539e201c85SYohann //------------------------------------------------------------------------------
154*6b92dc4bSJeremy 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) {
1569e201c85SYohann   data.slice[data.t_id_x + data.t_id_y * T_1D] = *U;
1579e201c85SYohann   __syncthreads();
1589e201c85SYohann   *V = 0.0;
1599e201c85SYohann   if (data.t_id_x < P_1D && data.t_id_y < P_1D) {
1609e201c85SYohann     for (CeedInt i = 0; i < Q_1D; i++) {
1619e201c85SYohann       *V += B[data.t_id_x + i * P_1D] * data.slice[i + data.t_id_y * T_1D];  // Contract x direction
1629e201c85SYohann     }
1639e201c85SYohann   }
1649e201c85SYohann   __syncthreads();
1659e201c85SYohann }
1669e201c85SYohann 
1679e201c85SYohann //------------------------------------------------------------------------------
1689e201c85SYohann // 2D transpose tensor contract and add x
1699e201c85SYohann //------------------------------------------------------------------------------
170*6b92dc4bSJeremy 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) {
1729e201c85SYohann   data.slice[data.t_id_x + data.t_id_y * T_1D] = *U;
1739e201c85SYohann   __syncthreads();
1749e201c85SYohann   if (data.t_id_x < P_1D && data.t_id_y < P_1D) {
1759e201c85SYohann     for (CeedInt i = 0; i < Q_1D; i++) {
1769e201c85SYohann       *V += B[data.t_id_x + i * P_1D] * data.slice[i + data.t_id_y * T_1D];  // Contract x direction
1779e201c85SYohann     }
1789e201c85SYohann   }
1799e201c85SYohann   __syncthreads();
1809e201c85SYohann }
1819e201c85SYohann 
1829e201c85SYohann //------------------------------------------------------------------------------
1839e201c85SYohann // 2D interpolate to quadrature points
1849e201c85SYohann //------------------------------------------------------------------------------
185*6b92dc4bSJeremy 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++) {
189*6b92dc4bSJeremy L Thompson     ContractX2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp], c_B, r_t);
190*6b92dc4bSJeremy 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 //------------------------------------------------------------------------------
197*6b92dc4bSJeremy 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++) {
202*6b92dc4bSJeremy L Thompson     ContractTransposeY2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp], c_B, r_t);
203*6b92dc4bSJeremy 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 //------------------------------------------------------------------------------
210*6b92dc4bSJeremy 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++) {
215*6b92dc4bSJeremy L Thompson     ContractX2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp], c_G, r_t);
216*6b92dc4bSJeremy L Thompson     ContractY2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t, c_B, &r_V[comp + 0 * NUM_COMP]);
217*6b92dc4bSJeremy L Thompson     ContractX2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp], c_B, r_t);
218*6b92dc4bSJeremy 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 //------------------------------------------------------------------------------
225*6b92dc4bSJeremy 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++) {
230*6b92dc4bSJeremy L Thompson     ContractTransposeY2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp + 0 * NUM_COMP], c_B, r_t);
231*6b92dc4bSJeremy L Thompson     ContractTransposeX2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t, c_G, &r_V[comp]);
232*6b92dc4bSJeremy L Thompson     ContractTransposeY2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp + 1 * NUM_COMP], c_G, r_t);
233*6b92dc4bSJeremy L Thompson     ContractTransposeAddX2d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t, c_B, &r_V[comp]);
2349e201c85SYohann   }
2359e201c85SYohann }
2369e201c85SYohann 
2379e201c85SYohann //------------------------------------------------------------------------------
2389e201c85SYohann // 2D quadrature weights
2399e201c85SYohann //------------------------------------------------------------------------------
2409e201c85SYohann template <int Q_1D>
2419e201c85SYohann inline __device__ void WeightTensor2d(SharedData_Hip &data, const CeedScalar *__restrict__ q_weight_1d, CeedScalar *w) {
2422b730f8bSJeremy 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;
2439e201c85SYohann }
2449e201c85SYohann 
2459e201c85SYohann //------------------------------------------------------------------------------
2469e201c85SYohann // 3D
2479e201c85SYohann //------------------------------------------------------------------------------
2489e201c85SYohann 
2499e201c85SYohann //------------------------------------------------------------------------------
2509e201c85SYohann // 3D tensor contract x
2519e201c85SYohann //------------------------------------------------------------------------------
252*6b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
2539e201c85SYohann inline __device__ void ContractX3d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
2549e201c85SYohann   CeedScalar r_B[P_1D];
2559e201c85SYohann   for (CeedInt i = 0; i < P_1D; i++) {
2569e201c85SYohann     r_B[i] = B[i + data.t_id_x * P_1D];
2579e201c85SYohann   }
2589e201c85SYohann 
2599e201c85SYohann   for (CeedInt k = 0; k < P_1D; k++) {
2609e201c85SYohann     data.slice[data.t_id_x + data.t_id_y * T_1D] = U[k];
2619e201c85SYohann     __syncthreads();
2629e201c85SYohann     V[k] = 0.0;
2639e201c85SYohann     if (data.t_id_x < Q_1D && data.t_id_y < P_1D) {
2649e201c85SYohann       for (CeedInt i = 0; i < P_1D; i++) {
2659e201c85SYohann         V[k] += r_B[i] * data.slice[i + data.t_id_y * T_1D];  // Contract x direction
2669e201c85SYohann       }
2679e201c85SYohann     }
2689e201c85SYohann     __syncthreads();
2699e201c85SYohann   }
2709e201c85SYohann }
2719e201c85SYohann 
2729e201c85SYohann //------------------------------------------------------------------------------
2739e201c85SYohann // 3D tensor contract y
2749e201c85SYohann //------------------------------------------------------------------------------
275*6b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
2769e201c85SYohann inline __device__ void ContractY3d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
2779e201c85SYohann   CeedScalar r_B[P_1D];
2789e201c85SYohann   for (CeedInt i = 0; i < P_1D; i++) {
2799e201c85SYohann     r_B[i] = B[i + data.t_id_y * P_1D];
2809e201c85SYohann   }
2819e201c85SYohann 
2829e201c85SYohann   for (CeedInt k = 0; k < P_1D; k++) {
2839e201c85SYohann     data.slice[data.t_id_x + data.t_id_y * T_1D] = U[k];
2849e201c85SYohann     __syncthreads();
2859e201c85SYohann     V[k] = 0.0;
2869e201c85SYohann     if (data.t_id_x < Q_1D && data.t_id_y < Q_1D) {
2879e201c85SYohann       for (CeedInt i = 0; i < P_1D; i++) {
2889e201c85SYohann         V[k] += r_B[i] * data.slice[data.t_id_x + i * T_1D];  // Contract y direction
2899e201c85SYohann       }
2909e201c85SYohann     }
2919e201c85SYohann     __syncthreads();
2929e201c85SYohann   }
2939e201c85SYohann }
2949e201c85SYohann 
2959e201c85SYohann //------------------------------------------------------------------------------
2969e201c85SYohann // 3D tensor contract z
2979e201c85SYohann //------------------------------------------------------------------------------
298*6b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
2999e201c85SYohann inline __device__ void ContractZ3d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
3009e201c85SYohann   for (CeedInt k = 0; k < Q_1D; k++) {
3019e201c85SYohann     V[k] = 0.0;
3029e201c85SYohann     if (data.t_id_x < Q_1D && data.t_id_y < Q_1D) {
3039e201c85SYohann       for (CeedInt i = 0; i < P_1D; i++) {
3049e201c85SYohann         V[k] += B[i + k * P_1D] * U[i];  // Contract z direction
3059e201c85SYohann       }
3069e201c85SYohann     }
3079e201c85SYohann   }
3089e201c85SYohann }
3099e201c85SYohann 
3109e201c85SYohann //------------------------------------------------------------------------------
3119e201c85SYohann // 3D transpose tensor contract z
3129e201c85SYohann //------------------------------------------------------------------------------
313*6b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
3149e201c85SYohann inline __device__ void ContractTransposeZ3d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
3159e201c85SYohann   for (CeedInt k = 0; k < P_1D; k++) {
3169e201c85SYohann     V[k] = 0.0;
3179e201c85SYohann     if (data.t_id_x < Q_1D && data.t_id_y < Q_1D) {
3189e201c85SYohann       for (CeedInt i = 0; i < Q_1D; i++) {
3199e201c85SYohann         V[k] += B[k + i * P_1D] * U[i];  // Contract z direction
3209e201c85SYohann       }
3219e201c85SYohann     }
3229e201c85SYohann   }
3239e201c85SYohann }
3249e201c85SYohann 
3259e201c85SYohann //------------------------------------------------------------------------------
3269e201c85SYohann // 3D transpose tensor contract y
3279e201c85SYohann //------------------------------------------------------------------------------
328*6b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
3299e201c85SYohann inline __device__ void ContractTransposeY3d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
3309e201c85SYohann   CeedScalar r_B[Q_1D];
3319e201c85SYohann   for (CeedInt i = 0; i < Q_1D; i++) {
3329e201c85SYohann     r_B[i] = B[data.t_id_y + i * P_1D];
3339e201c85SYohann   }
3349e201c85SYohann 
3359e201c85SYohann   for (CeedInt k = 0; k < P_1D; k++) {
3369e201c85SYohann     data.slice[data.t_id_x + data.t_id_y * T_1D] = U[k];
3379e201c85SYohann     __syncthreads();
3389e201c85SYohann     V[k] = 0.0;
3399e201c85SYohann     if (data.t_id_x < Q_1D && data.t_id_y < P_1D) {
3409e201c85SYohann       for (CeedInt i = 0; i < Q_1D; i++) {
3419e201c85SYohann         V[k] += r_B[i] * data.slice[data.t_id_x + i * T_1D];  // Contract y direction
3429e201c85SYohann       }
3439e201c85SYohann     }
3449e201c85SYohann     __syncthreads();
3459e201c85SYohann   }
3469e201c85SYohann }
3479e201c85SYohann 
3489e201c85SYohann //------------------------------------------------------------------------------
3499e201c85SYohann // 3D transpose tensor contract y
3509e201c85SYohann //------------------------------------------------------------------------------
351*6b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
3529e201c85SYohann inline __device__ void ContractTransposeAddY3d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
3539e201c85SYohann   CeedScalar r_B[Q_1D];
3549e201c85SYohann   for (CeedInt i = 0; i < Q_1D; i++) {
3559e201c85SYohann     r_B[i] = B[data.t_id_y + i * P_1D];
3569e201c85SYohann   }
3579e201c85SYohann 
3589e201c85SYohann   for (CeedInt k = 0; k < P_1D; k++) {
3599e201c85SYohann     data.slice[data.t_id_x + data.t_id_y * T_1D] = U[k];
3609e201c85SYohann     __syncthreads();
3619e201c85SYohann     if (data.t_id_x < Q_1D && data.t_id_y < P_1D) {
3629e201c85SYohann       for (CeedInt i = 0; i < Q_1D; i++) {
3639e201c85SYohann         V[k] += r_B[i] * data.slice[data.t_id_x + i * T_1D];  // Contract y direction
3649e201c85SYohann       }
3659e201c85SYohann     }
3669e201c85SYohann     __syncthreads();
3679e201c85SYohann   }
3689e201c85SYohann }
3699e201c85SYohann 
3709e201c85SYohann //------------------------------------------------------------------------------
3719e201c85SYohann // 3D transpose tensor contract x
3729e201c85SYohann //------------------------------------------------------------------------------
373*6b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
3749e201c85SYohann inline __device__ void ContractTransposeX3d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
3759e201c85SYohann   CeedScalar r_B[Q_1D];
3769e201c85SYohann   for (CeedInt i = 0; i < Q_1D; i++) {
3779e201c85SYohann     r_B[i] = B[data.t_id_x + i * P_1D];
3789e201c85SYohann   }
3799e201c85SYohann 
3809e201c85SYohann   for (CeedInt k = 0; k < P_1D; k++) {
3819e201c85SYohann     data.slice[data.t_id_x + data.t_id_y * T_1D] = U[k];
3829e201c85SYohann     __syncthreads();
3839e201c85SYohann     V[k] = 0.0;
3849e201c85SYohann     if (data.t_id_x < P_1D && data.t_id_y < P_1D) {
3859e201c85SYohann       for (CeedInt i = 0; i < Q_1D; i++) {
3869e201c85SYohann         V[k] += r_B[i] * data.slice[i + data.t_id_y * T_1D];  // Contract x direction
3879e201c85SYohann       }
3889e201c85SYohann     }
3899e201c85SYohann     __syncthreads();
3909e201c85SYohann   }
3919e201c85SYohann }
3929e201c85SYohann 
3939e201c85SYohann //------------------------------------------------------------------------------
3949e201c85SYohann // 3D transpose tensor contract add x
3959e201c85SYohann //------------------------------------------------------------------------------
396*6b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
3979e201c85SYohann inline __device__ void ContractTransposeAddX3d(SharedData_Hip &data, const CeedScalar *U, const CeedScalar *B, CeedScalar *V) {
3989e201c85SYohann   CeedScalar r_B[Q_1D];
3999e201c85SYohann   for (CeedInt i = 0; i < Q_1D; i++) {
4009e201c85SYohann     r_B[i] = B[data.t_id_x + i * P_1D];
4019e201c85SYohann   }
4029e201c85SYohann 
4039e201c85SYohann   for (CeedInt k = 0; k < P_1D; k++) {
4049e201c85SYohann     data.slice[data.t_id_x + data.t_id_y * T_1D] = U[k];
4059e201c85SYohann     __syncthreads();
4069e201c85SYohann     if (data.t_id_x < P_1D && data.t_id_y < P_1D) {
4079e201c85SYohann       for (CeedInt i = 0; i < Q_1D; i++) {
4089e201c85SYohann         V[k] += r_B[i] * data.slice[i + data.t_id_y * T_1D];  // Contract x direction
4099e201c85SYohann       }
4109e201c85SYohann     }
4119e201c85SYohann     __syncthreads();
4129e201c85SYohann   }
4139e201c85SYohann }
4149e201c85SYohann 
4159e201c85SYohann //------------------------------------------------------------------------------
4169e201c85SYohann // 3D interpolate to quadrature points
4179e201c85SYohann //------------------------------------------------------------------------------
418*6b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
4199e201c85SYohann inline __device__ void InterpTensor3d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, CeedScalar *__restrict__ r_V) {
4209e201c85SYohann   CeedScalar r_t1[T_1D];
4219e201c85SYohann   CeedScalar r_t2[T_1D];
4229e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
423*6b92dc4bSJeremy L Thompson     ContractX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * P_1D], c_B, r_t1);
424*6b92dc4bSJeremy L Thompson     ContractY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_B, r_t2);
425*6b92dc4bSJeremy L Thompson     ContractZ3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_B, &r_V[comp * Q_1D]);
4269e201c85SYohann   }
4279e201c85SYohann }
4289e201c85SYohann 
4299e201c85SYohann //------------------------------------------------------------------------------
4309e201c85SYohann // 3D interpolate transpose
4319e201c85SYohann //------------------------------------------------------------------------------
432*6b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
4332b730f8bSJeremy L Thompson inline __device__ void InterpTransposeTensor3d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B,
4342b730f8bSJeremy L Thompson                                                CeedScalar *__restrict__ r_V) {
4359e201c85SYohann   CeedScalar r_t1[T_1D];
4369e201c85SYohann   CeedScalar r_t2[T_1D];
4379e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
438*6b92dc4bSJeremy L Thompson     ContractTransposeZ3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * Q_1D], c_B, r_t1);
439*6b92dc4bSJeremy L Thompson     ContractTransposeY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_B, r_t2);
440*6b92dc4bSJeremy L Thompson     ContractTransposeX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_B, &r_V[comp * P_1D]);
4419e201c85SYohann   }
4429e201c85SYohann }
4439e201c85SYohann 
4449e201c85SYohann //------------------------------------------------------------------------------
4459e201c85SYohann // 3D derivatives at quadrature points
4469e201c85SYohann //------------------------------------------------------------------------------
447*6b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
4482b730f8bSJeremy L Thompson inline __device__ void GradTensor3d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, const CeedScalar *c_G,
4492b730f8bSJeremy L Thompson                                     CeedScalar *__restrict__ r_V) {
4509e201c85SYohann   CeedScalar r_t1[T_1D];
4519e201c85SYohann   CeedScalar r_t2[T_1D];
4529e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
453*6b92dc4bSJeremy L Thompson     ContractX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * P_1D], c_G, r_t1);
454*6b92dc4bSJeremy L Thompson     ContractY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_B, r_t2);
455*6b92dc4bSJeremy 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]);
456*6b92dc4bSJeremy L Thompson     ContractX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * P_1D], c_B, r_t1);
457*6b92dc4bSJeremy L Thompson     ContractY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_G, r_t2);
458*6b92dc4bSJeremy 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]);
459*6b92dc4bSJeremy L Thompson     ContractX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * P_1D], c_B, r_t1);
460*6b92dc4bSJeremy L Thompson     ContractY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_B, r_t2);
461*6b92dc4bSJeremy 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]);
4629e201c85SYohann   }
4639e201c85SYohann }
4649e201c85SYohann 
4659e201c85SYohann //------------------------------------------------------------------------------
4669e201c85SYohann // 3D derivatives transpose
4679e201c85SYohann //------------------------------------------------------------------------------
468*6b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
4692b730f8bSJeremy L Thompson inline __device__ void GradTransposeTensor3d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, const CeedScalar *c_G,
4702b730f8bSJeremy L Thompson                                              CeedScalar *__restrict__ r_V) {
4719e201c85SYohann   CeedScalar r_t1[T_1D];
4729e201c85SYohann   CeedScalar r_t2[T_1D];
4739e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
474*6b92dc4bSJeremy 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);
475*6b92dc4bSJeremy L Thompson     ContractTransposeY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_B, r_t2);
476*6b92dc4bSJeremy L Thompson     ContractTransposeX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_G, &r_V[comp * P_1D]);
477*6b92dc4bSJeremy 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);
478*6b92dc4bSJeremy L Thompson     ContractTransposeY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_G, r_t2);
479*6b92dc4bSJeremy L Thompson     ContractTransposeAddX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_B, &r_V[comp * P_1D]);
480*6b92dc4bSJeremy 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);
481*6b92dc4bSJeremy L Thompson     ContractTransposeY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_B, r_t2);
482*6b92dc4bSJeremy L Thompson     ContractTransposeAddX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_B, &r_V[comp * P_1D]);
4839e201c85SYohann   }
4849e201c85SYohann }
4859e201c85SYohann 
4869e201c85SYohann //------------------------------------------------------------------------------
4879e201c85SYohann // 3D derivatives at quadrature points
4889e201c85SYohann //------------------------------------------------------------------------------
489*6b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
4902b730f8bSJeremy L Thompson inline __device__ void GradTensorCollocated3d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B, const CeedScalar *c_G,
4912b730f8bSJeremy L Thompson                                               CeedScalar *__restrict__ r_V) {
4929e201c85SYohann   CeedScalar r_t1[T_1D];
4939e201c85SYohann   CeedScalar r_t2[T_1D];
4949e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
495*6b92dc4bSJeremy L Thompson     ContractX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, &r_U[comp * P_1D], c_B, r_t1);
496*6b92dc4bSJeremy L Thompson     ContractY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_B, r_t2);
497*6b92dc4bSJeremy L Thompson     ContractZ3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_B, r_t1);
498*6b92dc4bSJeremy 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]);
499*6b92dc4bSJeremy 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]);
500*6b92dc4bSJeremy 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]);
5019e201c85SYohann   }
5029e201c85SYohann }
5039e201c85SYohann 
5049e201c85SYohann //------------------------------------------------------------------------------
5059e201c85SYohann // 3D derivatives transpose
5069e201c85SYohann //------------------------------------------------------------------------------
507*6b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int Q_1D, int T_1D>
5082b730f8bSJeremy L Thompson inline __device__ void GradTransposeTensorCollocated3d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_B,
5092b730f8bSJeremy L Thompson                                                        const CeedScalar *c_G, CeedScalar *__restrict__ r_V) {
5109e201c85SYohann   CeedScalar r_t1[T_1D];
5119e201c85SYohann   CeedScalar r_t2[T_1D];
5129e201c85SYohann   for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
513*6b92dc4bSJeremy 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);
514*6b92dc4bSJeremy 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);
515*6b92dc4bSJeremy 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);
516*6b92dc4bSJeremy L Thompson     ContractTransposeZ3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_B, r_t1);
517*6b92dc4bSJeremy L Thompson     ContractTransposeY3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t1, c_B, r_t2);
518*6b92dc4bSJeremy L Thompson     ContractTransposeX3d<NUM_COMP, P_1D, Q_1D, T_1D>(data, r_t2, c_B, &r_V[comp * P_1D]);
5199e201c85SYohann   }
5209e201c85SYohann }
5219e201c85SYohann 
5229e201c85SYohann //------------------------------------------------------------------------------
5239e201c85SYohann // 3D quadrature weights
5249e201c85SYohann //------------------------------------------------------------------------------
5259e201c85SYohann template <int Q_1D>
5269e201c85SYohann inline __device__ void WeightTensor3d(SharedData_Hip &data, const CeedScalar *__restrict__ q_weight_1d, CeedScalar *w) {
5279e201c85SYohann   const bool       quad = (data.t_id_x < Q_1D && data.t_id_y < Q_1D);
5289e201c85SYohann   const CeedScalar pw   = quad ? q_weight_1d[data.t_id_x] * q_weight_1d[data.t_id_y] : 0.0;
5299e201c85SYohann   for (CeedInt q = 0; q < Q_1D; q++) {
5309e201c85SYohann     w[q] = quad ? pw * q_weight_1d[q] : 0.0;
5319e201c85SYohann   }
5329e201c85SYohann }
533