xref: /libCEED/include/ceed/jit-source/magma/magma-basis-interp-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 interpolation 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 
179e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
18f80f4a74SSebastian Grimberg // interp basis action (3D)
193c1e2affSSebastian Grimberg template <typename T, int DIM_U, int DIM_V, int NUM_COMP, int P, int Q, int rU_SIZE, int rV_SIZE>
207132caa0SSebastian 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,
217132caa0SSebastian Grimberg                                                          T rTmp[Q], T *swork) {
22f80f4a74SSebastian Grimberg   // Assumptions
233c1e2affSSebastian Grimberg   // 1. 1D threads of size max(P,Q)^2
243c1e2affSSebastian Grimberg   // 2. input:  rU[DIM_U x NUM_COMP x rU_SIZE] in registers (per thread)
253c1e2affSSebastian Grimberg   // 3. output: rV[DIM_V x NUM_COMP x rV_SIZE] in registers (per thread)
26f80f4a74SSebastian Grimberg   // 4. Three products per component
273c1e2affSSebastian Grimberg   //  4.1 Batch P^2 of (1xP) matrices times (PxQ) matrix => Batch P^2 of (1xQ) matrices
283c1e2affSSebastian Grimberg   //  4.2 Batch P   of (QxP) matrices times (PxQ) matrix => Batch P   of (QxQ) matrices
293c1e2affSSebastian Grimberg   //  4.3 Batch 1   of (Q^2xP_) matrix times (PxQ) matrix => (Q^2xQ_) matrix
30f80f4a74SSebastian Grimberg   // 5. Each thread computes one row of the output of each product
31f80f4a74SSebastian Grimberg   // 6. Sync is recommended before and after the call
32f80f4a74SSebastian Grimberg 
333c1e2affSSebastian Grimberg   for (int comp = 0; comp < NUM_COMP; comp++) {
343c1e2affSSebastian Grimberg     // Batch P^2 of (1xP) matrices [reg] times (PxQ) matrix [shmem] => Batch P^2 of (1xQ) matrices [shmem]
353c1e2affSSebastian Grimberg     if (tx < (P * P)) {
36f80f4a74SSebastian Grimberg       const int batchid = tx;
37f80f4a74SSebastian Grimberg       const int sld     = 1;
383c1e2affSSebastian Grimberg       T        *sTmp    = swork + batchid * (1 * Q);
393c1e2affSSebastian Grimberg       for (int j = 0; j < Q; j++) {
40f80f4a74SSebastian Grimberg         rTmp[0] = 0.0;
413c1e2affSSebastian Grimberg         for (int i = 0; i < P; i++) {
423c1e2affSSebastian Grimberg           rTmp[0] += rU[0][comp][i] * sT(i, j);
43f80f4a74SSebastian Grimberg         }
44f80f4a74SSebastian Grimberg         sTmp(0, j, sld) = rTmp[0];
45f80f4a74SSebastian Grimberg       }
463c1e2affSSebastian Grimberg     }  // end of: if (tx < P*P)
47f80f4a74SSebastian Grimberg     __syncthreads();
48f80f4a74SSebastian Grimberg 
493c1e2affSSebastian Grimberg     // Batch P of (QxP) matrices [shmem] times (PxQ) matrix [shmem] => Batch P of (QxQ) matrices [reg]
503c1e2affSSebastian Grimberg     if (tx < (P * Q)) {
513c1e2affSSebastian Grimberg       const int batchid = tx / Q;
523c1e2affSSebastian Grimberg       const int tx_     = tx % Q;
533c1e2affSSebastian Grimberg       const int sld     = Q;
543c1e2affSSebastian Grimberg       T        *sTmp    = swork + batchid * (Q * P);  // sTmp is input
553c1e2affSSebastian Grimberg       for (int j = 0; j < Q; j++) {
56f80f4a74SSebastian Grimberg         rTmp[j] = 0.0;
573c1e2affSSebastian Grimberg         for (int i = 0; i < P; i++) {
58f80f4a74SSebastian Grimberg           rTmp[j] += sTmp(tx_, i, sld) * sT(i, j);
59f80f4a74SSebastian Grimberg         }
60f80f4a74SSebastian Grimberg       }
61f80f4a74SSebastian Grimberg     }
62f80f4a74SSebastian Grimberg     __syncthreads();
63f80f4a74SSebastian Grimberg 
643c1e2affSSebastian Grimberg     // write rTmp[] into shmem as batch P of QxQ matrices
653c1e2affSSebastian Grimberg     if (tx < (P * Q)) {
663c1e2affSSebastian Grimberg       const int batchid = tx / Q;
673c1e2affSSebastian Grimberg       const int tx_     = tx % Q;
683c1e2affSSebastian Grimberg       const int sld     = Q;
693c1e2affSSebastian Grimberg       T        *sTmp    = swork + batchid * (Q * Q);
703c1e2affSSebastian Grimberg       for (int j = 0; j < Q; j++) {
71f80f4a74SSebastian Grimberg         sTmp(tx_, j, sld) = rTmp[j];
72f80f4a74SSebastian Grimberg       }
73f80f4a74SSebastian Grimberg     }
74f80f4a74SSebastian Grimberg     __syncthreads();
75f80f4a74SSebastian Grimberg 
763c1e2affSSebastian Grimberg     // Batch 1 of (Q^2xP_) matrices [shmem] times (PxQ) matrix [shmem] => Batch 1 of (Q^2xQ_) matrices [reg]
773c1e2affSSebastian Grimberg     if (tx < (Q * Q)) {
783c1e2affSSebastian Grimberg       // No need to declare batchid = (tx  / Q^2) = always zero
793c1e2affSSebastian Grimberg       // No need to declare tx_     = (tx_ % Q^2) = always tx
803c1e2affSSebastian Grimberg       const int sld  = Q * Q;
81f80f4a74SSebastian Grimberg       T        *sTmp = swork;
823c1e2affSSebastian Grimberg       for (int j = 0; j < Q; j++) {
83f80f4a74SSebastian Grimberg         rTmp[0] = 0.0;
843c1e2affSSebastian Grimberg         for (int i = 0; i < P; i++) {
85f80f4a74SSebastian Grimberg           rTmp[0] += sTmp(tx, i, sld) * sT(i, j);
86f80f4a74SSebastian Grimberg         }
873c1e2affSSebastian Grimberg         rV[0][comp][j] += rTmp[0];
88f80f4a74SSebastian Grimberg       }
89f80f4a74SSebastian Grimberg     }
90f80f4a74SSebastian Grimberg     __syncthreads();
91f80f4a74SSebastian Grimberg   }
92f80f4a74SSebastian Grimberg }
93f80f4a74SSebastian Grimberg 
949e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
953c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q *BASIS_MAX_P_Q, MAGMA_MAXTHREADS_3D)) __global__
96f80f4a74SSebastian Grimberg     void magma_interpn_3d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV,
97f80f4a74SSebastian Grimberg                                  const int cstrdV, const int nelem) {
98f80f4a74SSebastian Grimberg   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
99f80f4a74SSebastian Grimberg 
100f80f4a74SSebastian Grimberg   const int tx      = threadIdx.x;
101f80f4a74SSebastian Grimberg   const int ty      = threadIdx.y;
102f80f4a74SSebastian Grimberg   const int elem_id = (blockIdx.x * blockDim.y) + ty;
103f80f4a74SSebastian Grimberg 
104f80f4a74SSebastian Grimberg   if (elem_id >= nelem) return;
105f80f4a74SSebastian Grimberg 
1063c1e2affSSebastian Grimberg   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
1073c1e2affSSebastian Grimberg   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
1083c1e2affSSebastian Grimberg   CeedScalar rTmp[BASIS_Q]                  = {0.0};
109f80f4a74SSebastian Grimberg 
110f80f4a74SSebastian Grimberg   // shift global memory pointers by elem stride
111f80f4a74SSebastian Grimberg   dU += elem_id * estrdU;
112f80f4a74SSebastian Grimberg   dV += elem_id * estrdV;
113f80f4a74SSebastian Grimberg 
114f80f4a74SSebastian Grimberg   // assign shared memory pointers
1153c1e2affSSebastian Grimberg   CeedScalar *sT   = (CeedScalar *)shared_data;
1163c1e2affSSebastian Grimberg   CeedScalar *sTmp = sT + BASIS_P * BASIS_Q;
1173c1e2affSSebastian Grimberg   sTmp += ty * (max(BASIS_P * BASIS_P * BASIS_MAX_P_Q, BASIS_P * BASIS_Q * BASIS_Q));
118f80f4a74SSebastian Grimberg 
119f80f4a74SSebastian Grimberg   // read T
120f80f4a74SSebastian Grimberg   if (ty == 0) {
1219e0c01faSSebastian Grimberg     read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dT, sT);
122f80f4a74SSebastian Grimberg   }
123f80f4a74SSebastian Grimberg 
1243c1e2affSSebastian Grimberg   // read U (idim = 0 for dU, i_DIM = 0 for rU, u_dimstride is always 0)
1259e0c01faSSebastian Grimberg   read_U_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dU, cstrdU, rU, sTmp, tx);
126f80f4a74SSebastian Grimberg   // there is a sync at the end of this function
127f80f4a74SSebastian Grimberg 
1287132caa0SSebastian 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);
129f80f4a74SSebastian Grimberg   __syncthreads();
130f80f4a74SSebastian Grimberg 
131f80f4a74SSebastian Grimberg   // write V
1329e0c01faSSebastian Grimberg   write_V_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV, cstrdV, rV, tx);
133f80f4a74SSebastian Grimberg }
134f80f4a74SSebastian Grimberg 
1359e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
1363c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q *BASIS_MAX_P_Q, MAGMA_MAXTHREADS_3D)) __global__
137f80f4a74SSebastian Grimberg     void magma_interpt_3d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV,
138f80f4a74SSebastian Grimberg                                  const int cstrdV, const int nelem) {
139f80f4a74SSebastian Grimberg   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
140f80f4a74SSebastian Grimberg 
141f80f4a74SSebastian Grimberg   const int tx      = threadIdx.x;
142f80f4a74SSebastian Grimberg   const int ty      = threadIdx.y;
143f80f4a74SSebastian Grimberg   const int elem_id = (blockIdx.x * blockDim.y) + ty;
144f80f4a74SSebastian Grimberg 
145f80f4a74SSebastian Grimberg   if (elem_id >= nelem) return;
146f80f4a74SSebastian Grimberg 
1473c1e2affSSebastian Grimberg   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
1483c1e2affSSebastian Grimberg   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
1493c1e2affSSebastian Grimberg   CeedScalar rTmp[BASIS_P]                  = {0.0};
150f80f4a74SSebastian Grimberg 
151f80f4a74SSebastian Grimberg   // shift global memory pointers by elem stride
152f80f4a74SSebastian Grimberg   dU += elem_id * estrdU;
153f80f4a74SSebastian Grimberg   dV += elem_id * estrdV;
154f80f4a74SSebastian Grimberg 
155f80f4a74SSebastian Grimberg   // assign shared memory pointers
1563c1e2affSSebastian Grimberg   CeedScalar *sT   = (CeedScalar *)shared_data;
1573c1e2affSSebastian Grimberg   CeedScalar *sTmp = sT + BASIS_Q * BASIS_P;
1583c1e2affSSebastian Grimberg   sTmp += ty * (max(BASIS_Q * BASIS_Q * BASIS_MAX_P_Q, BASIS_Q * BASIS_P * BASIS_P));
159f80f4a74SSebastian Grimberg 
160f80f4a74SSebastian Grimberg   // read T
161f80f4a74SSebastian Grimberg   if (ty == 0) {
1629e0c01faSSebastian Grimberg     read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dT, sT);
163f80f4a74SSebastian Grimberg   }
164f80f4a74SSebastian Grimberg 
1653c1e2affSSebastian Grimberg   // read U (idim = 0 for dU, i_DIM = 0 for rU, u_dimstride is always 0)
1669e0c01faSSebastian Grimberg   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU, cstrdU, rU, sTmp, tx);
167f80f4a74SSebastian Grimberg   // there is a sync at the end of this function
168f80f4a74SSebastian Grimberg 
1697132caa0SSebastian 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);
170f80f4a74SSebastian Grimberg   __syncthreads();
171f80f4a74SSebastian Grimberg 
172f80f4a74SSebastian Grimberg   // write V
1739e0c01faSSebastian Grimberg   write_V_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV, cstrdV, rV, tx);
174f80f4a74SSebastian Grimberg }
175*db2becc9SJeremy L Thompson 
176*db2becc9SJeremy L Thompson ////////////////////////////////////////////////////////////////////////////////
177*db2becc9SJeremy L Thompson extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q *BASIS_MAX_P_Q, MAGMA_MAXTHREADS_3D)) __global__
178*db2becc9SJeremy L Thompson     void magma_interpta_3d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV,
179*db2becc9SJeremy L Thompson                                   const int cstrdV, const int nelem) {
180*db2becc9SJeremy L Thompson   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
181*db2becc9SJeremy L Thompson 
182*db2becc9SJeremy L Thompson   const int tx      = threadIdx.x;
183*db2becc9SJeremy L Thompson   const int ty      = threadIdx.y;
184*db2becc9SJeremy L Thompson   const int elem_id = (blockIdx.x * blockDim.y) + ty;
185*db2becc9SJeremy L Thompson 
186*db2becc9SJeremy L Thompson   if (elem_id >= nelem) return;
187*db2becc9SJeremy L Thompson 
188*db2becc9SJeremy L Thompson   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
189*db2becc9SJeremy L Thompson   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
190*db2becc9SJeremy L Thompson   CeedScalar rTmp[BASIS_P]                  = {0.0};
191*db2becc9SJeremy L Thompson 
192*db2becc9SJeremy L Thompson   // shift global memory pointers by elem stride
193*db2becc9SJeremy L Thompson   dU += elem_id * estrdU;
194*db2becc9SJeremy L Thompson   dV += elem_id * estrdV;
195*db2becc9SJeremy L Thompson 
196*db2becc9SJeremy L Thompson   // assign shared memory pointers
197*db2becc9SJeremy L Thompson   CeedScalar *sT   = (CeedScalar *)shared_data;
198*db2becc9SJeremy L Thompson   CeedScalar *sTmp = sT + BASIS_Q * BASIS_P;
199*db2becc9SJeremy L Thompson   sTmp += ty * (max(BASIS_Q * BASIS_Q * BASIS_MAX_P_Q, BASIS_Q * BASIS_P * BASIS_P));
200*db2becc9SJeremy L Thompson 
201*db2becc9SJeremy L Thompson   // read T
202*db2becc9SJeremy L Thompson   if (ty == 0) {
203*db2becc9SJeremy L Thompson     read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dT, sT);
204*db2becc9SJeremy L Thompson   }
205*db2becc9SJeremy L Thompson 
206*db2becc9SJeremy L Thompson   // read U (idim = 0 for dU, i_DIM = 0 for rU, u_dimstride is always 0)
207*db2becc9SJeremy L Thompson   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU, cstrdU, rU, sTmp, tx);
208*db2becc9SJeremy L Thompson   // there is a sync at the end of this function
209*db2becc9SJeremy L Thompson 
210*db2becc9SJeremy 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);
211*db2becc9SJeremy L Thompson   __syncthreads();
212*db2becc9SJeremy L Thompson 
213*db2becc9SJeremy L Thompson   // sum into V
214*db2becc9SJeremy L Thompson   sum_V_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV, cstrdV, rV, tx);
215*db2becc9SJeremy L Thompson }
216