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