Advanced Multi-Physics (AMP)
On-Line Documentation
MechanicsNonlinearUpdatedLagrangianElement.h
Go to the documentation of this file.
1#ifndef included_AMP_MechanicsNonlinearUpdatedLagrangianElement
2#define included_AMP_MechanicsNonlinearUpdatedLagrangianElement
3
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"
8
9#include <memory>
10#include <vector>
11
12
13namespace AMP::Operator {
14
15
21{
22public:
25 std::shared_ptr<const ElementOperationParameters> params )
26 : MechanicsElement( params ), d_elementOutputVector( nullptr )
27 {
28 d_JxW = &( d_fe->get_JxW() );
29
30 d_dphi = &( d_fe->get_dphi() );
31
32 d_phi = &( d_fe->get_phi() );
33
34 d_xyz = &( d_fe->get_xyz() );
35
37
38 d_gaussPtCnt = 0;
39 }
40
43
51 void setElementVectors( const std::vector<std::vector<double>> &elementInputVectors,
52 const std::vector<std::vector<double>> &elementInputVectors_pre,
53 std::vector<double> &elementOutputVector )
54 {
55 d_elementInputVectors = elementInputVectors;
56 d_elementInputVectors_pre = elementInputVectors_pre;
57 d_elementOutputVector = &( elementOutputVector );
58 }
59
63 void apply() override
64 {
67 } else {
69 }
70 }
71
84 void initMaterialModel( const std::vector<double> &initTempVector );
85
103 const std::vector<std::vector<double>> &elementInputVectors,
104 const std::vector<std::vector<double>> &elementInputVectors_pre );
105
114 void printStressAndStrain( FILE *fp, const std::vector<std::vector<double>> &inputVec );
115
125 void computeStressAndStrain( const std::vector<std::vector<double>> &inputVec,
126 std::vector<double> &stressVec,
127 std::vector<double> &strainVec );
128
132 void computeDeformationGradient( const std::vector<std::vector<libMesh::RealGradient>> &dphi,
133 const std::vector<libMesh::Point> &xyz,
134 unsigned int num_nodes,
135 unsigned int qp,
136 double F[3][3] );
137
141 void initializeReferenceXYZ( std::vector<double> &elementRefXYZ );
142
146 void assignReferenceXYZ( const std::vector<double> &elementRefXYZ )
147 {
148 d_elementRefXYZ = elementRefXYZ;
149 }
150
152
154
156 {
157 if ( ( d_useJaumannRate == false ) && ( d_useFlanaganTaylorElem == true ) ) {
160 }
161 }
162
163protected:
165
168
170 const std::vector<std::vector<double>> &,
171 double[3][3],
172 double[3][3] );
173
176
177 void
179
184
189
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>
204 std::vector<std::vector<double>> d_elementInputVectors;
208 std::vector<std::vector<double>>
212 std::vector<double> *d_elementOutputVector;
215
216 std::vector<double> d_elementRefXYZ;
217
218 std::vector<double> d_leftStretchV_n;
219
220 std::vector<double> d_leftStretchV_np1;
221
222 std::vector<double> d_rotationR_n;
223
224 std::vector<double> d_rotationR_np1;
225
227
228private:
229};
230
231
232} // namespace AMP::Operator
233
234#endif
std::shared_ptr< libMesh::FEBase > d_fe
void materialModelNonlinearGaussPointOperation(MechanicsNonlinearElement::MaterialUpdateType, const std::vector< std::vector< double > > &, double[3][3], double[3][3])
void materialModelPreNonlinearElementOperation(MechanicsNonlinearElement::MaterialUpdateType)
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)
void materialModelPostNonlinearGaussPointOperation(MechanicsNonlinearElement::MaterialUpdateType)
void printStressAndStrain(FILE *fp, const std::vector< std::vector< double > > &inputVec)
void materialModelPreNonlinearGaussPointOperation(MechanicsNonlinearElement::MaterialUpdateType)
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)
void materialModelPostNonlinearElementOperation(MechanicsNonlinearElement::MaterialUpdateType)
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])
MechanicsNonlinearUpdatedLagrangianElement(std::shared_ptr< const ElementOperationParameters > params)
Constructor.
void computeStressAndStrain(const std::vector< std::vector< double > > &inputVec, std::vector< double > &stressVec, std::vector< double > &strainVec)



Advanced Multi-Physics (AMP)
Oak Ridge National Laboratory
Idaho National Laboratory
Los Alamos National Laboratory
This page automatically produced from the
source code by doxygen
Last updated: Tue Mar 10 2026 13:06:41.
Comments on this page