GEOSX
Classes | Public Types | Public Member Functions | Protected Attributes | List of all members
geosx::SolidMechanicsLagrangianFEMKernels::ImplicitNewmark< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE > Class Template Reference

Implements kernels for solving the equations of motion using an implicit Newmark's method.. More...

#include <SolidMechanicsSmallStrainImplicitNewmarkKernel.hpp>

Inheritance diagram for geosx::SolidMechanicsLagrangianFEMKernels::ImplicitNewmark< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >:
Inheritance graph
[legend]

Classes

class  StackVariables
 Kernel variables allocated on the stack. More...
 

Public Types

using Base = QuasiStatic< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >
 Alias for the base class;.
 
- Public Types inherited from geosx::SolidMechanicsLagrangianFEMKernels::QuasiStatic< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >
using Base = finiteElement::ImplicitKernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, 3, 3 >
 Alias for the base class;.
 
- Public Types inherited from geosx::finiteElement::ImplicitKernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, 3, 3 >
using Base = KernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, NUM_DOF_PER_TEST_SP, NUM_DOF_PER_TRIAL_SP >
 Alias for the base class. (i.e. geosx::finiteElement::KernelBase)
 

Public Member Functions

 ImplicitNewmark (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 (&inputGravityVector)[3], real64 const inputNewmarkGamma, real64 const inputNewmarkBeta, real64 const inputMassDamping, real64 const inputStiffnessDamping, real64 const inputDt)
 Constructor. More...
 
GEOSX_HOST_DEVICE GEOSX_FORCE_INLINE void setup (localIndex const k, StackVariables &stack) const
 Copy global values from primary field to a local stack array. More...
 
GEOSX_DEVICE GEOSX_FORCE_INLINE void quadraturePointKernel (localIndex const k, localIndex const q, StackVariables &stack) const
 Performs a state update at a quadrature point. More...
 
GEOSX_HOST_DEVICE GEOSX_FORCE_INLINE real64 complete (localIndex const k, StackVariables &stack) const
 Performs the complete phase for the kernel. More...
 
- Public Member Functions inherited from geosx::SolidMechanicsLagrangianFEMKernels::QuasiStatic< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >
 QuasiStatic (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 (&inputGravityVector)[3])
 Constructor. More...
 
GEOSX_HOST_DEVICE GEOSX_FORCE_INLINE void setup (localIndex const k, StackVariables &stack) const
 Copy global values from primary field to a local stack array. More...
 
template<typename STRESS_MODIFIER = NoOpFunctors>
GEOSX_HOST_DEVICE GEOSX_FORCE_INLINE void quadraturePointKernel (localIndex const k, localIndex const q, StackVariables &stack, STRESS_MODIFIER &&stressModifier=NoOpFunctors{}) const
 Performs a state update at a quadrature point. More...
 
GEOSX_HOST_DEVICE GEOSX_FORCE_INLINE real64 complete (localIndex const k, StackVariables &stack) const
 Performs the complete phase for the kernel. More...
 
- Public Member Functions inherited from geosx::finiteElement::ImplicitKernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, 3, 3 >
 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)
 Constructor. More...
 
GEOSX_HOST_DEVICE GEOSX_FORCE_INLINE void setup (localIndex const k, StackVariables &stack) const
 Performs the setup phase for the kernel. More...
 
- Public Member Functions inherited from geosx::finiteElement::KernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, NUM_DOF_PER_TEST_SP, NUM_DOF_PER_TRIAL_SP >
 KernelBase (SUBREGION_TYPE const &elementSubRegion, FE_TYPE const &finiteElementSpace, CONSTITUTIVE_TYPE &inputConstitutiveType)
 Constructor. More...
 
GEOSX_HOST_DEVICE GEOSX_FORCE_INLINE void setup (localIndex const k, StackVariables &stack) const
 Performs the setup phase for the kernel. More...
 
GEOSX_HOST_DEVICE GEOSX_FORCE_INLINE void quadraturePointKernel (localIndex const k, localIndex const q, StackVariables &stack) const
 Performs a state update at a quadrature point. More...
 
GEOSX_HOST_DEVICE GEOSX_FORCE_INLINE real64 complete (localIndex const k, StackVariables &stack) const
 Performs the complete phase for the kernel. More...
 

Protected Attributes

arrayView2d< real64 const, nodes::TOTAL_DISPLACEMENT_USD > const m_vtilde
 The rank-global velocity predictor.
 
arrayView2d< real64 const, nodes::INCR_DISPLACEMENT_USD > const m_uhattilde
 The rank-global incremental displacement predictor.
 
real64 const m_newmarkGamma
 The Gamma parameter for Newmark's method.
 
real64 const m_newmarkBeta
 The Beta parameter for Newmark's method.
 
real64 const m_massDamping
 The mass damping coefficient.
 
