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 1D 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 (2D) 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_2d_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, T *swork) { 22 // Assumptions 23 // 1. 1D threads of size max(P,Q) 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. Two products per component 27 // 4.1 Batch P of (1xP) matrices times (PxQ) matrix => Batch P of (1xQ) matrices 28 // 4.2 Batch 1 of (QxP) matrix times (PxQ) matrix => (QxQ) matrix 29 // 5. Each thread computes one row of the output of each product 30 // 6. Sync is recommended before and after the call 31 32 for (int comp = 0; comp < NUM_COMP; comp++) { 33 // 1st product -- Batch P of (1xP) matrices [reg] x (PxQ) [shmem] => Batch P of (1xQ) matrices 34 // the batch output P x (1xQ) is written on the fly to shmem 35 if (tx < 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; 41 for (int i = 0; i < P; i++) { 42 rTmp += rU[0][comp][i] * sT(i, j); 43 } 44 sTmp(0, j, sld) = rTmp; 45 } 46 } // end of: if (tx < P) 47 __syncthreads(); 48 49 // 2nd product -- Batch 1 of a (QxP) matrix [shmem] x (PxQ) [shmem] => (QxQ) matrix [reg] 50 if (tx < Q) { 51 const int batchid = 0; 52 const int sld = Q; 53 T *sTmp = swork + batchid * (Q * P); 54 for (int j = 0; j < Q; j++) { 55 rTmp = 0.0; 56 for (int i = 0; i < P; i++) { 57 rTmp += sTmp(tx, i, sld) * sT(i, j); 58 } 59 rV[0][comp][j] += rTmp; 60 } 61 } 62 __syncthreads(); 63 } 64 } 65 66 //////////////////////////////////////////////////////////////////////////////// 67 extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__ 68 void magma_interpn_2d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV, 69 const int cstrdV, const int nelem) { 70 MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 71 72 const int tx = threadIdx.x; 73 const int ty = threadIdx.y; 74 const int elem_id = (blockIdx.x * blockDim.y) + ty; 75 76 if (elem_id >= nelem) return; 77 78 CeedScalar rU[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // for a non-fused operator BASIS_DIM is always 1 79 CeedScalar rV[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // for a non-fused operator BASIS_DIM is always 1 80 CeedScalar rTmp = 0.0; 81 82 // shift global memory pointers by elem stride 83 dU += elem_id * estrdU; 84 dV += elem_id * estrdV; 85 86 // assign shared memory pointers 87 CeedScalar *sT = (CeedScalar *)shared_data; 88 CeedScalar *sTmp = sT + BASIS_P * BASIS_Q; 89 sTmp += ty * (BASIS_P * BASIS_MAX_P_Q); 90 91 // read T 92 if (ty == 0) { 93 read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dT, sT); 94 } 95 96 // read U -- there is a sync at the end of this function 97 read_U_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dU, cstrdU, rU, sTmp, tx); 98 99 // no sync needed here -- read_U_2d already syncs at the end 100 magma_interp_2d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_P, BASIS_Q, BASIS_P, BASIS_Q>(sT, rU, rV, tx, rTmp, sTmp); 101 __syncthreads(); 102 103 // write V 104 write_V_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV, cstrdV, rV, tx); 105 } 106 107 //////////////////////////////////////////////////////////////////////////////// 108 extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__ 109 void magma_interpt_2d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV, 110 const int cstrdV, const int nelem) { 111 MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 112 113 const int tx = threadIdx.x; 114 const int ty = threadIdx.y; 115 const int elem_id = (blockIdx.x * blockDim.y) + ty; 116 117 if (elem_id >= nelem) return; 118 119 CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // for a non-fused operator BASIS_DIM is always 1 120 CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // for a non-fused operator BASIS_DIM is always 1 121 CeedScalar rTmp = 0.0; 122 123 // shift global memory pointers by elem stride 124 dU += elem_id * estrdU; 125 dV += elem_id * estrdV; 126 127 // assign shared memory pointers 128 CeedScalar *sT = (CeedScalar *)shared_data; 129 CeedScalar *sTmp = sT + BASIS_Q * BASIS_P; 130 sTmp += ty * (BASIS_Q * BASIS_MAX_P_Q); 131 132 // read T 133 if (ty == 0) { 134 read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dT, sT); 135 } 136 137 // read U -- there is a sync at the end of this function 138 read_U_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU, cstrdU, rU, sTmp, tx); 139 140 // no sync needed here -- read_U_2d already syncs at the end 141 magma_interp_2d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P>(sT, rU, rV, tx, rTmp, sTmp); 142 __syncthreads(); 143 144 // write V 145 write_V_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV, cstrdV, rV, tx); 146 } 147 148 //////////////////////////////////////////////////////////////////////////////// 149 extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__ 150 void magma_interpta_2d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV, 151 const int cstrdV, const int nelem) { 152 MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 153 154 const int tx = threadIdx.x; 155 const int ty = threadIdx.y; 156 const int elem_id = (blockIdx.x * blockDim.y) + ty; 157 158 if (elem_id >= nelem) return; 159 160 CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // for a non-fused operator BASIS_DIM is always 1 161 CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // for a non-fused operator BASIS_DIM is always 1 162 CeedScalar rTmp = 0.0; 163 164 // shift global memory pointers by elem stride 165 dU += elem_id * estrdU; 166 dV += elem_id * estrdV; 167 168 // assign shared memory pointers 169 CeedScalar *sT = (CeedScalar *)shared_data; 170 CeedScalar *sTmp = sT + BASIS_Q * BASIS_P; 171 sTmp += ty * (BASIS_Q * BASIS_MAX_P_Q); 172 173 // read T 174 if (ty == 0) { 175 read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dT, sT); 176 } 177 178 // read U -- there is a sync at the end of this function 179 read_U_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU, cstrdU, rU, sTmp, tx); 180 181 // no sync needed here -- read_U_2d already syncs at the end 182 magma_interp_2d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P>(sT, rU, rV, tx, rTmp, sTmp); 183 __syncthreads(); 184 185 // sum into V 186 sum_V_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV, cstrdV, rV, tx); 187 } 188