Advanced Multi-Physics (AMP)
On-Line Documentation
MechanicsLinearElement.h
Go to the documentation of this file.
1#ifndef included_AMP_MechanicsLinearElement
2#define included_AMP_MechanicsLinearElement
3
4#include "AMP/operators/mechanics/MechanicsElement.h"
5
6#include <memory>
7#include <vector>
8
9
10namespace AMP::Operator {
11
19{
20public:
22 explicit MechanicsLinearElement( std::shared_ptr<const ElementOperationParameters> params )
23 : MechanicsElement( params ), d_elementStiffnessMatrix( nullptr )
24 {
25 d_JxW = &( d_fe->get_JxW() );
26 d_dphi = &( d_fe->get_dphi() );
27 d_xyz = &( d_fe->get_xyz() );
28 AMP_INSIST( ( d_useJaumannRate == false ),
29 "Jaumann rate with small strain does not make any sense." );
30 }
31
34
40 void setElementStiffnessMatrix( std::vector<std::vector<double>> &elementStiffnessMatrix )
41 {
42 d_elementStiffnessMatrix = &( elementStiffnessMatrix );
43 }
44
48 void apply() override
49 {
52 } else {
54 }
55 }
56
64 void printStressAndStrain( FILE *fp, const std::vector<double> &dispVec );
65
74 void computeStressAndStrain( const std::vector<double> &dispVec,
75 std::vector<double> &stressVec,
76 std::vector<double> &strainVec );
77
78protected:
83
88
89 const std::vector<libMesh::Real> *d_JxW;
92 const std::vector<std::vector<libMesh::RealGradient>>
96 const std::vector<libMesh::Point>
99 std::vector<std::vector<double>> *d_elementStiffnessMatrix;
101private:
102};
103} // namespace AMP::Operator
104
105#endif
std::shared_ptr< libMesh::FEBase > d_fe
void setElementStiffnessMatrix(std::vector< std::vector< double > > &elementStiffnessMatrix)
const std::vector< libMesh::Real > * d_JxW
const std::vector< std::vector< libMesh::RealGradient > > * d_dphi
std::vector< std::vector< double > > * d_elementStiffnessMatrix
MechanicsLinearElement(std::shared_ptr< const ElementOperationParameters > params)
Constructor.
void computeStressAndStrain(const std::vector< double > &dispVec, std::vector< double > &stressVec, std::vector< double > &strainVec)
const std::vector< libMesh::Point > * d_xyz
void printStressAndStrain(FILE *fp, const std::vector< double > &dispVec)
#define AMP_INSIST(EXP, MSG)
Insist error.



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