xref: /libCEED/include/ceed/jit-source/magma/magma-basis-weight-1d.h (revision d4cc18453651bd0f94c1a2e078b2646a92dafdcc)
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 weight in 1D
103c1e2affSSebastian Grimberg #include "magma-common-tensor.h"
113c1e2affSSebastian Grimberg 
129e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
13f80f4a74SSebastian Grimberg // weight basis action -- 1D
143c1e2affSSebastian Grimberg template <typename T, int Q>
magma_weight_1d_device(const T * sTweight,T * sV,const int tx)153c1e2affSSebastian Grimberg static __device__ __inline__ void magma_weight_1d_device(const T *sTweight, T *sV, const int tx) {
16f80f4a74SSebastian Grimberg   // Assumptions
173c1e2affSSebastian Grimberg   // 1. 1D thread configuration of size Q
189d15e85bSSebastian Grimberg   // 2. The output sV is in shared memory -- size Q
193c1e2affSSebastian Grimberg   if (tx < Q) {
20f80f4a74SSebastian Grimberg     sV[tx] = sTweight[tx];
21f80f4a74SSebastian Grimberg   }
22f80f4a74SSebastian Grimberg }
23f80f4a74SSebastian Grimberg 
249e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
__launch_bounds__(MAGMA_BASIS_BOUNDS (BASIS_Q,MAGMA_MAXTHREADS_1D))253c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_Q, MAGMA_MAXTHREADS_1D)) __global__
26f80f4a74SSebastian Grimberg     void magma_weight_1d_kernel(const CeedScalar *dqweight1d, CeedScalar *dV, const int v_stride, const int nelem) {
27f80f4a74SSebastian Grimberg   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
28f80f4a74SSebastian Grimberg 
29f80f4a74SSebastian Grimberg   const int tx      = threadIdx.x;
30f80f4a74SSebastian Grimberg   const int ty      = threadIdx.y;
31f80f4a74SSebastian Grimberg   const int elem_id = (blockIdx.x * blockDim.y) + ty;
32f80f4a74SSebastian Grimberg 
33f80f4a74SSebastian Grimberg   if (elem_id >= nelem) return;
34f80f4a74SSebastian Grimberg 
35f80f4a74SSebastian Grimberg   // global memory pointers
36f80f4a74SSebastian Grimberg   dV += elem_id * v_stride;
37f80f4a74SSebastian Grimberg 
38f80f4a74SSebastian Grimberg   // shared memory pointers
39f80f4a74SSebastian Grimberg   CeedScalar *sTweight = (CeedScalar *)shared_data;
403c1e2affSSebastian Grimberg   CeedScalar *sV       = sTweight + BASIS_Q;
413c1e2affSSebastian Grimberg   sV += ty * BASIS_Q;
42f80f4a74SSebastian Grimberg 
43f80f4a74SSebastian Grimberg   // read dqweight_1d
443c1e2affSSebastian Grimberg   if (ty == 0 && tx < BASIS_Q) {
45f80f4a74SSebastian Grimberg     sTweight[tx] = dqweight1d[tx];
46f80f4a74SSebastian Grimberg   }
47f80f4a74SSebastian Grimberg 
48f80f4a74SSebastian Grimberg   __syncthreads();
493c1e2affSSebastian Grimberg   magma_weight_1d_device<CeedScalar, BASIS_Q>(sTweight, sV, tx);
50f80f4a74SSebastian Grimberg   __syncthreads();
51f80f4a74SSebastian Grimberg 
52f80f4a74SSebastian Grimberg   // write V
53f80f4a74SSebastian Grimberg   dV[tx] = sV[tx];
54f80f4a74SSebastian Grimberg }
55