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 2D 10 #ifndef CEED_MAGMA_BASIS_WEIGHT_2D_H 11 #define CEED_MAGMA_BASIS_WEIGHT_2D_H 12 13 #include "magma-common-tensor.h" 14 15 ////////////////////////////////////////////////////////////////////////////////////////// 16 // weight basis action -- 2D 17 template <typename T, int DIM, int NUM_COMP, int Q, int i_DIM, int i_COMP> 18 static __device__ __inline__ void magma_weight_2d_device(const T *sTweight, T rV[DIM][NUM_COMP][Q], const int tx) { 19 // Assumptions 20 // 1. 1D thread configuration of size Q 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) { 27 // x sTweight[j] for first update 28 // x sTweight[tx] for second update 29 for (int j = 0; j < Q; j++) { 30 rV[i_DIM][i_COMP][j] = sTweight[j] * sTweight[tx]; 31 } 32 } 33 } 34 35 ////////////////////////////////////////////////////////////////////////////////////////// 36 extern "C" __launch_bounds__(MAGMA_BASIS_BOUNDS(BASIS_Q, MAGMA_MAXTHREADS_2D)) __global__ 37 void magma_weight_2d_kernel(const CeedScalar *dqweight1d, CeedScalar *dV, const int v_stride, const int nelem) { 38 MAGMA_DEVICE_SHARED(CeedScalar, shared_data) 39 40 const int tx = threadIdx.x; 41 const int ty = threadIdx.y; 42 const int elem_id = (blockIdx.x * blockDim.y) + ty; 43 44 if (elem_id >= nelem) return; 45 46 CeedScalar rV[1][1][BASIS_Q]; // allocate with BASIS_DIM=BASIS_NUM_COMP=1, but sizes may differ for a fused operator 47 // global memory pointers 48 dV += elem_id * v_stride; 49 50 // shared memory pointers 51 CeedScalar *sTweight = (CeedScalar *)shared_data; 52 53 // read dqweight_1d 54 if (ty == 0 && tx < BASIS_Q) { 55 sTweight[tx] = dqweight1d[tx]; 56 } 57 58 __syncthreads(); 59 magma_weight_2d_device<CeedScalar, 1, 1, BASIS_Q, 0, 0>(sTweight, rV, tx); 60 61 // write V 62 if (tx < BASIS_Q) { 63 for (int j = 0; j < BASIS_Q; j++) { 64 dV[j * BASIS_Q + tx] = rV[0][0][j]; 65 } 66 } 67 } 68 69 #endif // CEED_MAGMA_BASIS_WEIGHT_2D_H 70