Advanced Multi-Physics (AMP)
On-Line Documentation
VonMises_IsotropicKinematicHardening.h
Go to the documentation of this file.
1#ifndef included_AMP_VonMises_IsotropicKinematicHardening
2#define included_AMP_VonMises_IsotropicKinematicHardening
3
4#include "AMP/operators/mechanics/MechanicsMaterialModel.h"
5
6#include <memory>
7#include <vector>
8
9
10namespace AMP::Operator {
11
13{
14public:
16 std::shared_ptr<MechanicsMaterialModelParameters> );
17
19
20 void getConstitutiveMatrix( double *& ) override;
21
22 void getInternalStress( const std::vector<std::vector<double>> &, double *& ) override;
23
24 void preLinearAssembly() override { d_gaussPtCnt = 0; }
25
27
28 void preNonlinearInit( bool, bool ) override;
29
30 void nonlinearInitGaussPointOperation( double ) override;
31
32 void preNonlinearAssembly() override
33 {
35
36 d_gaussPtCnt = 0;
37 }
38
40
41 void postNonlinearAssembly() override;
42
43 void preNonlinearReset() override { d_gaussPtCnt = 0; }
44
46
47 void nonlinearResetGaussPointOperation( const std::vector<std::vector<double>> & ) override;
48
49 void globalReset() override;
50
51 void postNonlinearReset() override;
52
53 void preNonlinearJacobian() override { d_gaussPtCnt = 0; }
54
56
57 void nonlinearJacobianGaussPointOperation( const std::vector<std::vector<double>> & ) override;
58
59protected:
60 void radialReturn( const double *stra_np1,
61 double *stre_np1,
62 double *back_stress_np1,
63 double *ystre_np1,
64 double *eph_bar_plas_np1,
65 const std::vector<std::vector<double>> &strain );
66
67 void computeEvalv( const std::vector<std::vector<double>> &strain );
68
70
72
74
75 std::vector<double> d_E;
76
77 std::vector<double> d_Nu;
78
79 double default_E;
80
81 double default_Nu;
82
83 double d_delta;
84
85 double d_K_0;
86
87 double d_K_inf;
88
89 double d_H;
90
91 double d_beta;
92
93 double d_Ep;
94
95 double d_Kin;
96
97 double d_Sig_0;
98
99 unsigned int d_gaussPtCnt;
100
101 unsigned int
104 unsigned int Plastic_Gauss_Point;
108
109 std::vector<double> d_EquilibriumStress;
110
111 std::vector<double> d_EquilibriumBackStress;
112
113 std::vector<double> d_EquilibriumStrain;
114
115 std::vector<double> d_EquilibriumYieldStress;
116
118
119 std::vector<double> d_Lambda;
120
121 std::vector<int> d_ElPl;
122
123 std::vector<double> d_tmp1Stress;
124
125 std::vector<double> d_tmp1BackStress;
126
127 std::vector<double> d_tmp1Strain;
128
129 std::vector<double> d_tmp1YieldStress;
130
131 std::vector<double> d_tmp1EffectivePlasticStrain;
132
133 std::vector<double> d_tmp2Stress;
134
135 std::vector<double> d_tmp2BackStress;
136
137 std::vector<double> d_tmp2YieldStress;
138
139 std::vector<double> d_tmp2EffectivePlasticStrain;
140
142
144
146
148
149private:
150};
151} // namespace AMP::Operator
152
153#endif
void getInternalStress(const std::vector< std::vector< double > > &, double *&) override
VonMises_IsotropicKinematicHardening(std::shared_ptr< MechanicsMaterialModelParameters >)
void radialReturn(const double *stra_np1, double *stre_np1, double *back_stress_np1, double *ystre_np1, double *eph_bar_plas_np1, const std::vector< std::vector< double > > &strain)
void nonlinearResetGaussPointOperation(const std::vector< std::vector< double > > &) override
void nonlinearJacobianGaussPointOperation(const std::vector< std::vector< double > > &) override
void computeEvalv(const std::vector< std::vector< double > > &strain)



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