xref: /libCEED/include/ceed/jit-source/magma/magma-basis-weight-3d.h (revision ca94c3ddc8f82b7d93a79f9e4812e99b8be840ff)
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 tensor basis weight in 3D
10 #ifndef CEED_MAGMA_BASIS_WEIGHT_3D_H
11 #define CEED_MAGMA_BASIS_WEIGHT_3D_H
12 
13 #include "magma-common-tensor.h"
14 
15 ////////////////////////////////////////////////////////////////////////////////
16 // weight basis action -- 3D
17 template <typename T, int DIM, int NUM_COMP, int Q, int i_DIM, int i_COMP>
18 static __device__ __inline__ void magma_weight_3d_device(const T *sTweight, T rV[DIM][NUM_COMP][Q], const int tx) {
19   // Assumptions
20   // 1. 1D thread configuration of size Q^2
21   // 2. rV[][][] matches the storage used in other actions (interp, grad, ... etc)
22   // 3. i_DIM and i_COMP specify which indexes to use in rV,
23   //    since the output per thread is a register array of size Q
24   // 4. Sync is recommended after the call (to make sure sTweight can be overwritten)
25 
26   if (tx < Q * Q) {
27     // x sTweight[j]    for first update
28     // x sTweight[tx%Q] for second update
29     // x sTweight[tx/Q] for third update
30     for (int j = 0; j < Q; j++) {
31       rV[i_DIM][i_COMP][j] = sTweight[j] * sTweight[tx % Q] * sTweight[tx / Q];
32     }
33   }
34 }
35 
36 ////////////////////////////////////////////////////////////////////////////////
37 extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_Q *BASIS_Q, MAGMA_MAXTHREADS_3D)) __global__
38     void magma_weight_3d_kernel(const CeedScalar *dqweight1d, CeedScalar *dV, const int v_stride, const int nelem) {
39   MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
40 
41   const int tx      = threadIdx.x;
42   const int ty      = threadIdx.y;
43   const int elem_id = (blockIdx.x * blockDim.y) + ty;
44 
45   if (elem_id >= nelem) return;
46 
47   CeedScalar rV[1][1][BASIS_Q];  // allocate with BASIS_DIM=BASIS_NUM_COMP=1, but sizes may differ for a fused operator
48   // global memory pointers
49   dV += elem_id * v_stride;
50 
51   // shared memory pointers
52   CeedScalar *sTweight = (CeedScalar *)shared_data;
53 
54   // read dqweight_1d
55   if (tx < BASIS_Q) {
56     sTweight[tx] = dqweight1d[tx];
57   }
58   __syncthreads();
59 
60   magma_weight_3d_device<CeedScalar, 1, 1, BASIS_Q, 0, 0>(sTweight, rV, tx);
61 
62   // write V
63   if (tx < (BASIS_Q * BASIS_Q)) {
64     for (int j = 0; j < BASIS_Q; j++) {
65       dV[j * (BASIS_Q * BASIS_Q) + tx] = rV[0][0][j];
66     }
67   }
68 }
69 
70 #endif  // CEED_MAGMA_BASIS_WEIGHT_3D_H
71