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 3D 103c1e2affSSebastian Grimberg #ifndef CEED_MAGMA_BASIS_GRAD_3D_H 113c1e2affSSebastian Grimberg #define CEED_MAGMA_BASIS_GRAD_3D_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 #define sTmp2(i, j, ldw) sTmp2[(j) * (ldw) + (i)] 19f80f4a74SSebastian Grimberg 20*9e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 21f80f4a74SSebastian Grimberg // grad basis action (3D) 22f80f4a74SSebastian Grimberg // This function is called three times at a higher level for 3D 233c1e2affSSebastian Grimberg // DIM_U -- for the size of rU[DIM_U * NUM_COMP * MAX_P_Q] 243c1e2affSSebastian Grimberg // DIM_V -- for the size of rV[DIM_V * NUM_COMP * MAX_P_Q] 253c1e2affSSebastian Grimberg // i_DIM -- the index of the outermost loop over dimensions in grad 263c1e2affSSebastian Grimberg // i_DIM_U -- which dim index of rU is accessed (always 0 for notrans, 0, 1, or 2 for trans) 273c1e2affSSebastian Grimberg // i_DIM_V -- which dim index of rV is accessed (0, 1, or 2 for notrans, always 0 for trans) 28f80f4a74SSebastian Grimberg // the scalar beta is used to specify whether to accumulate to rV, or overwrite it 293c1e2affSSebastian 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> 303c1e2affSSebastian Grimberg static __device__ __inline__ void magma_grad_3d_device(const T *sTinterp, const T *sTgrad, T rU[DIM_U][NUM_COMP][rU_SIZE], 313c1e2affSSebastian Grimberg T rV[DIM_V][NUM_COMP][rV_SIZE], T beta, const int tx, T rTmp, T *swork) { 32f80f4a74SSebastian Grimberg // Assumptions 333c1e2affSSebastian Grimberg // 0. This device routine applies grad for one dim only (i_DIM), so it should be thrice for 3D 343c1e2affSSebastian Grimberg // 1. 1D threads of size max(P,Q)^2 353c1e2affSSebastian Grimberg // 2. input: rU[DIM_U x NUM_COMP x rU_SIZE] in registers (per thread) 363c1e2affSSebastian Grimberg // 3. output: rV[DIM_V x NUM_COMP x rV_SIZE] in registers (per thread) 37f80f4a74SSebastian Grimberg // 4. Three products per each (dim,component) pair 383c1e2affSSebastian Grimberg // 4.1 Batch P^2 of (1xP) matrices times (PxQ) matrix => Batch P^2 of (1xQ) matrices 393c1e2affSSebastian Grimberg // 4.2 Batch P of (QxP) matrices times (PxQ) matrix => Batch P of (QxQ) matrices 403c1e2affSSebastian Grimberg // 4.3 Batch 1 of (Q^2xP_) matrix times (PxQ) matrix => (Q^2xQ_) matrix 41f80f4a74SSebastian Grimberg // 6. Each thread computes one row of the output of each product 42f80f4a74SSebastian Grimberg // 7. Sync is recommended before and after the call 43f80f4a74SSebastian Grimberg 44f80f4a74SSebastian Grimberg T *sW1 = swork; 453c1e2affSSebastian Grimberg T *sW2 = sW1 + P * P * Q; 463c1e2affSSebastian Grimberg for (int comp = 0; comp < NUM_COMP; comp++) { 473c1e2affSSebastian Grimberg // Batch P^2 of (1xP) matrices [reg] times (PxQ) matrix [shmem] => Batch P^2 of (1xQ) matrices [shmem] 483c1e2affSSebastian Grimberg if (tx < (P * P)) { 49f80f4a74SSebastian Grimberg const int batchid = tx; 50f80f4a74SSebastian Grimberg const int sld = 1; 513c1e2affSSebastian Grimberg const T *sT = (i_DIM == 0) ? sTgrad : sTinterp; 523c1e2affSSebastian Grimberg T *sTmp = sW1 + batchid * (1 * Q); 533c1e2affSSebastian Grimberg for (int j = 0; j < Q; j++) { 54f80f4a74SSebastian Grimberg rTmp = 0.0; 553c1e2affSSebastian Grimberg for (int i = 0; i < P; i++) { 563c1e2affSSebastian Grimberg rTmp += rU[i_DIM_U][comp][i] * sT(i, j); 57f80f4a74SSebastian Grimberg } 58f80f4a74SSebastian Grimberg sTmp(0, j, sld) = rTmp; 59f80f4a74SSebastian Grimberg } 603c1e2affSSebastian Grimberg } // end of: if (tx < P*P) 61f80f4a74SSebastian Grimberg __syncthreads(); 62f80f4a74SSebastian Grimberg 633c1e2affSSebastian Grimberg // Batch P of (QxP) matrices [shmem] times (PxQ) matrix [shmem] => Batch P of (QxQ) matrices [reg] 643c1e2affSSebastian Grimberg if (tx < (P * Q)) { 653c1e2affSSebastian Grimberg const int batchid = tx / Q; 663c1e2affSSebastian Grimberg const int tx_ = tx % Q; 673c1e2affSSebastian Grimberg const int sld = Q; 683c1e2affSSebastian Grimberg const T *sT = (i_DIM == 1) ? sTgrad : sTinterp; 693c1e2affSSebastian Grimberg T *sTmp = sW1 + batchid * (Q * P); // sTmp is input 703c1e2affSSebastian Grimberg T *sTmp2 = sW2 + batchid * (Q * Q); // sTmp2 is output 713c1e2affSSebastian Grimberg for (int j = 0; j < Q; j++) { 72f80f4a74SSebastian Grimberg rTmp = 0.0; 733c1e2affSSebastian Grimberg for (int i = 0; i < P; i++) { 74f80f4a74SSebastian Grimberg rTmp += sTmp(tx_, i, sld) * sT(i, j); 75f80f4a74SSebastian Grimberg } 76f80f4a74SSebastian Grimberg sTmp2(tx_, j, sld) = rTmp; 77f80f4a74SSebastian Grimberg } 78f80f4a74SSebastian Grimberg } 79f80f4a74SSebastian Grimberg __syncthreads(); 80f80f4a74SSebastian Grimberg 813c1e2affSSebastian Grimberg // Batch 1 of (Q^2xP_) matrices [shmem] times (PxQ) matrix [shmem] => Batch 1 of (Q^2xQ_) matrices [reg] 823c1e2affSSebastian Grimberg if (tx < (Q * Q)) { 833c1e2affSSebastian Grimberg // No need to declare batchid = (tx / Q^2) = always zero 843c1e2affSSebastian Grimberg // No need to declare tx_ = (tx_ % Q^2) = always tx 853c1e2affSSebastian Grimberg const int sld = Q * Q; 863c1e2affSSebastian Grimberg const T *sT = (i_DIM == 2) ? sTgrad : sTinterp; 87f80f4a74SSebastian Grimberg T *sTmp = sW2; // sTmp is input 883c1e2affSSebastian Grimberg for (int j = 0; j < Q; j++) { 89f80f4a74SSebastian Grimberg rTmp = 0.0; 903c1e2affSSebastian Grimberg for (int i = 0; i < P; i++) { 91f80f4a74SSebastian Grimberg rTmp += sTmp(tx, i, sld) * sT(i, j); 92f80f4a74SSebastian Grimberg } 933c1e2affSSebastian Grimberg rV[i_DIM_V][comp][j] *= beta; 943c1e2affSSebastian Grimberg rV[i_DIM_V][comp][j] += rTmp; 95f80f4a74SSebastian Grimberg } 96f80f4a74SSebastian Grimberg } 97f80f4a74SSebastian Grimberg __syncthreads(); 983c1e2affSSebastian Grimberg } // loop over NUM_COMP 99f80f4a74SSebastian Grimberg } 100f80f4a74SSebastian Grimberg 101*9e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 1023c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q *BASIS_MAX_P_Q, MAGMA_MAXTHREADS_3D)) __global__ 103f80f4a74SSebastian Grimberg void magma_gradn_3d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU, 104f80f4a74SSebastian Grimberg const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) { 105f80f4a74SSebastian Grimberg MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 106f80f4a74SSebastian Grimberg 107f80f4a74SSebastian Grimberg const int tx = threadIdx.x; 108f80f4a74SSebastian Grimberg const int ty = threadIdx.y; 109f80f4a74SSebastian Grimberg const int elem_id = (blockIdx.x * blockDim.y) + ty; 110f80f4a74SSebastian Grimberg magma_trans_t transT = MagmaNoTrans; 111f80f4a74SSebastian Grimberg 112f80f4a74SSebastian Grimberg if (elem_id >= nelem) return; 113f80f4a74SSebastian Grimberg 1143c1e2affSSebastian Grimberg CeedScalar rU[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // here DIM_U = 1, but might be different for a fused operator 1153c1e2affSSebastian Grimberg CeedScalar rV[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // here DIM_V = 1, but might be different for a fused operator 116f80f4a74SSebastian Grimberg CeedScalar rTmp = 0.0; 117f80f4a74SSebastian Grimberg 118f80f4a74SSebastian Grimberg // shift global memory pointers by elem stride 119f80f4a74SSebastian Grimberg dU += elem_id * estrdU; 120f80f4a74SSebastian Grimberg dV += elem_id * estrdV; 121f80f4a74SSebastian Grimberg 122f80f4a74SSebastian Grimberg // assign shared memory pointers 1233c1e2affSSebastian Grimberg CeedScalar *sTinterp = (CeedScalar *)shared_data; 1243c1e2affSSebastian Grimberg CeedScalar *sTgrad = sTinterp + BASIS_P * BASIS_Q; 1253c1e2affSSebastian Grimberg CeedScalar *sTmp = sTgrad + BASIS_P * BASIS_Q; 1263c1e2affSSebastian Grimberg sTmp += ty * (max(BASIS_P * BASIS_P * BASIS_P, (BASIS_P * BASIS_P * BASIS_Q) + (BASIS_P * BASIS_Q * BASIS_Q))); 127f80f4a74SSebastian Grimberg 128f80f4a74SSebastian Grimberg // read T 129f80f4a74SSebastian Grimberg if (ty == 0) { 130*9e0c01faSSebastian Grimberg read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dinterp1d, sTinterp); 131*9e0c01faSSebastian Grimberg read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dgrad1d, sTgrad); 132f80f4a74SSebastian Grimberg } 133f80f4a74SSebastian Grimberg __syncthreads(); 134f80f4a74SSebastian Grimberg 135f80f4a74SSebastian Grimberg // No need to read V ( required only in transposed grad ) 136f80f4a74SSebastian Grimberg const CeedScalar beta = 0.0; 137f80f4a74SSebastian Grimberg 1383c1e2affSSebastian Grimberg /* read U (idim = 0 for dU, i_DIM = 0 for rU) -- 139f80f4a74SSebastian Grimberg there is a sync at the end of this function */ 140*9e0c01faSSebastian Grimberg read_U_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx); 141f80f4a74SSebastian Grimberg 1423c1e2affSSebastian Grimberg /* first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) -- 143f80f4a74SSebastian Grimberg output from rV[0][][] into dV (idim = 0) */ 1443c1e2affSSebastian Grimberg magma_grad_3d_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); 145f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_3d_device */ 146*9e0c01faSSebastian Grimberg write_V_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (0 * dstrdV), cstrdV, rV, tx); 147f80f4a74SSebastian Grimberg 1483c1e2affSSebastian Grimberg /* second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) -- 149f80f4a74SSebastian Grimberg output from rV[0][][] into dV (idim = 1) */ 1503c1e2affSSebastian Grimberg magma_grad_3d_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); 151f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_3d_device */ 152*9e0c01faSSebastian Grimberg write_V_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (1 * dstrdV), cstrdV, rV, tx); 153f80f4a74SSebastian Grimberg 1543c1e2affSSebastian Grimberg /* third call (i_DIM = 2, i_DIM_U = 0, i_DIM_V = 0) -- 155f80f4a74SSebastian Grimberg output from rV[0][][] into dV (idim = 2) */ 1563c1e2affSSebastian Grimberg magma_grad_3d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_P, BASIS_Q, BASIS_P, BASIS_Q, 2, 0, 0>(sTinterp, sTgrad, rU, rV, beta, tx, rTmp, sTmp); 157f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_3d_device */ 158*9e0c01faSSebastian Grimberg write_V_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV + (2 * dstrdV), cstrdV, rV, tx); 159f80f4a74SSebastian Grimberg } 160f80f4a74SSebastian Grimberg 161*9e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 1623c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q *BASIS_MAX_P_Q, MAGMA_MAXTHREADS_3D)) __global__ 163f80f4a74SSebastian Grimberg void magma_gradt_3d_kernel(const CeedScalar *dinterp1d, const CeedScalar *dgrad1d, const CeedScalar *dU, const int estrdU, const int cstrdU, 164f80f4a74SSebastian Grimberg const int dstrdU, CeedScalar *dV, const int estrdV, const int cstrdV, const int dstrdV, const int nelem) { 165f80f4a74SSebastian Grimberg MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 166f80f4a74SSebastian Grimberg 167f80f4a74SSebastian Grimberg const int tx = threadIdx.x; 168f80f4a74SSebastian Grimberg const int ty = threadIdx.y; 169f80f4a74SSebastian Grimberg const int elem_id = (blockIdx.x * blockDim.y) + ty; 170f80f4a74SSebastian Grimberg magma_trans_t transT = MagmaTrans; 171f80f4a74SSebastian Grimberg 172f80f4a74SSebastian Grimberg if (elem_id >= nelem) return; 173f80f4a74SSebastian Grimberg 1743c1e2affSSebastian Grimberg CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // here DIM_U = 1, but might be different for a fused operator 1753c1e2affSSebastian Grimberg CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // here DIM_V = 1, but might be different for a fused operator 176f80f4a74SSebastian Grimberg CeedScalar rTmp = 0.0; 177f80f4a74SSebastian Grimberg 178f80f4a74SSebastian Grimberg // shift global memory pointers by elem stride 179f80f4a74SSebastian Grimberg dU += elem_id * estrdU; 180f80f4a74SSebastian Grimberg dV += elem_id * estrdV; 181f80f4a74SSebastian Grimberg 182f80f4a74SSebastian Grimberg // assign shared memory pointers 1833c1e2affSSebastian Grimberg CeedScalar *sTinterp = (CeedScalar *)shared_data; 1843c1e2affSSebastian Grimberg CeedScalar *sTgrad = sTinterp + BASIS_Q * BASIS_P; 1853c1e2affSSebastian Grimberg CeedScalar *sTmp = sTgrad + BASIS_Q * BASIS_P; 1863c1e2affSSebastian Grimberg sTmp += ty * (max(BASIS_Q * BASIS_Q * BASIS_Q, (BASIS_Q * BASIS_Q * BASIS_P) + (BASIS_Q * BASIS_P * BASIS_P))); 187f80f4a74SSebastian Grimberg 188f80f4a74SSebastian Grimberg // read T 189f80f4a74SSebastian Grimberg if (ty == 0) { 190*9e0c01faSSebastian Grimberg read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dinterp1d, sTinterp); 191*9e0c01faSSebastian Grimberg read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dgrad1d, sTgrad); 192f80f4a74SSebastian Grimberg } 193f80f4a74SSebastian Grimberg __syncthreads(); 194f80f4a74SSebastian Grimberg 195f80f4a74SSebastian Grimberg // read V (since this is transposed mode) 196f80f4a74SSebastian Grimberg const CeedScalar beta = 1.0; 197*9e0c01faSSebastian Grimberg read_V_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV + (0 * dstrdV), cstrdV, rV, tx); 198f80f4a74SSebastian Grimberg 1993c1e2affSSebastian Grimberg /* read U (idim = 0 for dU, i_DIM = 0 for rU) -- 200f80f4a74SSebastian Grimberg there is a sync at the end of this function */ 201*9e0c01faSSebastian Grimberg read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (0 * dstrdU), cstrdU, rU, sTmp, tx); 2023c1e2affSSebastian Grimberg /* then first call (i_DIM = 0, i_DIM_U = 0, i_DIM_V = 0) */ 2033c1e2affSSebastian Grimberg magma_grad_3d_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); 204f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_3d_device */ 205f80f4a74SSebastian Grimberg 2063c1e2affSSebastian Grimberg /* read U (idim = 1 for dU, i_DIM = 0 for rU) -- 207f80f4a74SSebastian Grimberg there is a sync at the end of this function */ 208*9e0c01faSSebastian Grimberg read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (1 * dstrdU), cstrdU, rU, sTmp, tx); 2093c1e2affSSebastian Grimberg /* then second call (i_DIM = 1, i_DIM_U = 0, i_DIM_V = 0) */ 2103c1e2affSSebastian Grimberg magma_grad_3d_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); 211f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_3d_device */ 212f80f4a74SSebastian Grimberg 2133c1e2affSSebastian Grimberg /* read U (idim = 2 for dU, i_DIM = 0 for rU) -- 214f80f4a74SSebastian Grimberg there is a sync at the end of this function */ 215*9e0c01faSSebastian Grimberg read_U_3d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU + (2 * dstrdU), cstrdU, rU, sTmp, tx); 2163c1e2affSSebastian Grimberg /* then third call (i_DIM = 2, i_DIM_U = 0, i_DIM_V = 0) */ 2173c1e2affSSebastian Grimberg magma_grad_3d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P, 2, 0, 0>(sTinterp, sTgrad, rU, rV, beta, tx, rTmp, sTmp); 218f80f4a74SSebastian Grimberg /* there is a sync at the end of magma_grad_3d_device */ 219f80f4a74SSebastian Grimberg 220f80f4a74SSebastian Grimberg // write V 221*9e0c01faSSebastian Grimberg write_V_3d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV + (0 * dstrdV), cstrdV, rV, tx); 222f80f4a74SSebastian Grimberg } 2233c1e2affSSebastian Grimberg 2243c1e2affSSebastian Grimberg #endif // CEED_MAGMA_BASIS_GRAD_3D_H 225