Implements kernels for solving the equations of motion using an implicit Newmark's method..
More...
|
| 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 |
|
| 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) |
|
| 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...
|
|
| 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...
|
|
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
-
nodeManager | Reference to the NodeManager object. |
edgeManager | Reference to the EdgeManager object. |
faceManager | Reference to the FaceManager object. |
targetRegionIndex | Index of the region the subregion belongs to. |
inputDofNumber | The dof number for the primary field. |
rankOffset | dof index offset of current rank |
inputMatrix | Reference to the Jacobian matrix. |
inputRhs | Reference to the RHS vector. |
inputDt | The timestep for the physics update. |
Constructor.
- Parameters
-
elementSubRegion | Reference to the SUBREGION_TYPE(class template parameter) object. |
inputConstitutiveType | The constitutive object. |
finiteElementSpace | Placeholder for the finite element space object, which currently doesn't do much. |
inputGravityVector | The 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 47 of file ImplicitSmallStrainNewmark.hpp.
template<typename SUBREGION_TYPE , typename CONSTITUTIVE_TYPE , typename FE_TYPE >
template<typename POLICY , typename KERNEL_TYPE >
Kernel Launcher.
- Template Parameters
-
POLICY | The RAJA policy to use for the launch. |
NUM_QUADRATURE_POINTS | The number of quadrature points per element. |
KERNEL_TYPE | The type of Kernel to execute. |
- Parameters
-
numElems | The number of elements to process in this launch. |
kernelComponent | The 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 188 of file ImplicitSmallStrainNewmark_impl.hpp.
template<typename SUBREGION_TYPE , typename CONSTITUTIVE_TYPE , typename FE_TYPE >
Performs a state update at a quadrature point.
- Template Parameters
-
STACK_VARIABLE_TYPE | The type of StackVariable that holds the stack variables. This is most likely a defined in a type that derives from KernelBase. |
- Parameters
-
k | The element index. |
q | The quadrature point index. |
stack | The 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.