xref: /libCEED/include/ceed/jit-source/magma/magma-basis-interp-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 interpolation 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 
169e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
17f80f4a74SSebastian Grimberg // interp basis action (3D)
183c1e2affSSebastian Grimberg template <typename T, int DIM_U, int DIM_V, int NUM_COMP, int P, int Q, int rU_SIZE, int rV_SIZE>
magma_interp_3d_device(const T * sT,T rU[DIM_U][NUM_COMP][rU_SIZE],T rV[DIM_V][NUM_COMP][rV_SIZE],const int tx,T rTmp[Q],T * swork)197132caa0SSebastian Grimberg static __device__ __inline__ void magma_interp_3d_device(const T *sT, T rU[DIM_U][NUM_COMP][rU_SIZE], T rV[DIM_V][NUM_COMP][rV_SIZE], const int tx,
207132caa0SSebastian Grimberg                                                          T rTmp[Q], T *swork) {
21f80f4a74SSebastian Grimberg   // Assumptions
223c1e2affSSebastian Grimberg   // 1. 1D threads of size max(P,Q)^2
233c1e2affSSebastian Grimberg   // 2. input:  rU[DIM_U x NUM_COMP x rU_SIZE] in registers (per thread)
243c1e2affSSebastian Grimberg   // 3. output: rV[DIM_V x NUM_COMP x rV_SIZE] in registers (per thread)
25f80f4a74SSebastian Grimberg   // 4. Three products per component
263c1e2affSSebastian Grimberg   //  4.1 Batch P^2 of (1xP) matrices times (PxQ) matrix => Batch P^2 of (1xQ) matrices
273c1e2affSSebastian Grimberg   //  4.2 Batch P   of (QxP) matrices times (PxQ) matrix => Batch P   of (QxQ) matrices
283c1e2affSSebastian Grimberg   //  4.3 Batch 1   of (Q^2xP_) matrix times (PxQ) matrix => (Q^2xQ_) matrix
29f80f4a74SSebastian Grimberg   // 5. Each thread computes one row of the output of each product
30f80f4a74SSebastian Grimberg   // 6. Sync is recommended before and after the call
31f80f4a74SSebastian Grimberg 
323c1e2affSSebastian Grimberg   for (int comp = 0; comp < NUM_COMP; comp++) {
333c1e2affSSebastian Grimberg     // Batch P^2 of (1xP) matrices [reg] times (PxQ) matrix [shmem] => Batch P^2 of (1xQ) matrices [shmem]
343c1e2affSSebastian Grimberg     if (tx < (P * P)) {
35f80f4a74SSebastian Grimberg       const int batchid = tx;
36f80f4a74SSebastian Grimberg       const int sld     = 1;
373c1e2affSSebastian Grimberg       T        *sTmp    = swork + batchid * (1 * Q);
383c1e2affSSebastian Grimberg       for (int j = 0; j < Q; j++) {
39f80f4a74SSebastian Grimberg         rTmp[0] = 0.0;
403c1e2affSSebastian Grimberg         for (int i = 0; i < P; i++) {
413c1e2affSSebastian Grimberg           rTmp[0] += rU[0][comp][i] * sT(i, j);
42f80f4a74SSebastian Grimberg         }
43f80f4a74SSebastian Grimberg         sTmp(0, j, sld) = rTmp[0];
44f80f4a74SSebastian Grimberg       }
453c1e2affSSebastian Grimberg     }  // end of: if (tx < P*P)
46f80f4a74SSebastian Grimberg     __syncthreads();
47f80f4a74SSebastian Grimberg 
483c1e2affSSebastian Grimberg     // Batch P of (QxP) matrices [shmem] times (PxQ) matrix [shmem] => Batch P of (QxQ) matrices [reg]
493c1e2affSSebastian Grimberg     if (tx < (P * Q)) {
503c1e2affSSebastian Grimberg       const int batchid = tx / Q;
513c1e2affSSebastian Grimberg       const int tx_     = tx % Q;
523c1e2affSSebastian Grimberg       const int sld     = Q;
533c1e2affSSebastian Grimberg       T        *sTmp    = swork + batchid * (Q * P);  // sTmp is input
543c1e2affSSebastian Grimberg       for (int j = 0; j < Q; j++) {
55f80f4a74SSebastian Grimberg         rTmp[j] = 0.0;
563c1e2affSSebastian Grimberg         for (int i = 0; i < P; i++) {
57f80f4a74SSebastian Grimberg           rTmp[j] += sTmp(tx_, i, sld) * sT(i, j);
58f80f4a74SSebastian Grimberg         }
59f80f4a74SSebastian Grimberg       }
60f80f4a74SSebastian Grimberg     }
61f80f4a74SSebastian Grimberg     __syncthreads();
62f80f4a74SSebastian Grimberg 
633c1e2affSSebastian Grimberg     // write rTmp[] into shmem as batch P of QxQ matrices
643c1e2affSSebastian Grimberg     if (tx < (P * Q)) {
653c1e2affSSebastian Grimberg       const int batchid = tx / Q;
663c1e2affSSebastian Grimberg       const int tx_     = tx % Q;
673c1e2affSSebastian Grimberg       const int sld     = Q;
683c1e2affSSebastian Grimberg       T        *sTmp    = swork + batchid * (Q * Q);
693c1e2affSSebastian Grimberg       for (int j = 0; j < Q; j++) {
70f80f4a74SSebastian Grimberg         sTmp(tx_, j, sld) = rTmp[j];
71f80f4a74SSebastian Grimberg       }
72f80f4a74SSebastian Grimberg     }
73f80f4a74SSebastian Grimberg     __syncthreads();
74f80f4a74SSebastian Grimberg 
753c1e2affSSebastian Grimberg     // Batch 1 of (Q^2xP_) matrices [shmem] times (PxQ) matrix [shmem] => Batch 1 of (Q^2xQ_) matrices [reg]
763c1e2affSSebastian Grimberg     if (tx < (Q * Q)) {
773c1e2affSSebastian Grimberg       // No need to declare batchid = (tx  / Q^2) = always zero
783c1e2affSSebastian Grimberg       // No need to declare tx_     = (tx_ % Q^2) = always tx
793c1e2affSSebastian Grimberg       const int sld  = Q * Q;
80f80f4a74SSebastian Grimberg       T        *sTmp = swork;
813c1e2affSSebastian Grimberg       for (int j = 0; j < Q; j++) {
82f80f4a74SSebastian Grimberg         rTmp[0] = 0.0;
833c1e2affSSebastian Grimberg         for (int i = 0; i < P; i++) {
84f80f4a74SSebastian Grimberg           rTmp[0] += sTmp(tx, i, sld) * sT(i, j);
85f80f4a74SSebastian Grimberg         }
863c1e2affSSebastian Grimberg         rV[0][comp][j] += rTmp[0];
87f80f4a74SSebastian Grimberg       }
88f80f4a74SSebastian Grimberg     }
89f80f4a74SSebastian Grimberg     __syncthreads();
90f80f4a74SSebastian Grimberg   }
91f80f4a74SSebastian Grimberg }
92f80f4a74SSebastian Grimberg 
939e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
__launch_bounds__(MAGMA_BASIS_BOUNDS (BASIS_MAX_P_Q * BASIS_MAX_P_Q,MAGMA_MAXTHREADS_3D))943c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q *BASIS_MAX_P_Q, MAGMA_MAXTHREADS_3D)) __global__
95f80f4a74SSebastian Grimberg     void magma_interpn_3d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV,
96f80f4a74SSebastian Grimberg                                  const int cstrdV, const int nelem) {
97f80f4a74SSebastian Grimberg   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
98f80f4a74SSebastian Grimberg 
99f80f4a74SSebastian Grimberg   const int tx      = threadIdx.x;
100f80f4a74SSebastian Grimberg   const int ty      = threadIdx.y;
101f80f4a74SSebastian Grimberg   const int elem_id = (blockIdx.x * blockDim.y) + ty;
102f80f4a74SSebastian Grimberg 
103f80f4a74SSebastian Grimberg   if (elem_id >= nelem) return;
104f80f4a74SSebastian Grimberg 
1053c1e2affSSebastian Grimberg   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
1063c1e2affSSebastian Grimberg   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
1073c1e2affSSebastian Grimberg   CeedScalar rTmp[BASIS_Q]                  = {0.0};
108f80f4a74SSebastian Grimberg 
109f80f4a74SSebastian Grimberg   // shift global memory pointers by elem stride
110f80f4a74SSebastian Grimberg   dU += elem_id * estrdU;
111f80f4a74SSebastian Grimberg   dV += elem_id * estrdV;
112f80f4a74SSebastian Grimberg 
113f80f4a74SSebastian Grimberg   // assign shared memory pointers
1143c1e2affSSebastian Grimberg   CeedScalar *sT   = (CeedScalar *)shared_data;
1153c1e2affSSebastian Grimberg   CeedScalar *sTmp = sT + BASIS_P * BASIS_Q;
1163c1e2affSSebastian Grimberg   sTmp += ty * (max(BASIS_P * BASIS_P * BASIS_MAX_P_Q, BASIS_P * BASIS_Q * BASIS_Q));
117f80f4a74SSebastian Grimberg 
118f80f4a74SSebastian Grimberg   // read T
119f80f4a74SSebastian Grimberg   if (ty == 0) {
1209e0c01faSSebastian Grimberg     read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dT, sT);
121f80f4a74SSebastian Grimberg   }
122f80f4a74SSebastian Grimberg 
1233c1e2affSSebastian Grimberg   // read U (idim = 0 for dU, i_DIM = 0 for rU, u_dimstride is always 0)
1249e0c01faSSebastian Grimberg   read_U_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dU, cstrdU, rU, sTmp, tx);
125f80f4a74SSebastian Grimberg   // there is a sync at the end of this function
126f80f4a74SSebastian Grimberg 
1277132caa0SSebastian Grimberg   magma_interp_3d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_P, BASIS_Q, BASIS_P, BASIS_Q>(sT, rU, rV, tx, rTmp, sTmp);
128f80f4a74SSebastian Grimberg   __syncthreads();
129f80f4a74SSebastian Grimberg 
130f80f4a74SSebastian Grimberg   // write V
1319e0c01faSSebastian Grimberg   write_V_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV, cstrdV, rV, tx);
132f80f4a74SSebastian Grimberg }
133f80f4a74SSebastian Grimberg 
1349e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
__launch_bounds__(MAGMA_BASIS_BOUNDS (BASIS_MAX_P_Q * BASIS_MAX_P_Q,MAGMA_MAXTHREADS_3D))1353c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q *BASIS_MAX_P_Q, MAGMA_MAXTHREADS_3D)) __global__
136f80f4a74SSebastian Grimberg     void magma_interpt_3d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV,
137f80f4a74SSebastian Grimberg                                  const int cstrdV, const int nelem) {
138f80f4a74SSebastian Grimberg   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
139f80f4a74SSebastian Grimberg 
140f80f4a74SSebastian Grimberg   const int tx      = threadIdx.x;
141f80f4a74SSebastian Grimberg   const int ty      = threadIdx.y;
142f80f4a74SSebastian Grimberg   const int elem_id = (blockIdx.x * blockDim.y) + ty;
143f80f4a74SSebastian Grimberg 
144f80f4a74SSebastian Grimberg   if (elem_id >= nelem) return;
145f80f4a74SSebastian Grimberg 
1463c1e2affSSebastian Grimberg   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
1473c1e2affSSebastian Grimberg   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
1483c1e2affSSebastian Grimberg   CeedScalar rTmp[BASIS_P]                  = {0.0};
149f80f4a74SSebastian Grimberg 
150f80f4a74SSebastian Grimberg   // shift global memory pointers by elem stride
151f80f4a74SSebastian Grimberg   dU += elem_id * estrdU;
152f80f4a74SSebastian Grimberg   dV += elem_id * estrdV;
153f80f4a74SSebastian Grimberg 
154f80f4a74SSebastian Grimberg   // assign shared memory pointers
1553c1e2affSSebastian Grimberg   CeedScalar *sT   = (CeedScalar *)shared_data;
1563c1e2affSSebastian Grimberg   CeedScalar *sTmp = sT + BASIS_Q * BASIS_P;
1573c1e2affSSebastian Grimberg   sTmp += ty * (max(BASIS_Q * BASIS_Q * BASIS_MAX_P_Q, BASIS_Q * BASIS_P * BASIS_P));
158f80f4a74SSebastian Grimberg 
159f80f4a74SSebastian Grimberg   // read T
160f80f4a74SSebastian Grimberg   if (ty == 0) {
1619e0c01faSSebastian Grimberg     read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dT, sT);
162f80f4a74SSebastian Grimberg   }
163f80f4a74SSebastian Grimberg 
1643c1e2affSSebastian Grimberg   // read U (idim = 0 for dU, i_DIM = 0 for rU, u_dimstride is always 0)
1659e0c01faSSebastian Grimberg   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU, cstrdU, rU, sTmp, tx);
166f80f4a74SSebastian Grimberg   // there is a sync at the end of this function
167f80f4a74SSebastian Grimberg 
1687132caa0SSebastian Grimberg   magma_interp_3d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P>(sT, rU, rV, tx, rTmp, sTmp);
169f80f4a74SSebastian Grimberg   __syncthreads();
170f80f4a74SSebastian Grimberg 
171f80f4a74SSebastian Grimberg   // write V
1729e0c01faSSebastian Grimberg   write_V_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV, cstrdV, rV, tx);
173f80f4a74SSebastian Grimberg }
174db2becc9SJeremy L Thompson 
175db2becc9SJeremy L Thompson ////////////////////////////////////////////////////////////////////////////////
__launch_bounds__(MAGMA_BASIS_BOUNDS (BASIS_MAX_P_Q * BASIS_MAX_P_Q,MAGMA_MAXTHREADS_3D))176db2becc9SJeremy L Thompson extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q *BASIS_MAX_P_Q, MAGMA_MAXTHREADS_3D)) __global__
177db2becc9SJeremy L Thompson     void magma_interpta_3d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV,
178db2becc9SJeremy L Thompson                                   const int cstrdV, const int nelem) {
179db2becc9SJeremy L Thompson   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
180db2becc9SJeremy L Thompson 
181db2becc9SJeremy L Thompson   const int tx      = threadIdx.x;
182db2becc9SJeremy L Thompson   const int ty      = threadIdx.y;
183db2becc9SJeremy L Thompson   const int elem_id = (blockIdx.x * blockDim.y) + ty;
184db2becc9SJeremy L Thompson 
185db2becc9SJeremy L Thompson   if (elem_id >= nelem) return;
186db2becc9SJeremy L Thompson 
187db2becc9SJeremy L Thompson   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
188db2becc9SJeremy L Thompson   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
189db2becc9SJeremy L Thompson   CeedScalar rTmp[BASIS_P]                  = {0.0};
190db2becc9SJeremy L Thompson 
191db2becc9SJeremy L Thompson   // shift global memory pointers by elem stride
192db2becc9SJeremy L Thompson   dU += elem_id * estrdU;
193db2becc9SJeremy L Thompson   dV += elem_id * estrdV;
194db2becc9SJeremy L Thompson 
195db2becc9SJeremy L Thompson   // assign shared memory pointers
196db2becc9SJeremy L Thompson   CeedScalar *sT   = (CeedScalar *)shared_data;
197db2becc9SJeremy L Thompson   CeedScalar *sTmp = sT + BASIS_Q * BASIS_P;
198db2becc9SJeremy L Thompson   sTmp += ty * (max(BASIS_Q * BASIS_Q * BASIS_MAX_P_Q, BASIS_Q * BASIS_P * BASIS_P));
199db2becc9SJeremy L Thompson 
200db2becc9SJeremy L Thompson   // read T
201db2becc9SJeremy L Thompson   if (ty == 0) {
202db2becc9SJeremy L Thompson     read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dT, sT);
203db2becc9SJeremy L Thompson   }
204db2becc9SJeremy L Thompson 
205db2becc9SJeremy L Thompson   // read U (idim = 0 for dU, i_DIM = 0 for rU, u_dimstride is always 0)
206db2becc9SJeremy L Thompson   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU, cstrdU, rU, sTmp, tx);
207db2becc9SJeremy L Thompson   // there is a sync at the end of this function
208db2becc9SJeremy L Thompson 
209db2becc9SJeremy L Thompson   magma_interp_3d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P>(sT, rU, rV, tx, rTmp, sTmp);
210db2becc9SJeremy L Thompson   __syncthreads();
211db2becc9SJeremy L Thompson 
212db2becc9SJeremy L Thompson   // sum into V
213db2becc9SJeremy L Thompson   sum_V_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV, cstrdV, rV, tx);
214db2becc9SJeremy L Thompson }
215