xref: /libCEED/include/ceed/jit-source/magma/magma-basis-interp-3d.h (revision 392d940c90dcc80d2c8cd0455511b2e57938a796)
1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Internal header for MAGMA tensor basis interpolation in 3D
10 
11 #include "magma-common-tensor.h"
12 
13 // macros to abstract access of shared memory and reg. file
14 #define sT(i, j) sT[(j) * P + (i)]
15 #define sTmp(i, j, ldw) sTmp[(j) * (ldw) + (i)]
16 
17 ////////////////////////////////////////////////////////////////////////////////
18 // interp basis action (3D)
19 template <typename T, int DIM_U, int DIM_V, int NUM_COMP, int P, int Q, int rU_SIZE, int rV_SIZE>
20 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,
21                                                          T rTmp[Q], T *swork) {
22   // Assumptions
23   // 1. 1D threads of size max(P,Q)^2
24   // 2. input:  rU[DIM_U x NUM_COMP x rU_SIZE] in registers (per thread)
25   // 3. output: rV[DIM_V x NUM_COMP x rV_SIZE] in registers (per thread)
26   // 4. Three products per component
27   //  4.1 Batch P^2 of (1xP) matrices times (PxQ) matrix => Batch P^2 of (1xQ) matrices
28   //  4.2 Batch P   of (QxP) matrices times (PxQ) matrix => Batch P   of (QxQ) matrices
29   //  4.3 Batch 1   of (Q^2xP_) matrix times (PxQ) matrix => (Q^2xQ_) matrix
30   // 5. Each thread computes one row of the output of each product
31   // 6. Sync is recommended before and after the call
32 
33   for (int comp = 0; comp < NUM_COMP; comp++) {
34     // Batch P^2 of (1xP) matrices [reg] times (PxQ) matrix [shmem] => Batch P^2 of (1xQ) matrices [shmem]
35     if (tx < (P * P)) {
36       const int batchid = tx;
37       const int sld     = 1;
38       T        *sTmp    = swork + batchid * (1 * Q);
39       for (int j = 0; j < Q; j++) {
40         rTmp[0] = 0.0;
41         for (int i = 0; i < P; i++) {
42           rTmp[0] += rU[0][comp][i] * sT(i, j);
43         }
44         sTmp(0, j, sld) = rTmp[0];
45       }
46     }  // end of: if (tx < P*P)
47     __syncthreads();
48 
49     // Batch P of (QxP) matrices [shmem] times (PxQ) matrix [shmem] => Batch P of (QxQ) matrices [reg]
50     if (tx < (P * Q)) {
51       const int batchid = tx / Q;
52       const int tx_     = tx % Q;
53       const int sld     = Q;
54       T        *sTmp    = swork + batchid * (Q * P);  // sTmp is input
55       for (int j = 0; j < Q; j++) {
56         rTmp[j] = 0.0;
57         for (int i = 0; i < P; i++) {
58           rTmp[j] += sTmp(tx_, i, sld) * sT(i, j);
59         }
60       }
61     }
62     __syncthreads();
63 
64     // write rTmp[] into shmem as batch P of QxQ matrices
65     if (tx < (P * Q)) {
66       const int batchid = tx / Q;
67       const int tx_     = tx % Q;
68       const int sld     = Q;
69       T        *sTmp    = swork + batchid * (Q * Q);
70       for (int j = 0; j < Q; j++) {
71         sTmp(tx_, j, sld) = rTmp[j];
72       }
73     }
74     __syncthreads();
75 
76     // Batch 1 of (Q^2xP_) matrices [shmem] times (PxQ) matrix [shmem] => Batch 1 of (Q^2xQ_) matrices [reg]
77     if (tx < (Q * Q)) {
78       // No need to declare batchid = (tx  / Q^2) = always zero
79       // No need to declare tx_     = (tx_ % Q^2) = always tx
80       const int sld  = Q * Q;
81       T        *sTmp = swork;
82       for (int j = 0; j < Q; j++) {
83         rTmp[0] = 0.0;
84         for (int i = 0; i < P; i++) {
85           rTmp[0] += sTmp(tx, i, sld) * sT(i, j);
86         }
87         rV[0][comp][j] += rTmp[0];
88       }
89     }
90     __syncthreads();
91   }
92 }
93 
94 ////////////////////////////////////////////////////////////////////////////////
95 extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q *BASIS_MAX_P_Q, MAGMA_MAXTHREADS_3D)) __global__
96     void magma_interpn_3d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV,
97                                  const int cstrdV, const int nelem) {
98   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
99 
100   const int tx      = threadIdx.x;
101   const int ty      = threadIdx.y;
102   const int elem_id = (blockIdx.x * blockDim.y) + ty;
103 
104   if (elem_id >= nelem) return;
105 
106   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
107   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
108   CeedScalar rTmp[BASIS_Q]                  = {0.0};
109 
110   // shift global memory pointers by elem stride
111   dU += elem_id * estrdU;
112   dV += elem_id * estrdV;
113 
114   // assign shared memory pointers
115   CeedScalar *sT   = (CeedScalar *)shared_data;
116   CeedScalar *sTmp = sT + BASIS_P * BASIS_Q;
117   sTmp += ty * (max(BASIS_P * BASIS_P * BASIS_MAX_P_Q, BASIS_P * BASIS_Q * BASIS_Q));
118 
119   // read T
120   if (ty == 0) {
121     read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dT, sT);
122   }
123 
124   // read U (idim = 0 for dU, i_DIM = 0 for rU, u_dimstride is always 0)
125   read_U_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dU, cstrdU, rU, sTmp, tx);
126   // there is a sync at the end of this function
127 
128   magma_interp_3d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_P, BASIS_Q, BASIS_P, BASIS_Q>(sT, rU, rV, tx, rTmp, sTmp);
129   __syncthreads();
130 
131   // write V
132   write_V_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV, cstrdV, rV, tx);
133 }
134 
135 ////////////////////////////////////////////////////////////////////////////////
136 extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q *BASIS_MAX_P_Q, MAGMA_MAXTHREADS_3D)) __global__
137     void magma_interpt_3d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV,
138                                  const int cstrdV, const int nelem) {
139   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
140 
141   const int tx      = threadIdx.x;
142   const int ty      = threadIdx.y;
143   const int elem_id = (blockIdx.x * blockDim.y) + ty;
144 
145   if (elem_id >= nelem) return;
146 
147   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
148   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
149   CeedScalar rTmp[BASIS_P]                  = {0.0};
150 
151   // shift global memory pointers by elem stride
152   dU += elem_id * estrdU;
153   dV += elem_id * estrdV;
154 
155   // assign shared memory pointers
156   CeedScalar *sT   = (CeedScalar *)shared_data;
157   CeedScalar *sTmp = sT + BASIS_Q * BASIS_P;
158   sTmp += ty * (max(BASIS_Q * BASIS_Q * BASIS_MAX_P_Q, BASIS_Q * BASIS_P * BASIS_P));
159 
160   // read T
161   if (ty == 0) {
162     read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dT, sT);
163   }
164 
165   // read U (idim = 0 for dU, i_DIM = 0 for rU, u_dimstride is always 0)
166   read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU, cstrdU, rU, sTmp, tx);
167   // there is a sync at the end of this function
168 
169   magma_interp_3d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P>(sT, rU, rV, tx, rTmp, sTmp);
170   __syncthreads();
171 
172   // write V
173   write_V_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV, cstrdV, rV, tx);
174 }
175