1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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 #include "magma-common-tensor.h" 113c1e2affSSebastian Grimberg 12f80f4a74SSebastian Grimberg // macros to abstract access of shared memory and reg. file 133c1e2affSSebastian Grimberg #define sT(i, j) sT[(j) * P + (i)] 14f80f4a74SSebastian Grimberg #define sTmp(i, j, ldw) sTmp[(j) * (ldw) + (i)] 15f80f4a74SSebastian Grimberg 169e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 177132caa0SSebastian Grimberg // Helper function to add or set into V 187132caa0SSebastian Grimberg template <typename T, bool Add> 197132caa0SSebastian Grimberg struct magma_grad_2d_device_accumulate; 207132caa0SSebastian Grimberg 217132caa0SSebastian Grimberg template <typename T> 227132caa0SSebastian Grimberg struct magma_grad_2d_device_accumulate<T, true> { 237132caa0SSebastian Grimberg static __device__ __inline__ void op(T &rV, const T &rTmp) { rV += rTmp; } 247132caa0SSebastian Grimberg }; 257132caa0SSebastian Grimberg 267132caa0SSebastian Grimberg template <typename T> 277132caa0SSebastian Grimberg struct magma_grad_2d_device_accumulate<T, false> { 287132caa0SSebastian Grimberg static __device__ __inline__ void op(T &rV, const T &rTmp) { rV = rTmp; } 297132caa0SSebastian Grimberg }; 307132caa0SSebastian Grimberg 317132caa0SSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 32f80f4a74SSebastian Grimberg // grad basis action (2D) 33f80f4a74SSebastian Grimberg // This function is called two times at a higher level for 2D 343c1e2affSSebastian Grimberg // DIM_U -- for the size of rU[DIM_U * NUM_COMP * MAX_P_Q] 353c1e2affSSebastian Grimberg // DIM_V -- for the size of rV[DIM_V * NUM_COMP * MAX_P_Q] 363c1e2affSSebastian Grimberg // i_DIM -- the index of the outermost loop over dimensions in grad 373c1e2affSSebastian Grimberg // i_DIM_U -- which dim index of rU is accessed (always 0 for notrans, 0 or 1 for trans) 383c1e2affSSebastian Grimberg // i_DIM_V -- which dim index of rV is accessed (0 or 1 for notrans, always 0 for trans) 397132caa0SSebastian 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> 403c1e2affSSebastian Grimberg static __device__ __inline__ void magma_grad_2d_device(const T *sTinterp, const T *sTgrad, T rU[DIM_U][NUM_COMP][rU_SIZE], 417132caa0SSebastian Grimberg T rV[DIM_V][NUM_COMP][rV_SIZE], const int tx, T rTmp, T *swork) { 42f80f4a74SSebastian Grimberg // Assumptions 433c1e2affSSebastian Grimberg // 0. This device routine applies grad for one dim only (i_DIM), so it should be called twice for 2D 443c1e2affSSebastian Grimberg // 1. 1D threads of size max(P,Q) 453c1e2affSSebastian Grimberg // 2. input: rU[DIM_U x NUM_COMP x P] in registers (per thread) 463c1e2affSSebastian Grimberg // 3. output: rV[DIM_V x NUM_COMP x Q] in registers (per thread) 47f80f4a74SSebastian Grimberg // 4. Two products per each (dim,component) pair 483c1e2affSSebastian Grimberg // 4.1 Batch P of (1xP) matrices times (PxQ) matrix => Batch P of (1xQ) matrices 493c1e2affSSebastian Grimberg // 4.2 Batch 1 of (QxP) matrix times (PxQ) matrix => (QxQ) matrix 50f80f4a74SSebastian Grimberg // 6. Each thread computes one row of the output of each product 51f80f4a74SSebastian Grimberg // 7. Sync is recommended before and after the call 52f80f4a74SSebastian Grimberg 533c1e2affSSebastian Grimberg for (int comp = 0; comp < NUM_COMP; comp++) { 543c1e2affSSebastian Grimberg // 1st product -- Batch P of (1xP) matrices [reg] x (PxQ) [shmem] => Batch P of (1xQ) matrices 553c1e2affSSebastian Grimberg // the batch output P x (1xQ) is written on the fly to shmem 563c1e2affSSebastian Grimberg if (tx < P) { 57f80f4a74SSebastian Grimberg const int batchid = tx; 58f80f4a74SSebastian Grimberg const int sld = 1; 593c1e2affSSebastian Grimberg const T *sT = (i_DIM == 0) ? sTgrad : sTinterp; 603c1e2affSSebastian Grimberg T *sTmp = swork + batchid * (1 * Q); 613c1e2affSSebastian Grimberg for (int j = 0; j < Q; j++) { 62f80f4a74SSebastian Grimberg rTmp = 0.0; 633c1e2affSSebastian Grimberg for (int i = 0; i < P; i++) { 643c1e2affSSebastian Grimberg rTmp += rU[i_DIM_U][comp][i] * sT(i, j); 65f80f4a74SSebastian Grimberg } 66f80f4a74SSebastian Grimberg sTmp(0, j, sld) = rTmp; 67f80f4a74SSebastian Grimberg } 683c1e2affSSebastian Grimberg } // end of: if (tx < P) 69f80f4a74SSebastian Grimberg __syncthreads(); 70f80f4a74SSebastian Grimberg 713c1e2affSSebastian Grimberg // 2nd product -- Batch 1 of a (QxP) matrix [shmem] x (PxQ) [shmem] => (QxQ) matrix [reg] 723c1e2affSSebastian Grimberg if (tx < Q) { 73f80f4a74SSebastian Grimberg const int batchid = 0; 743c1e2affSSebastian Grimberg const int sld = Q; 753c1e2affSSebastian Grimberg const T *sT = (i_DIM == 1) ? sTgrad : sTinterp; 763c1e2affSSebastian Grimberg T *sTmp = swork + batchid * (Q * P); 773c1e2affSSebastian Grimberg for (int j = 0; j < Q; j++) { 78f80f4a74SSebastian Grimberg rTmp = 0.0; 793c1e2affSSebastian Grimberg for (int i = 0; i < P; i++) { 80f80f4a74SSebastian Grimberg rTmp += sTmp(tx, i, sld) * sT(i, j); 81f80f4a74SSebastian Grimberg } 827132caa0SSebastian Grimberg magma_grad_2d_device_accumulate<T, ADD>::op(rV[i_DIM_V][comp][j], rTmp); 83f80f4a74SSebastian Grimberg } 84f80f4a74SSebastian Grimberg } 85f80f4a74SSebastian Grimberg __syncthreads(); 863c1e2affSSebastian Grimberg } // loop over NUM_COMP 87f80f4a74SSebastian Grimberg } 88f80f4a74SSebastian Grimberg 899e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 903c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__ 91f80f4a74SSebastian Grimberg void magma_gradn_2d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU, 92f80f4a74SSebastian Grimberg const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) { 93f80f4a74SSebastian Grimberg MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 94f80f4a74SSebastian Grimberg 95f80f4a74SSebastian Grimberg const int tx = threadIdx.x; 96f80f4a74SSebastian Grimberg const int ty = threadIdx.y; 97f80f4a74SSebastian Grimberg const int elem_id = (blockIdx.x * blockDim.y) + ty; 98f80f4a74SSebastian Grimberg 99f80f4a74SSebastian Grimberg if (elem_id >= nelem) return; 100f80f4a74SSebastian Grimberg 1013c1e2affSSebastian Grimberg CeedScalar rU[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // here DIM_U = 1, but might be different for a fused operator 1023c1e2affSSebastian Grimberg CeedScalar rV[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // here DIM_V = 1, but might be different for a fused operator 103f80f4a74SSebastian Grimberg CeedScalar rTmp = 0.0; 104f80f4a74SSebastian Grimberg 105f80f4a74SSebastian Grimberg // shift global memory pointers by elem stride 106f80f4a74SSebastian Grimberg dU += elem_id * estrdU; 107f80f4a74SSebastian Grimberg dV += elem_id * estrdV; 108f80f4a74SSebastian Grimberg 109f80f4a74SSebastian Grimberg // assign shared memory pointers 1103c1e2affSSebastian Grimberg CeedScalar *sTinterp = (CeedScalar *)shared_data; 1113c1e2affSSebastian Grimberg CeedScalar *sTgrad = sTinterp + BASIS_P * BASIS_Q; 1123c1e2affSSebastian Grimberg CeedScalar *sTmp = sTgrad + BASIS_P * BASIS_Q; 1133c1e2affSSebastian Grimberg sTmp += ty * (BASIS_P * BASIS_MAX_P_Q); 114f80f4a74SSebastian Grimberg 115f80f4a74SSebastian Grimberg // read T 116f80f4a74SSebastian Grimberg if (ty == 0) { 1179e0c01faSSebastian Grimberg read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dinterp1d, sTinterp); 1189e0c01faSSebastian Grimberg read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dgrad1d, sTgrad); 119f80f4a74SSebastian Grimberg } 120f80f4a74SSebastian Grimberg 1213c1e2affSSebastian Grimberg /* read U (idim = 0 for dU, i_DIM = 0 for rU) -- 122f80f4a74SSebastian Grimberg there is a sync at the end of this function */ 1239e0c01faSSebastian Grimberg read_U_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx); 124f80f4a74SSebastian Grimberg 1253c1e2affSSebastian Grimberg /* first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) -- 126f80f4a74SSebastian Grimberg output from rV[0][][] into dV (idim = 0) */ 1277132caa0SSebastian 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, 1287132caa0SSebastian Grimberg sTmp); 129f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_2d_device */ 1309e0c01faSSebastian Grimberg write_V_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (0 * dstrdV), cstrdV, rV, tx); 131f80f4a74SSebastian Grimberg 1323c1e2affSSebastian Grimberg /* second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) -- 133f80f4a74SSebastian Grimberg output from rV[0][][] into dV (idim = 1) */ 1347132caa0SSebastian 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, 1357132caa0SSebastian Grimberg sTmp); 136f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_2d_device */ 1379e0c01faSSebastian Grimberg write_V_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (1 * dstrdV), cstrdV, rV, tx); 138f80f4a74SSebastian Grimberg } 139f80f4a74SSebastian Grimberg 1409e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 1413c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__ 142f80f4a74SSebastian Grimberg void magma_gradt_2d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU, 143f80f4a74SSebastian Grimberg const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) { 144f80f4a74SSebastian Grimberg MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 145f80f4a74SSebastian Grimberg 146f80f4a74SSebastian Grimberg const int tx = threadIdx.x; 147f80f4a74SSebastian Grimberg const int ty = threadIdx.y; 148f80f4a74SSebastian Grimberg const int elem_id = (blockIdx.x * blockDim.y) + ty; 149f80f4a74SSebastian Grimberg 150f80f4a74SSebastian Grimberg if (elem_id >= nelem) return; 151f80f4a74SSebastian Grimberg 1523c1e2affSSebastian Grimberg CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // here DIM_U = 1, but might be different for a fused operator 1533c1e2affSSebastian Grimberg CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // here DIM_V = 1, but might be different for a fused operator 154f80f4a74SSebastian Grimberg CeedScalar rTmp = 0.0; 155f80f4a74SSebastian Grimberg 156f80f4a74SSebastian Grimberg // shift global memory pointers by elem stride 157f80f4a74SSebastian Grimberg dU += elem_id * estrdU; 158f80f4a74SSebastian Grimberg dV += elem_id * estrdV; 159f80f4a74SSebastian Grimberg 160f80f4a74SSebastian Grimberg // assign shared memory pointers 1613c1e2affSSebastian Grimberg CeedScalar *sTinterp = (CeedScalar *)shared_data; 1623c1e2affSSebastian Grimberg CeedScalar *sTgrad = sTinterp + BASIS_Q * BASIS_P; 1633c1e2affSSebastian Grimberg CeedScalar *sTmp = sTgrad + BASIS_Q * BASIS_P; 1643c1e2affSSebastian Grimberg sTmp += ty * (BASIS_Q * BASIS_MAX_P_Q); 165f80f4a74SSebastian Grimberg 166f80f4a74SSebastian Grimberg // read T 167f80f4a74SSebastian Grimberg if (ty == 0) { 1689e0c01faSSebastian Grimberg read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dinterp1d, sTinterp); 1699e0c01faSSebastian Grimberg read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dgrad1d, sTgrad); 170f80f4a74SSebastian Grimberg } 171f80f4a74SSebastian Grimberg __syncthreads(); 172f80f4a74SSebastian Grimberg 1733c1e2affSSebastian Grimberg /* read U (idim = 0 for dU, i_DIM = 0 for rU) -- 174f80f4a74SSebastian Grimberg there is a sync at the end of this function */ 1759e0c01faSSebastian Grimberg read_U_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx); 1763c1e2affSSebastian Grimberg /* first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) */ 1777132caa0SSebastian 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); 178f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_2d_device */ 179f80f4a74SSebastian Grimberg 1803c1e2affSSebastian Grimberg /* read U (idim = 1 for dU, i_DIM = 0 for rU) -- 181f80f4a74SSebastian Grimberg there is a sync at the end of this function */ 1829e0c01faSSebastian Grimberg read_U_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (1 * dstrdU), cstrdU, rU, sTmp, tx); 1833c1e2affSSebastian Grimberg /* second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) */ 1847132caa0SSebastian 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); 185f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_2d_device */ 186f80f4a74SSebastian Grimberg 187f80f4a74SSebastian Grimberg // write V 1889e0c01faSSebastian Grimberg write_V_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV + (0 * dstrdV), cstrdV, rV, tx); 189f80f4a74SSebastian Grimberg } 190db2becc9SJeremy L Thompson 191db2becc9SJeremy L Thompson //////////////////////////////////////////////////////////////////////////////// 192db2becc9SJeremy L Thompson extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__ 193db2becc9SJeremy L Thompson void magma_gradta_2d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU, 194db2becc9SJeremy L Thompson const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) { 195db2becc9SJeremy L Thompson MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 196db2becc9SJeremy L Thompson 197db2becc9SJeremy L Thompson const int tx = threadIdx.x; 198db2becc9SJeremy L Thompson const int ty = threadIdx.y; 199db2becc9SJeremy L Thompson const int elem_id = (blockIdx.x * blockDim.y) + ty; 200db2becc9SJeremy L Thompson 201db2becc9SJeremy L Thompson if (elem_id >= nelem) return; 202db2becc9SJeremy L Thompson 203db2becc9SJeremy L Thompson CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // here DIM_U = 1, but might be different for a fused operator 204db2becc9SJeremy L Thompson CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // here DIM_V = 1, but might be different for a fused operator 205db2becc9SJeremy L Thompson CeedScalar rTmp = 0.0; 206db2becc9SJeremy L Thompson 207db2becc9SJeremy L Thompson // shift global memory pointers by elem stride 208db2becc9SJeremy L Thompson dU += elem_id * estrdU; 209db2becc9SJeremy L Thompson dV += elem_id * estrdV; 210db2becc9SJeremy L Thompson 211db2becc9SJeremy L Thompson // assign shared memory pointers 212db2becc9SJeremy L Thompson CeedScalar *sTinterp = (CeedScalar *)shared_data; 213db2becc9SJeremy L Thompson CeedScalar *sTgrad = sTinterp + BASIS_Q * BASIS_P; 214db2becc9SJeremy L Thompson CeedScalar *sTmp = sTgrad + BASIS_Q * BASIS_P; 215db2becc9SJeremy L Thompson sTmp += ty * (BASIS_Q * BASIS_MAX_P_Q); 216db2becc9SJeremy L Thompson 217db2becc9SJeremy L Thompson // read T 218db2becc9SJeremy L Thompson if (ty == 0) { 219db2becc9SJeremy L Thompson read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dinterp1d, sTinterp); 220db2becc9SJeremy L Thompson read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dgrad1d, sTgrad); 221db2becc9SJeremy L Thompson } 222db2becc9SJeremy L Thompson __syncthreads(); 223db2becc9SJeremy L Thompson 224db2becc9SJeremy L Thompson /* read U (idim = 0 for dU, i_DIM = 0 for rU) -- 225db2becc9SJeremy L Thompson there is a sync at the end of this function */ 226db2becc9SJeremy L Thompson read_U_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx); 227db2becc9SJeremy L Thompson /* first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) */ 228db2becc9SJeremy L Thompson 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); 229db2becc9SJeremy L Thompson /* there is a sync at the end of magma_grad_2d_device */ 230db2becc9SJeremy L Thompson 231db2becc9SJeremy L Thompson /* read U (idim = 1 for dU, i_DIM = 0 for rU) -- 232db2becc9SJeremy L Thompson there is a sync at the end of this function */ 233db2becc9SJeremy L Thompson read_U_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (1 * dstrdU), cstrdU, rU, sTmp, tx); 234db2becc9SJeremy L Thompson /* second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) */ 235db2becc9SJeremy L Thompson 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); 236db2becc9SJeremy L Thompson /* there is a sync at the end of magma_grad_2d_device */ 237db2becc9SJeremy L Thompson 238db2becc9SJeremy L Thompson // sum into V 239db2becc9SJeremy L Thompson sum_V_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV + (0 * dstrdV), cstrdV, rV, tx); 240db2becc9SJeremy L Thompson } 241