1*5aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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 199e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 207132caa0SSebastian Grimberg // Helper function to add or set into V 217132caa0SSebastian Grimberg template <typename T, bool Add> 227132caa0SSebastian Grimberg struct magma_grad_2d_device_accumulate; 237132caa0SSebastian Grimberg 247132caa0SSebastian Grimberg template <typename T> 257132caa0SSebastian Grimberg struct magma_grad_2d_device_accumulate<T, true> { 267132caa0SSebastian Grimberg static __device__ __inline__ void op(T &rV, const T &rTmp) { rV += rTmp; } 277132caa0SSebastian Grimberg }; 287132caa0SSebastian Grimberg 297132caa0SSebastian Grimberg template <typename T> 307132caa0SSebastian Grimberg struct magma_grad_2d_device_accumulate<T, false> { 317132caa0SSebastian Grimberg static __device__ __inline__ void op(T &rV, const T &rTmp) { rV = rTmp; } 327132caa0SSebastian Grimberg }; 337132caa0SSebastian Grimberg 347132caa0SSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 35f80f4a74SSebastian Grimberg // grad basis action (2D) 36f80f4a74SSebastian Grimberg // This function is called two times at a higher level for 2D 373c1e2affSSebastian Grimberg // DIM_U -- for the size of rU[DIM_U * NUM_COMP * MAX_P_Q] 383c1e2affSSebastian Grimberg // DIM_V -- for the size of rV[DIM_V * NUM_COMP * MAX_P_Q] 393c1e2affSSebastian Grimberg // i_DIM -- the index of the outermost loop over dimensions in grad 403c1e2affSSebastian Grimberg // i_DIM_U -- which dim index of rU is accessed (always 0 for notrans, 0 or 1 for trans) 413c1e2affSSebastian Grimberg // i_DIM_V -- which dim index of rV is accessed (0 or 1 for notrans, always 0 for trans) 427132caa0SSebastian 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, bool ADD> 433c1e2affSSebastian Grimberg static __device__ __inline__ void magma_grad_2d_device(const T *sTinterp, const T *sTgrad, T rU[DIM_U][NUM_COMP][rU_SIZE], 447132caa0SSebastian Grimberg T rV[DIM_V][NUM_COMP][rV_SIZE], const int tx, T rTmp, T *swork) { 45f80f4a74SSebastian Grimberg // Assumptions 463c1e2affSSebastian Grimberg // 0. This device routine applies grad for one dim only (i_DIM), so it should be called twice for 2D 473c1e2affSSebastian Grimberg // 1. 1D threads of size max(P,Q) 483c1e2affSSebastian Grimberg // 2. input: rU[DIM_U x NUM_COMP x P] in registers (per thread) 493c1e2affSSebastian Grimberg // 3. output: rV[DIM_V x NUM_COMP x Q] in registers (per thread) 50f80f4a74SSebastian Grimberg // 4. Two products per each (dim,component) pair 513c1e2affSSebastian Grimberg // 4.1 Batch P of (1xP) matrices times (PxQ) matrix => Batch P of (1xQ) matrices 523c1e2affSSebastian Grimberg // 4.2 Batch 1 of (QxP) matrix times (PxQ) matrix => (QxQ) matrix 53f80f4a74SSebastian Grimberg // 6. Each thread computes one row of the output of each product 54f80f4a74SSebastian Grimberg // 7. Sync is recommended before and after the call 55f80f4a74SSebastian Grimberg 563c1e2affSSebastian Grimberg for (int comp = 0; comp < NUM_COMP; comp++) { 573c1e2affSSebastian Grimberg // 1st product -- Batch P of (1xP) matrices [reg] x (PxQ) [shmem] => Batch P of (1xQ) matrices 583c1e2affSSebastian Grimberg // the batch output P x (1xQ) is written on the fly to shmem 593c1e2affSSebastian Grimberg if (tx < P) { 60f80f4a74SSebastian Grimberg const int batchid = tx; 61f80f4a74SSebastian Grimberg const int sld = 1; 623c1e2affSSebastian Grimberg const T *sT = (i_DIM == 0) ? sTgrad : sTinterp; 633c1e2affSSebastian Grimberg T *sTmp = swork + batchid * (1 * Q); 643c1e2affSSebastian Grimberg for (int j = 0; j < Q; j++) { 65f80f4a74SSebastian Grimberg rTmp = 0.0; 663c1e2affSSebastian Grimberg for (int i = 0; i < P; i++) { 673c1e2affSSebastian Grimberg rTmp += rU[i_DIM_U][comp][i] * sT(i, j); 68f80f4a74SSebastian Grimberg } 69f80f4a74SSebastian Grimberg sTmp(0, j, sld) = rTmp; 70f80f4a74SSebastian Grimberg } 713c1e2affSSebastian Grimberg } // end of: if (tx < P) 72f80f4a74SSebastian Grimberg __syncthreads(); 73f80f4a74SSebastian Grimberg 743c1e2affSSebastian Grimberg // 2nd product -- Batch 1 of a (QxP) matrix [shmem] x (PxQ) [shmem] => (QxQ) matrix [reg] 753c1e2affSSebastian Grimberg if (tx < Q) { 76f80f4a74SSebastian Grimberg const int batchid = 0; 773c1e2affSSebastian Grimberg const int sld = Q; 783c1e2affSSebastian Grimberg const T *sT = (i_DIM == 1) ? sTgrad : sTinterp; 793c1e2affSSebastian Grimberg T *sTmp = swork + batchid * (Q * P); 803c1e2affSSebastian Grimberg for (int j = 0; j < Q; j++) { 81f80f4a74SSebastian Grimberg rTmp = 0.0; 823c1e2affSSebastian Grimberg for (int i = 0; i < P; i++) { 83f80f4a74SSebastian Grimberg rTmp += sTmp(tx, i, sld) * sT(i, j); 84f80f4a74SSebastian Grimberg } 857132caa0SSebastian Grimberg magma_grad_2d_device_accumulate<T, ADD>::op(rV[i_DIM_V][comp][j], rTmp); 86f80f4a74SSebastian Grimberg } 87f80f4a74SSebastian Grimberg } 88f80f4a74SSebastian Grimberg __syncthreads(); 893c1e2affSSebastian Grimberg } // loop over NUM_COMP 90f80f4a74SSebastian Grimberg } 91f80f4a74SSebastian Grimberg 929e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 933c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__ 94f80f4a74SSebastian Grimberg void magma_gradn_2d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU, 95f80f4a74SSebastian Grimberg const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) { 96f80f4a74SSebastian Grimberg MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 97f80f4a74SSebastian Grimberg 98f80f4a74SSebastian Grimberg const int tx = threadIdx.x; 99f80f4a74SSebastian Grimberg const int ty = threadIdx.y; 100f80f4a74SSebastian Grimberg const int elem_id = (blockIdx.x * blockDim.y) + ty; 101f80f4a74SSebastian Grimberg 102f80f4a74SSebastian Grimberg if (elem_id >= nelem) return; 103f80f4a74SSebastian Grimberg 1043c1e2affSSebastian Grimberg CeedScalar rU[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // here DIM_U = 1, but might be different for a fused operator 1053c1e2affSSebastian Grimberg CeedScalar rV[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // here DIM_V = 1, but might be different for a fused operator 106f80f4a74SSebastian Grimberg CeedScalar rTmp = 0.0; 107f80f4a74SSebastian Grimberg 108f80f4a74SSebastian Grimberg // shift global memory pointers by elem stride 109f80f4a74SSebastian Grimberg dU += elem_id * estrdU; 110f80f4a74SSebastian Grimberg dV += elem_id * estrdV; 111f80f4a74SSebastian Grimberg 112f80f4a74SSebastian Grimberg // assign shared memory pointers 1133c1e2affSSebastian Grimberg CeedScalar *sTinterp = (CeedScalar *)shared_data; 1143c1e2affSSebastian Grimberg CeedScalar *sTgrad = sTinterp + BASIS_P * BASIS_Q; 1153c1e2affSSebastian Grimberg CeedScalar *sTmp = sTgrad + BASIS_P * BASIS_Q; 1163c1e2affSSebastian Grimberg sTmp += ty * (BASIS_P * BASIS_MAX_P_Q); 117f80f4a74SSebastian Grimberg 118f80f4a74SSebastian Grimberg // read T 119f80f4a74SSebastian Grimberg if (ty == 0) { 1209e0c01faSSebastian Grimberg read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dinterp1d, sTinterp); 1219e0c01faSSebastian Grimberg read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dgrad1d, sTgrad); 122f80f4a74SSebastian Grimberg } 123f80f4a74SSebastian Grimberg 1243c1e2affSSebastian Grimberg /* read U (idim = 0 for dU, i_DIM = 0 for rU) -- 125f80f4a74SSebastian Grimberg there is a sync at the end of this function */ 1269e0c01faSSebastian Grimberg read_U_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx); 127f80f4a74SSebastian Grimberg 1283c1e2affSSebastian Grimberg /* first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) -- 129f80f4a74SSebastian Grimberg output from rV[0][][] into dV (idim = 0) */ 1307132caa0SSebastian Grimberg 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, 1317132caa0SSebastian Grimberg sTmp); 132f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_2d_device */ 1339e0c01faSSebastian Grimberg write_V_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (0 * dstrdV), cstrdV, rV, tx); 134f80f4a74SSebastian Grimberg 1353c1e2affSSebastian Grimberg /* second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) -- 136f80f4a74SSebastian Grimberg output from rV[0][][] into dV (idim = 1) */ 1377132caa0SSebastian Grimberg 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, 1387132caa0SSebastian Grimberg sTmp); 139f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_2d_device */ 1409e0c01faSSebastian Grimberg write_V_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (1 * dstrdV), cstrdV, rV, tx); 141f80f4a74SSebastian Grimberg } 142f80f4a74SSebastian Grimberg 1439e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 1443c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__ 145f80f4a74SSebastian Grimberg void magma_gradt_2d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU, 146f80f4a74SSebastian Grimberg const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) { 147f80f4a74SSebastian Grimberg MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 148f80f4a74SSebastian Grimberg 149f80f4a74SSebastian Grimberg const int tx = threadIdx.x; 150f80f4a74SSebastian Grimberg const int ty = threadIdx.y; 151f80f4a74SSebastian Grimberg const int elem_id = (blockIdx.x * blockDim.y) + ty; 152f80f4a74SSebastian Grimberg 153f80f4a74SSebastian Grimberg if (elem_id >= nelem) return; 154f80f4a74SSebastian Grimberg 1553c1e2affSSebastian Grimberg CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // here DIM_U = 1, but might be different for a fused operator 1563c1e2affSSebastian Grimberg CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // here DIM_V = 1, but might be different for a fused operator 157f80f4a74SSebastian Grimberg CeedScalar rTmp = 0.0; 158f80f4a74SSebastian Grimberg 159f80f4a74SSebastian Grimberg // shift global memory pointers by elem stride 160f80f4a74SSebastian Grimberg dU += elem_id * estrdU; 161f80f4a74SSebastian Grimberg dV += elem_id * estrdV; 162f80f4a74SSebastian Grimberg 163f80f4a74SSebastian Grimberg // assign shared memory pointers 1643c1e2affSSebastian Grimberg CeedScalar *sTinterp = (CeedScalar *)shared_data; 1653c1e2affSSebastian Grimberg CeedScalar *sTgrad = sTinterp + BASIS_Q * BASIS_P; 1663c1e2affSSebastian Grimberg CeedScalar *sTmp = sTgrad + BASIS_Q * BASIS_P; 1673c1e2affSSebastian Grimberg sTmp += ty * (BASIS_Q * BASIS_MAX_P_Q); 168f80f4a74SSebastian Grimberg 169f80f4a74SSebastian Grimberg // read T 170f80f4a74SSebastian Grimberg if (ty == 0) { 1719e0c01faSSebastian Grimberg read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dinterp1d, sTinterp); 1729e0c01faSSebastian Grimberg read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dgrad1d, sTgrad); 173f80f4a74SSebastian Grimberg } 174f80f4a74SSebastian Grimberg __syncthreads(); 175f80f4a74SSebastian Grimberg 1763c1e2affSSebastian Grimberg /* read U (idim = 0 for dU, i_DIM = 0 for rU) -- 177f80f4a74SSebastian Grimberg there is a sync at the end of this function */ 1789e0c01faSSebastian Grimberg read_U_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx); 1793c1e2affSSebastian Grimberg /* first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) */ 1807132caa0SSebastian Grimberg 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); 181f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_2d_device */ 182f80f4a74SSebastian Grimberg 1833c1e2affSSebastian Grimberg /* read U (idim = 1 for dU, i_DIM = 0 for rU) -- 184f80f4a74SSebastian Grimberg there is a sync at the end of this function */ 1859e0c01faSSebastian Grimberg read_U_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (1 * dstrdU), cstrdU, rU, sTmp, tx); 1863c1e2affSSebastian Grimberg /* second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) */ 1877132caa0SSebastian Grimberg 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); 188f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_2d_device */ 189f80f4a74SSebastian Grimberg 190f80f4a74SSebastian Grimberg // write V 1919e0c01faSSebastian Grimberg write_V_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV + (0 * dstrdV), cstrdV, rV, tx); 192f80f4a74SSebastian Grimberg } 1933c1e2affSSebastian Grimberg 1943c1e2affSSebastian Grimberg #endif // CEED_MAGMA_BASIS_GRAD_2D_H 195