xref: /libCEED/include/ceed/jit-source/magma/magma-basis-interp-1d.h (revision 5aed82e4fa97acf4ba24a7f10a35f5303a6798e0)
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_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 
189e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
19f80f4a74SSebastian Grimberg // interp basis action (1D)
203c1e2affSSebastian Grimberg template <typename T, int DIM, int NUM_COMP, int P, int Q>
217132caa0SSebastian Grimberg static __device__ __inline__ void magma_interp_1d_device(const T *sT, 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 
313c1e2affSSebastian Grimberg   if (tx < Q) {
323c1e2affSSebastian Grimberg     for (int comp = 0; comp < NUM_COMP; comp++) {
337132caa0SSebastian Grimberg       T rv = 0.0;
343c1e2affSSebastian Grimberg       for (int i = 0; i < P; i++) {
353c1e2affSSebastian Grimberg         rv += sU[comp][i] * sT(i, tx);  // sT[tx * P + i];
36f80f4a74SSebastian Grimberg       }
373c1e2affSSebastian Grimberg       sV[comp][tx] = rv;
38f80f4a74SSebastian Grimberg     }
39f80f4a74SSebastian Grimberg   }
40f80f4a74SSebastian Grimberg }
41f80f4a74SSebastian Grimberg 
429e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
433c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_1D)) __global__
44f80f4a74SSebastian Grimberg     void magma_interpn_1d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV,
45f80f4a74SSebastian Grimberg                                  const int cstrdV, const int nelem) {
46f80f4a74SSebastian Grimberg   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
47f80f4a74SSebastian Grimberg 
48f80f4a74SSebastian Grimberg   const int tx      = threadIdx.x;
49f80f4a74SSebastian Grimberg   const int ty      = threadIdx.y;
50f80f4a74SSebastian Grimberg   const int elem_id = (blockIdx.x * blockDim.y) + ty;
51f80f4a74SSebastian Grimberg 
52f80f4a74SSebastian Grimberg   if (elem_id >= nelem) return;
53f80f4a74SSebastian Grimberg 
543c1e2affSSebastian Grimberg   CeedScalar *sU[BASIS_NUM_COMP];
553c1e2affSSebastian Grimberg   CeedScalar *sV[BASIS_NUM_COMP];
56f80f4a74SSebastian Grimberg 
57f80f4a74SSebastian Grimberg   // shift global memory pointers by elem stride
58f80f4a74SSebastian Grimberg   dU += elem_id * estrdU;
59f80f4a74SSebastian Grimberg   dV += elem_id * estrdV;
60f80f4a74SSebastian Grimberg 
61f80f4a74SSebastian Grimberg   // assign shared memory pointers
623c1e2affSSebastian Grimberg   CeedScalar *sT = (CeedScalar *)shared_data;
633c1e2affSSebastian Grimberg   CeedScalar *sW = sT + BASIS_P * BASIS_Q;
643c1e2affSSebastian Grimberg   sU[0]          = sW + ty * BASIS_NUM_COMP * (BASIS_P + BASIS_Q);
653c1e2affSSebastian Grimberg   sV[0]          = sU[0] + (BASIS_NUM_COMP * 1 * BASIS_P);
663c1e2affSSebastian Grimberg   for (int comp = 1; comp < BASIS_NUM_COMP; comp++) {
673c1e2affSSebastian Grimberg     sU[comp] = sU[comp - 1] + (1 * BASIS_P);
683c1e2affSSebastian Grimberg     sV[comp] = sV[comp - 1] + (1 * BASIS_Q);
69f80f4a74SSebastian Grimberg   }
70f80f4a74SSebastian Grimberg 
71f80f4a74SSebastian Grimberg   // read T
72f80f4a74SSebastian Grimberg   if (ty == 0) {
739e0c01faSSebastian Grimberg     read_T_notrans_gm2sm<BASIS_P, BASIS_Q>(tx, dT, sT);
74f80f4a74SSebastian Grimberg   }
75f80f4a74SSebastian Grimberg 
76f80f4a74SSebastian Grimberg   // read U
773c1e2affSSebastian Grimberg   read_1d<CeedScalar, BASIS_P, BASIS_NUM_COMP>(dU, cstrdU, sU, tx);
78f80f4a74SSebastian Grimberg 
79f80f4a74SSebastian Grimberg   __syncthreads();
807132caa0SSebastian Grimberg   magma_interp_1d_device<CeedScalar, BASIS_DIM, BASIS_NUM_COMP, BASIS_P, BASIS_Q>(sT, sU, sV, tx);
81f80f4a74SSebastian Grimberg   __syncthreads();
82f80f4a74SSebastian Grimberg 
83f80f4a74SSebastian Grimberg   // write V
843c1e2affSSebastian Grimberg   write_1d<CeedScalar, BASIS_Q, BASIS_NUM_COMP>(sV, dV, cstrdV, tx);
85f80f4a74SSebastian Grimberg }
86f80f4a74SSebastian Grimberg 
879e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
883c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_MAX_P_Q, MAGMA_MAXTHREADS_1D)) __global__
89f80f4a74SSebastian Grimberg     void magma_interpt_1d_kernel(const CeedScalar *dT, const CeedScalar *dU, const int estrdU, const int cstrdU, CeedScalar *dV, const int estrdV,
90f80f4a74SSebastian Grimberg                                  const int cstrdV, const int nelem) {
91f80f4a74SSebastian Grimberg   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
92f80f4a74SSebastian Grimberg 
93f80f4a74SSebastian Grimberg   const int tx      = threadIdx.x;
94f80f4a74SSebastian Grimberg   const int ty      = threadIdx.y;
95f80f4a74SSebastian Grimberg   const int elem_id = (blockIdx.x * blockDim.y) + ty;
96f80f4a74SSebastian Grimberg 
97f80f4a74SSebastian Grimberg   if (elem_id >= nelem) return;
98f80f4a74SSebastian Grimberg 
993c1e2affSSebastian Grimberg   CeedScalar *sU[BASIS_NUM_COMP];
1003c1e2affSSebastian Grimberg   CeedScalar *sV[BASIS_NUM_COMP];
101f80f4a74SSebastian Grimberg 
102f80f4a74SSebastian Grimberg   // shift global memory pointers by elem stride
103f80f4a74SSebastian Grimberg   dU += elem_id * estrdU;
104f80f4a74SSebastian Grimberg   dV += elem_id * estrdV;
105f80f4a74SSebastian Grimberg 
106f80f4a74SSebastian Grimberg   // assign shared memory pointers
1073c1e2affSSebastian Grimberg   CeedScalar *sT = (CeedScalar *)shared_data;
1083c1e2affSSebastian Grimberg   CeedScalar *sW = sT + BASIS_Q * BASIS_P;
1093c1e2affSSebastian Grimberg   sU[0]          = sW + ty * BASIS_NUM_COMP * (BASIS_Q + BASIS_P);
1103c1e2affSSebastian Grimberg   sV[0]          = sU[0] + (BASIS_NUM_COMP * 1 * BASIS_Q);
1113c1e2affSSebastian Grimberg   for (int comp = 1; comp < BASIS_NUM_COMP; comp++) {
1123c1e2affSSebastian Grimberg     sU[comp] = sU[comp - 1] + (1 * BASIS_Q);
1133c1e2affSSebastian Grimberg     sV[comp] = sV[comp - 1] + (1 * BASIS_P);
114f80f4a74SSebastian Grimberg   }
115f80f4a74SSebastian Grimberg 
116f80f4a74SSebastian Grimberg   // read T
117f80f4a74SSebastian Grimberg   if (ty == 0) {
1189e0c01faSSebastian Grimberg     read_T_trans_gm2sm<BASIS_Q, BASIS_P>(tx, dT, sT);
119f80f4a74SSebastian Grimberg   }
120f80f4a74SSebastian Grimberg 
121f80f4a74SSebastian Grimberg   // read U
1223c1e2affSSebastian Grimberg   read_1d<CeedScalar, BASIS_Q, BASIS_NUM_COMP>(dU, cstrdU, sU, tx);
123f80f4a74SSebastian Grimberg 
124f80f4a74SSebastian Grimberg   __syncthreads();
1257132caa0SSebastian Grimberg   magma_interp_1d_device<CeedScalar, BASIS_DIM, BASIS_NUM_COMP, BASIS_Q, BASIS_P>(sT, sU, sV, tx);
126f80f4a74SSebastian Grimberg   __syncthreads();
127f80f4a74SSebastian Grimberg 
128f80f4a74SSebastian Grimberg   // write V
1293c1e2affSSebastian Grimberg   write_1d<CeedScalar, BASIS_P, BASIS_NUM_COMP>(sV, dV, cstrdV, tx);
130f80f4a74SSebastian Grimberg }
1313c1e2affSSebastian Grimberg 
1323c1e2affSSebastian Grimberg #endif  // CEED_MAGMA_BASIS_INTERP_1D_H
133