xref: /libCEED/include/ceed/jit-source/magma/magma-basis-weight-nontensor.h (revision 9d15e85b4f78ffb2d2860753c87a3b1789cc3bb6)
1940a72f1SSebastian Grimberg // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2940a72f1SSebastian Grimberg // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3940a72f1SSebastian Grimberg //
4940a72f1SSebastian Grimberg // SPDX-License-Identifier: BSD-2-Clause
5940a72f1SSebastian Grimberg //
6940a72f1SSebastian Grimberg // This file is part of CEED:  http://github.com/ceed
7940a72f1SSebastian Grimberg 
8940a72f1SSebastian Grimberg /// @file
9940a72f1SSebastian Grimberg /// Internal header for MAGMA non-tensor basis weight
10940a72f1SSebastian Grimberg #ifndef CEED_MAGMA_BASIS_WEIGHT_NONTENSOR_H
11940a72f1SSebastian Grimberg #define CEED_MAGMA_BASIS_WEIGHT_NONTENSOR_H
12940a72f1SSebastian Grimberg 
13940a72f1SSebastian Grimberg #include "magma-common-nontensor.h"
14940a72f1SSebastian Grimberg 
15940a72f1SSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
16940a72f1SSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_Q, MAGMA_MAXTHREADS_1D)) __global__
17*9d15e85bSSebastian Grimberg     void magma_weight_nontensor(int n, const CeedScalar *__restrict__ dqweight, CeedScalar *__restrict__ dV) {
18940a72f1SSebastian Grimberg   MAGMA_DEVICE_SHARED(CeedScalar, shared_data);
19940a72f1SSebastian Grimberg 
20940a72f1SSebastian Grimberg   const int tx = threadIdx.x;
21940a72f1SSebastian Grimberg   const int ty = threadIdx.y;
22940a72f1SSebastian Grimberg   const int id = blockIdx.x * blockDim.y + ty;
23940a72f1SSebastian Grimberg 
24940a72f1SSebastian Grimberg   // terminate threads with no work
25940a72f1SSebastian Grimberg   if (id >= n) return;
26940a72f1SSebastian Grimberg 
27*9d15e85bSSebastian Grimberg   dV += id * BASIS_Q;
28940a72f1SSebastian Grimberg 
29940a72f1SSebastian Grimberg   // shared memory pointers
30940a72f1SSebastian Grimberg   CeedScalar *sqweight = (CeedScalar *)shared_data;
31940a72f1SSebastian Grimberg   CeedScalar *sV       = sqweight + BASIS_Q;
32940a72f1SSebastian Grimberg   sV += ty * BASIS_Q;
33940a72f1SSebastian Grimberg 
34940a72f1SSebastian Grimberg   // read qweight
35940a72f1SSebastian Grimberg   if (ty == 0 && tx < BASIS_Q) {
36940a72f1SSebastian Grimberg     sqweight[tx] = dqweight[tx];
37940a72f1SSebastian Grimberg   }
38940a72f1SSebastian Grimberg   __syncthreads();
39940a72f1SSebastian Grimberg 
40940a72f1SSebastian Grimberg   if (tx < BASIS_Q) {
41940a72f1SSebastian Grimberg     sV[tx] = sqweight[tx];
42940a72f1SSebastian Grimberg   }
43940a72f1SSebastian Grimberg 
44940a72f1SSebastian Grimberg   // write V
45940a72f1SSebastian Grimberg   dV[tx] = sV[tx];
46940a72f1SSebastian Grimberg }
47940a72f1SSebastian Grimberg 
48940a72f1SSebastian Grimberg #endif  // CEED_MAGMA_BASIS_WEIGHT_NONTENSOR_H
49