Advanced Multi-Physics (AMP)
On-Line Documentation
MechanicsMaterialModel.h
Go to the documentation of this file.
1#ifndef included_AMP_MechanicsMaterialModel
2#define included_AMP_MechanicsMaterialModel
3
4#include <cstring>
5
6#include "AMP/materials/Material.h"
7#include "AMP/operators/ElementPhysicsModel.h"
8#include "AMP/operators/mechanics/UpdatedLagrangianUtils.h"
9#include <memory>
10
11namespace AMP::Operator {
12
14
25{
26public:
32 explicit MechanicsMaterialModel( std::shared_ptr<MechanicsMaterialModelParameters> params )
33 : ElementPhysicsModel( params )
34 {
36 params->d_db->getWithDefault<bool>( "USE_MATERIALS_LIBRARY", false );
37
39 params->d_db->getWithDefault<bool>( "USE_UPDATED_LAGRANGIAN", false );
40
41 d_useJaumannRate = params->d_db->getWithDefault<bool>( "USE_JAUMANN_RATE", false );
42
44 params->d_db->getWithDefault<bool>( "USE_CONTINUUM_TANGENT", false );
45
46 d_checkCladOrPellet = params->d_db->getWithDefault<bool>( "IS_IT_CLAD", false );
47
48 if ( d_useMaterialsLibrary == true ) {
49 AMP_INSIST( ( params->d_db->keyExists( "Material" ) ), "Key ''Material'' is missing!" );
50 std::string matname = params->d_db->getString( "Material" );
52 }
53
54 d_currentTime = 0;
56 }
57
62
67 virtual void getConstitutiveMatrix( double *& ) {}
68
73 virtual void getConstitutiveMatrixUpdatedLagrangian( double[6][6], double[3][3] ) {}
74
79 virtual void getStressForUpdatedLagrangian( double[6] ) {}
80
86 virtual void getInternalStress( const std::vector<std::vector<double>> &, double *& ) {}
87
89 const std::vector<std::vector<double>> &, double *&, double[3][3], double[3][3], double )
90 {
91 }
92
96 virtual void getEffectiveStress( double *& ) {}
97
102 virtual void getEquivalentStrain( double *& ) {}
103
109 virtual void getExternalStress( double *& ) {}
110
111 /*
112 A bunch of virtual functions have been implemented to allow the
113 programmer to initialize, reset or update any of the variables within
114 the material model from outside the gauss point calculation.
115 Most of these functions are not implemented. The are added whenever it
116 is required.
117 */
118
119 // For LinearOperator's Reset/Init
120
121 virtual void preLinearAssembly() {}
122
123 virtual void postLinearAssembly() {}
124
126
128
130
132
133 // For NonlinearOperator's Init
134
135 virtual void preNonlinearInit( bool, bool ) {}
136
137 virtual void postNonlinearInit() {}
138
140
142
144
146
152 virtual void nonlinearInitGaussPointOperation( double ) {}
153
154 // For NonlinearOperator's Apply
155
156 virtual void preNonlinearAssembly() {}
157
158 virtual void postNonlinearAssembly() {}
159
161
163
165
167
168 // For NonlinearOperator's Reset
169
170 virtual void globalReset() {}
171
172 virtual void preNonlinearReset() {}
173
174 virtual void postNonlinearReset() {}
175
177
179
181
183
190 virtual void nonlinearResetGaussPointOperation( const std::vector<std::vector<double>> & ) {}
191
192 virtual void nonlinearResetGaussPointOperation_UL( const std::vector<std::vector<double>> &,
193 double[3][3],
194 double[3][3] )
195 {
196 }
197
198 // For NonlinearOperator's GetJacobianParameters
199
200 virtual void preNonlinearJacobian() {}
201
202 virtual void postNonlinearJacobian() {}
203
205
207
209
211
212 virtual void nonlinearJacobianGaussPointOperation( const std::vector<std::vector<double>> & ) {}
213
214 virtual void nonlinearJacobianGaussPointOperation_UL( const std::vector<std::vector<double>> &,
215 double[3][3],
216 double[3][3] )
217 {
218 }
219
220 void updateTime( double currTime )
221 {
223 d_currentTime = currTime;
224 }
225
226 std::shared_ptr<AMP::Materials::Material> getMaterial() { return d_material; }
227
228protected:
245
246 std::shared_ptr<AMP::Materials::Material>
249private:
250};
251} // namespace AMP::Operator
252
253#endif
virtual void nonlinearResetGaussPointOperation(const std::vector< std::vector< double > > &)
MechanicsMaterialModel(std::shared_ptr< MechanicsMaterialModelParameters > params)
virtual void getConstitutiveMatrixUpdatedLagrangian(double[6][6], double[3][3])
virtual void nonlinearJacobianGaussPointOperation_UL(const std::vector< std::vector< double > > &, double[3][3], double[3][3])
std::shared_ptr< AMP::Materials::Material > d_material
virtual void getStressForUpdatedLagrangian(double[6])
virtual void nonlinearResetGaussPointOperation_UL(const std::vector< std::vector< double > > &, double[3][3], double[3][3])
virtual void getInternalStress(const std::vector< std::vector< double > > &, double *&)
virtual void nonlinearJacobianGaussPointOperation(const std::vector< std::vector< double > > &)
std::shared_ptr< AMP::Materials::Material > getMaterial()
virtual void getInternalStress_UL(const std::vector< std::vector< double > > &, double *&, double[3][3], double[3][3], double)
#define AMP_INSIST(EXP, MSG)
Insist error.
std::unique_ptr< Material > getMaterial(const std::string &name)
Get a material.
ElementPhysicsModelParameters MechanicsMaterialModelParameters



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