xref: /libCEED/rust/libceed-sys/c-src/include/ceed/jit-source/magma/magma-basis-weight-1d.h (revision 9d15e85b4f78ffb2d2860753c87a3b1789cc3bb6)
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 weight in 1D
103c1e2affSSebastian Grimberg #ifndef CEED_MAGMA_BASIS_WEIGHT_1D_H
113c1e2affSSebastian Grimberg #define CEED_MAGMA_BASIS_WEIGHT_1D_H
123c1e2affSSebastian Grimberg 
133c1e2affSSebastian Grimberg #include "magma-common-tensor.h"
143c1e2affSSebastian Grimberg 
15f80f4a74SSebastian Grimberg //////////////////////////////////////////////////////////////////////////////////////////
16f80f4a74SSebastian Grimberg // weight basis action -- 1D
173c1e2affSSebastian Grimberg template <typename T, int Q>
183c1e2affSSebastian Grimberg static __device__ __inline__ void magma_weight_1d_device(const T *sTweight, T *sV, const int tx) {
19f80f4a74SSebastian Grimberg   // Assumptions
203c1e2affSSebastian Grimberg   // 1. 1D thread configuration of size Q
21*9d15e85bSSebastian Grimberg   // 2. The output sV is in shared memory -- size Q
223c1e2affSSebastian Grimberg   if (tx < Q) {
23f80f4a74SSebastian Grimberg     sV[tx] = sTweight[tx];
24f80f4a74SSebastian Grimberg   }
25f80f4a74SSebastian Grimberg }
26f80f4a74SSebastian Grimberg 
27f80f4a74SSebastian Grimberg //////////////////////////////////////////////////////////////////////////////////////////
283c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_Q, MAGMA_MAXTHREADS_1D)) __global__
29f80f4a74SSebastian Grimberg     void magma_weight_1d_kernel(const CeedScalar *dqweight1d, CeedScalar *dV, const int v_stride, const int nelem) {
30f80f4a74SSebastian Grimberg   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
31f80f4a74SSebastian Grimberg 
32f80f4a74SSebastian Grimberg   const int tx      = threadIdx.x;
33f80f4a74SSebastian Grimberg   const int ty      = threadIdx.y;
34f80f4a74SSebastian Grimberg   const int elem_id = (blockIdx.x * blockDim.y) + ty;
35f80f4a74SSebastian Grimberg 
36f80f4a74SSebastian Grimberg   if (elem_id >= nelem) return;
37f80f4a74SSebastian Grimberg 
38f80f4a74SSebastian Grimberg   // global memory pointers
39f80f4a74SSebastian Grimberg   dV += elem_id * v_stride;
40f80f4a74SSebastian Grimberg 
41f80f4a74SSebastian Grimberg   // shared memory pointers
42f80f4a74SSebastian Grimberg   CeedScalar *sTweight = (CeedScalar *)shared_data;
433c1e2affSSebastian Grimberg   CeedScalar *sV       = sTweight + BASIS_Q;
443c1e2affSSebastian Grimberg   sV += ty * BASIS_Q;
45f80f4a74SSebastian Grimberg 
46f80f4a74SSebastian Grimberg   // read dqweight_1d
473c1e2affSSebastian Grimberg   if (ty == 0 && tx < BASIS_Q) {
48f80f4a74SSebastian Grimberg     sTweight[tx] = dqweight1d[tx];
49f80f4a74SSebastian Grimberg   }
50f80f4a74SSebastian Grimberg 
51f80f4a74SSebastian Grimberg   __syncthreads();
523c1e2affSSebastian Grimberg   magma_weight_1d_device<CeedScalar, BASIS_Q>(sTweight, sV, tx);
53f80f4a74SSebastian Grimberg   __syncthreads();
54f80f4a74SSebastian Grimberg 
55f80f4a74SSebastian Grimberg   // write V
56f80f4a74SSebastian Grimberg   dV[tx] = sV[tx];
57f80f4a74SSebastian Grimberg }
583c1e2affSSebastian Grimberg 
593c1e2affSSebastian Grimberg #endif  // CEED_MAGMA_BASIS_WEIGHT_1D_H
60