21 #ifndef GEOSX_FINITEELEMENT_IMPLICITKERNELBASE_HPP_ 22 #define GEOSX_FINITEELEMENT_IMPLICITKERNELBASE_HPP_ 45 template<
typename SUBREGION_TYPE,
46 typename CONSTITUTIVE_TYPE,
48 int NUM_DOF_PER_TEST_SP,
49 int NUM_DOF_PER_TRIAL_SP >
54 NUM_DOF_PER_TRIAL_SP >
62 NUM_DOF_PER_TRIAL_SP >;
85 SUBREGION_TYPE
const & elementSubRegion,
86 FE_TYPE
const & finiteElementSpace,
87 CONSTITUTIVE_TYPE *
const inputConstitutiveType,
92 Base( elementSubRegion,
94 inputConstitutiveType ),
164 stack.localRowDofIndex[a*numDofPerTestSupportPoint+i] =
m_dofNumber[localNodeIndex]+i;
173 stack.localColDofIndex[a*numDofPerTrialSupportPoint+i] =
m_dofNumber[localNodeIndex]+i;
#define GEOSX_HOST_DEVICE
Marks a host-device function.
Define the base interface for implicit finite element kernels.
real64 localResidual[numRows]
C-array storage for the element local residual vector.
Kernel variables allocated on the stack.
long long int globalIndex
Global index type (for indexing objects across MPI partitions).
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data...
ImplicitKernelBase(NodeManager const &nodeManager, EdgeManager const &edgeManager, FaceManager const &faceManager, SUBREGION_TYPE const &elementSubRegion, FE_TYPE const &finiteElementSpace, CONSTITUTIVE_TYPE *const inputConstitutiveType, arrayView1d< globalIndex const > const &inputDofNumber, globalIndex const rankOffset, CRSMatrixView< real64, globalIndex const > const &inputMatrix, arrayView1d< real64 > const &inputRhs)
Constructor.
Kernel variables allocated on the stack.
GEOSX_HOST_DEVICE GEOSX_FORCE_INLINE void setup(localIndex const k, StackVariables &stack) const
Performs the setup phase for the kernel.
traits::ViewTypeConst< typename SUBREGION_TYPE::NodeMapType::base_type > const m_elemsToNodes
The element to nodes map.
This class serves to provide a "view" of a multidimensional array.
static constexpr int numRows
The number of rows in the element local jacobian matrix.
arrayView1d< globalIndex const > const m_dofNumber
The global degree of freedom number.
globalIndex const m_dofRankOffset
The global rank offset.
double real64
64-bit floating point type.
Define the base interface for finite element kernels.
static constexpr int numTrialSupportPointsPerElem
GEOSX_HOST_DEVICE StackVariables()
globalIndex localColDofIndex[numCols]
C-array storage for the element local column degrees of freedom.
static constexpr int numTestSupportPointsPerElem
static constexpr int numDofPerTestSupportPoint
This class provides an interface to ObjectManagerBase in order to manage edge data.
#define GEOSX_FORCE_INLINE
Marks a function or lambda for inlining.
static constexpr int numCols
The number of columns in the element local jacobian matrix.
#define GEOSX_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
static constexpr int numDofPerTrialSupportPoint
std::ptrdiff_t localIndex
Local index type (for indexing objects within an MPI partition).
CRSMatrixView< real64, globalIndex const > const m_matrix
The global Jacobian matrix.
The FaceManager class provides an interface to ObjectManagerBase in order to manage face data...
globalIndex localRowDofIndex[numRows]
C-array storage for the element local row degrees of freedom.
arrayView1d< real64 > const m_rhs
The global residaul vector.
real64 localJacobian[numRows][numCols]
C-array storage for the element local Jacobian matrix.