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 interpolation in 1D 103c1e2affSSebastian Grimberg #ifndef CEED_MAGMA_BASIS_INTERP_2D_H 113c1e2affSSebastian Grimberg #define CEED_MAGMA_BASIS_INTERP_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 //////////////////////////////////////////////////////////////////////////////// 20f80f4a74SSebastian Grimberg // interp basis action (2D) 213c1e2affSSebastian Grimberg template <typename T, int DIM_U, int DIM_V, int NUM_COMP, int P, int Q, int rU_SIZE, int rV_SIZE> 227132caa0SSebastian Grimberg 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, 237132caa0SSebastian Grimberg T rTmp, T *swork) { 24f80f4a74SSebastian Grimberg // Assumptions 253c1e2affSSebastian Grimberg // 1. 1D threads of size max(P,Q) 263c1e2affSSebastian Grimberg // 2. input: rU[DIM_U x NUM_COMP x rU_SIZE] in registers (per thread) 273c1e2affSSebastian Grimberg // 3. output: rV[DIM_V x NUM_COMP x rV_SIZE] in registers (per thread) 28f80f4a74SSebastian Grimberg // 4. Two products per component 293c1e2affSSebastian Grimberg // 4.1 Batch P of (1xP) matrices times (PxQ) matrix => Batch P of (1xQ) matrices 303c1e2affSSebastian Grimberg // 4.2 Batch 1 of (QxP) matrix times (PxQ) matrix => (QxQ) matrix 31f80f4a74SSebastian Grimberg // 5. Each thread computes one row of the output of each product 32f80f4a74SSebastian Grimberg // 6. Sync is recommended before and after the call 33f80f4a74SSebastian Grimberg 343c1e2affSSebastian Grimberg for (int comp = 0; comp < NUM_COMP; comp++) { 353c1e2affSSebastian Grimberg // 1st product -- Batch P of (1xP) matrices [reg] x (PxQ) [shmem] => Batch P of (1xQ) matrices 363c1e2affSSebastian Grimberg // the batch output P x (1xQ) is written on the fly to shmem 373c1e2affSSebastian Grimberg if (tx < P) { 38f80f4a74SSebastian Grimberg const int batchid = tx; 39f80f4a74SSebastian Grimberg const int sld = 1; 403c1e2affSSebastian Grimberg T *sTmp = swork + batchid * (1 * Q); 413c1e2affSSebastian Grimberg for (int j = 0; j < Q; j++) { 42f80f4a74SSebastian Grimberg rTmp = 0.0; 433c1e2affSSebastian Grimberg for (int i = 0; i < P; i++) { 443c1e2affSSebastian Grimberg rTmp += rU[0][comp][i] * sT(i, j); 45f80f4a74SSebastian Grimberg } 46f80f4a74SSebastian Grimberg sTmp(0, j, sld) = rTmp; 47f80f4a74SSebastian Grimberg } 483c1e2affSSebastian Grimberg } // end of: if (tx < P) 49f80f4a74SSebastian Grimberg __syncthreads(); 50f80f4a74SSebastian Grimberg 513c1e2affSSebastian Grimberg // 2nd product -- Batch 1 of a (QxP) matrix [shmem] x (PxQ) [shmem] => (QxQ) matrix [reg] 523c1e2affSSebastian Grimberg if (tx < Q) { 53f80f4a74SSebastian Grimberg const int batchid = 0; 543c1e2affSSebastian Grimberg const int sld = Q; 553c1e2affSSebastian Grimberg T *sTmp = swork + batchid * (Q * P); 563c1e2affSSebastian Grimberg for (int j = 0; j < Q; j++) { 57f80f4a74SSebastian Grimberg rTmp = 0.0; 583c1e2affSSebastian Grimberg for (int i = 0; i < P; i++) { 59f80f4a74SSebastian Grimberg rTmp += sTmp(tx, i, sld) * sT(i, j); 60f80f4a74SSebastian Grimberg } 613c1e2affSSebastian Grimberg rV[0][comp][j] += rTmp; 62f80f4a74SSebastian Grimberg } 63f80f4a74SSebastian Grimberg } 64f80f4a74SSebastian Grimberg __syncthreads(); 65f80f4a74SSebastian Grimberg } 66f80f4a74SSebastian Grimberg } 67f80f4a74SSebastian Grimberg 689e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 693c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__ 70f80f4a74SSebastian Grimberg void magma_interpn_2d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV, 71f80f4a74SSebastian Grimberg const int cstrdV, const int nelem) { 72f80f4a74SSebastian Grimberg MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 73f80f4a74SSebastian Grimberg 74f80f4a74SSebastian Grimberg const int tx = threadIdx.x; 75f80f4a74SSebastian Grimberg const int ty = threadIdx.y; 76f80f4a74SSebastian Grimberg const int elem_id = (blockIdx.x * blockDim.y) + ty; 77f80f4a74SSebastian Grimberg 78f80f4a74SSebastian Grimberg if (elem_id >= nelem) return; 79f80f4a74SSebastian Grimberg 803c1e2affSSebastian Grimberg CeedScalar rU[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // for a non-fused operator BASIS_DIM is always 1 813c1e2affSSebastian Grimberg CeedScalar rV[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // for a non-fused operator BASIS_DIM is always 1 82f80f4a74SSebastian Grimberg CeedScalar rTmp = 0.0; 83f80f4a74SSebastian Grimberg 84f80f4a74SSebastian Grimberg // shift global memory pointers by elem stride 85f80f4a74SSebastian Grimberg dU += elem_id * estrdU; 86f80f4a74SSebastian Grimberg dV += elem_id * estrdV; 87f80f4a74SSebastian Grimberg 88f80f4a74SSebastian Grimberg // assign shared memory pointers 893c1e2affSSebastian Grimberg CeedScalar *sT = (CeedScalar *)shared_data; 903c1e2affSSebastian Grimberg CeedScalar *sTmp = sT + BASIS_P * BASIS_Q; 913c1e2affSSebastian Grimberg sTmp += ty * (BASIS_P * BASIS_MAX_P_Q); 92f80f4a74SSebastian Grimberg 93f80f4a74SSebastian Grimberg // read T 94f80f4a74SSebastian Grimberg if (ty == 0) { 959e0c01faSSebastian Grimberg read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dT, sT); 96f80f4a74SSebastian Grimberg } 97f80f4a74SSebastian Grimberg 98f80f4a74SSebastian Grimberg // read U -- there is a sync at the end of this function 999e0c01faSSebastian Grimberg read_U_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dU, cstrdU, rU, sTmp, tx); 100f80f4a74SSebastian Grimberg 1019e0c01faSSebastian Grimberg // no sync needed here -- read_U_2d already syncs at the end 1027132caa0SSebastian Grimberg magma_interp_2d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_P, BASIS_Q, BASIS_P, BASIS_Q>(sT, rU, rV, tx, rTmp, sTmp); 103f80f4a74SSebastian Grimberg __syncthreads(); 104f80f4a74SSebastian Grimberg 105f80f4a74SSebastian Grimberg // write V 1069e0c01faSSebastian Grimberg write_V_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV, cstrdV, rV, tx); 107f80f4a74SSebastian Grimberg } 108f80f4a74SSebastian Grimberg 1099e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 1103c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__ 111f80f4a74SSebastian Grimberg void magma_interpt_2d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV, 112f80f4a74SSebastian Grimberg const int cstrdV, const int nelem) { 113f80f4a74SSebastian Grimberg MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 114f80f4a74SSebastian Grimberg 115f80f4a74SSebastian Grimberg const int tx = threadIdx.x; 116f80f4a74SSebastian Grimberg const int ty = threadIdx.y; 117f80f4a74SSebastian Grimberg const int elem_id = (blockIdx.x * blockDim.y) + ty; 118f80f4a74SSebastian Grimberg 119f80f4a74SSebastian Grimberg if (elem_id >= nelem) return; 120f80f4a74SSebastian Grimberg 1213c1e2affSSebastian Grimberg CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // for a non-fused operator BASIS_DIM is always 1 1223c1e2affSSebastian Grimberg CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // for a non-fused operator BASIS_DIM is always 1 123f80f4a74SSebastian Grimberg CeedScalar rTmp = 0.0; 124f80f4a74SSebastian Grimberg 125f80f4a74SSebastian Grimberg // shift global memory pointers by elem stride 126f80f4a74SSebastian Grimberg dU += elem_id * estrdU; 127f80f4a74SSebastian Grimberg dV += elem_id * estrdV; 128f80f4a74SSebastian Grimberg 129f80f4a74SSebastian Grimberg // assign shared memory pointers 1303c1e2affSSebastian Grimberg CeedScalar *sT = (CeedScalar *)shared_data; 1313c1e2affSSebastian Grimberg CeedScalar *sTmp = sT + BASIS_Q * BASIS_P; 1323c1e2affSSebastian Grimberg sTmp += ty * (BASIS_Q * BASIS_MAX_P_Q); 133f80f4a74SSebastian Grimberg 134f80f4a74SSebastian Grimberg // read T 135f80f4a74SSebastian Grimberg if (ty == 0) { 1369e0c01faSSebastian Grimberg read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dT, sT); 137f80f4a74SSebastian Grimberg } 138f80f4a74SSebastian Grimberg 139f80f4a74SSebastian Grimberg // read U -- there is a sync at the end of this function 1409e0c01faSSebastian Grimberg read_U_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU, cstrdU, rU, sTmp, tx); 141f80f4a74SSebastian Grimberg 1429e0c01faSSebastian Grimberg // no sync needed here -- read_U_2d already syncs at the end 1437132caa0SSebastian Grimberg magma_interp_2d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P>(sT, rU, rV, tx, rTmp, sTmp); 144f80f4a74SSebastian Grimberg __syncthreads(); 145f80f4a74SSebastian Grimberg 146f80f4a74SSebastian Grimberg // write V 1479e0c01faSSebastian Grimberg write_V_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV, cstrdV, rV, tx); 148f80f4a74SSebastian Grimberg } 1493c1e2affSSebastian Grimberg 1503c1e2affSSebastian Grimberg #endif // CEED_MAGMA_BASIS_INTERP_2D_H 151