xref: /libCEED/rust/libceed-sys/c-src/include/ceed/jit-source/magma/magma-basis-interp-2d.h (revision 3c1e2aff6d111c93ca8797996aaf987f66b08927)
1f80f4a74SSebastian Grimberg // Copyright (c) 2017-2022, 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 
8*3c1e2affSSebastian Grimberg /// @file
9*3c1e2affSSebastian Grimberg /// Internal header for MAGMA tensor basis interpolation in 1D
10*3c1e2affSSebastian Grimberg #ifndef CEED_MAGMA_BASIS_INTERP_2D_H
11*3c1e2affSSebastian Grimberg #define CEED_MAGMA_BASIS_INTERP_2D_H
12*3c1e2affSSebastian Grimberg 
13*3c1e2affSSebastian Grimberg #include "magma-common-tensor.h"
14*3c1e2affSSebastian Grimberg 
15f80f4a74SSebastian Grimberg // macros to abstract access of shared memory and reg. file
16*3c1e2affSSebastian Grimberg #define sT(i, j) sT[(j)*P + (i)]
17f80f4a74SSebastian Grimberg #define sTmp(i, j, ldw) sTmp[(j) * (ldw) + (i)]
18f80f4a74SSebastian Grimberg 
19f80f4a74SSebastian Grimberg //////////////////////////////////////////////////////////////////////////////////////////
20f80f4a74SSebastian Grimberg // interp basis action (2D)
21*3c1e2affSSebastian Grimberg template <typename T, int DIM_U, int DIM_V, int NUM_COMP, int P, int Q, int rU_SIZE, int rV_SIZE>
22*3c1e2affSSebastian Grimberg static __device__ __inline__ void magma_interp_2d_device(const T *sT, magma_trans_t transT, T rU[DIM_U][NUM_COMP][rU_SIZE],
23*3c1e2affSSebastian Grimberg                                                          T rV[DIM_V][NUM_COMP][rV_SIZE], const int tx, T rTmp, T *swork) {
24f80f4a74SSebastian Grimberg   // Assumptions
25*3c1e2affSSebastian Grimberg   // 1. 1D threads of size max(P,Q)
26*3c1e2affSSebastian Grimberg   // 2. input:  rU[DIM_U x NUM_COMP x rU_SIZE] in registers (per thread)
27*3c1e2affSSebastian Grimberg   // 3. output: rV[DIM_V x NUM_COMP x rV_SIZE] in registers (per thread)
28f80f4a74SSebastian Grimberg   // 4. Two products per component
29*3c1e2affSSebastian Grimberg   //  4.1 Batch P of (1xP) matrices times (PxQ) matrix => Batch P of (1xQ) matrices
30*3c1e2affSSebastian Grimberg   //  4.2 Batch 1 of (QxP) matrix   times (PxQ) matrix => (QxQ) matrix
31f80f4a74SSebastian Grimberg   // 5. Each thread computes one row of the output of each product
32f80f4a74SSebastian Grimberg   // 6. Sync is recommended before and after the call
33f80f4a74SSebastian Grimberg 
34*3c1e2affSSebastian Grimberg   for (int comp = 0; comp < NUM_COMP; comp++) {
35*3c1e2affSSebastian Grimberg     // 1st product -- Batch P of (1xP) matrices [reg] x (PxQ) [shmem] => Batch P of (1xQ) matrices
36*3c1e2affSSebastian Grimberg     // the batch output P x (1xQ) is written on the fly to shmem
37*3c1e2affSSebastian Grimberg     if (tx < P) {
38f80f4a74SSebastian Grimberg       const int batchid = tx;
39f80f4a74SSebastian Grimberg       const int sld     = 1;
40*3c1e2affSSebastian Grimberg       T        *sTmp    = swork + batchid * (1 * Q);
41*3c1e2affSSebastian Grimberg       for (int j = 0; j < Q; j++) {
42f80f4a74SSebastian Grimberg         rTmp = 0.0;
43*3c1e2affSSebastian Grimberg         for (int i = 0; i < P; i++) {
44*3c1e2affSSebastian Grimberg           rTmp += rU[0][comp][i] * sT(i, j);
45f80f4a74SSebastian Grimberg         }
46f80f4a74SSebastian Grimberg         sTmp(0, j, sld) = rTmp;
47f80f4a74SSebastian Grimberg       }
48*3c1e2affSSebastian Grimberg     }  // end of: if (tx < P)
49f80f4a74SSebastian Grimberg     __syncthreads();
50f80f4a74SSebastian Grimberg 
51*3c1e2affSSebastian Grimberg     // 2nd product -- Batch 1 of a (QxP) matrix [shmem] x (PxQ) [shmem] => (QxQ) matrix [reg]
52*3c1e2affSSebastian Grimberg     if (tx < Q) {
53f80f4a74SSebastian Grimberg       const int batchid = 0;
54*3c1e2affSSebastian Grimberg       const int sld     = Q;
55*3c1e2affSSebastian Grimberg       T        *sTmp    = swork + batchid * (Q * P);
56*3c1e2affSSebastian Grimberg       for (int j = 0; j < Q; j++) {
57f80f4a74SSebastian Grimberg         rTmp = 0.0;
58*3c1e2affSSebastian Grimberg         for (int i = 0; i < P; i++) {
59f80f4a74SSebastian Grimberg           rTmp += sTmp(tx, i, sld) * sT(i, j);
60f80f4a74SSebastian Grimberg         }
61*3c1e2affSSebastian Grimberg         rV[0][comp][j] += rTmp;
62f80f4a74SSebastian Grimberg       }
63f80f4a74SSebastian Grimberg     }
64f80f4a74SSebastian Grimberg     __syncthreads();
65f80f4a74SSebastian Grimberg   }
66f80f4a74SSebastian Grimberg }
67f80f4a74SSebastian Grimberg 
68f80f4a74SSebastian Grimberg //////////////////////////////////////////////////////////////////////////////////////////
69*3c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__
70f80f4a74SSebastian Grimberg     void magma_interpn_2d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV,
71f80f4a74SSebastian Grimberg                                  const int cstrdV, const int nelem) {
72f80f4a74SSebastian Grimberg   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
73f80f4a74SSebastian Grimberg 
74f80f4a74SSebastian Grimberg   const int     tx      = threadIdx.x;
75f80f4a74SSebastian Grimberg   const int     ty      = threadIdx.y;
76f80f4a74SSebastian Grimberg   const int     elem_id = (blockIdx.x * blockDim.y) + ty;
77f80f4a74SSebastian Grimberg   magma_trans_t transT  = MagmaNoTrans;
78f80f4a74SSebastian Grimberg 
79f80f4a74SSebastian Grimberg   if (elem_id >= nelem) return;
80f80f4a74SSebastian Grimberg 
81*3c1e2affSSebastian Grimberg   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
82*3c1e2affSSebastian Grimberg   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
83f80f4a74SSebastian Grimberg   CeedScalar rTmp                           = 0.0;
84f80f4a74SSebastian Grimberg 
85f80f4a74SSebastian Grimberg   // shift global memory pointers by elem stride
86f80f4a74SSebastian Grimberg   dU += elem_id * estrdU;
87f80f4a74SSebastian Grimberg   dV += elem_id * estrdV;
88f80f4a74SSebastian Grimberg 
89f80f4a74SSebastian Grimberg   // assign shared memory pointers
90*3c1e2affSSebastian Grimberg   CeedScalar *sT   = (CeedScalar *)shared_data;
91*3c1e2affSSebastian Grimberg   CeedScalar *sTmp = sT + BASIS_P * BASIS_Q;
92*3c1e2affSSebastian Grimberg   sTmp += ty * (BASIS_P * BASIS_MAX_P_Q);
93f80f4a74SSebastian Grimberg 
94f80f4a74SSebastian Grimberg   // read T
95f80f4a74SSebastian Grimberg   if (ty == 0) {
96*3c1e2affSSebastian Grimberg     dread_T_gm2sm<BASIS_P, BASIS_Q>(tx, transT, dT, sT);
97f80f4a74SSebastian Grimberg   }
98f80f4a74SSebastian Grimberg 
99f80f4a74SSebastian Grimberg   // read U -- there is a sync at the end of this function
100*3c1e2affSSebastian Grimberg   readU_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dU, cstrdU, rU, sTmp, tx);
101f80f4a74SSebastian Grimberg 
102f80f4a74SSebastian Grimberg   // no sync needed here -- readU_2d already syncs at the end
103*3c1e2affSSebastian Grimberg   magma_interp_2d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_P, BASIS_Q, BASIS_P, BASIS_Q>(sT, transT, rU, rV, tx, rTmp, sTmp);
104f80f4a74SSebastian Grimberg   __syncthreads();
105f80f4a74SSebastian Grimberg 
106f80f4a74SSebastian Grimberg   // write V
107*3c1e2affSSebastian Grimberg   writeV_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV, cstrdV, rV, tx);
108f80f4a74SSebastian Grimberg }
109f80f4a74SSebastian Grimberg 
110f80f4a74SSebastian Grimberg //////////////////////////////////////////////////////////////////////////////////////////
111*3c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__
112f80f4a74SSebastian Grimberg     void magma_interpt_2d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV,
113f80f4a74SSebastian Grimberg                                  const int cstrdV, const int nelem) {
114f80f4a74SSebastian Grimberg   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
115f80f4a74SSebastian Grimberg 
116f80f4a74SSebastian Grimberg   const int     tx      = threadIdx.x;
117f80f4a74SSebastian Grimberg   const int     ty      = threadIdx.y;
118f80f4a74SSebastian Grimberg   const int     elem_id = (blockIdx.x * blockDim.y) + ty;
119f80f4a74SSebastian Grimberg   magma_trans_t transT  = MagmaTrans;
120f80f4a74SSebastian Grimberg 
121f80f4a74SSebastian Grimberg   if (elem_id >= nelem) return;
122f80f4a74SSebastian Grimberg 
123*3c1e2affSSebastian Grimberg   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
124*3c1e2affSSebastian Grimberg   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // for a non-fused operator BASIS_DIM is always 1
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
132*3c1e2affSSebastian Grimberg   CeedScalar *sT   = (CeedScalar *)shared_data;
133*3c1e2affSSebastian Grimberg   CeedScalar *sTmp = sT + BASIS_Q * BASIS_P;
134*3c1e2affSSebastian Grimberg   sTmp += ty * (BASIS_Q * BASIS_MAX_P_Q);
135f80f4a74SSebastian Grimberg 
136f80f4a74SSebastian Grimberg   // read T
137f80f4a74SSebastian Grimberg   if (ty == 0) {
138*3c1e2affSSebastian Grimberg     dread_T_gm2sm<BASIS_Q, BASIS_P>(tx, transT, dT, sT);
139f80f4a74SSebastian Grimberg   }
140f80f4a74SSebastian Grimberg 
141f80f4a74SSebastian Grimberg   // read V
142*3c1e2affSSebastian Grimberg   readV_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV, cstrdV, rV, tx);
143f80f4a74SSebastian Grimberg 
144f80f4a74SSebastian Grimberg   // read U -- there is a sync at the end of this function
145*3c1e2affSSebastian Grimberg   readU_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU, cstrdU, rU, sTmp, tx);
146f80f4a74SSebastian Grimberg 
147f80f4a74SSebastian Grimberg   // no sync needed here -- readU_2d already syncs at the end
148*3c1e2affSSebastian Grimberg   magma_interp_2d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P>(sT, transT, rU, rV, tx, rTmp, sTmp);
149f80f4a74SSebastian Grimberg   __syncthreads();
150f80f4a74SSebastian Grimberg 
151f80f4a74SSebastian Grimberg   // write V
152*3c1e2affSSebastian Grimberg   writeV_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV, cstrdV, rV, tx);
153f80f4a74SSebastian Grimberg }
154*3c1e2affSSebastian Grimberg 
155*3c1e2affSSebastian Grimberg #endif  // CEED_MAGMA_BASIS_INTERP_2D_H
156