1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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 2D
103c1e2affSSebastian Grimberg #include "magma-common-tensor.h"
113c1e2affSSebastian Grimberg
129e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
13f80f4a74SSebastian Grimberg // weight basis action -- 2D
143c1e2affSSebastian Grimberg template <typename T, int DIM, int NUM_COMP, int Q, int i_DIM, int i_COMP>
magma_weight_2d_device(const T * sTweight,T rV[DIM][NUM_COMP][Q],const int tx)153c1e2affSSebastian Grimberg static __device__ __inline__ void magma_weight_2d_device(const T *sTweight, T rV[DIM][NUM_COMP][Q], const int tx) {
16f80f4a74SSebastian Grimberg // Assumptions
173c1e2affSSebastian Grimberg // 1. 1D thread configuration of size Q
18f80f4a74SSebastian Grimberg // 2. rV[][][] matches the storage used in other actions (interp, grad, ... etc)
193c1e2affSSebastian Grimberg // 3. i_DIM and i_COMP specify which indexes to use in rV,
203c1e2affSSebastian Grimberg // since the output per thread is a register array of size Q
21f80f4a74SSebastian Grimberg // 4. Sync is recommended after the call (to make sure sTweight can be overwritten)
22f80f4a74SSebastian Grimberg
233c1e2affSSebastian Grimberg if (tx < Q) {
24f80f4a74SSebastian Grimberg // x sTweight[j] for first update
25f80f4a74SSebastian Grimberg // x sTweight[tx] for second update
263c1e2affSSebastian Grimberg for (int j = 0; j < Q; j++) {
273c1e2affSSebastian Grimberg rV[i_DIM][i_COMP][j] = sTweight[j] * sTweight[tx];
28f80f4a74SSebastian Grimberg }
29f80f4a74SSebastian Grimberg }
30f80f4a74SSebastian Grimberg }
31f80f4a74SSebastian Grimberg
329e0c01faSSebastian Grimberg ////////////////////////////////////////////////////////////////////////////////
__launch_bounds__(MAGMA_BASIS_BOUNDS (BASIS_Q,MAGMA_MAXTHREADS_2D))333c1e2affSSebastian Grimberg extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_Q, MAGMA_MAXTHREADS_2D)) __global__
34f80f4a74SSebastian Grimberg void magma_weight_2d_kernel(const CeedScalar *dqweight1d, CeedScalar *dV, const int v_stride, const int nelem) {
35f80f4a74SSebastian Grimberg MAGMA_DEVICE_SHARED(CeedScalar, shared_data)
36f80f4a74SSebastian Grimberg
37f80f4a74SSebastian Grimberg const int tx = threadIdx.x;
38f80f4a74SSebastian Grimberg const int ty = threadIdx.y;
39f80f4a74SSebastian Grimberg const int elem_id = (blockIdx.x * blockDim.y) + ty;
40f80f4a74SSebastian Grimberg
41f80f4a74SSebastian Grimberg if (elem_id >= nelem) return;
42f80f4a74SSebastian Grimberg
433c1e2affSSebastian Grimberg CeedScalar rV[1][1][BASIS_Q]; // allocate with BASIS_DIM=BASIS_NUM_COMP=1, but sizes may differ for a fused operator
44f80f4a74SSebastian Grimberg // global memory pointers
45f80f4a74SSebastian Grimberg dV += elem_id * v_stride;
46f80f4a74SSebastian Grimberg
47f80f4a74SSebastian Grimberg // shared memory pointers
48f80f4a74SSebastian Grimberg CeedScalar *sTweight = (CeedScalar *)shared_data;
49f80f4a74SSebastian Grimberg
50f80f4a74SSebastian Grimberg // read dqweight_1d
513c1e2affSSebastian Grimberg if (ty == 0 && tx < BASIS_Q) {
52f80f4a74SSebastian Grimberg sTweight[tx] = dqweight1d[tx];
53f80f4a74SSebastian Grimberg }
54f80f4a74SSebastian Grimberg
55f80f4a74SSebastian Grimberg __syncthreads();
563c1e2affSSebastian Grimberg magma_weight_2d_device<CeedScalar, 1, 1, BASIS_Q, 0, 0>(sTweight, rV, tx);
57f80f4a74SSebastian Grimberg
58f80f4a74SSebastian Grimberg // write V
593c1e2affSSebastian Grimberg if (tx < BASIS_Q) {
603c1e2affSSebastian Grimberg for (int j = 0; j < BASIS_Q; j++) {
613c1e2affSSebastian Grimberg dV[j * BASIS_Q + tx] = rV[0][0][j];
62f80f4a74SSebastian Grimberg }
63f80f4a74SSebastian Grimberg }
64f80f4a74SSebastian Grimberg }
65