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 interpolation in 1D 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 //////////////////////////////////////////////////////////////////////////////// 17f80f4a74SSebastian Grimberg // interp basis action (2D) 183c1e2affSSebastian Grimberg template <typename T, int DIM_U, int DIM_V, int NUM_COMP, int P, int Q, int rU_SIZE, int rV_SIZE> 197132caa0SSebastian 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, 207132caa0SSebastian Grimberg T rTmp, T *swork) { 21f80f4a74SSebastian Grimberg // Assumptions 223c1e2affSSebastian Grimberg // 1. 1D threads of size max(P,Q) 233c1e2affSSebastian Grimberg // 2. input: rU[DIM_U x NUM_COMP x rU_SIZE] in registers (per thread) 243c1e2affSSebastian Grimberg // 3. output: rV[DIM_V x NUM_COMP x rV_SIZE] in registers (per thread) 25f80f4a74SSebastian Grimberg // 4. Two products per component 263c1e2affSSebastian Grimberg // 4.1 Batch P of (1xP) matrices times (PxQ) matrix => Batch P of (1xQ) matrices 273c1e2affSSebastian Grimberg // 4.2 Batch 1 of (QxP) matrix times (PxQ) matrix => (QxQ) matrix 28f80f4a74SSebastian Grimberg // 5. Each thread computes one row of the output of each product 29f80f4a74SSebastian Grimberg // 6. Sync is recommended before and after the call 30f80f4a74SSebastian Grimberg 313c1e2affSSebastian Grimberg for (int comp = 0; comp < NUM_COMP; comp++) { 323c1e2affSSebastian Grimberg // 1st product -- Batch P of (1xP) matrices [reg] x (PxQ) [shmem] => Batch P of (1xQ) matrices 333c1e2affSSebastian Grimberg // the batch output P x (1xQ) is written on the fly to shmem 343c1e2affSSebastian Grimberg if (tx < P) { 35f80f4a74SSebastian Grimberg const int batchid = tx; 36f80f4a74SSebastian Grimberg const int sld = 1; 373c1e2affSSebastian Grimberg T *sTmp = swork + batchid * (1 * Q); 383c1e2affSSebastian Grimberg for (int j = 0; j < Q; j++) { 39f80f4a74SSebastian Grimberg rTmp = 0.0; 403c1e2affSSebastian Grimberg for (int i = 0; i < P; i++) { 413c1e2affSSebastian Grimberg rTmp += rU[0][comp][i] * sT(i, j); 42f80f4a74SSebastian Grimberg } 43f80f4a74SSebastian Grimberg sTmp(0, j, sld) = rTmp; 44f80f4a74SSebastian Grimberg } 453c1e2affSSebastian Grimberg } // end of: if (tx < P) 46f80f4a74SSebastian Grimberg __syncthreads(); 47f80f4a74SSebastian Grimberg 483c1e2affSSebastian Grimberg // 2nd product -- Batch 1 of a (QxP) matrix [shmem] x (PxQ) [shmem] => (QxQ) matrix [reg] 493c1e2affSSebastian Grimberg if (tx < Q) { 50f80f4a74SSebastian Grimberg const int batchid = 0; 513c1e2affSSebastian Grimberg const int sld = Q; 523c1e2affSSebastian Grimberg T *sTmp = swork + batchid * (Q * P); 533c1e2affSSebastian Grimberg for (int j = 0; j < Q; j++) { 54f80f4a74SSebastian Grimberg rTmp = 0.0; 553c1e2affSSebastian Grimberg for (int i = 0; i < P; i++) { 56f80f4a74SSebastian Grimberg rTmp += sTmp(tx, i, sld) * sT(i, j); 57f80f4a74SSebastian Grimberg } 583c1e2affSSebastian Grimberg rV[0][comp][j] += rTmp; 59f80f4a74SSebastian Grimberg } 60f80f4a74SSebastian Grimberg } 61f80f4a74SSebastian Grimberg __syncthreads(); 62f80f4a74SSebastian Grimberg } 63f80f4a74SSebastian Grimberg } 64f80f4a74SSebastian Grimberg 659e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 663c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__ 67f80f4a74SSebastian Grimberg void magma_interpn_2d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV, 68f80f4a74SSebastian Grimberg const int cstrdV, const int nelem) { 69f80f4a74SSebastian Grimberg MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 70f80f4a74SSebastian Grimberg 71f80f4a74SSebastian Grimberg const int tx = threadIdx.x; 72f80f4a74SSebastian Grimberg const int ty = threadIdx.y; 73f80f4a74SSebastian Grimberg const int elem_id = (blockIdx.x * blockDim.y) + ty; 74f80f4a74SSebastian Grimberg 75f80f4a74SSebastian Grimberg if (elem_id >= nelem) return; 76f80f4a74SSebastian Grimberg 773c1e2affSSebastian Grimberg CeedScalar rU[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // for a non-fused operator BASIS_DIM is always 1 783c1e2affSSebastian Grimberg CeedScalar rV[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // for a non-fused operator BASIS_DIM is always 1 79f80f4a74SSebastian Grimberg CeedScalar rTmp = 0.0; 80f80f4a74SSebastian Grimberg 81f80f4a74SSebastian Grimberg // shift global memory pointers by elem stride 82f80f4a74SSebastian Grimberg dU += elem_id * estrdU; 83f80f4a74SSebastian Grimberg dV += elem_id * estrdV; 84f80f4a74SSebastian Grimberg 85f80f4a74SSebastian Grimberg // assign shared memory pointers 863c1e2affSSebastian Grimberg CeedScalar *sT = (CeedScalar *)shared_data; 873c1e2affSSebastian Grimberg CeedScalar *sTmp = sT + BASIS_P * BASIS_Q; 883c1e2affSSebastian Grimberg sTmp += ty * (BASIS_P * BASIS_MAX_P_Q); 89f80f4a74SSebastian Grimberg 90f80f4a74SSebastian Grimberg // read T 91f80f4a74SSebastian Grimberg if (ty == 0) { 929e0c01faSSebastian Grimberg read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dT, sT); 93f80f4a74SSebastian Grimberg } 94f80f4a74SSebastian Grimberg 95f80f4a74SSebastian Grimberg // read U -- there is a sync at the end of this function 969e0c01faSSebastian Grimberg read_U_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dU, cstrdU, rU, sTmp, tx); 97f80f4a74SSebastian Grimberg 989e0c01faSSebastian Grimberg // no sync needed here -- read_U_2d already syncs at the end 997132caa0SSebastian 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); 100f80f4a74SSebastian Grimberg __syncthreads(); 101f80f4a74SSebastian Grimberg 102f80f4a74SSebastian Grimberg // write V 1039e0c01faSSebastian Grimberg write_V_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dV, cstrdV, rV, tx); 104f80f4a74SSebastian Grimberg } 105f80f4a74SSebastian Grimberg 1069e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 1073c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__ 108f80f4a74SSebastian Grimberg void magma_interpt_2d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV, 109f80f4a74SSebastian Grimberg const int cstrdV, const int nelem) { 110f80f4a74SSebastian Grimberg MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 111f80f4a74SSebastian Grimberg 112f80f4a74SSebastian Grimberg const int tx = threadIdx.x; 113f80f4a74SSebastian Grimberg const int ty = threadIdx.y; 114f80f4a74SSebastian Grimberg const int elem_id = (blockIdx.x * blockDim.y) + ty; 115f80f4a74SSebastian Grimberg 116f80f4a74SSebastian Grimberg if (elem_id >= nelem) return; 117f80f4a74SSebastian Grimberg 1183c1e2affSSebastian Grimberg CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // for a non-fused operator BASIS_DIM is always 1 1193c1e2affSSebastian Grimberg CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // for a non-fused operator BASIS_DIM is always 1 120f80f4a74SSebastian Grimberg CeedScalar rTmp = 0.0; 121f80f4a74SSebastian Grimberg 122f80f4a74SSebastian Grimberg // shift global memory pointers by elem stride 123f80f4a74SSebastian Grimberg dU += elem_id * estrdU; 124f80f4a74SSebastian Grimberg dV += elem_id * estrdV; 125f80f4a74SSebastian Grimberg 126f80f4a74SSebastian Grimberg // assign shared memory pointers 1273c1e2affSSebastian Grimberg CeedScalar *sT = (CeedScalar *)shared_data; 1283c1e2affSSebastian Grimberg CeedScalar *sTmp = sT + BASIS_Q * BASIS_P; 1293c1e2affSSebastian Grimberg sTmp += ty * (BASIS_Q * BASIS_MAX_P_Q); 130f80f4a74SSebastian Grimberg 131f80f4a74SSebastian Grimberg // read T 132f80f4a74SSebastian Grimberg if (ty == 0) { 1339e0c01faSSebastian Grimberg read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dT, sT); 134f80f4a74SSebastian Grimberg } 135f80f4a74SSebastian Grimberg 136f80f4a74SSebastian Grimberg // read U -- there is a sync at the end of this function 1379e0c01faSSebastian Grimberg read_U_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU, cstrdU, rU, sTmp, tx); 138f80f4a74SSebastian Grimberg 1399e0c01faSSebastian Grimberg // no sync needed here -- read_U_2d already syncs at the end 1407132caa0SSebastian 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); 141f80f4a74SSebastian Grimberg __syncthreads(); 142f80f4a74SSebastian Grimberg 143f80f4a74SSebastian Grimberg // write V 1449e0c01faSSebastian Grimberg write_V_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV, cstrdV, rV, tx); 145f80f4a74SSebastian Grimberg } 146db2becc9SJeremy L Thompson 147db2becc9SJeremy L Thompson //////////////////////////////////////////////////////////////////////////////// 148db2becc9SJeremy L Thompson extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_2D)) __global__ 149db2becc9SJeremy L Thompson void magma_interpta_2d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV, 150db2becc9SJeremy L Thompson const int cstrdV, const int nelem) { 151db2becc9SJeremy L Thompson MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 152db2becc9SJeremy L Thompson 153db2becc9SJeremy L Thompson const int tx = threadIdx.x; 154db2becc9SJeremy L Thompson const int ty = threadIdx.y; 155db2becc9SJeremy L Thompson const int elem_id = (blockIdx.x * blockDim.y) + ty; 156db2becc9SJeremy L Thompson 157db2becc9SJeremy L Thompson if (elem_id >= nelem) return; 158db2becc9SJeremy L Thompson 159db2becc9SJeremy L Thompson CeedScalar rU[1][BASIS_NUM_COMP][BASIS_Q] = {0.0}; // for a non-fused operator BASIS_DIM is always 1 160db2becc9SJeremy L Thompson CeedScalar rV[1][BASIS_NUM_COMP][BASIS_P] = {0.0}; // for a non-fused operator BASIS_DIM is always 1 161db2becc9SJeremy L Thompson CeedScalar rTmp = 0.0; 162db2becc9SJeremy L Thompson 163db2becc9SJeremy L Thompson // shift global memory pointers by elem stride 164db2becc9SJeremy L Thompson dU += elem_id * estrdU; 165db2becc9SJeremy L Thompson dV += elem_id * estrdV; 166db2becc9SJeremy L Thompson 167db2becc9SJeremy L Thompson // assign shared memory pointers 168db2becc9SJeremy L Thompson CeedScalar *sT = (CeedScalar *)shared_data; 169db2becc9SJeremy L Thompson CeedScalar *sTmp = sT + BASIS_Q * BASIS_P; 170db2becc9SJeremy L Thompson sTmp += ty * (BASIS_Q * BASIS_MAX_P_Q); 171db2becc9SJeremy L Thompson 172db2becc9SJeremy L Thompson // read T 173db2becc9SJeremy L Thompson if (ty == 0) { 174db2becc9SJeremy L Thompson read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dT, sT); 175db2becc9SJeremy L Thompson } 176db2becc9SJeremy L Thompson 177db2becc9SJeremy L Thompson // read U -- there is a sync at the end of this function 178db2becc9SJeremy L Thompson read_U_2d<CeedScalar, BASIS_Q, 1, BASIS_NUM_COMP, BASIS_Q, 0>(dU, cstrdU, rU, sTmp, tx); 179db2becc9SJeremy L Thompson 180db2becc9SJeremy L Thompson // no sync needed here -- read_U_2d already syncs at the end 181db2becc9SJeremy L Thompson magma_interp_2d_device<CeedScalar, 1, 1, BASIS_NUM_COMP, BASIS_Q, BASIS_P, BASIS_Q, BASIS_P>(sT, rU, rV, tx, rTmp, sTmp); 182db2becc9SJeremy L Thompson __syncthreads(); 183db2becc9SJeremy L Thompson 184db2becc9SJeremy L Thompson // sum into V 185db2becc9SJeremy L Thompson sum_V_2d<CeedScalar, BASIS_P, 1, BASIS_NUM_COMP, BASIS_P, 0>(dV, cstrdV, rV, tx); 186db2becc9SJeremy L Thompson } 187