xref: /libCEED/include/ceed/jit-source/magma/magma-basis-grad-2d.h (revision ca94c3ddc8f82b7d93a79f9e4812e99b8be840ff)
1 // Copyright (c) 2017-2022, 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 gradient in 2D
10 #ifndef CEED_MAGMA_BASIS_GRAD_2D_H
11 #define CEED_MAGMA_BASIS_GRAD_2D_H
12 
13 #include "magma-common-tensor.h"
14 
15 // macros to abstract access of shared memory and reg. file
16 #define sT(i, j) sT[(j)*P + (i)]
17 #define sTmp(i, j, ldw) sTmp[(j) * (ldw) + (i)]
18 
19 ////////////////////////////////////////////////////////////////////////////////
20 // Helper function to add or set into V
21 template <typename T, bool Add>
22 struct magma_grad_2d_device_accumulate;
23 
24 template <typename T>
25 struct magma_grad_2d_device_accumulate<T, true> {
26   static __device__ __inline__ void op(T &rV, const T &rTmp) { rV += rTmp; }
27 };
28 
29 template <typename T>
30 struct magma_grad_2d_device_accumulate<T, false> {
31   static __device__ __inline__ void op(T &rV, const T &rTmp) { rV = rTmp; }
32 };
33 
34 ////////////////////////////////////////////////////////////////////////////////
35 // grad basis action (2D)
36 // This function is called two times at a higher level for 2D
37 // DIM_U   -- for the size of rU[DIM_U * NUM_COMP * MAX_P_Q]
38 // DIM_V   -- for the size of rV[DIM_V * NUM_COMP * MAX_P_Q]
39 // i_DIM   -- the index of the outermost loop over dimensions in grad
40 // i_DIM_U -- which dim index of rU is accessed (always 0 for notrans, 0 or 1 for trans)
41 // i_DIM_V -- which dim index of rV is accessed (0 or 1 for notrans, always 0 for trans)
42 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>
43 static __device__ __inline__ void magma_grad_2d_device(const T *sTinterp, const T *sTgrad, T rU[DIM_U][NUM_COMP][rU_SIZE],
44                                                        T rV[DIM_V][NUM_COMP][rV_SIZE], const int tx, T rTmp, T *swork) {
45   // Assumptions
46   // 0. This device routine applies grad for one dim only (i_DIM), so it should be called twice for 2D
47   // 1. 1D threads of size max(P,Q)
48   // 2. input:  rU[DIM_U x NUM_COMP x P] in registers (per thread)
49   // 3. output: rV[DIM_V x NUM_COMP x Q] in registers (per thread)
50   // 4. Two products per each (dim,component) pair
51   //  4.1 Batch P of (1xP) matrices times (PxQ) matrix => Batch P of (1xQ) matrices
52   //  4.2 Batch 1 of (QxP) matrix   times (PxQ) matrix => (QxQ) matrix
53   // 6. Each thread computes one row of the output of each product
54   // 7. Sync is recommended before and after the call
55 
56   for (int comp = 0; comp < NUM_COMP; comp++) {
57     // 1st product -- Batch P of (1xP) matrices [reg] x (PxQ) [shmem] => Batch P of (1xQ) matrices
58     // the batch output P x (1xQ) is written on the fly to shmem
59     if (tx < P) {
60       const int batchid = tx;
61       const int sld     = 1;
62       const T  *sT      = (i_DIM == 0) ? sTgrad : sTinterp;
63       T        *sTmp    = swork + batchid * (1 * Q);
64       for (int j = 0; j < Q; j++) {
65         rTmp = 0.0;
66         for (int i = 0; i < P; i++) {
67           rTmp += rU[i_DIM_U][comp][i] * sT(i, j);
68         }
69         sTmp(0, j, sld) = rTmp;
70       }
71     }  // end of: if (tx < P)
72     __syncthreads();
73 
74     // 2nd product -- Batch 1 of a (QxP) matrix [shmem] x (PxQ) [shmem] => (QxQ) matrix [reg]
75     if (tx < Q) {
76       const int batchid = 0;
77       const int sld     = Q;
78       const T  *sT      = (i_DIM == 1) ? sTgrad : sTinterp;
79       T        *sTmp    = swork + batchid * (Q * P);
80       for (int j = 0; j < Q; j++) {
81         rTmp = 0.0;
82         for (int i = 0; i < P; i++) {
83           rTmp += sTmp(tx, i, sld) * sT(i, j);
84         }
85         magma_grad_2d_device_accumulate<T, ADD>::op(rV[i_DIM_V][comp][j], rTmp);
86       }
87     }
88     __syncthreads();
89   }  // loop over NUM_COMP
90 }
91 
92 ////////////////////////////////////////////////////////////////////////////////
93 extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__
94     void magma_gradn_2d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU,
95                                const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) {
96   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
97 
98   const int tx      = threadIdx.x;
99   const int ty      = threadIdx.y;
100   const int elem_id = (blockIdx.x * blockDim.y) + ty;
101 
102   if (elem_id >= nelem) return;
103 
104   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // here DIM_U = 1, but might be different for a fused operator
105   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // here DIM_V = 1, but might be different for a fused operator
106   CeedScalar rTmp                           = 0.0;
107 
108   // shift global memory pointers by elem stride
109   dU += elem_id * estrdU;
110   dV += elem_id * estrdV;
111 
112   // assign shared memory pointers
113   CeedScalar *sTinterp = (CeedScalar *)shared_data;
114   CeedScalar *sTgrad   = sTinterp + BASIS_P * BASIS_Q;
115   CeedScalar *sTmp     = sTgrad + BASIS_P * BASIS_Q;
116   sTmp += ty * (BASIS_P * BASIS_MAX_P_Q);
117 
118   // read T
119   if (ty == 0) {
120     read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dinterp1d, sTinterp);
121     read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dgrad1d, sTgrad);
122   }
123 
124   /* read U (idim = 0 for dU, i_DIM = 0 for rU) --
125      there is a sync at the end of this function */
126   read_U_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx);
127 
128   /* first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) --
129      output from rV[0][][] into dV (idim = 0) */
130   magma_grad_2d_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,
131                                                                                                              sTmp);
132   /* there is a sync at the end of magma_grad_2d_device */
133   write_V_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (0 * dstrdV), cstrdV, rV, tx);
134 
135   /* second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) --
136   output from rV[0][][] into dV (idim = 1) */
137   magma_grad_2d_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,
138                                                                                                              sTmp);
139   /* there is a sync at the end of magma_grad_2d_device */
140   write_V_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (1 * dstrdV), cstrdV, rV, tx);
141 }
142 
143 ////////////////////////////////////////////////////////////////////////////////
144 extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__
145     void magma_gradt_2d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU,
146                                const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) {
147   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
148 
149   const int tx      = threadIdx.x;
150   const int ty      = threadIdx.y;
151   const int elem_id = (blockIdx.x * blockDim.y) + ty;
152 
153   if (elem_id >= nelem) return;
154 
155   CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0};  // here DIM_U = 1, but might be different for a fused operator
156   CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0};  // here DIM_V = 1, but might be different for a fused operator
157   CeedScalar rTmp                           = 0.0;
158 
159   // shift global memory pointers by elem stride
160   dU += elem_id * estrdU;
161   dV += elem_id * estrdV;
162 
163   // assign shared memory pointers
164   CeedScalar *sTinterp = (CeedScalar *)shared_data;
165   CeedScalar *sTgrad   = sTinterp + BASIS_Q * BASIS_P;
166   CeedScalar *sTmp     = sTgrad + BASIS_Q * BASIS_P;
167   sTmp += ty * (BASIS_Q * BASIS_MAX_P_Q);
168 
169   // read T
170   if (ty == 0) {
171     read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dinterp1d, sTinterp);
172     read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dgrad1d, sTgrad);
173   }
174   __syncthreads();
175 
176   /* read U (idim = 0 for dU, i_DIM = 0 for rU) --
177      there is a sync at the end of this function */
178   read_U_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx);
179   /* first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) */
180   magma_grad_2d_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);
181   /* there is a sync at the end of magma_grad_2d_device */
182 
183   /* read U (idim = 1 for dU, i_DIM = 0 for rU) --
184      there is a sync at the end of this function */
185   read_U_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (1 * dstrdU), cstrdU, rU, sTmp, tx);
186   /* second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) */
187   magma_grad_2d_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);
188   /* there is a sync at the end of magma_grad_2d_device */
189 
190   // write V
191   write_V_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV + (0 * dstrdV), cstrdV, rV, tx);
192 }
193 
194 #endif  // CEED_MAGMA_BASIS_GRAD_2D_H
195