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 83c1e2affSSebastian Grimberg /// @file 93c1e2affSSebastian Grimberg /// Internal header for MAGMA tensor basis gradient in 2D 103c1e2affSSebastian Grimberg #ifndef CEED_MAGMA_BASIS_GRAD_2D_H 113c1e2affSSebastian Grimberg #define CEED_MAGMA_BASIS_GRAD_2D_H 123c1e2affSSebastian Grimberg 133c1e2affSSebastian Grimberg #include "magma-common-tensor.h" 143c1e2affSSebastian Grimberg 15f80f4a74SSebastian Grimberg // macros to abstract access of shared memory and reg. file 163c1e2affSSebastian Grimberg #define sT(i, j) sT[(j)*P + (i)] 17f80f4a74SSebastian Grimberg #define sTmp(i, j, ldw) sTmp[(j) * (ldw) + (i)] 18f80f4a74SSebastian Grimberg 19*9e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 20f80f4a74SSebastian Grimberg // grad basis action (2D) 21f80f4a74SSebastian Grimberg // This function is called two times at a higher level for 2D 223c1e2affSSebastian Grimberg // DIM_U -- for the size of rU[DIM_U * NUM_COMP * MAX_P_Q] 233c1e2affSSebastian Grimberg // DIM_V -- for the size of rV[DIM_V * NUM_COMP * MAX_P_Q] 243c1e2affSSebastian Grimberg // i_DIM -- the index of the outermost loop over dimensions in grad 253c1e2affSSebastian Grimberg // i_DIM_U -- which dim index of rU is accessed (always 0 for notrans, 0 or 1 for trans) 263c1e2affSSebastian Grimberg // i_DIM_V -- which dim index of rV is accessed (0 or 1 for notrans, always 0 for trans) 27f80f4a74SSebastian Grimberg // the scalar beta is used to specify whether to accumulate to rV, or overwrite it 283c1e2affSSebastian Grimberg 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> 293c1e2affSSebastian Grimberg static __device__ __inline__ void magma_grad_2d_device(const T *sTinterp, const T *sTgrad, T rU[DIM_U][NUM_COMP][rU_SIZE], 303c1e2affSSebastian Grimberg T rV[DIM_V][NUM_COMP][rV_SIZE], T beta, const int tx, T rTmp, T *swork) { 31f80f4a74SSebastian Grimberg // Assumptions 323c1e2affSSebastian Grimberg // 0. This device routine applies grad for one dim only (i_DIM), so it should be called twice for 2D 333c1e2affSSebastian Grimberg // 1. 1D threads of size max(P,Q) 343c1e2affSSebastian Grimberg // 2. input: rU[DIM_U x NUM_COMP x P] in registers (per thread) 353c1e2affSSebastian Grimberg // 3. output: rV[DIM_V x NUM_COMP x Q] in registers (per thread) 36f80f4a74SSebastian Grimberg // 4. Two products per each (dim,component) pair 373c1e2affSSebastian Grimberg // 4.1 Batch P of (1xP) matrices times (PxQ) matrix => Batch P of (1xQ) matrices 383c1e2affSSebastian Grimberg // 4.2 Batch 1 of (QxP) matrix times (PxQ) matrix => (QxQ) matrix 39f80f4a74SSebastian Grimberg // 6. Each thread computes one row of the output of each product 40f80f4a74SSebastian Grimberg // 7. Sync is recommended before and after the call 41f80f4a74SSebastian Grimberg 423c1e2affSSebastian Grimberg for (int comp = 0; comp < NUM_COMP; comp++) { 433c1e2affSSebastian Grimberg // 1st product -- Batch P of (1xP) matrices [reg] x (PxQ) [shmem] => Batch P of (1xQ) matrices 443c1e2affSSebastian Grimberg // the batch output P x (1xQ) is written on the fly to shmem 453c1e2affSSebastian Grimberg if (tx < P) { 46f80f4a74SSebastian Grimberg const int batchid = tx; 47f80f4a74SSebastian Grimberg const int sld = 1; 483c1e2affSSebastian Grimberg const T *sT = (i_DIM == 0) ? sTgrad : sTinterp; 493c1e2affSSebastian Grimberg T *sTmp = swork + batchid * (1 * Q); 503c1e2affSSebastian Grimberg for (int j = 0; j < Q; j++) { 51f80f4a74SSebastian Grimberg rTmp = 0.0; 523c1e2affSSebastian Grimberg for (int i = 0; i < P; i++) { 533c1e2affSSebastian Grimberg rTmp += rU[i_DIM_U][comp][i] * sT(i, j); 54f80f4a74SSebastian Grimberg } 55f80f4a74SSebastian Grimberg sTmp(0, j, sld) = rTmp; 56f80f4a74SSebastian Grimberg } 573c1e2affSSebastian Grimberg } // end of: if (tx < P) 58f80f4a74SSebastian Grimberg __syncthreads(); 59f80f4a74SSebastian Grimberg 603c1e2affSSebastian Grimberg // 2nd product -- Batch 1 of a (QxP) matrix [shmem] x (PxQ) [shmem] => (QxQ) matrix [reg] 613c1e2affSSebastian Grimberg if (tx < Q) { 62f80f4a74SSebastian Grimberg const int batchid = 0; 633c1e2affSSebastian Grimberg const int sld = Q; 643c1e2affSSebastian Grimberg const T *sT = (i_DIM == 1) ? sTgrad : sTinterp; 653c1e2affSSebastian Grimberg T *sTmp = swork + batchid * (Q * P); 663c1e2affSSebastian Grimberg for (int j = 0; j < Q; j++) { 67f80f4a74SSebastian Grimberg rTmp = 0.0; 683c1e2affSSebastian Grimberg for (int i = 0; i < P; i++) { 69f80f4a74SSebastian Grimberg rTmp += sTmp(tx, i, sld) * sT(i, j); 70f80f4a74SSebastian Grimberg } 713c1e2affSSebastian Grimberg rV[i_DIM_V][comp][j] *= beta; 723c1e2affSSebastian Grimberg rV[i_DIM_V][comp][j] += rTmp; 73f80f4a74SSebastian Grimberg } 74f80f4a74SSebastian Grimberg } 75f80f4a74SSebastian Grimberg __syncthreads(); 763c1e2affSSebastian Grimberg } // loop over NUM_COMP 77f80f4a74SSebastian Grimberg } 78f80f4a74SSebastian Grimberg 79*9e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 803c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__ 81f80f4a74SSebastian Grimberg void magma_gradn_2d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU, 82f80f4a74SSebastian Grimberg const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) { 83f80f4a74SSebastian Grimberg MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 84f80f4a74SSebastian Grimberg 85f80f4a74SSebastian Grimberg const int tx = threadIdx.x; 86f80f4a74SSebastian Grimberg const int ty = threadIdx.y; 87f80f4a74SSebastian Grimberg const int elem_id = (blockIdx.x * blockDim.y) + ty; 88f80f4a74SSebastian Grimberg magma_trans_t transT = MagmaNoTrans; 89f80f4a74SSebastian Grimberg 90f80f4a74SSebastian Grimberg if (elem_id >= nelem) return; 91f80f4a74SSebastian Grimberg 923c1e2affSSebastian Grimberg CeedScalar rU[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // here DIM_U = 1, but might be different for a fused operator 933c1e2affSSebastian Grimberg CeedScalar rV[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // here DIM_V = 1, but might be different for a fused operator 94f80f4a74SSebastian Grimberg CeedScalar rTmp = 0.0; 95f80f4a74SSebastian Grimberg 96f80f4a74SSebastian Grimberg // shift global memory pointers by elem stride 97f80f4a74SSebastian Grimberg dU += elem_id * estrdU; 98f80f4a74SSebastian Grimberg dV += elem_id * estrdV; 99f80f4a74SSebastian Grimberg 100f80f4a74SSebastian Grimberg // assign shared memory pointers 1013c1e2affSSebastian Grimberg CeedScalar *sTinterp = (CeedScalar *)shared_data; 1023c1e2affSSebastian Grimberg CeedScalar *sTgrad = sTinterp + BASIS_P * BASIS_Q; 1033c1e2affSSebastian Grimberg CeedScalar *sTmp = sTgrad + BASIS_P * BASIS_Q; 1043c1e2affSSebastian Grimberg sTmp += ty * (BASIS_P * BASIS_MAX_P_Q); 105f80f4a74SSebastian Grimberg 106f80f4a74SSebastian Grimberg // read T 107f80f4a74SSebastian Grimberg if (ty == 0) { 108*9e0c01faSSebastian Grimberg read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dinterp1d, sTinterp); 109*9e0c01faSSebastian Grimberg read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dgrad1d, sTgrad); 110f80f4a74SSebastian Grimberg } 111f80f4a74SSebastian Grimberg 112f80f4a74SSebastian Grimberg // No need to read V ( required only in transposed grad ) 113f80f4a74SSebastian Grimberg const CeedScalar beta = 0.0; 114f80f4a74SSebastian Grimberg 1153c1e2affSSebastian Grimberg /* read U (idim = 0 for dU, i_DIM = 0 for rU) -- 116f80f4a74SSebastian Grimberg there is a sync at the end of this function */ 117*9e0c01faSSebastian Grimberg read_U_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx); 118f80f4a74SSebastian Grimberg 1193c1e2affSSebastian Grimberg /* first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) -- 120f80f4a74SSebastian Grimberg output from rV[0][][] into dV (idim = 0) */ 1213c1e2affSSebastian Grimberg magma_grad_2d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_P, BASIS_Q, BASIS_P, BASIS_Q, 0, 0, 0>(sTinterp, sTgrad, rU, rV, beta, tx, rTmp, sTmp); 122f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_2d_device */ 123*9e0c01faSSebastian Grimberg write_V_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (0 * dstrdV), cstrdV, rV, tx); 124f80f4a74SSebastian Grimberg 1253c1e2affSSebastian Grimberg /* second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) -- 126f80f4a74SSebastian Grimberg output from rV[0][][] into dV (idim = 1) */ 1273c1e2affSSebastian Grimberg magma_grad_2d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_P, BASIS_Q, BASIS_P, BASIS_Q, 1, 0, 0>(sTinterp, sTgrad, rU, rV, beta, tx, rTmp, sTmp); 128f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_2d_device */ 129*9e0c01faSSebastian Grimberg write_V_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (1 * dstrdV), cstrdV, rV, tx); 130f80f4a74SSebastian Grimberg } 131f80f4a74SSebastian Grimberg 132*9e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 1333c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__ 134f80f4a74SSebastian Grimberg void magma_gradt_2d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU, 135f80f4a74SSebastian Grimberg const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) { 136f80f4a74SSebastian Grimberg MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 137f80f4a74SSebastian Grimberg 138f80f4a74SSebastian Grimberg const int tx = threadIdx.x; 139f80f4a74SSebastian Grimberg const int ty = threadIdx.y; 140f80f4a74SSebastian Grimberg const int elem_id = (blockIdx.x * blockDim.y) + ty; 141f80f4a74SSebastian Grimberg magma_trans_t transT = MagmaTrans; 142f80f4a74SSebastian Grimberg 143f80f4a74SSebastian Grimberg if (elem_id >= nelem) return; 144f80f4a74SSebastian Grimberg 1453c1e2affSSebastian Grimberg CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // here DIM_U = 1, but might be different for a fused operator 1463c1e2affSSebastian Grimberg CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // here DIM_V = 1, but might be different for a fused operator 147f80f4a74SSebastian Grimberg CeedScalar rTmp = 0.0; 148f80f4a74SSebastian Grimberg 149f80f4a74SSebastian Grimberg // shift global memory pointers by elem stride 150f80f4a74SSebastian Grimberg dU += elem_id * estrdU; 151f80f4a74SSebastian Grimberg dV += elem_id * estrdV; 152f80f4a74SSebastian Grimberg 153f80f4a74SSebastian Grimberg // assign shared memory pointers 1543c1e2affSSebastian Grimberg CeedScalar *sTinterp = (CeedScalar *)shared_data; 1553c1e2affSSebastian Grimberg CeedScalar *sTgrad = sTinterp + BASIS_Q * BASIS_P; 1563c1e2affSSebastian Grimberg CeedScalar *sTmp = sTgrad + BASIS_Q * BASIS_P; 1573c1e2affSSebastian Grimberg sTmp += ty * (BASIS_Q * BASIS_MAX_P_Q); 158f80f4a74SSebastian Grimberg 159f80f4a74SSebastian Grimberg // read T 160f80f4a74SSebastian Grimberg if (ty == 0) { 161*9e0c01faSSebastian Grimberg read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dinterp1d, sTinterp); 162*9e0c01faSSebastian Grimberg read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dgrad1d, sTgrad); 163f80f4a74SSebastian Grimberg } 164f80f4a74SSebastian Grimberg __syncthreads(); 165f80f4a74SSebastian Grimberg 166f80f4a74SSebastian Grimberg /* read V (since this is transposed mode -- 1673c1e2affSSebastian Grimberg idim = 0 for dV, i_DIM = 0 for rV) */ 168f80f4a74SSebastian Grimberg const CeedScalar beta = 1.0; 169*9e0c01faSSebastian Grimberg read_V_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV + (0 * dstrdV), cstrdV, rV, tx); 170f80f4a74SSebastian Grimberg 1713c1e2affSSebastian Grimberg /* read U (idim = 0 for dU, i_DIM = 0 for rU) -- 172f80f4a74SSebastian Grimberg there is a sync at the end of this function */ 173*9e0c01faSSebastian Grimberg read_U_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx); 1743c1e2affSSebastian Grimberg /* first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) */ 1753c1e2affSSebastian Grimberg magma_grad_2d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P, 0, 0, 0>(sTinterp, sTgrad, rU, rV, beta, tx, rTmp, sTmp); 176f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_2d_device */ 177f80f4a74SSebastian Grimberg 1783c1e2affSSebastian Grimberg /* read U (idim = 1 for dU, i_DIM = 0 for rU) -- 179f80f4a74SSebastian Grimberg there is a sync at the end of this function */ 180*9e0c01faSSebastian Grimberg read_U_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (1 * dstrdU), cstrdU, rU, sTmp, tx); 1813c1e2affSSebastian Grimberg /* second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) */ 1823c1e2affSSebastian Grimberg magma_grad_2d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P, 1, 0, 0>(sTinterp, sTgrad, rU, rV, beta, tx, rTmp, sTmp); 183f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_2d_device */ 184f80f4a74SSebastian Grimberg 185f80f4a74SSebastian Grimberg // write V 186*9e0c01faSSebastian Grimberg write_V_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV + (0 * dstrdV), cstrdV, rV, tx); 187f80f4a74SSebastian Grimberg } 1883c1e2affSSebastian Grimberg 1893c1e2affSSebastian Grimberg #endif // CEED_MAGMA_BASIS_GRAD_2D_H 190