Advanced Multi-Physics (AMP)
On-Line Documentation
MechanicsLinearUpdatedLagrangianElement.h
Go to the documentation of this file.
1#ifndef included_AMP_MechanicsLinearUpdatedLagrangianElement
2#define included_AMP_MechanicsLinearUpdatedLagrangianElement
3
4#include "AMP/operators/mechanics/MechanicsConstants.h"
5#include "AMP/operators/mechanics/MechanicsElement.h"
6#include "AMP/operators/mechanics/UpdatedLagrangianUtils.h"
7
8#include <memory>
9#include <vector>
10
11
12namespace AMP::Operator {
13
14
22{
23public:
26 std::shared_ptr<const ElementOperationParameters> params )
27 : MechanicsElement( params ), d_elementStiffnessMatrix( nullptr )
28 {
29 d_JxW = &( d_fe->get_JxW() );
30 d_dphi = &( d_fe->get_dphi() );
31 d_xyz = &( d_fe->get_xyz() );
33 }
34
37
43 void setElementStiffnessMatrix( std::vector<std::vector<double>> &elementStiffnessMatrix )
44 {
45 d_elementStiffnessMatrix = &( elementStiffnessMatrix );
46 }
47
51 void apply() override
52 {
55 } else {
57 }
58 }
59
67 void printStressAndStrain( FILE *fp, const std::vector<double> &dispVec );
68
77 void computeStressAndStrain( const std::vector<double> &dispVec,
78 std::vector<double> &stressVec,
79 std::vector<double> &strainVec );
80
86 void setElementVectors( const std::vector<std::vector<double>> &elementInputVectors )
87 {
88 d_elementInputVectors = elementInputVectors;
89 }
90
100 void initializeReferenceXYZ( std::vector<double> &elementRefXYZ );
101
106 void assignReferenceXYZ( const std::vector<double> &elementRefXYZ )
107 {
108 d_elementRefXYZ = elementRefXYZ;
109 }
110
111protected:
116
121
122 const std::vector<libMesh::Real> *d_JxW;
125 const std::vector<std::vector<libMesh::RealGradient>>
129 const std::vector<libMesh::Point>
132 std::vector<std::vector<double>> *d_elementStiffnessMatrix;
134 std::vector<std::vector<double>> d_elementInputVectors;
139
140 std::vector<double> d_elementRefXYZ;
141
142private:
143};
144
145} // namespace AMP::Operator
146
147#endif
std::shared_ptr< libMesh::FEBase > d_fe
void printStressAndStrain(FILE *fp, const std::vector< double > &dispVec)
const std::vector< std::vector< libMesh::RealGradient > > * d_dphi
void setElementStiffnessMatrix(std::vector< std::vector< double > > &elementStiffnessMatrix)
void initializeReferenceXYZ(std::vector< double > &elementRefXYZ)
void computeStressAndStrain(const std::vector< double > &dispVec, std::vector< double > &stressVec, std::vector< double > &strainVec)
void setElementVectors(const std::vector< std::vector< double > > &elementInputVectors)
MechanicsLinearUpdatedLagrangianElement(std::shared_ptr< const ElementOperationParameters > params)
Constructor.



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