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 interpolation in 1D 103c1e2affSSebastian Grimberg #ifndef CEED_MAGMA_BASIS_INTERP_1D_H 113c1e2affSSebastian Grimberg #define CEED_MAGMA_BASIS_INTERP_1D_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 18*9e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 19f80f4a74SSebastian Grimberg // interp basis action (1D) 203c1e2affSSebastian Grimberg template <typename T, int DIM, int NUM_COMP, int P, int Q> 213c1e2affSSebastian Grimberg static __device__ __inline__ void magma_interp_1d_device(const T *sT, magma_trans_t transT, T *sU[NUM_COMP], T *sV[NUM_COMP], const int tx) { 22f80f4a74SSebastian Grimberg // Assumptions 233c1e2affSSebastian Grimberg // 1. 1D threads of size max(P,Q) 243c1e2affSSebastian Grimberg // 2. sU[i] is 1xP: in shared memory 253c1e2affSSebastian Grimberg // 3. sV[i] is 1xQ: in shared memory 263c1e2affSSebastian Grimberg // 4. P_roduct per component is one row (1xP) times T matrix (PxQ) => one row (1xQ) 27f80f4a74SSebastian Grimberg // 5. Each thread computes one entry in sV[i] 28f80f4a74SSebastian Grimberg // 6. Must sync before and after call 29f80f4a74SSebastian Grimberg // 7. Note that the layout for U and V is different from 2D/3D problem 30f80f4a74SSebastian Grimberg 31f80f4a74SSebastian Grimberg T rv; 323c1e2affSSebastian Grimberg if (tx < Q) { 333c1e2affSSebastian Grimberg for (int comp = 0; comp < NUM_COMP; comp++) { 343c1e2affSSebastian Grimberg rv = (transT == MagmaTrans) ? sV[comp][tx] : 0.0; 353c1e2affSSebastian Grimberg for (int i = 0; i < P; i++) { 363c1e2affSSebastian Grimberg rv += sU[comp][i] * sT(i, tx); // sT[tx * P + i]; 37f80f4a74SSebastian Grimberg } 383c1e2affSSebastian Grimberg sV[comp][tx] = rv; 39f80f4a74SSebastian Grimberg } 40f80f4a74SSebastian Grimberg } 41f80f4a74SSebastian Grimberg } 42f80f4a74SSebastian Grimberg 43*9e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 443c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_1D)) __global__ 45f80f4a74SSebastian Grimberg void magma_interpn_1d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV, 46f80f4a74SSebastian Grimberg const int cstrdV, const int nelem) { 47f80f4a74SSebastian Grimberg MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 48f80f4a74SSebastian Grimberg 49f80f4a74SSebastian Grimberg const int tx = threadIdx.x; 50f80f4a74SSebastian Grimberg const int ty = threadIdx.y; 51f80f4a74SSebastian Grimberg const int elem_id = (blockIdx.x * blockDim.y) + ty; 52f80f4a74SSebastian Grimberg magma_trans_t transT = MagmaNoTrans; 53f80f4a74SSebastian Grimberg 54f80f4a74SSebastian Grimberg if (elem_id >= nelem) return; 55f80f4a74SSebastian Grimberg 563c1e2affSSebastian Grimberg CeedScalar *sU[BASIS_NUM_COMP]; 573c1e2affSSebastian Grimberg CeedScalar *sV[BASIS_NUM_COMP]; 58f80f4a74SSebastian Grimberg 59f80f4a74SSebastian Grimberg // shift global memory pointers by elem stride 60f80f4a74SSebastian Grimberg dU += elem_id * estrdU; 61f80f4a74SSebastian Grimberg dV += elem_id * estrdV; 62f80f4a74SSebastian Grimberg 63f80f4a74SSebastian Grimberg // assign shared memory pointers 643c1e2affSSebastian Grimberg CeedScalar *sT = (CeedScalar *)shared_data; 653c1e2affSSebastian Grimberg CeedScalar *sW = sT + BASIS_P * BASIS_Q; 663c1e2affSSebastian Grimberg sU[0] = sW + ty * BASIS_NUM_COMP * (BASIS_P + BASIS_Q); 673c1e2affSSebastian Grimberg sV[0] = sU[0] + (BASIS_NUM_COMP * 1 * BASIS_P); 683c1e2affSSebastian Grimberg for (int comp = 1; comp < BASIS_NUM_COMP; comp++) { 693c1e2affSSebastian Grimberg sU[comp] = sU[comp - 1] + (1 * BASIS_P); 703c1e2affSSebastian Grimberg sV[comp] = sV[comp - 1] + (1 * BASIS_Q); 71f80f4a74SSebastian Grimberg } 72f80f4a74SSebastian Grimberg 73f80f4a74SSebastian Grimberg // read T 74f80f4a74SSebastian Grimberg if (ty == 0) { 75*9e0c01faSSebastian Grimberg read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dT, sT); 76f80f4a74SSebastian Grimberg } 77f80f4a74SSebastian Grimberg 78f80f4a74SSebastian Grimberg // read U 793c1e2affSSebastian Grimberg read_1d<CeedScalar, BASIS_P, BASIS_NUM_COMP>(dU, cstrdU, sU, tx); 80f80f4a74SSebastian Grimberg 81f80f4a74SSebastian Grimberg __syncthreads(); 823c1e2affSSebastian Grimberg magma_interp_1d_device<CeedScalar, BASIS_DIM, BASIS_NUM_COMP, BASIS_P, BASIS_Q>(sT, transT, sU, sV, tx); 83f80f4a74SSebastian Grimberg __syncthreads(); 84f80f4a74SSebastian Grimberg 85f80f4a74SSebastian Grimberg // write V 863c1e2affSSebastian Grimberg write_1d<CeedScalar, BASIS_Q, BASIS_NUM_COMP>(sV, dV, cstrdV, tx); 87f80f4a74SSebastian Grimberg } 88f80f4a74SSebastian Grimberg 89*9e0c01faSSebastian Grimberg //////////////////////////////////////////////////////////////////////////////// 903c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_1D)) __global__ 91f80f4a74SSebastian Grimberg void magma_interpt_1d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV, 92f80f4a74SSebastian Grimberg const int cstrdV, 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 magma_trans_t transT = MagmaTrans; 99f80f4a74SSebastian Grimberg 100f80f4a74SSebastian Grimberg if (elem_id >= nelem) return; 101f80f4a74SSebastian Grimberg 1023c1e2affSSebastian Grimberg CeedScalar *sU[BASIS_NUM_COMP]; 1033c1e2affSSebastian Grimberg CeedScalar *sV[BASIS_NUM_COMP]; 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 *sT = (CeedScalar *)shared_data; 1113c1e2affSSebastian Grimberg CeedScalar *sW = sT + BASIS_Q * BASIS_P; 1123c1e2affSSebastian Grimberg sU[0] = sW + ty * BASIS_NUM_COMP * (BASIS_Q + BASIS_P); 1133c1e2affSSebastian Grimberg sV[0] = sU[0] + (BASIS_NUM_COMP * 1 * BASIS_Q); 1143c1e2affSSebastian Grimberg for (int comp = 1; comp < BASIS_NUM_COMP; comp++) { 1153c1e2affSSebastian Grimberg sU[comp] = sU[comp - 1] + (1 * BASIS_Q); 1163c1e2affSSebastian Grimberg sV[comp] = sV[comp - 1] + (1 * BASIS_P); 117f80f4a74SSebastian Grimberg } 118f80f4a74SSebastian Grimberg 119f80f4a74SSebastian Grimberg // read T 120f80f4a74SSebastian Grimberg if (ty == 0) { 121*9e0c01faSSebastian Grimberg read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dT, sT); 122f80f4a74SSebastian Grimberg } 123f80f4a74SSebastian Grimberg 124f80f4a74SSebastian Grimberg // read U 1253c1e2affSSebastian Grimberg read_1d<CeedScalar, BASIS_Q, BASIS_NUM_COMP>(dU, cstrdU, sU, tx); 126f80f4a74SSebastian Grimberg 127f80f4a74SSebastian Grimberg // read V 1283c1e2affSSebastian Grimberg read_1d<CeedScalar, BASIS_P, BASIS_NUM_COMP>(dV, cstrdV, sV, tx); 129f80f4a74SSebastian Grimberg 130f80f4a74SSebastian Grimberg __syncthreads(); 1313c1e2affSSebastian Grimberg magma_interp_1d_device<CeedScalar, BASIS_DIM, BASIS_NUM_COMP, BASIS_Q, BASIS_P>(sT, transT, sU, sV, tx); 132f80f4a74SSebastian Grimberg __syncthreads(); 133f80f4a74SSebastian Grimberg 134f80f4a74SSebastian Grimberg // write V 1353c1e2affSSebastian Grimberg write_1d<CeedScalar, BASIS_P, BASIS_NUM_COMP>(sV, dV, cstrdV, tx); 136f80f4a74SSebastian Grimberg } 1373c1e2affSSebastian Grimberg 1383c1e2affSSebastian Grimberg #endif // CEED_MAGMA_BASIS_INTERP_1D_H 139