20 #ifndef GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_KERNELS_IMPLCITSMALLSTRAINQUASISTATIC_HPP_
21 #define GEOS_PHYSICSSOLVERS_SOLIDMECHANICS_KERNELS_IMPLCITSMALLSTRAINQUASISTATIC_HPP_
29 namespace solidMechanicsLagrangianFEMKernels
56 template<
typename SUBREGION_TYPE,
57 typename CONSTITUTIVE_TYPE,
99 SUBREGION_TYPE
const & elementSubRegion,
100 FE_TYPE
const & finiteElementSpace,
101 CONSTITUTIVE_TYPE & inputConstitutiveType,
107 real64 const (&inputGravityVector)[3] );
131 #if !defined(CALC_FEM_SHAPE_IN_KERNEL)
205 template<
typename STRESS_MODIFIER = NoOpFunctors >
210 STRESS_MODIFIER && stressModifier =
NoOpFunctors{} )
const;
218 StackVariables & stack )
const;
223 template<
typename POLICY,
224 typename KERNEL_TYPE >
227 KERNEL_TYPE
const & kernelComponent );
#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.
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.
CONSTITUTIVE_TYPE::KernelWrapper const m_constitutiveUpdate
Used to forward arguments to a class that implements the KernelBase interface.
Implements kernels for solving quasi-static equilibrium.
arrayView2d< real64 const, nodes::INCR_DISPLACEMENT_USD > const m_uhat
The rank-global incremental displacement array.
static real64 kernelLaunch(localIndex const numElems, KERNEL_TYPE const &kernelComponent)
Kernel Launcher.
GEOS_HOST_DEVICE void quadraturePointKernel(localIndex const k, localIndex const q, StackVariables &stack, STRESS_MODIFIER &&stressModifier=NoOpFunctors{}) const
Performs a state update at a quadrature point.
arrayView2d< real64 const, nodes::TOTAL_DISPLACEMENT_USD > const m_disp
The rank-global displacement array.
static constexpr int numNodesPerElem
arrayView2d< real64 const > const m_density
The rank global density.
ImplicitSmallStrainQuasiStatic(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, real64 const (&inputGravityVector)[3])
Constructor.
real64 const m_gravityVector[3]
The gravity vector.
GEOS_HOST_DEVICE real64 complete(localIndex const k, StackVariables &stack) const
Performs the complete phase for the kernel.
GEOS_HOST_DEVICE real64 computeStabilizationScaling(localIndex const k) const
Get a parameter representative of the stiffness, used as physical scaling for the stabilization matri...
GEOS_HOST_DEVICE void setup(localIndex const k, StackVariables &stack) const
Copy global values from primary field to a local stack array.
arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const m_X
The array containing the nodal position array.
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).
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Kernel variables allocated on the stack.
Internal struct to provide no-op defaults used in the inclusion of lambda functions into kernel compo...
constexpr GEOS_HOST_DEVICE void operator()(localIndex const a, localIndex const b)
operator() no-op used for adding an additional dynamics term inside the jacobian assembly loop.
Kernel variables allocated on the stack.
real64 uhat_local[numNodesPerElem][numDofPerTrialSupportPoint]
Stack storage for the element local nodal incremental displacement.
real64 u_local[numNodesPerElem][numDofPerTrialSupportPoint]
Stack storage for the element local nodal displacement.
GEOS_HOST_DEVICE StackVariables()
Constructor.
real64 constitutiveStiffness[6][6]
Stack storage for the constitutive stiffness at a quadrature point.