GEOSX
Classes | Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Static Public Attributes | Protected Attributes | List of all members
geos::solidMechanicsLagrangianFEMKernels::ImplicitSmallStrainNewmark< 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 <ImplicitSmallStrainNewmark.hpp>

Inheritance diagram for geos::solidMechanicsLagrangianFEMKernels::ImplicitSmallStrainNewmark< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >:
Inheritance graph
[legend]

Classes

class  StackVariables
 

Public Types

using Base = ImplicitSmallStrainQuasiStatic< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >
 Alias for the base class;.
 
- Public Types inherited from geos::solidMechanicsLagrangianFEMKernels::ImplicitSmallStrainQuasiStatic< 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 geos::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. geos::finiteElement::KernelBase)
 

Public Member Functions

 ImplicitSmallStrainNewmark (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], real64 const inputNewmarkGamma, real64 const inputNewmarkBeta, real64 const inputMassDamping, real64 const inputStiffnessDamping)
 Constructor. More...
 
GEOS_HOST_DEVICE void setup (localIndex const k, StackVariables &stack) const
 Copy global values from primary field to a local stack array. More...
 
GEOS_HOST_DEVICE void quadraturePointKernel (localIndex const k, localIndex const q, StackVariables &stack) const
 Performs a state update at a quadrature point. More...
 
GEOS_HOST_DEVICE real64 complete (localIndex const k, StackVariables &stack) const
 
- Public Member Functions inherited from geos::solidMechanicsLagrangianFEMKernels::ImplicitSmallStrainQuasiStatic< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >
 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. More...
 
GEOS_HOST_DEVICE 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>
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. More...
 
GEOS_HOST_DEVICE real64 complete (localIndex const k, StackVariables &stack) const
 Performs the complete phase for the kernel. More...
 
template<typename STRESS_MODIFIER >
GEOS_HOST_DEVICE GEOS_FORCE_INLINE void quadraturePointKernel (localIndex const k, localIndex const q, StackVariables &stack, STRESS_MODIFIER &&stressModifier) const
 
template<typename POLICY , typename KERNEL_TYPE >
GEOS_FORCE_INLINE real64 kernelLaunch (localIndex const numElems, KERNEL_TYPE const &kernelComponent)
 
- Public Member Functions inherited from geos::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, real64 const inputDt)
 Constructor. More...
 
GEOS_HOST_DEVICE void setup (localIndex const k, StackVariables &stack) const
 Performs the setup phase for the kernel. More...
 
- Public Member Functions inherited from geos::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...
 
GEOS_HOST_DEVICE void setup (localIndex const k, StackVariables &stack) const
 Performs the setup phase for the kernel. More...
 
GEOS_HOST_DEVICE GEOS_FORCE_INLINE void quadraturePointKernel (localIndex const k, localIndex const q, StackVariables &stack) const
 Performs a state update at a quadrature point. More...
 
GEOS_HOST_DEVICE real64 complete (localIndex const k, StackVariables &stack) const
 Performs the complete phase for the kernel. More...
 

Static Public Member Functions

template<typename POLICY , typename KERNEL_TYPE >
static real64 kernelLaunch (localIndex const numElems, KERNEL_TYPE const &kernelComponent)
 Kernel Launcher. More...
 
- Static Public Member Functions inherited from geos::solidMechanicsLagrangianFEMKernels::ImplicitSmallStrainQuasiStatic< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >
template<typename POLICY , typename KERNEL_TYPE >
static real64 kernelLaunch (localIndex const numElems, KERNEL_TYPE const &kernelComponent)
 Kernel Launcher. More...
 
- Static Public Member Functions inherited from geos::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...
 

Public Attributes

const arrayView2d< real64 const, nodes::TOTAL_DISPLACEMENT_USD > m_disp
 The rank-global displacement array.
 
const arrayView2d< real64 const, nodes::INCR_DISPLACEMENT_USD > m_uhat
 The rank-global incremental displacement array.
 
const arrayView2d< real64 const > m_density
 The rank global density.
 
const arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > m_X
 The array containing the nodal position array.
 
const real64 m_gravityVector [3]
 The gravity vector.
 
- Public Attributes inherited from geos::finiteElement::ImplicitKernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, 3, 3 >
const traits::ViewTypeConst< typename SUBREGION_TYPE::NodeMapType::base_type > m_elemsToNodes
 The element to nodes map.
 
const FE_TYPE & m_finiteElementSpace
 

Static Public Attributes

static constexpr int numNodesPerElem
 
- Static Public Attributes inherited from geos::solidMechanicsLagrangianFEMKernels::ImplicitSmallStrainQuasiStatic< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >
static constexpr int numNodesPerElem = Base::maxNumTestSupportPointsPerElem
 
- Static Public Attributes inherited from geos::finiteElement::ImplicitKernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, 3, 3 >
static constexpr int maxNumTestSupportPointsPerElem
 
