xref: /libCEED/rust/libceed-sys/c-src/include/ceed/jit-source/magma/magma-basis-grad-3d.h (revision db2becc9f302fe8eb3a32ace50ce3f3a5d42e6c4)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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 
113c1e2affSSebastian Grimberg #include "magma-common-tensor.h"
123c1e2affSSebastian Grimberg 
13f80f4a74SSebastian Grimberg // macros to abstract access of shared memory and reg. file
143c1e2affSSebastian Grimberg #define sT(i, j) sT[(j) * P + (i)]
15f80f4a74SSebastian Grimberg #define sTmp(i, j, ldw) sTmp[(j) * (ldw) + (i)]
16f80f4a74SSebastian Grimberg #define sTmp2(i, j, ldw) sTmp2[(j) * (ldw) + (i)]
17f80f4a74SSebastian Grimberg 
189e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
197132caa0SSebastian Grimberg // Helper function to add or set into V
207132caa0SSebastian Grimberg template <typename T, bool Add>
217132caa0SSebastian Grimberg struct magma_grad_3d_device_accumulate;
227132caa0SSebastian Grimberg 
237132caa0SSebastian Grimberg template <typename T>
247132caa0SSebastian Grimberg struct magma_grad_3d_device_accumulate<T, true> {
257132caa0SSebastian Grimberg   static __device__ __inline__ void op(T &rV, const T &rTmp) { rV += rTmp; }
267132caa0SSebastian Grimberg };
277132caa0SSebastian Grimberg 
287132caa0SSebastian Grimberg template <typename T>
297132caa0SSebastian Grimberg struct magma_grad_3d_device_accumulate<T, false> {
307132caa0SSebastian Grimberg   static __device__ __inline__ void op(T &rV, const T &rTmp) { rV = rTmp; }
317132caa0SSebastian Grimberg };
327132caa0SSebastian Grimberg 
337132caa0SSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
34f80f4a74SSebastian Grimberg // grad basis action (3D)
35f80f4a74SSebastian Grimberg // This function is called three times at a higher level for 3D
363c1e2affSSebastian Grimberg // DIM_U   -- for the size of rU[DIM_U * NUM_COMP * MAX_P_Q]
373c1e2affSSebastian Grimberg // DIM_V   -- for the size of rV[DIM_V * NUM_COMP * MAX_P_Q]
383c1e2affSSebastian Grimberg // i_DIM   -- the index of the outermost loop over dimensions in grad
393c1e2affSSebastian Grimberg // i_DIM_U -- which dim index of rU is accessed (always 0 for notrans, 0, 1, or 2 for trans)
403c1e2affSSebastian Grimberg // i_DIM_V -- which dim index of rV is accessed (0, 1, or 2 for notrans, always 0 for trans)
417132caa0SSebastian 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>
423c1e2affSSebastian Grimberg static __device__ __inline__ void magma_grad_3d_device(const T *sTinterp, const T *sTgrad, T rU[DIM_U][NUM_COMP][rU_SIZE],
437132caa0SSebastian Grimberg                                                        T rV[DIM_V][NUM_COMP][rV_SIZE], const int tx, T rTmp, T *swork) {
44f80f4a74SSebastian Grimberg   // Assumptions
453c1e2affSSebastian Grimberg   // 0. This device routine applies grad for one dim only (i_DIM), so it should be thrice for 3D
463c1e2affSSebastian Grimberg   // 1. 1D threads of size max(P,Q)^2
473c1e2affSSebastian Grimberg   // 2. input:  rU[DIM_U x NUM_COMP x rU_SIZE] in registers (per thread)
483c1e2affSSebastian Grimberg   // 3. output: rV[DIM_V x NUM_COMP x rV_SIZE] in registers (per thread)
49f80f4a74SSebastian Grimberg   // 4. Three products per each (dim,component) pair
503c1e2affSSebastian Grimberg   //  4.1 Batch P^2 of (1xP) matrices times (PxQ) matrix => Batch P^2 of (1xQ) matrices
513c1e2affSSebastian Grimberg   //  4.2 Batch P   of (QxP) matrices times (PxQ) matrix => Batch P   of (QxQ) matrices
523c1e2affSSebastian Grimberg   //  4.3 Batch 1   of (Q^2xP_) matrix times (PxQ) matrix => (Q^2xQ_) matrix
53f80f4a74SSebastian Grimberg   // 6. Each thread computes one row of the output of each product
54f80f4a74SSebastian Grimberg   // 7. Sync is recommended before and after the call
55f80f4a74SSebastian Grimberg 
56f80f4a74SSebastian Grimberg   T *sW1 = swork;
573c1e2affSSebastian Grimberg   T *sW2 = sW1 + P * P * Q;
583c1e2affSSebastian Grimberg   for (int comp = 0; comp < NUM_COMP; comp++) {
593c1e2affSSebastian Grimberg     // Batch P^2 of (1xP) matrices [reg] times (PxQ) matrix [shmem] => Batch P^2 of (1xQ) matrices [shmem]
603c1e2affSSebastian Grimberg     if (tx < (P * P)) {
61f80f4a74SSebastian Grimberg       const int batchid = tx;
62f80f4a74SSebastian Grimberg       const int sld     = 1;
633c1e2affSSebastian Grimberg       const T  *sT      = (i_DIM == 0) ? sTgrad : sTinterp;
643c1e2affSSebastian Grimberg       T        *sTmp    = sW1 + batchid * (1 * Q);
653c1e2affSSebastian Grimberg       for (int j = 0; j < Q; j++) {
66f80f4a74SSebastian Grimberg         rTmp = 0.0;
673c1e2affSSebastian Grimberg         for (int i = 0; i < P; i++) {
683c1e2affSSebastian Grimberg           rTmp += rU[i_DIM_U][comp][i] * sT(i, j);
69f80f4a74SSebastian Grimberg         }
70f80f4a74SSebastian Grimberg         sTmp(0, j, sld) = rTmp;
71f80f4a74SSebastian Grimberg       }
723c1e2affSSebastian Grimberg     }  // end of: if (tx < P*P)
73f80f4a74SSebastian Grimberg     __syncthreads();
74f80f4a74SSebastian Grimberg 
753c1e2affSSebastian Grimberg     // Batch P of (QxP) matrices [shmem] times (PxQ) matrix [shmem] => Batch P of (QxQ) matrices [reg]
763c1e2affSSebastian Grimberg     if (tx < (P * Q)) {
773c1e2affSSebastian Grimberg       const int batchid = tx / Q;
783c1e2affSSebastian Grimberg       const int tx_     = tx % Q;
793c1e2affSSebastian Grimberg       const int sld     = Q;
803c1e2affSSebastian Grimberg       const T  *sT      = (i_DIM == 1) ? sTgrad : sTinterp;
813c1e2affSSebastian Grimberg       T        *sTmp    = sW1 + batchid * (Q * P);  // sTmp is input
823c1e2affSSebastian Grimberg       T        *sTmp2   = sW2 + batchid * (Q * Q);  // sTmp2 is output
833c1e2affSSebastian Grimberg       for (int j = 0; j < Q; j++) {
84f80f4a74SSebastian Grimberg         rTmp = 0.0;
853c1e2affSSebastian Grimberg         for (int i = 0; i < P; i++) {
86f80f4a74SSebastian Grimberg           rTmp += sTmp(tx_, i, sld) * sT(i, j);
87f80f4a74SSebastian Grimberg         }
88f80f4a74SSebastian Grimberg         sTmp2(tx_, j, sld) = rTmp;
89f80f4a74SSebastian Grimberg       }
90f80f4a74SSebastian Grimberg     }
91f80f4a74SSebastian Grimberg     __syncthreads();
92f80f4a74SSebastian Grimberg 
933c1e2affSSebastian Grimberg     // Batch 1 of (Q^2xP_) matrices [shmem] times (PxQ) matrix [shmem] => Batch 1 of (Q^2xQ_) matrices [reg]
943c1e2affSSebastian Grimberg     if (tx < (Q * Q)) {
953c1e2affSSebastian Grimberg       // No need to declare batchid = (tx  / Q^2) = always zero
963c1e2affSSebastian Grimberg       // No need to declare tx_     = (tx_ % Q^2) = always tx
973c1e2affSSebastian Grimberg       const int sld  = Q * Q;
983c1e2affSSebastian Grimberg       const T  *sT   = (i_DIM == 2) ? sTgrad : sTinterp;
99f80f4a74SSebastian Grimberg       T        *sTmp = sW2;  // sTmp is input
1003c1e2affSSebastian Grimberg       for (int j = 0; j < Q; j++) {
101f80f4a74SSebastian Grimberg         rTmp = 0.0;
1023c1e2affSSebastian Grimberg         for (int i = 0; i < P; i++) {
103f80f4a74SSebastian Grimberg           rTmp += sTmp(tx, i, sld) * sT(i, j);
104f80f4a74SSebastian Grimberg         }
1057132caa0SSebastian Grimberg         magma_grad_3d_device_accumulate<T, ADD>::op(rV[i_DIM_V][comp][j], rTmp);
106f80f4a74SSebastian Grimberg       }
107f80f4a74SSebastian Grimberg     }
108f80f4a74SSebastian Grimberg     __syncthreads();
1093c1e2affSSebastian Grimberg   }  // loop over NUM_COMP
110f80f4a74SSebastian Grimberg }
111f80f4a74SSebastian Grimberg 
1129e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
1133c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q *BASIS_MAX_P_Q, MAGMA_MAXTHREADS_3D)) __global__
114f80f4a74SSebastian Grimberg     void magma_gradn_3d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU,
115f80f4a74SSebastian Grimberg                                const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) {
116f80f4a74SSebastian Grimberg   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
117f80f4a74SSebastian Grimberg 
118f80f4a74SSebastian Grimberg   const int tx      = threadIdx.x;
119f80f4a74SSebastian Grimberg   const int ty      = threadIdx.y;
120f80f4a74SSebastian Grimberg   const int elem_id = (blockIdx.x * blockDim.y) + ty;
121f80f4a74SSebastian Grimberg 
122f80f4a74SSebastian Grimberg   if (elem_id >= nelem) return;
123f80f4a74SSebastian Grimberg 
1243c1e2affSSebastian Grimberg   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // here DIM_U = 1, but might be different for a fused operator
1253c1e2affSSebastian Grimberg   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // here DIM_V = 1, but might be different for a fused operator
126f80f4a74SSebastian Grimberg   CeedScalar rTmp                           = 0.0;
127f80f4a74SSebastian Grimberg 
128f80f4a74SSebastian Grimberg   // shift global memory pointers by elem stride
129f80f4a74SSebastian Grimberg   dU += elem_id * estrdU;
130f80f4a74SSebastian Grimberg   dV += elem_id * estrdV;
131f80f4a74SSebastian Grimberg 
132f80f4a74SSebastian Grimberg   // assign shared memory pointers
1333c1e2affSSebastian Grimberg   CeedScalar *sTinterp = (CeedScalar *)shared_data;
1343c1e2affSSebastian Grimberg   CeedScalar *sTgrad   = sTinterp + BASIS_P * BASIS_Q;
1353c1e2affSSebastian Grimberg   CeedScalar *sTmp     = sTgrad + BASIS_P * BASIS_Q;
1363c1e2affSSebastian Grimberg   sTmp += ty * (max(BASIS_P * BASIS_P * BASIS_P, (BASIS_P * BASIS_P * BASIS_Q) + (BASIS_P * BASIS_Q * BASIS_Q)));
137f80f4a74SSebastian Grimberg 
138f80f4a74SSebastian Grimberg   // read T
139f80f4a74SSebastian Grimberg   if (ty == 0) {
1409e0c01faSSebastian Grimberg     read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dinterp1d, sTinterp);
1419e0c01faSSebastian Grimberg     read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dgrad1d, sTgrad);
142f80f4a74SSebastian Grimberg   }
143f80f4a74SSebastian Grimberg   __syncthreads();
144f80f4a74SSebastian Grimberg 
1453c1e2affSSebastian Grimberg   /* read U (idim = 0 for dU, i_DIM = 0 for rU) --
146f80f4a74SSebastian Grimberg      there is a sync at the end of this function */
1479e0c01faSSebastian Grimberg   read_U_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx);
148f80f4a74SSebastian Grimberg 
1493c1e2affSSebastian Grimberg   /* first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) --
150f80f4a74SSebastian Grimberg      output from rV[0][][] into dV (idim = 0) */
1517132caa0SSebastian 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,
1527132caa0SSebastian Grimberg                                                                                                              sTmp);
153f80f4a74SSebastian Grimberg   /* there is a sync at the end of magma_grad_3d_device */
1549e0c01faSSebastian Grimberg   write_V_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (0 * dstrdV), cstrdV, rV, tx);
155f80f4a74SSebastian Grimberg 
1563c1e2affSSebastian Grimberg   /* second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) --
157f80f4a74SSebastian Grimberg      output from rV[0][][] into dV (idim = 1) */
1587132caa0SSebastian 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,
1597132caa0SSebastian Grimberg                                                                                                              sTmp);
160f80f4a74SSebastian Grimberg   /* there is a sync at the end of magma_grad_3d_device */
1619e0c01faSSebastian Grimberg   write_V_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (1 * dstrdV), cstrdV, rV, tx);
162f80f4a74SSebastian Grimberg 
1633c1e2affSSebastian Grimberg   /* third call (i_DIM = 2, i_DIM_U = 0, i_DIM_V = 0) --
164f80f4a74SSebastian Grimberg      output from rV[0][][] into dV (idim = 2) */
1657132caa0SSebastian 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,
1667132caa0SSebastian Grimberg                                                                                                              sTmp);
167f80f4a74SSebastian Grimberg   /* there is a sync at the end of magma_grad_3d_device */
1689e0c01faSSebastian Grimberg   write_V_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (2 * dstrdV), cstrdV, rV, tx);
169f80f4a74SSebastian Grimberg }
170f80f4a74SSebastian Grimberg 
1719e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
1723c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q *BASIS_MAX_P_Q, MAGMA_MAXTHREADS_3D)) __global__
173f80f4a74SSebastian Grimberg     void magma_gradt_3d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU,
174f80f4a74SSebastian Grimberg                                const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) {
175f80f4a74SSebastian Grimberg   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
176f80f4a74SSebastian Grimberg 
177f80f4a74SSebastian Grimberg   const int tx      = threadIdx.x;
178f80f4a74SSebastian Grimberg   const int ty      = threadIdx.y;
179f80f4a74SSebastian Grimberg   const int elem_id = (blockIdx.x * blockDim.y) + ty;
180f80f4a74SSebastian Grimberg 
181f80f4a74SSebastian Grimberg   if (elem_id >= nelem) return;
182f80f4a74SSebastian Grimberg 
1833c1e2affSSebastian Grimberg   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // here DIM_U = 1, but might be different for a fused operator
1843c1e2affSSebastian Grimberg   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // here DIM_V = 1, but might be different for a fused operator
185f80f4a74SSebastian Grimberg   CeedScalar rTmp                           = 0.0;
186f80f4a74SSebastian Grimberg 
187f80f4a74SSebastian Grimberg   // shift global memory pointers by elem stride
188f80f4a74SSebastian Grimberg   dU += elem_id * estrdU;
189f80f4a74SSebastian Grimberg   dV += elem_id * estrdV;
190f80f4a74SSebastian Grimberg 
191f80f4a74SSebastian Grimberg   // assign shared memory pointers
1923c1e2affSSebastian Grimberg   CeedScalar *sTinterp = (CeedScalar *)shared_data;
1933c1e2affSSebastian Grimberg   CeedScalar *sTgrad   = sTinterp + BASIS_Q * BASIS_P;
1943c1e2affSSebastian Grimberg   CeedScalar *sTmp     = sTgrad + BASIS_Q * BASIS_P;
1953c1e2affSSebastian Grimberg   sTmp += ty * (max(BASIS_Q * BASIS_Q * BASIS_Q, (BASIS_Q * BASIS_Q * BASIS_P) + (BASIS_Q * BASIS_P * BASIS_P)));
196f80f4a74SSebastian Grimberg 
197f80f4a74SSebastian Grimberg   // read T
198f80f4a74SSebastian Grimberg   if (ty == 0) {
1999e0c01faSSebastian Grimberg     read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dinterp1d, sTinterp);
2009e0c01faSSebastian Grimberg     read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dgrad1d, sTgrad);
201f80f4a74SSebastian Grimberg   }
202f80f4a74SSebastian Grimberg   __syncthreads();
203f80f4a74SSebastian Grimberg 
2043c1e2affSSebastian Grimberg   /* read U (idim = 0 for dU, i_DIM = 0 for rU) --
205f80f4a74SSebastian Grimberg      there is a sync at the end of this function */
2069e0c01faSSebastian Grimberg   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx);
2073c1e2affSSebastian Grimberg   /* then first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) */
2087132caa0SSebastian 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);
209f80f4a74SSebastian Grimberg   /* there is a sync at the end of magma_grad_3d_device */
210f80f4a74SSebastian Grimberg 
2113c1e2affSSebastian Grimberg   /* read U (idim = 1 for dU, i_DIM = 0 for rU) --
212f80f4a74SSebastian Grimberg      there is a sync at the end of this function */
2139e0c01faSSebastian Grimberg   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (1 * dstrdU), cstrdU, rU, sTmp, tx);
2143c1e2affSSebastian Grimberg   /* then second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) */
2157132caa0SSebastian 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);
216f80f4a74SSebastian Grimberg   /* there is a sync at the end of magma_grad_3d_device */
217f80f4a74SSebastian Grimberg 
2183c1e2affSSebastian Grimberg   /* read U (idim = 2 for dU, i_DIM = 0 for rU) --
219f80f4a74SSebastian Grimberg      there is a sync at the end of this function */
2209e0c01faSSebastian Grimberg   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (2 * dstrdU), cstrdU, rU, sTmp, tx);
2213c1e2affSSebastian Grimberg   /* then third call (i_DIM = 2, i_DIM_U = 0, i_DIM_V = 0) */
2227132caa0SSebastian 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);
223f80f4a74SSebastian Grimberg   /* there is a sync at the end of magma_grad_3d_device */
224f80f4a74SSebastian Grimberg 
225f80f4a74SSebastian Grimberg   // write V
2269e0c01faSSebastian Grimberg   write_V_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV + (0 * dstrdV), cstrdV, rV, tx);
227f80f4a74SSebastian Grimberg }
228*db2becc9SJeremy L Thompson 
229*db2becc9SJeremy L Thompson ////////////////////////////////////////////////////////////////////////////////
230*db2becc9SJeremy L Thompson extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q *BASIS_MAX_P_Q, MAGMA_MAXTHREADS_3D)) __global__
231*db2becc9SJeremy L Thompson     void magma_gradta_3d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU,
232*db2becc9SJeremy L Thompson                                 const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) {
233*db2becc9SJeremy L Thompson   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
234*db2becc9SJeremy L Thompson 
235*db2becc9SJeremy L Thompson   const int tx      = threadIdx.x;
236*db2becc9SJeremy L Thompson   const int ty      = threadIdx.y;
237*db2becc9SJeremy L Thompson   const int elem_id = (blockIdx.x * blockDim.y) + ty;
238*db2becc9SJeremy L Thompson 
239*db2becc9SJeremy L Thompson   if (elem_id >= nelem) return;
240*db2becc9SJeremy L Thompson 
241*db2becc9SJeremy L Thompson   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // here DIM_U = 1, but might be different for a fused operator
242*db2becc9SJeremy L Thompson   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // here DIM_V = 1, but might be different for a fused operator
243*db2becc9SJeremy L Thompson   CeedScalar rTmp                           = 0.0;
244*db2becc9SJeremy L Thompson 
245*db2becc9SJeremy L Thompson   // shift global memory pointers by elem stride
246*db2becc9SJeremy L Thompson   dU += elem_id * estrdU;
247*db2becc9SJeremy L Thompson   dV += elem_id * estrdV;
248*db2becc9SJeremy L Thompson 
249*db2becc9SJeremy L Thompson   // assign shared memory pointers
250*db2becc9SJeremy L Thompson   CeedScalar *sTinterp = (CeedScalar *)shared_data;
251*db2becc9SJeremy L Thompson   CeedScalar *sTgrad   = sTinterp + BASIS_Q * BASIS_P;
252*db2becc9SJeremy L Thompson   CeedScalar *sTmp     = sTgrad + BASIS_Q * BASIS_P;
253*db2becc9SJeremy L Thompson   sTmp += ty * (max(BASIS_Q * BASIS_Q * BASIS_Q, (BASIS_Q * BASIS_Q * BASIS_P) + (BASIS_Q * BASIS_P * BASIS_P)));
254*db2becc9SJeremy L Thompson 
255*db2becc9SJeremy L Thompson   // read T
256*db2becc9SJeremy L Thompson   if (ty == 0) {
257*db2becc9SJeremy L Thompson     read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dinterp1d, sTinterp);
258*db2becc9SJeremy L Thompson     read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dgrad1d, sTgrad);
259*db2becc9SJeremy L Thompson   }
260*db2becc9SJeremy L Thompson   __syncthreads();
261*db2becc9SJeremy L Thompson 
262*db2becc9SJeremy L Thompson   /* read U (idim = 0 for dU, i_DIM = 0 for rU) --
263*db2becc9SJeremy L Thompson      there is a sync at the end of this function */
264*db2becc9SJeremy L Thompson   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx);
265*db2becc9SJeremy L Thompson   /* then first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) */
266*db2becc9SJeremy 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);
267*db2becc9SJeremy L Thompson   /* there is a sync at the end of magma_grad_3d_device */
268*db2becc9SJeremy L Thompson 
269*db2becc9SJeremy L Thompson   /* read U (idim = 1 for dU, i_DIM = 0 for rU) --
270*db2becc9SJeremy L Thompson      there is a sync at the end of this function */
271*db2becc9SJeremy L Thompson   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (1 * dstrdU), cstrdU, rU, sTmp, tx);
272*db2becc9SJeremy L Thompson   /* then second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) */
273*db2becc9SJeremy 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);
274*db2becc9SJeremy L Thompson   /* there is a sync at the end of magma_grad_3d_device */
275*db2becc9SJeremy L Thompson 
276*db2becc9SJeremy L Thompson   /* read U (idim = 2 for dU, i_DIM = 0 for rU) --
277*db2becc9SJeremy L Thompson      there is a sync at the end of this function */
278*db2becc9SJeremy L Thompson   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (2 * dstrdU), cstrdU, rU, sTmp, tx);
279*db2becc9SJeremy L Thompson   /* then third call (i_DIM = 2, i_DIM_U = 0, i_DIM_V = 0) */
280*db2becc9SJeremy 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);
281*db2becc9SJeremy L Thompson   /* there is a sync at the end of magma_grad_3d_device */
282*db2becc9SJeremy L Thompson 
283*db2becc9SJeremy L Thompson   // sum into V
284*db2becc9SJeremy L Thompson   sum_V_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV + (0 * dstrdV), cstrdV, rV, tx);
285*db2becc9SJeremy L Thompson }
286