real64 const m_stiffnessDamping
 The stiffness damping coefficient.
 
real64 const m_dt
 The timestep for the update.
 
- Protected Attributes inherited from geosx::SolidMechanicsLagrangianFEMKernels::QuasiStatic< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >
arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const m_X
 The array containing the nodal position array.
 
arrayView2d< real64 const, nodes::TOTAL_DISPLACEMENT_USD > const m_disp
 The rank-global displacement array.
 
arrayView2d< real64 const, nodes::INCR_DISPLACEMENT_USD > const m_uhat
 The rank-global incremental displacement array.
 
real64 const m_gravityVector [3]
 The gravity vector.
 
arrayView2d< real64 const > const m_density
 The rank global density.
 
- Protected Attributes inherited from geosx::finiteElement::ImplicitKernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, 3, 3 >
arrayView1d< globalIndex const > const m_dofNumber
 The global degree of freedom number.
 
globalIndex const m_dofRankOffset
 The global rank offset.
 
CRSMatrixView< real64, globalIndex const > const m_matrix
 The global Jacobian matrix.
 
arrayView1d< real64 > const m_rhs
 The global residaul vector.
 
- Protected Attributes inherited from geosx::finiteElement::KernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, NUM_DOF_PER_TEST_SP, NUM_DOF_PER_TRIAL_SP >
traits::ViewTypeConst< typename SUBREGION_TYPE::NodeMapType::base_type > const m_elemsToNodes
 The element to nodes map.
 
arrayView1d< integer const > const m_elemGhostRank
 The element ghost rank array.
 
CONSTITUTIVE_TYPE::KernelWrapper const m_constitutiveUpdate
 
FE_TYPE const & m_finiteElementSpace
 

Additional Inherited Members

- Static Public Member Functions inherited from geosx::finiteElement::KernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, NUM_DOF_PER_TEST_SP, NUM_DOF_PER_TRIAL_SP >
template<typename POLICY , typename KERNEL_TYPE >
static real64 kernelLaunch (localIndex const numElems, KERNEL_TYPE const &kernelComponent)
 Kernel Launcher. More...
 
- Static Public Attributes inherited from geosx::SolidMechanicsLagrangianFEMKernels::QuasiStatic< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >
static constexpr int numNodesPerElem = Base::numTestSupportPointsPerElem
 
- Static Public Attributes inherited from geosx::finiteElement::KernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, NUM_DOF_PER_TEST_SP, NUM_DOF_PER_TRIAL_SP >
static constexpr int numTestSupportPointsPerElem = FE_TYPE::numNodes
 
static constexpr int numTrialSupportPointsPerElem = FE_TYPE::numNodes
 
static constexpr int numDofPerTestSupportPoint = NUM_DOF_PER_TEST_SP
 
static constexpr int numDofPerTrialSupportPoint = NUM_DOF_PER_TRIAL_SP
 
static constexpr int numQuadraturePointsPerElem = FE_TYPE::numQuadraturePoints
 Compile time value for the number of quadrature points per element.
 

Detailed Description

template<typename SUBREGION_TYPE, typename CONSTITUTIVE_TYPE, typename FE_TYPE>
class geosx::SolidMechanicsLagrangianFEMKernels::ImplicitNewmark< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >

Implements kernels for solving the equations of motion using an implicit Newmark's method..

Constructor. Constructor.

Parameters
nodeManagerReference to the NodeManager object.
edgeManagerReference to the EdgeManager object.
faceManagerReference to the FaceManager object.
targetRegionIndexIndex of the region the subregion belongs to.
inputDofNumberThe dof number for the primary field.
rankOffsetdof index offset of current rank
inputMatrixReference to the Jacobian matrix.
inputRhsReference to the RHS vector. Constructor.
elementSubRegionReference to the SUBREGION_TYPE(class template parameter) object.
inputConstitutiveTypeThe constitutive object.
finiteElementSpacePlaceholder for the finite element space object, which currently doesn't do much.
inputGravityVectorThe gravity vector.

Implicit Newmark Description

Implements the KernelBase interface functions required for solving the equations of motion using with an Implicit Newmark's Method with one of the "finite element kernel application" functions such as geosx::finiteElement::RegionBasedKernelApplication.

Definition at line 46 of file SolidMechanicsSmallStrainImplicitNewmarkKernel.hpp.

Constructor & Destructor Documentation

◆ ImplicitNewmark()

template<typename SUBREGION_TYPE , typename CONSTITUTIVE_TYPE , typename FE_TYPE >
geosx::SolidMechanicsLagrangianFEMKernels::ImplicitNewmark< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >::ImplicitNewmark ( 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 (&)  inputGravityVector[3],
real64 const  inputNewmarkGamma,
real64 const  inputNewmarkBeta,
real64 const  inputMassDamping,
real64 const  inputStiffnessDamping,
real64 const  inputDt 
)
inline

