xref: /libCEED/include/ceed/jit-source/magma/magma-basis-weight-nontensor.h (revision d4cc18453651bd0f94c1a2e078b2646a92dafdcc)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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 #include "magma-common-nontensor.h"
11940a72f1SSebastian Grimberg 
12940a72f1SSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
__launch_bounds__(MAGMA_BASIS_BOUNDS (BASIS_Q,MAGMA_MAXTHREADS_1D))13940a72f1SSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_Q, MAGMA_MAXTHREADS_1D)) __global__
149d15e85bSSebastian Grimberg     void magma_weight_nontensor(int n, const CeedScalar *__restrict__ dqweight, CeedScalar *__restrict__ dV) {
15940a72f1SSebastian Grimberg   MAGMA_DEVICE_SHARED(CeedScalar, shared_data);
16940a72f1SSebastian Grimberg 
17940a72f1SSebastian Grimberg   const int tx = threadIdx.x;
18940a72f1SSebastian Grimberg   const int ty = threadIdx.y;
19940a72f1SSebastian Grimberg   const int id = blockIdx.x * blockDim.y + ty;
20940a72f1SSebastian Grimberg 
21940a72f1SSebastian Grimberg   // terminate threads with no work
22940a72f1SSebastian Grimberg   if (id >= n) return;
23940a72f1SSebastian Grimberg 
249d15e85bSSebastian Grimberg   dV += id * BASIS_Q;
25940a72f1SSebastian Grimberg 
26940a72f1SSebastian Grimberg   // shared memory pointers
27940a72f1SSebastian Grimberg   CeedScalar *sqweight = (CeedScalar *)shared_data;
28940a72f1SSebastian Grimberg   CeedScalar *sV       = sqweight + BASIS_Q;
29940a72f1SSebastian Grimberg   sV += ty * BASIS_Q;
30940a72f1SSebastian Grimberg 
31940a72f1SSebastian Grimberg   // read qweight
32940a72f1SSebastian Grimberg   if (ty == 0 && tx < BASIS_Q) {
33940a72f1SSebastian Grimberg     sqweight[tx] = dqweight[tx];
34940a72f1SSebastian Grimberg   }
35940a72f1SSebastian Grimberg   __syncthreads();
36940a72f1SSebastian Grimberg 
37940a72f1SSebastian Grimberg   if (tx < BASIS_Q) {
38940a72f1SSebastian Grimberg     sV[tx] = sqweight[tx];
39940a72f1SSebastian Grimberg   }
40940a72f1SSebastian Grimberg 
41940a72f1SSebastian Grimberg   // write V
42940a72f1SSebastian Grimberg   dV[tx] = sV[tx];
43940a72f1SSebastian Grimberg }
44