1#ifndef included_AMP_MechanicsNonlinearUpdatedLagrangianElement
2#define included_AMP_MechanicsNonlinearUpdatedLagrangianElement
4#include "AMP/operators/mechanics/MechanicsConstants.h"
5#include "AMP/operators/mechanics/MechanicsElement.h"
6#include "AMP/operators/mechanics/MechanicsNonlinearElement.h"
7#include "AMP/operators/mechanics/UpdatedLagrangianUtils.h"
25 std::shared_ptr<const ElementOperationParameters> params )
52 const std::vector<std::vector<double>> &elementInputVectors_pre,
53 std::vector<double> &elementOutputVector )
103 const std::vector<std::vector<double>> &elementInputVectors,
104 const std::vector<std::vector<double>> &elementInputVectors_pre );
126 std::vector<double> &stressVec,
127 std::vector<double> &strainVec );
133 const std::vector<libMesh::Point> &xyz,
134 unsigned int num_nodes,
170 const std::vector<std::vector<double>> &,
190 const std::vector<libMesh::Real> *
d_JxW;
193 const std::vector<std::vector<libMesh::RealGradient>>
197 const std::vector<std::vector<libMesh::Real>>
201 const std::vector<libMesh::Point>
208 std::vector<std::vector<double>>
bool d_useReducedIntegration
std::shared_ptr< libMesh::FEBase > d_fe
bool d_useFlanaganTaylorElem
std::vector< double > d_leftStretchV_np1
void zeroOutGaussPointCount()
void materialModelNonlinearGaussPointOperation(MechanicsNonlinearElement::MaterialUpdateType, const std::vector< std::vector< double > > &, double[3][3], double[3][3])
void preNonlinearElementInit()
void materialModelPreNonlinearElementOperation(MechanicsNonlinearElement::MaterialUpdateType)
const std::vector< libMesh::Real > * d_JxW
void initializeReferenceXYZ(std::vector< double > &elementRefXYZ)
void updateMaterialModel(MechanicsNonlinearElement::MaterialUpdateType type, const std::vector< std::vector< double > > &elementInputVectors, const std::vector< std::vector< double > > &elementInputVectors_pre)
std::vector< double > d_leftStretchV_n
void materialModelPostNonlinearGaussPointOperation(MechanicsNonlinearElement::MaterialUpdateType)
void printStressAndStrain(FILE *fp, const std::vector< std::vector< double > > &inputVec)
std::vector< std::vector< double > > d_elementInputVectors
std::vector< double > d_rotationR_np1
const std::vector< std::vector< libMesh::Real > > * d_phi
void materialModelPreNonlinearGaussPointOperation(MechanicsNonlinearElement::MaterialUpdateType)
const std::vector< libMesh::Point > * d_xyz
std::vector< double > d_elementRefXYZ
void setElementVectors(const std::vector< std::vector< double > > &elementInputVectors, const std::vector< std::vector< double > > &elementInputVectors_pre, std::vector< double > &elementOutputVector)
void initMaterialModel(const std::vector< double > &initTempVector)
bool d_onePointShearIntegration
void materialModelPostNonlinearElementOperation(MechanicsNonlinearElement::MaterialUpdateType)
virtual ~MechanicsNonlinearUpdatedLagrangianElement()
Destructor.
const std::vector< std::vector< libMesh::RealGradient > > * d_dphi
void computeDeformationGradient(const std::vector< std::vector< libMesh::RealGradient > > &dphi, const std::vector< libMesh::Point > &xyz, unsigned int num_nodes, unsigned int qp, double F[3][3])
std::vector< double > d_rotationR_n
std::vector< std::vector< double > > d_elementInputVectors_pre
void assignReferenceXYZ(const std::vector< double > &elementRefXYZ)
MechanicsNonlinearUpdatedLagrangianElement(std::shared_ptr< const ElementOperationParameters > params)
Constructor.
std::vector< double > * d_elementOutputVector
void computeStressAndStrain(const std::vector< std::vector< double > > &inputVec, std::vector< double > &stressVec, std::vector< double > &strainVec)