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