Constructor.

Constructor. Constructor.

Parameters
nodeManagerReference to the NodeManager object.
edgeManagerReference to the EdgeManager object.
faceManagerReference to the FaceManager object.
targetRegionIndexIndex of the region the subregion belongs to.
inputDofNumberThe dof number for the primary field.
rankOffsetdof index offset of current rank
inputMatrixReference to the Jacobian matrix.
inputRhsReference to the RHS vector. Constructor.
elementSubRegionReference to the SUBREGION_TYPE(class template parameter) object.
inputConstitutiveTypeThe constitutive object.
finiteElementSpacePlaceholder for the finite element space object, which currently doesn't do much.
inputGravityVectorThe gravity vector.
inputNewmarkGammaThe Gamma parameter of the Newmark method.
inputNewmarkBetaThe Beta parameter for the Newmark method.
inputMassDampingThe mass damping coefficient.
inputStiffnessDampingThe stiffness damping coefficient.
inputDtThe timestep for the physics update.

Definition at line 82 of file SolidMechanicsSmallStrainImplicitNewmarkKernel.hpp.

Member Function Documentation

◆ complete()

template<typename SUBREGION_TYPE , typename CONSTITUTIVE_TYPE , typename FE_TYPE >
GEOSX_HOST_DEVICE GEOSX_FORCE_INLINE real64 geosx::SolidMechanicsLagrangianFEMKernels::ImplicitNewmark< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >::complete ( localIndex const  k,
StackVariables stack 
) const
inline

Performs the complete phase for the kernel.

Template Parameters
STACK_VARIABLE_TYPEThe type of StackVariable that holds the stack variables. This is most likely a defined in a type that derives from KernelBase.
Parameters
kThe element index.
stackThe StackVariable object that hold the stack variables.
Returns
The maximum contribution to the residual.

KernelBase::complete() Description

The operations typically found in complete are the mapping of the local Jacobian and Residual into the global Jacobian and Residual.

The ImplicitNewmark implementation adds residual and jacobian contributions from stiffness based damping.

Definition at line 230 of file SolidMechanicsSmallStrainImplicitNewmarkKernel.hpp.

◆ quadraturePointKernel()

template<typename SUBREGION_TYPE , typename CONSTITUTIVE_TYPE , typename FE_TYPE >
GEOSX_DEVICE GEOSX_FORCE_INLINE void geosx::SolidMechanicsLagrangianFEMKernels::ImplicitNewmark< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >::quadraturePointKernel ( localIndex const  k,
localIndex const  q,
StackVariables stack 
) const
inline

Performs a state update at a quadrature point.

Template Parameters
STACK_VARIABLE_TYPEThe type of StackVariable that holds the stack variables. This is most likely a defined in a type that derives from KernelBase.
Parameters
kThe element index.
qThe quadrature point index.
stackThe StackVariable object that hold the stack variables.

KernelBase::quadraturePointKernel() Description

The operations found here are the mapping from the support points to the quadrature point, calculation of gradients, etc. From this data the state of the constitutive model is updated if required by the physics package.

The ImplcitNewmark kernel adds the calculation of the inertia damping, jacobian and residual contributions.

Definition at line 188 of file SolidMechanicsSmallStrainImplicitNewmarkKernel.hpp.

◆ setup()

template<typename SUBREGION_TYPE , typename CONSTITUTIVE_TYPE , typename FE_TYPE >
GEOSX_HOST_DEVICE GEOSX_FORCE_INLINE void geosx::SolidMechanicsLagrangianFEMKernels::ImplicitNewmark< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >::setup ( localIndex const  k,
StackVariables stack 
) const
inline

Copy global values from primary field to a local stack array.

Copy global values from primary field to a local stack array. Performs the setup phase for the kernel.

Template Parameters
STACK_VARIABLE_TYPEThe type of StackVariable that holds the stack variables. This is most likely a defined in a type that derives from KernelBase.
Parameters
kThe element index.
stackThe StackVariable object that hold the stack variables.

KernelBase::setup() Description

The operations typically found in setup are thing such as the collection of global data into local stack storage.

ImplicitKernelBase::setup() Description

In this implementation, the element local Row and Column DOF stack arrays are filled for when we fill the global matrix and rhs.

Note
This seems like a waste of register space. We should do this in complete() unless we actually need these dof somewhere else in the kernel.

For the QuasiStatic implementation, global values from the displacement, incremental displacement, and degree of freedom numbers are placed into element local stack storage.

For the ImplicitNewmark implementation, global values from the velocity predictor, and the incremental displacement predictor are placed into element local stack storage.

Definition at line 165 of file SolidMechanicsSmallStrainImplicitNewmarkKernel.hpp.


The documentation for this class was generated from the following file: