xref: /libCEED/include/ceed/jit-source/magma/magma-basis-weight-nontensor.h (revision 49a40c8a2d720db341b0b117b89656b473cbebfb)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Internal header for MAGMA non-tensor basis weight
10 #ifndef CEED_MAGMA_BASIS_WEIGHT_NONTENSOR_H
11 #define CEED_MAGMA_BASIS_WEIGHT_NONTENSOR_H
12 
13 #include "magma-common-nontensor.h"
14 
15 ////////////////////////////////////////////////////////////////////////////////
16 extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_Q, MAGMA_MAXTHREADS_1D)) __global__
17     void magma_weight_nontensor(int n, const CeedScalar *__restrict__ dqweight, CeedScalar *__restrict__ dV) {
18   MAGMA_DEVICE_SHARED(CeedScalar, shared_data);
19 
20   const int tx = threadIdx.x;
21   const int ty = threadIdx.y;
22   const int id = blockIdx.x * blockDim.y + ty;
23 
24   // terminate threads with no work
25   if (id >= n) return;
26 
27   dV += id * BASIS_Q;
28 
29   // shared memory pointers
30   CeedScalar *sqweight = (CeedScalar *)shared_data;
31   CeedScalar *sV       = sqweight + BASIS_Q;
32   sV += ty * BASIS_Q;
33 
34   // read qweight
35   if (ty == 0 && tx < BASIS_Q) {
36     sqweight[tx] = dqweight[tx];
37   }
38   __syncthreads();
39 
40   if (tx < BASIS_Q) {
41     sV[tx] = sqweight[tx];
42   }
43 
44   // write V
45   dV[tx] = sV[tx];
46 }
47 
48 #endif  // CEED_MAGMA_BASIS_WEIGHT_NONTENSOR_H
49