#include <MechanicsLinearElement.h>

Public Member Functions | |
| void | apply () override |
| void | computeStressAndStrain (const std::vector< double > &dispVec, std::vector< double > &stressVec, std::vector< double > &strainVec) |
| unsigned int | getNumberOfGaussPoints () |
| void | initializeForCurrentElement (const libMesh::Elem *elem, std::shared_ptr< MechanicsMaterialModel > materialModel) |
| MechanicsLinearElement (std::shared_ptr< const ElementOperationParameters > params) | |
| Constructor. | |
| void | printStressAndStrain (FILE *fp, const std::vector< double > &dispVec) |
| void | setElementStiffnessMatrix (std::vector< std::vector< double > > &elementStiffnessMatrix) |
| virtual | ~MechanicsLinearElement () |
| Destructor. | |
Protected Member Functions | |
| void | apply_Normal () |
| void | apply_Reduced () |
Protected Attributes | |
| const std::vector< std::vector< libMesh::RealGradient > > * | d_dphi |
| const libMesh::Elem * | d_elem |
| std::vector< std::vector< double > > * | d_elementStiffnessMatrix |
| std::shared_ptr< libMesh::FEBase > | d_fe |
| std::shared_ptr< libMesh::FEType > | d_feType |
| int | d_iDebugPrintInfoLevel |
| const std::vector< libMesh::Real > * | d_JxW |
| std::shared_ptr< MechanicsMaterialModel > | d_materialModel |
| std::shared_ptr< libMesh::QBase > | d_qrule |
| bool | d_useFlanaganTaylorElem |
| bool | d_useJaumannRate |
| bool | d_useReducedIntegration |
| const std::vector< libMesh::Point > * | d_xyz |
A class for representing the element level computation performed within a linear finite element operator for modelling solid mechanics. The linear operator could either be a linear elasticity operator or it could be the jacobian of a nonlinear elasticity/elasto-plasticity operator.
Definition at line 18 of file MechanicsLinearElement.h.
|
inlineexplicit |
Constructor.
Definition at line 22 of file MechanicsLinearElement.h.
References AMP_INSIST, d_dphi, AMP::Operator::MechanicsElement::d_fe, d_JxW, AMP::Operator::MechanicsElement::d_useJaumannRate, and d_xyz.
|
inlinevirtual |
Destructor.
Definition at line 33 of file MechanicsLinearElement.h.
|
inlineoverridevirtual |
Element stiffness matrix computation.
Implements AMP::Operator::ElementOperation.
Definition at line 48 of file MechanicsLinearElement.h.
References apply_Normal(), apply_Reduced(), and AMP::Operator::MechanicsElement::d_useReducedIntegration.
|
protected |
Element stiffness matrix computation using normal integration scheme.
Referenced by apply().
|
protected |
Element stiffness matrix computation using reduced integration scheme.
Referenced by apply().
| void AMP::Operator::MechanicsLinearElement::computeStressAndStrain | ( | const std::vector< double > & | dispVec, |
| std::vector< double > & | stressVec, | ||
| std::vector< double > & | strainVec | ||
| ) |
Computes the stress and strain values at the Gauss points within the current element The 6 components of stress and strain at each Gauss point are arranged in the order: xx, yy, zz, yz, xz and xy.
| [in] | dispVec | Displacements at the nodes of the current element. |
| [out] | stressVec | Stresses at the Gauss points of the current element. |
| [out] | strainVec | Strains at the Gauss points of the current element. |
|
inlineinherited |
Definition at line 72 of file MechanicsElement.h.
References AMP::Operator::MechanicsElement::d_qrule.
|
inlineinherited |
This function is used by the Linear and Nonlinear mechanics FEOperators to pass the current element and material model to this class during the finite element assembly operation.
| [in] | elem | Pointer to the current element within a finite element assembly. |
| [in] | materialModel | Shared pointer to the mechanics material model used in the current element. |
Definition at line 62 of file MechanicsElement.h.
References AMP::Operator::MechanicsElement::d_elem, and AMP::Operator::MechanicsElement::d_materialModel.
| void AMP::Operator::MechanicsLinearElement::printStressAndStrain | ( | FILE * | fp, |
| const std::vector< double > & | dispVec | ||
| ) |
Writes the stess and strain values at the Gauss points within the current element to the file. The 6 components of stress and strain at each Gauss point are arranged in the order: xx, yy, zz, yz, xz and xy.
| [in] | fp | File pointer |
| [in] | dispVec | Displacements at the nodes of the current element. |
|
inline |
This function is used by MechanicsLinearFEOperator to pass the address of the element stiffness matrix to this class.
| [in] | elementStiffnessMatrix | Element stiffness matrix |
Definition at line 40 of file MechanicsLinearElement.h.
References d_elementStiffnessMatrix.
|
protected |
Spatial Derivatives of the shape functions at the Gauss points in the current element.
Definition at line 93 of file MechanicsLinearElement.h.
Referenced by MechanicsLinearElement().
|
protectedinherited |
Pointer to the current element within the finite element assembly.
Definition at line 90 of file MechanicsElement.h.
Referenced by AMP::Operator::MechanicsElement::initializeForCurrentElement().
|
protected |
Element stiffness matrix.
Definition at line 99 of file MechanicsLinearElement.h.
Referenced by setElementStiffnessMatrix().
|
protectedinherited |
Finite element shape functions.
Definition at line 84 of file MechanicsElement.h.
Referenced by MechanicsLinearElement(), AMP::Operator::MechanicsLinearUpdatedLagrangianElement::MechanicsLinearUpdatedLagrangianElement(), AMP::Operator::MechanicsNonlinearElement::MechanicsNonlinearElement(), and AMP::Operator::MechanicsNonlinearUpdatedLagrangianElement::MechanicsNonlinearUpdatedLagrangianElement().
|
protectedinherited |
Type of polynomial used for the finite element shape functions. This includes both the polynomial order: First order/Second order etc. and polynomial family: Lagrange/Hierarchic/Hermite etc.
Definition at line 78 of file MechanicsElement.h.
|
protectedinherited |
< Inside Green-Naghdi stress-rate whether to use Flanagan Taylor stress-srate or not.
Definition at line 105 of file MechanicsElement.h.
|
protected |
Product of the determinant of Jacobian and the quadrature weight at the Gauss points in the current element.
Definition at line 89 of file MechanicsLinearElement.h.
Referenced by MechanicsLinearElement().
|
protectedinherited |
Shared pointer to the mechanics material model used in the current element.
Definition at line 93 of file MechanicsElement.h.
Referenced by AMP::Operator::MechanicsElement::initializeForCurrentElement().
|
protectedinherited |
Quadtrature rule used for numerical integration.
Definition at line 87 of file MechanicsElement.h.
Referenced by AMP::Operator::MechanicsElement::getNumberOfGaussPoints().
|
protectedinherited |
Definition at line 101 of file MechanicsElement.h.
Referenced by AMP::Operator::MechanicsNonlinearUpdatedLagrangianElement::resetElementInfo().
|
protectedinherited |
A flag that checks whether to use Jaumann Rate in Updated Lagrangian formulation or not.
Definition at line 97 of file MechanicsElement.h.
Referenced by MechanicsLinearElement(), AMP::Operator::MechanicsNonlinearElement::MechanicsNonlinearElement(), and AMP::Operator::MechanicsNonlinearUpdatedLagrangianElement::resetElementInfo().
|
protectedinherited |
A flag that is true if reduced integration scheme is used and false otherwise.
Definition at line 75 of file MechanicsElement.h.
Referenced by apply(), AMP::Operator::MechanicsLinearUpdatedLagrangianElement::apply(), AMP::Operator::MechanicsNonlinearElement::apply(), and AMP::Operator::MechanicsNonlinearUpdatedLagrangianElement::apply().
|
protected |
Locations of the Gauss points in the current element.
Definition at line 97 of file MechanicsLinearElement.h.
Referenced by MechanicsLinearElement().
|
Advanced Multi-Physics (AMP) Oak Ridge National Laboratory Idaho National Laboratory Los Alamos National Laboratory |
This page automatically produced from the source code by Last updated: Tue Mar 10 2026 13:06:44. Comments on this page |