static constexpr int maxNumTrialSupportPointsPerElem
 
static constexpr int numDofPerTestSupportPoint
 
static constexpr int numDofPerTrialSupportPoint
 
- Static Public Attributes inherited from geos::finiteElement::KernelBase< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE, NUM_DOF_PER_TEST_SP, NUM_DOF_PER_TRIAL_SP >
static constexpr int maxNumTestSupportPointsPerElem = FE_TYPE::maxSupportPoints
 
static constexpr int maxNumTrialSupportPointsPerElem = FE_TYPE::maxSupportPoints
 
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.
 

Protected Attributes

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

Additional Inherited Members

- Protected Member Functions inherited from geos::solidMechanicsLagrangianFEMKernels::ImplicitSmallStrainQuasiStatic< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >
GEOS_HOST_DEVICE real64 computeStabilizationScaling (localIndex const k) const
 Get a parameter representative of the stiffness, used as physical scaling for the stabilization matrix. More...
 

Detailed Description

template<typename SUBREGION_TYPE, typename CONSTITUTIVE_TYPE, typename FE_TYPE>
class geos::solidMechanicsLagrangianFEMKernels::ImplicitSmallStrainNewmark< 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.
inputDtThe timestep for the physics update.

Constructor.

Parameters
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 geos::finiteElement::RegionBasedKernelApplication.

Definition at line 46 of file ImplicitSmallStrainNewmark.hpp.

Constructor & Destructor Documentation

◆ ImplicitSmallStrainNewmark()

template<typename SUBREGION_TYPE , typename CONSTITUTIVE_TYPE , typename FE_TYPE >
geos::solidMechanicsLagrangianFEMKernels::ImplicitSmallStrainNewmark< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >::ImplicitSmallStrainNewmark ( 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],
real64 const  inputNewmarkGamma,
real64 const  inputNewmarkBeta,
real64 const  inputMassDamping,
real64 const  inputStiffnessDamping 
)

Constructor.

Parameters
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 35 of file ImplicitSmallStrainNewmark_impl.hpp.

Member Function Documentation

◆ complete()

template<typename SUBREGION_TYPE , typename CONSTITUTIVE_TYPE , typename FE_TYPE >
GEOS_HOST_DEVICE real64 geos::solidMechanicsLagrangianFEMKernels::ImplicitSmallStrainNewmark< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >::complete ( localIndex const  k,
StackVariables stack 
) const
inline

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

Definition at line 143 of file ImplicitSmallStrainNewmark_impl.hpp.

◆ kernelLaunch()

template<typename SUBREGION_TYPE , typename CONSTITUTIVE_TYPE , typename FE_TYPE >
template<typename POLICY , typename KERNEL_TYPE >
real64 geos::solidMechanicsLagrangianFEMKernels::ImplicitSmallStrainNewmark< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >::kernelLaunch ( localIndex const  numElems,
KERNEL_TYPE const &  kernelComponent 
)
static

Kernel Launcher.

Template Parameters
POLICYThe RAJA policy to use for the launch.
NUM_QUADRATURE_POINTSThe number of quadrature points per element.
KERNEL_TYPEThe type of Kernel to execute.
Parameters
numElemsThe number of elements to process in this launch.
kernelComponentThe instantiation of KERNEL_TYPE to execute.
Returns
The maximum residual contribution.

This is a generic launching function for all of the finite element kernels that follow the interface set by KernelBase.

Definition at line 187 of file ImplicitSmallStrainNewmark_impl.hpp.

◆ quadraturePointKernel()

template<typename SUBREGION_TYPE , typename CONSTITUTIVE_TYPE , typename FE_TYPE >
GEOS_HOST_DEVICE void geos::solidMechanicsLagrangianFEMKernels::ImplicitSmallStrainNewmark< 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 103 of file ImplicitSmallStrainNewmark_impl.hpp.

◆ setup()

template<typename SUBREGION_TYPE , typename CONSTITUTIVE_TYPE , typename FE_TYPE >
GEOS_HOST_DEVICE void geos::solidMechanicsLagrangianFEMKernels::ImplicitSmallStrainNewmark< 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.

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 82 of file ImplicitSmallStrainNewmark_impl.hpp.

Member Data Documentation

◆ numNodesPerElem

template<typename SUBREGION_TYPE , typename CONSTITUTIVE_TYPE , typename FE_TYPE >
constexpr int geos::solidMechanicsLagrangianFEMKernels::ImplicitSmallStrainQuasiStatic< SUBREGION_TYPE, CONSTITUTIVE_TYPE, FE_TYPE >::numNodesPerElem
staticconstexpr

Maximum number of nodes per element, which is equal to the maxNumTestSupportPointPerElem and maxNumTrialSupportPointPerElem by definition. When the FE_TYPE is not a Virtual Element, this will be the actual number of nodes per element.

Definition at line 76 of file ImplicitSmallStrainQuasiStatic.hpp.


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