xref: /libCEED/include/ceed/jit-source/magma/magma-basis-grad-3d.h (revision d4cc18453651bd0f94c1a2e078b2646a92dafdcc)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
2f80f4a74SSebastian Grimberg // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3f80f4a74SSebastian Grimberg //
4f80f4a74SSebastian Grimberg // SPDX-License-Identifier: BSD-2-Clause
5f80f4a74SSebastian Grimberg //
6f80f4a74SSebastian Grimberg // This file is part of CEED:  http://github.com/ceed
7f80f4a74SSebastian Grimberg 
83c1e2affSSebastian Grimberg /// @file
93c1e2affSSebastian Grimberg /// Internal header for MAGMA tensor basis gradient in 3D
103c1e2affSSebastian Grimberg #include "magma-common-tensor.h"
113c1e2affSSebastian Grimberg 
12f80f4a74SSebastian Grimberg // macros to abstract access of shared memory and reg. file
133c1e2affSSebastian Grimberg #define sT(i, j) sT[(j) * P + (i)]
14f80f4a74SSebastian Grimberg #define sTmp(i, j, ldw) sTmp[(j) * (ldw) + (i)]
15f80f4a74SSebastian Grimberg #define sTmp2(i, j, ldw) sTmp2[(j) * (ldw) + (i)]
16f80f4a74SSebastian Grimberg 
179e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
187132caa0SSebastian Grimberg // Helper function to add or set into V
197132caa0SSebastian Grimberg template <typename T, bool Add>
207132caa0SSebastian Grimberg struct magma_grad_3d_device_accumulate;
217132caa0SSebastian Grimberg 
227132caa0SSebastian Grimberg template <typename T>
237132caa0SSebastian Grimberg struct magma_grad_3d_device_accumulate<T, true> {
247132caa0SSebastian Grimberg   static __device__ __inline__ void op(T &rV, const T &rTmp) { rV += rTmp; }
257132caa0SSebastian Grimberg };
267132caa0SSebastian Grimberg 
277132caa0SSebastian Grimberg template <typename T>
287132caa0SSebastian Grimberg struct magma_grad_3d_device_accumulate<T, false> {
297132caa0SSebastian Grimberg   static __device__ __inline__ void op(T &rV, const T &rTmp) { rV = rTmp; }
307132caa0SSebastian Grimberg };
317132caa0SSebastian Grimberg 
327132caa0SSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
33f80f4a74SSebastian Grimberg // grad basis action (3D)
34f80f4a74SSebastian Grimberg // This function is called three times at a higher level for 3D
353c1e2affSSebastian Grimberg // DIM_U   -- for the size of rU[DIM_U * NUM_COMP * MAX_P_Q]
363c1e2affSSebastian Grimberg // DIM_V   -- for the size of rV[DIM_V * NUM_COMP * MAX_P_Q]
373c1e2affSSebastian Grimberg // i_DIM   -- the index of the outermost loop over dimensions in grad
383c1e2affSSebastian Grimberg // i_DIM_U -- which dim index of rU is accessed (always 0 for notrans, 0, 1, or 2 for trans)
393c1e2affSSebastian Grimberg // i_DIM_V -- which dim index of rV is accessed (0, 1, or 2 for notrans, always 0 for trans)
407132caa0SSebastian Grimberg template <typename T, int DIM_U, int DIM_V, int NUM_COMP, int P, int Q, int rU_SIZE, int rV_SIZE, int i_DIM, int i_DIM_U, int i_DIM_V, bool ADD>
413c1e2affSSebastian Grimberg static __device__ __inline__ void magma_grad_3d_device(const T *sTinterp, const T *sTgrad, T rU[DIM_U][NUM_COMP][rU_SIZE],
427132caa0SSebastian Grimberg                                                        T rV[DIM_V][NUM_COMP][rV_SIZE], const int tx, T rTmp, T *swork) {
43f80f4a74SSebastian Grimberg   // Assumptions
443c1e2affSSebastian Grimberg   // 0. This device routine applies grad for one dim only (i_DIM), so it should be thrice for 3D
453c1e2affSSebastian Grimberg   // 1. 1D threads of size max(P,Q)^2
463c1e2affSSebastian Grimberg   // 2. input:  rU[DIM_U x NUM_COMP x rU_SIZE] in registers (per thread)
473c1e2affSSebastian Grimberg   // 3. output: rV[DIM_V x NUM_COMP x rV_SIZE] in registers (per thread)
48f80f4a74SSebastian Grimberg   // 4. Three products per each (dim,component) pair
493c1e2affSSebastian Grimberg   //  4.1 Batch P^2 of (1xP) matrices times (PxQ) matrix => Batch P^2 of (1xQ) matrices
503c1e2affSSebastian Grimberg   //  4.2 Batch P   of (QxP) matrices times (PxQ) matrix => Batch P   of (QxQ) matrices
513c1e2affSSebastian Grimberg   //  4.3 Batch 1   of (Q^2xP_) matrix times (PxQ) matrix => (Q^2xQ_) matrix
52f80f4a74SSebastian Grimberg   // 6. Each thread computes one row of the output of each product
53f80f4a74SSebastian Grimberg   // 7. Sync is recommended before and after the call
54f80f4a74SSebastian Grimberg 
55f80f4a74SSebastian Grimberg   T *sW1 = swork;
563c1e2affSSebastian Grimberg   T *sW2 = sW1 + P * P * Q;
573c1e2affSSebastian Grimberg   for (int comp = 0; comp < NUM_COMP; comp++) {
583c1e2affSSebastian Grimberg     // Batch P^2 of (1xP) matrices [reg] times (PxQ) matrix [shmem] => Batch P^2 of (1xQ) matrices [shmem]
593c1e2affSSebastian Grimberg     if (tx < (P * P)) {
60f80f4a74SSebastian Grimberg       const int batchid = tx;
61f80f4a74SSebastian Grimberg       const int sld     = 1;
623c1e2affSSebastian Grimberg       const T  *sT      = (i_DIM == 0) ? sTgrad : sTinterp;
633c1e2affSSebastian Grimberg       T        *sTmp    = sW1 + batchid * (1 * Q);
643c1e2affSSebastian Grimberg       for (int j = 0; j < Q; j++) {
65f80f4a74SSebastian Grimberg         rTmp = 0.0;
663c1e2affSSebastian Grimberg         for (int i = 0; i < P; i++) {
673c1e2affSSebastian Grimberg           rTmp += rU[i_DIM_U][comp][i] * sT(i, j);
68f80f4a74SSebastian Grimberg         }
69f80f4a74SSebastian Grimberg         sTmp(0, j, sld) = rTmp;
70f80f4a74SSebastian Grimberg       }
713c1e2affSSebastian Grimberg     }  // end of: if (tx < P*P)
72f80f4a74SSebastian Grimberg     __syncthreads();
73f80f4a74SSebastian Grimberg 
743c1e2affSSebastian Grimberg     // Batch P of (QxP) matrices [shmem] times (PxQ) matrix [shmem] => Batch P of (QxQ) matrices [reg]
753c1e2affSSebastian Grimberg     if (tx < (P * Q)) {
763c1e2affSSebastian Grimberg       const int batchid = tx / Q;
773c1e2affSSebastian Grimberg       const int tx_     = tx % Q;
783c1e2affSSebastian Grimberg       const int sld     = Q;
793c1e2affSSebastian Grimberg       const T  *sT      = (i_DIM == 1) ? sTgrad : sTinterp;
803c1e2affSSebastian Grimberg       T        *sTmp    = sW1 + batchid * (Q * P);  // sTmp is input
813c1e2affSSebastian Grimberg       T        *sTmp2   = sW2 + batchid * (Q * Q);  // sTmp2 is output
823c1e2affSSebastian Grimberg       for (int j = 0; j < Q; j++) {
83f80f4a74SSebastian Grimberg         rTmp = 0.0;
843c1e2affSSebastian Grimberg         for (int i = 0; i < P; i++) {
85f80f4a74SSebastian Grimberg           rTmp += sTmp(tx_, i, sld) * sT(i, j);
86f80f4a74SSebastian Grimberg         }
87f80f4a74SSebastian Grimberg         sTmp2(tx_, j, sld) = rTmp;
88f80f4a74SSebastian Grimberg       }
89f80f4a74SSebastian Grimberg     }
90f80f4a74SSebastian Grimberg     __syncthreads();
91f80f4a74SSebastian Grimberg 
923c1e2affSSebastian Grimberg     // Batch 1 of (Q^2xP_) matrices [shmem] times (PxQ) matrix [shmem] => Batch 1 of (Q^2xQ_) matrices [reg]
933c1e2affSSebastian Grimberg     if (tx < (Q * Q)) {
943c1e2affSSebastian Grimberg       // No need to declare batchid = (tx  / Q^2) = always zero
953c1e2affSSebastian Grimberg       // No need to declare tx_     = (tx_ % Q^2) = always tx
963c1e2affSSebastian Grimberg       const int sld  = Q * Q;
973c1e2affSSebastian Grimberg       const T  *sT   = (i_DIM == 2) ? sTgrad : sTinterp;
98f80f4a74SSebastian Grimberg       T        *sTmp = sW2;  // sTmp is input
993c1e2affSSebastian Grimberg       for (int j = 0; j < Q; j++) {
100f80f4a74SSebastian Grimberg         rTmp = 0.0;
1013c1e2affSSebastian Grimberg         for (int i = 0; i < P; i++) {
102f80f4a74SSebastian Grimberg           rTmp += sTmp(tx, i, sld) * sT(i, j);
103f80f4a74SSebastian Grimberg         }
1047132caa0SSebastian Grimberg         magma_grad_3d_device_accumulate<T, ADD>::op(rV[i_DIM_V][comp][j], rTmp);
105f80f4a74SSebastian Grimberg       }
106f80f4a74SSebastian Grimberg     }
107f80f4a74SSebastian Grimberg     __syncthreads();
1083c1e2affSSebastian Grimberg   }  // loop over NUM_COMP
109f80f4a74SSebastian Grimberg }
110f80f4a74SSebastian Grimberg 
1119e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
1123c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q *BASIS_MAX_P_Q, MAGMA_MAXTHREADS_3D)) __global__
113f80f4a74SSebastian Grimberg     void magma_gradn_3d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU,
114f80f4a74SSebastian Grimberg                                const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) {
115f80f4a74SSebastian Grimberg   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
116f80f4a74SSebastian Grimberg 
117f80f4a74SSebastian Grimberg   const int tx      = threadIdx.x;
118f80f4a74SSebastian Grimberg   const int ty      = threadIdx.y;
119f80f4a74SSebastian Grimberg   const int elem_id = (blockIdx.x * blockDim.y) + ty;
120f80f4a74SSebastian Grimberg 
121f80f4a74SSebastian Grimberg   if (elem_id >= nelem) return;
122f80f4a74SSebastian Grimberg 
1233c1e2affSSebastian Grimberg   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // here DIM_U = 1, but might be different for a fused operator
1243c1e2affSSebastian Grimberg   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // here DIM_V = 1, but might be different for a fused operator
125f80f4a74SSebastian Grimberg   CeedScalar rTmp                           = 0.0;
126f80f4a74SSebastian Grimberg 
127f80f4a74SSebastian Grimberg   // shift global memory pointers by elem stride
128f80f4a74SSebastian Grimberg   dU += elem_id * estrdU;
129f80f4a74SSebastian Grimberg   dV += elem_id * estrdV;
130f80f4a74SSebastian Grimberg 
131f80f4a74SSebastian Grimberg   // assign shared memory pointers
1323c1e2affSSebastian Grimberg   CeedScalar *sTinterp = (CeedScalar *)shared_data;
1333c1e2affSSebastian Grimberg   CeedScalar *sTgrad   = sTinterp + BASIS_P * BASIS_Q;
1343c1e2affSSebastian Grimberg   CeedScalar *sTmp     = sTgrad + BASIS_P * BASIS_Q;
1353c1e2affSSebastian Grimberg   sTmp += ty * (max(BASIS_P * BASIS_P * BASIS_P, (BASIS_P * BASIS_P * BASIS_Q) + (BASIS_P * BASIS_Q * BASIS_Q)));
136f80f4a74SSebastian Grimberg 
137f80f4a74SSebastian Grimberg   // read T
138f80f4a74SSebastian Grimberg   if (ty == 0) {
1399e0c01faSSebastian Grimberg     read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dinterp1d, sTinterp);
1409e0c01faSSebastian Grimberg     read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dgrad1d, sTgrad);
141f80f4a74SSebastian Grimberg   }
142f80f4a74SSebastian Grimberg   __syncthreads();
143f80f4a74SSebastian Grimberg 
1443c1e2affSSebastian Grimberg   /* read U (idim = 0 for dU, i_DIM = 0 for rU) --
145f80f4a74SSebastian Grimberg      there is a sync at the end of this function */
1469e0c01faSSebastian Grimberg   read_U_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx);
147f80f4a74SSebastian Grimberg 
1483c1e2affSSebastian Grimberg   /* first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) --
149f80f4a74SSebastian Grimberg      output from rV[0][][] into dV (idim = 0) */
1507132caa0SSebastian Grimberg   magma_grad_3d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_P, BASIS_Q, BASIS_P, BASIS_Q, 0, 0, 0, false>(sTinterp, sTgrad, rU, rV, tx, rTmp,
1517132caa0SSebastian Grimberg                                                                                                              sTmp);
152f80f4a74SSebastian Grimberg   /* there is a sync at the end of magma_grad_3d_device */
1539e0c01faSSebastian Grimberg   write_V_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (0 * dstrdV), cstrdV, rV, tx);
154f80f4a74SSebastian Grimberg 
1553c1e2affSSebastian Grimberg   /* second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) --
156f80f4a74SSebastian Grimberg      output from rV[0][][] into dV (idim = 1) */
1577132caa0SSebastian Grimberg   magma_grad_3d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_P, BASIS_Q, BASIS_P, BASIS_Q, 1, 0, 0, false>(sTinterp, sTgrad, rU, rV, tx, rTmp,
1587132caa0SSebastian Grimberg                                                                                                              sTmp);
159f80f4a74SSebastian Grimberg   /* there is a sync at the end of magma_grad_3d_device */
1609e0c01faSSebastian Grimberg   write_V_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (1 * dstrdV), cstrdV, rV, tx);
161f80f4a74SSebastian Grimberg 
1623c1e2affSSebastian Grimberg   /* third call (i_DIM = 2, i_DIM_U = 0, i_DIM_V = 0) --
163f80f4a74SSebastian Grimberg      output from rV[0][][] into dV (idim = 2) */
1647132caa0SSebastian Grimberg   magma_grad_3d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_P, BASIS_Q, BASIS_P, BASIS_Q, 2, 0, 0, false>(sTinterp, sTgrad, rU, rV, tx, rTmp,
1657132caa0SSebastian Grimberg                                                                                                              sTmp);
166f80f4a74SSebastian Grimberg   /* there is a sync at the end of magma_grad_3d_device */
1679e0c01faSSebastian Grimberg   write_V_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (2 * dstrdV), cstrdV, rV, tx);
168f80f4a74SSebastian Grimberg }
169f80f4a74SSebastian Grimberg 
1709e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
1713c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q *BASIS_MAX_P_Q, MAGMA_MAXTHREADS_3D)) __global__
172f80f4a74SSebastian Grimberg     void magma_gradt_3d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU,
173f80f4a74SSebastian Grimberg                                const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) {
174f80f4a74SSebastian Grimberg   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
175f80f4a74SSebastian Grimberg 
176f80f4a74SSebastian Grimberg   const int tx      = threadIdx.x;
177f80f4a74SSebastian Grimberg   const int ty      = threadIdx.y;
178f80f4a74SSebastian Grimberg   const int elem_id = (blockIdx.x * blockDim.y) + ty;
179f80f4a74SSebastian Grimberg 
180f80f4a74SSebastian Grimberg   if (elem_id >= nelem) return;
181f80f4a74SSebastian Grimberg 
1823c1e2affSSebastian Grimberg   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // here DIM_U = 1, but might be different for a fused operator
1833c1e2affSSebastian Grimberg   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // here DIM_V = 1, but might be different for a fused operator
184f80f4a74SSebastian Grimberg   CeedScalar rTmp                           = 0.0;
185f80f4a74SSebastian Grimberg 
186f80f4a74SSebastian Grimberg   // shift global memory pointers by elem stride
187f80f4a74SSebastian Grimberg   dU += elem_id * estrdU;
188f80f4a74SSebastian Grimberg   dV += elem_id * estrdV;
189f80f4a74SSebastian Grimberg 
190f80f4a74SSebastian Grimberg   // assign shared memory pointers
1913c1e2affSSebastian Grimberg   CeedScalar *sTinterp = (CeedScalar *)shared_data;
1923c1e2affSSebastian Grimberg   CeedScalar *sTgrad   = sTinterp + BASIS_Q * BASIS_P;
1933c1e2affSSebastian Grimberg   CeedScalar *sTmp     = sTgrad + BASIS_Q * BASIS_P;
1943c1e2affSSebastian Grimberg   sTmp += ty * (max(BASIS_Q * BASIS_Q * BASIS_Q, (BASIS_Q * BASIS_Q * BASIS_P) + (BASIS_Q * BASIS_P * BASIS_P)));
195f80f4a74SSebastian Grimberg 
196f80f4a74SSebastian Grimberg   // read T
197f80f4a74SSebastian Grimberg   if (ty == 0) {
1989e0c01faSSebastian Grimberg     read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dinterp1d, sTinterp);
1999e0c01faSSebastian Grimberg     read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dgrad1d, sTgrad);
200f80f4a74SSebastian Grimberg   }
201f80f4a74SSebastian Grimberg   __syncthreads();
202f80f4a74SSebastian Grimberg 
2033c1e2affSSebastian Grimberg   /* read U (idim = 0 for dU, i_DIM = 0 for rU) --
204f80f4a74SSebastian Grimberg      there is a sync at the end of this function */
2059e0c01faSSebastian Grimberg   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx);
2063c1e2affSSebastian Grimberg   /* then first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) */
2077132caa0SSebastian Grimberg   magma_grad_3d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P, 0, 0, 0, true>(sTinterp, sTgrad, rU, rV, tx, rTmp, sTmp);
208f80f4a74SSebastian Grimberg   /* there is a sync at the end of magma_grad_3d_device */
209f80f4a74SSebastian Grimberg 
2103c1e2affSSebastian Grimberg   /* read U (idim = 1 for dU, i_DIM = 0 for rU) --
211f80f4a74SSebastian Grimberg      there is a sync at the end of this function */
2129e0c01faSSebastian Grimberg   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (1 * dstrdU), cstrdU, rU, sTmp, tx);
2133c1e2affSSebastian Grimberg   /* then second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) */
2147132caa0SSebastian Grimberg   magma_grad_3d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P, 1, 0, 0, true>(sTinterp, sTgrad, rU, rV, tx, rTmp, sTmp);
215f80f4a74SSebastian Grimberg   /* there is a sync at the end of magma_grad_3d_device */
216f80f4a74SSebastian Grimberg 
2173c1e2affSSebastian Grimberg   /* read U (idim = 2 for dU, i_DIM = 0 for rU) --
218f80f4a74SSebastian Grimberg      there is a sync at the end of this function */
2199e0c01faSSebastian Grimberg   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (2 * dstrdU), cstrdU, rU, sTmp, tx);
2203c1e2affSSebastian Grimberg   /* then third call (i_DIM = 2, i_DIM_U = 0, i_DIM_V = 0) */
2217132caa0SSebastian Grimberg   magma_grad_3d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P, 2, 0, 0, true>(sTinterp, sTgrad, rU, rV, tx, rTmp, sTmp);
222f80f4a74SSebastian Grimberg   /* there is a sync at the end of magma_grad_3d_device */
223f80f4a74SSebastian Grimberg 
224f80f4a74SSebastian Grimberg   // write V
2259e0c01faSSebastian Grimberg   write_V_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV + (0 * dstrdV), cstrdV, rV, tx);
226f80f4a74SSebastian Grimberg }
227db2becc9SJeremy L Thompson 
228db2becc9SJeremy L Thompson ////////////////////////////////////////////////////////////////////////////////
229db2becc9SJeremy L Thompson extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q *BASIS_MAX_P_Q, MAGMA_MAXTHREADS_3D)) __global__
230db2becc9SJeremy L Thompson     void magma_gradta_3d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU,
231db2becc9SJeremy L Thompson                                 const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) {
232db2becc9SJeremy L Thompson   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
233db2becc9SJeremy L Thompson 
234db2becc9SJeremy L Thompson   const int tx      = threadIdx.x;
235db2becc9SJeremy L Thompson   const int ty      = threadIdx.y;
236db2becc9SJeremy L Thompson   const int elem_id = (blockIdx.x * blockDim.y) + ty;
237db2becc9SJeremy L Thompson 
238db2becc9SJeremy L Thompson   if (elem_id >= nelem) return;
239db2becc9SJeremy L Thompson 
240db2becc9SJeremy L Thompson   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // here DIM_U = 1, but might be different for a fused operator
241db2becc9SJeremy L Thompson   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // here DIM_V = 1, but might be different for a fused operator
242db2becc9SJeremy L Thompson   CeedScalar rTmp                           = 0.0;
243db2becc9SJeremy L Thompson 
244db2becc9SJeremy L Thompson   // shift global memory pointers by elem stride
245db2becc9SJeremy L Thompson   dU += elem_id * estrdU;
246db2becc9SJeremy L Thompson   dV += elem_id * estrdV;
247db2becc9SJeremy L Thompson 
248db2becc9SJeremy L Thompson   // assign shared memory pointers
249db2becc9SJeremy L Thompson   CeedScalar *sTinterp = (CeedScalar *)shared_data;
250db2becc9SJeremy L Thompson   CeedScalar *sTgrad   = sTinterp + BASIS_Q * BASIS_P;
251db2becc9SJeremy L Thompson   CeedScalar *sTmp     = sTgrad + BASIS_Q * BASIS_P;
252db2becc9SJeremy L Thompson   sTmp += ty * (max(BASIS_Q * BASIS_Q * BASIS_Q, (BASIS_Q * BASIS_Q * BASIS_P) + (BASIS_Q * BASIS_P * BASIS_P)));
253db2becc9SJeremy L Thompson 
254db2becc9SJeremy L Thompson   // read T
255db2becc9SJeremy L Thompson   if (ty == 0) {
256db2becc9SJeremy L Thompson     read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dinterp1d, sTinterp);
257db2becc9SJeremy L Thompson     read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dgrad1d, sTgrad);
258db2becc9SJeremy L Thompson   }
259db2becc9SJeremy L Thompson   __syncthreads();
260db2becc9SJeremy L Thompson 
261db2becc9SJeremy L Thompson   /* read U (idim = 0 for dU, i_DIM = 0 for rU) --
262db2becc9SJeremy L Thompson      there is a sync at the end of this function */
263db2becc9SJeremy L Thompson   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx);
264db2becc9SJeremy L Thompson   /* then first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) */
265db2becc9SJeremy L Thompson   magma_grad_3d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P, 0, 0, 0, true>(sTinterp, sTgrad, rU, rV, tx, rTmp, sTmp);
266db2becc9SJeremy L Thompson   /* there is a sync at the end of magma_grad_3d_device */
267db2becc9SJeremy L Thompson 
268db2becc9SJeremy L Thompson   /* read U (idim = 1 for dU, i_DIM = 0 for rU) --
269db2becc9SJeremy L Thompson      there is a sync at the end of this function */
270db2becc9SJeremy L Thompson   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (1 * dstrdU), cstrdU, rU, sTmp, tx);
271db2becc9SJeremy L Thompson   /* then second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) */
272db2becc9SJeremy L Thompson   magma_grad_3d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P, 1, 0, 0, true>(sTinterp, sTgrad, rU, rV, tx, rTmp, sTmp);
273db2becc9SJeremy L Thompson   /* there is a sync at the end of magma_grad_3d_device */
274db2becc9SJeremy L Thompson 
275db2becc9SJeremy L Thompson   /* read U (idim = 2 for dU, i_DIM = 0 for rU) --
276db2becc9SJeremy L Thompson      there is a sync at the end of this function */
277db2becc9SJeremy L Thompson   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (2 * dstrdU), cstrdU, rU, sTmp, tx);
278db2becc9SJeremy L Thompson   /* then third call (i_DIM = 2, i_DIM_U = 0, i_DIM_V = 0) */
279db2becc9SJeremy L Thompson   magma_grad_3d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P, 2, 0, 0, true>(sTinterp, sTgrad, rU, rV, tx, rTmp, sTmp);
280db2becc9SJeremy L Thompson   /* there is a sync at the end of magma_grad_3d_device */
281db2becc9SJeremy L Thompson 
282db2becc9SJeremy L Thompson   // sum into V
283db2becc9SJeremy L Thompson   sum_V_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV + (0 * dstrdV), cstrdV, rV, tx);
284db2becc9SJeremy L Thompson }
285