21 #ifndef GEOS_PHYSICSSOLVERS_SIMPLEPDE_LAPLACEFEMKERNELS_HPP_ 
   22 #define GEOS_PHYSICSSOLVERS_SIMPLEPDE_LAPLACEFEMKERNELS_HPP_ 
   24 #define GEOS_DISPATCH_VEM  
   55 template< 
typename SUBREGION_TYPE,
 
   56           typename CONSTITUTIVE_TYPE,
 
   95                     SUBREGION_TYPE 
const & elementSubRegion,
 
   96                     FE_TYPE 
const & finiteElementSpace,
 
   97                     CONSTITUTIVE_TYPE & inputConstitutiveType,
 
  103                     string const fieldName ):
 
  110           inputConstitutiveType,
 
  116     m_X( nodeManager.referencePosition() ),
 
  141 #if !defined(CALC_FEM_SHAPE_IN_KERNEL) 
  169     stack.numCols = stack.numRows;
 
  170     for( 
localIndex a = 0; a < stack.numRows; ++a )
 
  174 #if defined(CALC_FEM_SHAPE_IN_KERNEL) 
  175       for( 
int i=0; i<3; ++i )
 
  177         stack.
xLocal[ a ][ i ] = 
m_X[ localNodeIndex ][ i ];
 
  182       stack.localRowDofIndex[a] = 
m_dofNumber[localNodeIndex];
 
  183       stack.localColDofIndex[a] = 
m_dofNumber[localNodeIndex];
 
  186     addGradGradStabilizationMatrix< FE_TYPE, numDofPerTrialSupportPoint >( stack.feStack,
 
  187                                                                            stack.localJacobian );
 
  201                                                                            stack.feStack, dNdX );
 
  202     for( 
localIndex a = 0; a < stack.numRows; ++a )
 
  204       for( 
localIndex b = 0; b < stack.numCols; ++b )
 
  206         stack.localJacobian[ a ][ b ] += LvArray::tensorOps::AiBi< 3 >( dNdX[a], dNdX[b] ) * detJ;
 
  226     for( 
localIndex a = 0; a < stack.numRows; ++a )
 
  228       for( 
localIndex b = 0; b < stack.numCols; ++b )
 
  230         stack.localResidual[ a ] += stack.localJacobian[ a ][ b ] * stack.
primaryField_local[ b ];
 
  234     for( 
int a = 0; a < stack.numRows; ++a )
 
  237       if( dof < 0 || dof >= 
m_matrix.numRows() ) 
continue;
 
  238       m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
 
  239                                                                               stack.localColDofIndex,
 
  240                                                                               stack.localJacobian[ a ],
 
  243       RAJA::atomicAdd< parallelDeviceAtomic >( &
m_rhs[ dof ], stack.localResidual[ a ] );
 
  244       maxForce = fmax( maxForce, fabs( stack.localResidual[ a ] ) );
 
#define GEOS_HOST_DEVICE
Marks a host-device function.
 
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
 
This class provides an interface to ObjectManagerBase in order to manage edge data.
 
The FaceManager class provides an interface to ObjectManagerBase in order to manage face data.
 
Implements kernels for solving Laplace's equation.
 
GEOS_HOST_DEVICE real64 complete(localIndex const k, StackVariables &stack) const
Performs the complete phase for the kernel.
 
LaplaceFEMKernel(NodeManager const &nodeManager, EdgeManager const &edgeManager, FaceManager const &faceManager, localIndex const targetRegionIndex, SUBREGION_TYPE const &elementSubRegion, FE_TYPE const &finiteElementSpace, CONSTITUTIVE_TYPE &inputConstitutiveType, arrayView1d< globalIndex const > const inputDofNumber, globalIndex const rankOffset, CRSMatrixView< real64, globalIndex const > const inputMatrix, arrayView1d< real64 > const inputRhs, real64 const inputDt, string const fieldName)
Constructor.
 
GEOS_HOST_DEVICE void quadraturePointKernel(localIndex const k, localIndex const q, StackVariables &stack) const
Performs a state update at a quadrature point.
 
arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const m_X
The array containing the nodal position array.
 
GEOS_HOST_DEVICE void setup(localIndex const k, StackVariables &stack) const
Copy global values from primary field to a local stack array.
 
arrayView1d< real64 const > const m_primaryField
The global primary field array.
 
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
 
Define the base interface for implicit finite element kernels.
 
globalIndex const m_dofRankOffset
The global rank offset.
 
traits::ViewTypeConst< typename SUBREGION_TYPE::NodeMapType::base_type > const m_elemsToNodes
The element to nodes map.
 
FE_TYPE const m_finiteElementSpace
 
CRSMatrixView< real64, globalIndex const > const m_matrix
The global Jacobian matrix.
 
FE_TYPE::template MeshData< SUBREGION_TYPE > m_meshData
Data structure containing mesh data used to setup the finite element.
 
static constexpr int numDofPerTestSupportPoint
 
static constexpr int numDofPerTrialSupportPoint
 
arrayView1d< globalIndex const > const m_dofNumber
The global degree of freedom number.
 
real64 const m_dt
time increment
 
static constexpr int maxNumTestSupportPointsPerElem
 
arrayView1d< real64 > const m_rhs
The global residaul vector.
 
Used to forward arguments to a class that implements the KernelBase interface.
 
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
 
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
 
double real64
64-bit floating point type.
 
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
 
LvArray::CRSMatrixView< T, COL_INDEX, INDEX_TYPE const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
 
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
 
Array< T, 1 > array1d
Alias for 1D array.
 
Kernel variables allocated on the stack.
 
GEOS_HOST_DEVICE StackVariables()
Constructor.
 
real64 primaryField_local[maxNumTestSupportPointsPerElem]
C-array storage for the element local primary field variable.
 
Kernel variables allocated on the stack.