Advanced Multi-Physics (AMP)
On-Line Documentation
MechanicsNonlinearElement.h
Go to the documentation of this file.
1#ifndef included_AMP_MechanicsNonlinearElement
2#define included_AMP_MechanicsNonlinearElement
3
4#include "AMP/operators/mechanics/MechanicsConstants.h"
5#include "AMP/operators/mechanics/MechanicsElement.h"
6
7#include <memory>
8#include <vector>
9
10
11namespace AMP::Operator {
12
13
19{
20public:
25
27 explicit MechanicsNonlinearElement( std::shared_ptr<const ElementOperationParameters> params )
28 : MechanicsElement( params ), d_elementOutputVector( nullptr )
29 {
30 d_JxW = &( d_fe->get_JxW() );
31
32 d_dphi = &( d_fe->get_dphi() );
33
34 d_phi = &( d_fe->get_phi() );
35
36 d_xyz = &( d_fe->get_xyz() );
37
38 AMP_INSIST( ( d_useJaumannRate == false ),
39 "Jaumann rate with small strain does not make any sense." );
40 }
41
44
51 void setElementVectors( const std::vector<std::vector<double>> &elementInputVectors,
52 std::vector<double> &elementOutputVector )
53 {
54 d_elementInputVectors = elementInputVectors;
55 d_elementOutputVector = &( elementOutputVector );
56 }
57
61 void apply() override
62 {
65 } else {
67 }
68 }
69
82 void initMaterialModel( const std::vector<double> &initTempVector );
83
96 const std::vector<std::vector<double>> &elementInputVectors );
97
106 void printStressAndStrain( FILE *fp, const std::vector<std::vector<double>> &inputVec );
107
117 void computeStressAndStrain( const std::vector<std::vector<double>> &inputVec,
118 std::vector<double> &stressVec,
119 std::vector<double> &strainVec );
120
121protected:
123
125
127 const std::vector<std::vector<double>> & );
128
130
132
137
142
143 const std::vector<libMesh::Real> *d_JxW;
146 const std::vector<std::vector<libMesh::RealGradient>>
150 const std::vector<std::vector<libMesh::Real>>
154 const std::vector<libMesh::Point>
157 std::vector<std::vector<double>> d_elementInputVectors;
161 std::vector<double> *d_elementOutputVector;
163private:
164};
165
166
167} // namespace AMP::Operator
168
169#endif
std::shared_ptr< libMesh::FEBase > d_fe
std::vector< std::vector< double > > d_elementInputVectors
void computeStressAndStrain(const std::vector< std::vector< double > > &inputVec, std::vector< double > &stressVec, std::vector< double > &strainVec)
void materialModelPreNonlinearElementOperation(MaterialUpdateType)
void printStressAndStrain(FILE *fp, const std::vector< std::vector< double > > &inputVec)
void initMaterialModel(const std::vector< double > &initTempVector)
void updateMaterialModel(MaterialUpdateType type, const std::vector< std::vector< double > > &elementInputVectors)
const std::vector< std::vector< libMesh::Real > > * d_phi
const std::vector< std::vector< libMesh::RealGradient > > * d_dphi
void materialModelPostNonlinearElementOperation(MaterialUpdateType)
const std::vector< libMesh::Point > * d_xyz
void materialModelNonlinearGaussPointOperation(MaterialUpdateType, const std::vector< std::vector< double > > &)
void setElementVectors(const std::vector< std::vector< double > > &elementInputVectors, std::vector< double > &elementOutputVector)
void materialModelPostNonlinearGaussPointOperation(MaterialUpdateType)
MechanicsNonlinearElement(std::shared_ptr< const ElementOperationParameters > params)
Constructor.
const std::vector< libMesh::Real > * d_JxW
void materialModelPreNonlinearGaussPointOperation(MaterialUpdateType)
#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