22 #ifndef GEOS_FINITEELEMENT_IMPLICITKERNELBASE_HPP_
23 #define GEOS_FINITEELEMENT_IMPLICITKERNELBASE_HPP_
46 template<
typename SUBREGION_TYPE,
47 typename CONSTITUTIVE_TYPE,
49 int NUM_DOF_PER_TEST_SP,
50 int NUM_DOF_PER_TRIAL_SP >
55 NUM_DOF_PER_TRIAL_SP >
63 NUM_DOF_PER_TRIAL_SP >;
90 SUBREGION_TYPE
const & elementSubRegion,
91 FE_TYPE
const & finiteElementSpace,
92 CONSTITUTIVE_TYPE & inputConstitutiveType,
98 Base( elementSubRegion,
100 inputConstitutiveType ),
107 FiniteElementBase::initialize< FE_TYPE >( nodeManager,
138 for(
int ii = 0; ii < maxNumRows; ++ii )
140 for(
int jj = 0; jj < maxNumCols; ++jj )
142 localJacobian[ii][jj] = 0.0;
163 real64 localJacobian[maxNumRows][maxNumCols];
186 m_finiteElementSpace.template setup< FE_TYPE >( k, m_meshData, stack.
feStack );
187 localIndex numTestSupportPoints = m_finiteElementSpace.template numSupportPoints< FE_TYPE >( stack.
feStack );
188 localIndex numTrialSupportPoints = numTestSupportPoints;
189 stack.
numRows = numTestSupportPoints * numDofPerTestSupportPoint;
190 stack.
numCols = numTrialSupportPoints * numDofPerTrialSupportPoint;
191 for(
localIndex a=0; a<numTestSupportPoints; ++a )
193 localIndex const localNodeIndex = m_elemsToNodes[k][a];
194 for(
int i=0; i<numDofPerTestSupportPoint; ++i )
196 stack.
localRowDofIndex[a*numDofPerTestSupportPoint+i] = m_dofNumber[localNodeIndex]+i;
200 for(
localIndex a=0; a<numTrialSupportPoints; ++a )
202 localIndex const localNodeIndex = m_elemsToNodes[k][a];
203 for(
int i=0; i<numDofPerTrialSupportPoint; ++i )
205 stack.
localColDofIndex[a*numDofPerTrialSupportPoint+i] = m_dofNumber[localNodeIndex]+i;
225 typename FE_TYPE::template MeshData< SUBREGION_TYPE >
m_meshData;
#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.
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.
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
ImplicitKernelBase(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)
Constructor.
static constexpr int numDofPerTrialSupportPoint
GEOS_HOST_DEVICE void setup(localIndex const k, StackVariables &stack) const
Performs the setup phase for the kernel.
arrayView1d< globalIndex const > const m_dofNumber
The global degree of freedom number.
static constexpr int maxNumTrialSupportPointsPerElem
real64 const m_dt
time increment
static constexpr int maxNumTestSupportPointsPerElem
arrayView1d< real64 > const m_rhs
The global residaul vector.
Define the base interface for finite element kernels.
traits::ViewTypeConst< typename SUBREGION_TYPE::NodeMapType::base_type > const m_elemsToNodes
The element to nodes map.
FE_TYPE const m_finiteElementSpace
static constexpr int numDofPerTestSupportPoint
static constexpr int numDofPerTrialSupportPoint
static constexpr int maxNumTrialSupportPointsPerElem
static constexpr int maxNumTestSupportPointsPerElem
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix 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).
Kernel variables allocated on the stack.
localIndex numCols
The actual number of columns in the element local jacobian matrix (<= maxNumCols).
real64 localJacobian[maxNumRows][maxNumCols]
C-array storage for the element local Jacobian matrix.
GEOS_HOST_DEVICE StackVariables()
globalIndex localRowDofIndex[maxNumRows]
C-array storage for the element local row degrees of freedom.
globalIndex localColDofIndex[maxNumCols]
C-array storage for the element local column degrees of freedom.
FE_TYPE::StackVariables feStack
Stack variables needed for the underlying FEM type.
static constexpr int maxNumRows
The number of rows in the pre-allocated element local jacobian matrix (upper bound for numRows).
real64 localResidual[maxNumRows]
C-array storage for the element local residual vector.
static constexpr int maxNumCols
The number of columns in the pre-allocated element local jacobian matrix (upper bound for numCols).
localIndex numRows
The actual number of rows in the element local jacobian matrix (<= maxNumRows).
Kernel variables allocated on the stack.