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