www.mooseframework.org
TwoParameterPlasticityStressUpdate.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 
12 #include "Conversion.h" // for stringify
13 #include "libmesh/utility.h" // for Utility::pow
14 
16 
17 InputParameters
19 {
20  InputParameters params = MultiParameterPlasticityStressUpdate::validParams();
21  params.addClassDescription("Return-map and Jacobian algorithms for (P, Q) plastic models");
22  return params;
23 }
24 
26  const InputParameters & parameters, unsigned num_yf, unsigned num_intnl)
27  : MultiParameterPlasticityStressUpdate(parameters, _num_pq, num_yf, num_intnl),
28  _p_trial(0.0),
29  _q_trial(0.0),
30  _Epp(0.0),
31  _Eqq(0.0),
32  _dgaE_dpt(0.0),
33  _dgaE_dqt(0.0),
34  _dp_dpt(0.0),
35  _dq_dpt(0.0),
36  _dp_dqt(0.0),
37  _dq_dqt(0.0)
38 {
39 }
40 
41 void
42 TwoParameterPlasticityStressUpdate::yieldFunctionValuesV(const std::vector<Real> & stress_params,
43  const std::vector<Real> & intnl,
44  std::vector<Real> & yf) const
45 {
46  const Real p = stress_params[0];
47  const Real q = stress_params[1];
48  yieldFunctionValues(p, q, intnl, yf);
49 }
50 
51 void
53 {
54  setEppEqq(Eijkl, _Epp, _Eqq);
55  _Eij[0][0] = _Epp;
56  _Eij[1][0] = _Eij[0][1] = 0.0;
57  _Eij[1][1] = _Eqq;
58  _En = _Epp;
59  _Cij[0][0] = 1.0 / _Epp;
60  _Cij[1][0] = _Cij[0][1] = 0.0;
61  _Cij[1][1] = 1.0 / _Eqq;
62 }
63 
64 void
65 TwoParameterPlasticityStressUpdate::preReturnMapV(const std::vector<Real> & trial_stress_params,
66  const RankTwoTensor & stress_trial,
67  const std::vector<Real> & intnl_old,
68  const std::vector<Real> & yf,
69  const RankFourTensor & Eijkl)
70 {
71  const Real p_trial = trial_stress_params[0];
72  const Real q_trial = trial_stress_params[1];
73  preReturnMap(p_trial, q_trial, stress_trial, intnl_old, yf, Eijkl);
74 }
75 
76 void
78  Real /*q_trial*/,
79  const RankTwoTensor & /*stress_trial*/,
80  const std::vector<Real> & /*intnl_old*/,
81  const std::vector<Real> & /*yf*/,
82  const RankFourTensor & /*Eijkl*/)
83 {
84  return;
85 }
86 
87 void
89  Real q_trial,
90  const std::vector<Real> & intnl_old,
91  Real & p,
92  Real & q,
93  Real & gaE,
94  std::vector<Real> & intnl) const
95 {
96  p = p_trial;
97  q = q_trial;
98  gaE = 0.0;
99  std::copy(intnl_old.begin(), intnl_old.end(), intnl.begin());
100 }
101 
102 void
103 TwoParameterPlasticityStressUpdate::computeAllQV(const std::vector<Real> & stress_params,
104  const std::vector<Real> & intnl,
105  std::vector<yieldAndFlow> & all_q) const
106 {
107  const Real p = stress_params[0];
108  const Real q = stress_params[1];
109  computeAllQ(p, q, intnl, all_q);
110 }
111 
112 void
113 TwoParameterPlasticityStressUpdate::initializeVarsV(const std::vector<Real> & trial_stress_params,
114  const std::vector<Real> & intnl_old,
115  std::vector<Real> & stress_params,
116  Real & gaE,
117  std::vector<Real> & intnl) const
118 {
119  const Real p_trial = trial_stress_params[0];
120  const Real q_trial = trial_stress_params[1];
121  Real p;
122  Real q;
123  initializeVars(p_trial, q_trial, intnl_old, p, q, gaE, intnl);
124  stress_params[0] = p;
125  stress_params[1] = q;
126 }
127 
128 void
129 TwoParameterPlasticityStressUpdate::setIntnlValuesV(const std::vector<Real> & trial_stress_params,
130  const std::vector<Real> & current_stress_params,
131  const std::vector<Real> & intnl_old,
132  std::vector<Real> & intnl) const
133 {
134  const Real p_trial = trial_stress_params[0];
135  const Real q_trial = trial_stress_params[1];
136  const Real p = current_stress_params[0];
137  const Real q = current_stress_params[1];
138  setIntnlValues(p_trial, q_trial, p, q, intnl_old, intnl);
139 }
140 
141 void
143  const std::vector<Real> & trial_stress_params,
144  const std::vector<Real> & current_stress_params,
145  const std::vector<Real> & intnl_old,
146  std::vector<std::vector<Real>> & dintnl) const
147 {
148  const Real p_trial = trial_stress_params[0];
149  const Real q_trial = trial_stress_params[1];
150  const Real p = current_stress_params[0];
151  const Real q = current_stress_params[1];
152  setIntnlDerivatives(p_trial, q_trial, p, q, intnl_old, dintnl);
153 }
154 
155 void
157  std::vector<Real> & stress_params) const
158 {
159  Real p;
160  Real q;
161  computePQ(stress, p, q);
162  stress_params[0] = p;
163  stress_params[1] = q;
164 }
165 
166 void
168  const RankTwoTensor & stress_trial,
169  const std::vector<Real> & trial_stress_params,
170  const RankTwoTensor & stress,
171  const std::vector<Real> & stress_params,
172  Real gaE,
173  const yieldAndFlow & smoothed_q,
174  const RankFourTensor & Eijkl,
175  bool compute_full_tangent_operator,
176  const std::vector<std::vector<Real>> & dvar_dtrial,
177  RankFourTensor & cto)
178 {
179  const Real p_trial = trial_stress_params[0];
180  const Real q_trial = trial_stress_params[1];
181  const Real p = stress_params[0];
182  const Real q = stress_params[1];
183  _dp_dpt = dvar_dtrial[0][0];
184  _dp_dqt = dvar_dtrial[0][1];
185  _dq_dpt = dvar_dtrial[1][0];
186  _dq_dqt = dvar_dtrial[1][1];
187  _dgaE_dpt = dvar_dtrial[2][0];
188  _dgaE_dqt = dvar_dtrial[2][1];
189  consistentTangentOperator(stress_trial,
190  p_trial,
191  q_trial,
192  stress,
193  p,
194  q,
195  gaE,
196  smoothed_q,
197  Eijkl,
198  compute_full_tangent_operator,
199  cto);
200 }
201 
202 void
204  const RankTwoTensor & stress_trial,
205  Real /*p_trial*/,
206  Real /*q_trial*/,
207  const RankTwoTensor & stress,
208  Real /*p*/,
209  Real /*q*/,
210  Real gaE,
211  const yieldAndFlow & smoothed_q,
212  const RankFourTensor & elasticity_tensor,
213  bool compute_full_tangent_operator,
214  RankFourTensor & cto) const
215 {
216  cto = elasticity_tensor;
217  if (!compute_full_tangent_operator)
218  return;
219 
220  const RankTwoTensor dpdsig = dpdstress(stress);
221  const RankTwoTensor dpdsig_trial = dpdstress(stress_trial);
222  const RankTwoTensor dqdsig = dqdstress(stress);
223  const RankTwoTensor dqdsig_trial = dqdstress(stress_trial);
224 
225  const RankTwoTensor s1 = elasticity_tensor * ((1.0 / _Epp) * (1.0 - _dp_dpt) * dpdsig +
226  (1.0 / _Eqq) * (-_dq_dpt) * dqdsig);
227  const RankTwoTensor s2 = elasticity_tensor * ((1.0 / _Epp) * (-_dp_dqt) * dpdsig +
228  (1.0 / _Eqq) * (1.0 - _dq_dqt) * dqdsig);
229  const RankTwoTensor t1 = elasticity_tensor * dpdsig_trial;
230  const RankTwoTensor t2 = elasticity_tensor * dqdsig_trial;
231 
232  for (unsigned i = 0; i < _tensor_dimensionality; ++i)
233  for (unsigned j = 0; j < _tensor_dimensionality; ++j)
234  for (unsigned k = 0; k < _tensor_dimensionality; ++k)
235  for (unsigned l = 0; l < _tensor_dimensionality; ++l)
236  cto(i, j, k, l) -= s1(i, j) * t1(k, l) + s2(i, j) * t2(k, l);
237 
238  const RankFourTensor d2pdsig2 = d2pdstress2(stress);
239  const RankFourTensor d2qdsig2 = d2qdstress2(stress);
240 
241  const RankFourTensor Tijab = elasticity_tensor * (gaE / _Epp) *
242  (smoothed_q.dg[0] * d2pdsig2 + smoothed_q.dg[1] * d2qdsig2);
243 
244  RankFourTensor inv = RankFourTensor(RankFourTensor::initIdentityFour) + Tijab;
245  try
246  {
247  inv = inv.transposeMajor().invSymm();
248  }
249  catch (const MooseException & e)
250  {
251  // Cannot form the inverse, so probably at some degenerate place in stress space.
252  // Just return with the "best estimate" of the cto.
253  mooseWarning("TwoParameterPlasticityStressUpdate: Cannot invert 1+T in consistent tangent "
254  "operator computation at quadpoint ",
255  _qp,
256  " of element ",
257  _current_elem->id());
258  return;
259  }
260 
261  cto = (cto.transposeMajor() * inv).transposeMajor();
262 }
263 
264 void
266  const std::vector<Real> & stress_params,
267  Real gaE,
268  const std::vector<Real> & intnl,
269  const yieldAndFlow & smoothed_q,
270  const RankFourTensor & Eijkl,
271  RankTwoTensor & stress) const
272 {
273  const Real p_ok = stress_params[0];
274  const Real q_ok = stress_params[1];
275  setStressAfterReturn(stress_trial, p_ok, q_ok, gaE, intnl, smoothed_q, Eijkl, stress);
276 }
277 
278 void
280  const RankTwoTensor & /*stress_trial*/,
281  Real gaE,
282  const yieldAndFlow & smoothed_q,
283  const RankFourTensor & /*elasticity_tensor*/,
284  const RankTwoTensor & returned_stress,
285  RankTwoTensor & inelastic_strain_increment) const
286 {
287  inelastic_strain_increment = (gaE / _Epp) * (smoothed_q.dg[0] * dpdstress(returned_stress) +
288  smoothed_q.dg[1] * dqdstress(returned_stress));
289 }
290 
291 std::vector<RankTwoTensor>
293 {
294  std::vector<RankTwoTensor> dsp(_num_pq, RankTwoTensor());
295  dsp[0] = dpdstress(stress);
296  dsp[1] = dqdstress(stress);
297  return dsp;
298 }
299 
300 std::vector<RankFourTensor>
302 {
303  std::vector<RankFourTensor> d2(_num_pq, RankFourTensor());
304  d2[0] = d2pdstress2(stress);
305  d2[1] = d2qdstress2(stress);
306  return d2;
307 }
defineLegacyParams
defineLegacyParams(TwoParameterPlasticityStressUpdate)
TwoParameterPlasticityStressUpdate::d2qdstress2
virtual RankFourTensor d2qdstress2(const RankTwoTensor &stress) const =0
d2(q)/d(stress)/d(stress) Derived classes must override this
TwoParameterPlasticityStressUpdate::computeStressParams
virtual void computeStressParams(const RankTwoTensor &stress, std::vector< Real > &stress_params) const override
Computes stress_params, given stress.
Definition: TwoParameterPlasticityStressUpdate.C:156
MultiParameterPlasticityStressUpdate::yieldAndFlow
Struct designed to hold info about a single yield function and its derivatives, as well as the flow d...
Definition: MultiParameterPlasticityStressUpdate.h:214
MultiParameterPlasticityStressUpdate::yieldAndFlow::dg
std::vector< Real > dg
Definition: MultiParameterPlasticityStressUpdate.h:219
TwoParameterPlasticityStressUpdate::setIntnlValuesV
void setIntnlValuesV(const std::vector< Real > &trial_stress_params, const std::vector< Real > &current_stress_params, const std::vector< Real > &intnl_old, std::vector< Real > &intnl) const override
Sets the internal parameters based on the trial values of stress_params, their current values,...
Definition: TwoParameterPlasticityStressUpdate.C:129
MultiParameterPlasticityStressUpdate::_Cij
std::vector< std::vector< Real > > _Cij
_Cij[i, j] * _Eij[j, k] = 1 iff j == k
Definition: MultiParameterPlasticityStressUpdate.h:144
TwoParameterPlasticityStressUpdate::d2pdstress2
virtual RankFourTensor d2pdstress2(const RankTwoTensor &stress) const =0
d2(p)/d(stress)/d(stress) Derived classes must override this
TwoParameterPlasticityStressUpdate::yieldFunctionValues
virtual void yieldFunctionValues(Real p, Real q, const std::vector< Real > &intnl, std::vector< Real > &yf) const =0
Computes the values of the yield functions, given p, q and intnl parameters.
TwoParameterPlasticityStressUpdate::d2stress_param_dstress
virtual std::vector< RankFourTensor > d2stress_param_dstress(const RankTwoTensor &stress) const override
d2(stress_param[i])/d(stress)/d(stress) at given stress
Definition: TwoParameterPlasticityStressUpdate.C:301
TwoParameterPlasticityStressUpdate::setEppEqq
virtual void setEppEqq(const RankFourTensor &Eijkl, Real &Epp, Real &Eqq) const =0
Set Epp and Eqq based on the elasticity tensor Derived classes must override this.
TwoParameterPlasticityStressUpdate::setIntnlValues
virtual void setIntnlValues(Real p_trial, Real q_trial, Real p, Real q, const std::vector< Real > &intnl_old, std::vector< Real > &intnl) const =0
Sets the internal parameters based on the trial values of p and q, their current values,...
TwoParameterPlasticityStressUpdate::validParams
static InputParameters validParams()
Definition: TwoParameterPlasticityStressUpdate.C:18
TwoParameterPlasticityStressUpdate::setEffectiveElasticity
virtual void setEffectiveElasticity(const RankFourTensor &Eijkl) override
Sets _Eij and _En and _Cij.
Definition: TwoParameterPlasticityStressUpdate.C:52
TwoParameterPlasticityStressUpdate::_dp_dqt
Real _dp_dqt
derivative of Variable with respect to trial variable (used in consistent-tangent-operator calculatio...
Definition: TwoParameterPlasticityStressUpdate.h:63
TwoParameterPlasticityStressUpdate::_Epp
Real _Epp
elasticity tensor in p direction
Definition: TwoParameterPlasticityStressUpdate.h:49
TwoParameterPlasticityStressUpdate::dstress_param_dstress
virtual std::vector< RankTwoTensor > dstress_param_dstress(const RankTwoTensor &stress) const override
d(stress_param[i])/d(stress) at given stress
Definition: TwoParameterPlasticityStressUpdate.C:292
TwoParameterPlasticityStressUpdate::_dgaE_dqt
Real _dgaE_dqt
derivative of Variable with respect to trial variable (used in consistent-tangent-operator calculatio...
Definition: TwoParameterPlasticityStressUpdate.h:57
TwoParameterPlasticityStressUpdate::setIntnlDerivativesV
virtual void setIntnlDerivativesV(const std::vector< Real > &trial_stress_params, const std::vector< Real > &current_stress_params, const std::vector< Real > &intnl, std::vector< std::vector< Real >> &dintnl) const override
Sets the derivatives of internal parameters, based on the trial values of stress_params,...
Definition: TwoParameterPlasticityStressUpdate.C:142
MultiParameterPlasticityStressUpdate::validParams
static InputParameters validParams()
Definition: MultiParameterPlasticityStressUpdate.C:23
TwoParameterPlasticityStressUpdate::_dq_dqt
Real _dq_dqt
derivative of Variable with respect to trial variable (used in consistent-tangent-operator calculatio...
Definition: TwoParameterPlasticityStressUpdate.h:65
TwoParameterPlasticityStressUpdate::yieldFunctionValuesV
void yieldFunctionValuesV(const std::vector< Real > &stress_params, const std::vector< Real > &intnl, std::vector< Real > &yf) const override
Computes the values of the yield functions, given stress_params and intnl parameters.
Definition: TwoParameterPlasticityStressUpdate.C:42
MultiParameterPlasticityStressUpdate::_Eij
std::vector< std::vector< Real > > _Eij
E[i, j] in the system of equations to be solved.
Definition: MultiParameterPlasticityStressUpdate.h:138
TwoParameterPlasticityStressUpdate
TwoParameterPlasticityStressUpdate performs the return-map algorithm and associated stress updates fo...
Definition: TwoParameterPlasticityStressUpdate.h:29
TwoParameterPlasticityStressUpdate::TwoParameterPlasticityStressUpdate
TwoParameterPlasticityStressUpdate(const InputParameters &parameters, unsigned num_yf, unsigned num_intnl)
Definition: TwoParameterPlasticityStressUpdate.C:25
TwoParameterPlasticityStressUpdate::setStressAfterReturnV
void setStressAfterReturnV(const RankTwoTensor &stress_trial, const std::vector< Real > &stress_params, Real gaE, const std::vector< Real > &intnl, const yieldAndFlow &smoothed_q, const RankFourTensor &Eijkl, RankTwoTensor &stress) const override
Sets stress from the admissible parameters.
Definition: TwoParameterPlasticityStressUpdate.C:265
TwoParameterPlasticityStressUpdate::computePQ
virtual void computePQ(const RankTwoTensor &stress, Real &p, Real &q) const =0
Computes p and q, given stress.
TwoParameterPlasticityStressUpdate::initializeVarsV
void initializeVarsV(const std::vector< Real > &trial_stress_params, const std::vector< Real > &intnl_old, std::vector< Real > &stress_params, Real &gaE, std::vector< Real > &intnl) const override
Sets (stress_params, intnl) at "good guesses" of the solution to the Return-Map algorithm.
Definition: TwoParameterPlasticityStressUpdate.C:113
TwoParameterPlasticityStressUpdate::computeAllQV
void computeAllQV(const std::vector< Real > &stress_params, const std::vector< Real > &intnl, std::vector< yieldAndFlow > &all_q) const override
Completely fills all_q with correct values.
Definition: TwoParameterPlasticityStressUpdate.C:103
TwoParameterPlasticityStressUpdate::preReturnMapV
void preReturnMapV(const std::vector< Real > &trial_stress_params, const RankTwoTensor &stress_trial, const std::vector< Real > &intnl_old, const std::vector< Real > &yf, const RankFourTensor &Eijkl) override
Derived classes may employ this function to record stuff or do other computations prior to the return...
Definition: TwoParameterPlasticityStressUpdate.C:65
TwoParameterPlasticityStressUpdate::_dp_dpt
Real _dp_dpt
derivative of Variable with respect to trial variable (used in consistent-tangent-operator calculatio...
Definition: TwoParameterPlasticityStressUpdate.h:59
TwoParameterPlasticityStressUpdate::preReturnMap
virtual void preReturnMap(Real p_trial, Real q_trial, const RankTwoTensor &stress_trial, const std::vector< Real > &intnl_old, const std::vector< Real > &yf, const RankFourTensor &Eijkl)
Derived classes may employ this function to record stuff or do other computations prior to the return...
Definition: TwoParameterPlasticityStressUpdate.C:77
TwoParameterPlasticityStressUpdate::dpdstress
virtual RankTwoTensor dpdstress(const RankTwoTensor &stress) const =0
d(p)/d(stress) Derived classes must override this
TwoParameterPlasticityStressUpdate::_dgaE_dpt
Real _dgaE_dpt
derivative of Variable with respect to trial variable (used in consistent-tangent-operator calculatio...
Definition: TwoParameterPlasticityStressUpdate.h:55
RankTwoTensor
RankTwoTensorTempl< Real > RankTwoTensor
Definition: ACGrGrElasticDrivingForce.h:17
TwoParameterPlasticityStressUpdate::initializeVars
virtual void initializeVars(Real p_trial, Real q_trial, const std::vector< Real > &intnl_old, Real &p, Real &q, Real &gaE, std::vector< Real > &intnl) const
Sets (p, q, gaE, intnl) at "good guesses" of the solution to the Return-Map algorithm.
Definition: TwoParameterPlasticityStressUpdate.C:88
MultiParameterPlasticityStressUpdate::_En
Real _En
normalising factor
Definition: MultiParameterPlasticityStressUpdate.h:141
TwoParameterPlasticityStressUpdate::_num_pq
constexpr static int _num_pq
Number of variables = 2 = (p, q)
Definition: TwoParameterPlasticityStressUpdate.h:40
TwoParameterPlasticityStressUpdate::setIntnlDerivatives
virtual void setIntnlDerivatives(Real p_trial, Real q_trial, Real p, Real q, const std::vector< Real > &intnl, std::vector< std::vector< Real >> &dintnl) const =0
Sets the derivatives of internal parameters, based on the trial values of p and q,...
RankFourTensorTempl< Real >
TwoParameterPlasticityStressUpdate::setStressAfterReturn
virtual void setStressAfterReturn(const RankTwoTensor &stress_trial, Real p_ok, Real q_ok, Real gaE, const std::vector< Real > &intnl, const yieldAndFlow &smoothed_q, const RankFourTensor &Eijkl, RankTwoTensor &stress) const =0
Sets stress from the admissible parameters.
TwoParameterPlasticityStressUpdate::_dq_dpt
Real _dq_dpt
derivative of Variable with respect to trial variable (used in consistent-tangent-operator calculatio...
Definition: TwoParameterPlasticityStressUpdate.h:61
TwoParameterPlasticityStressUpdate::setInelasticStrainIncrementAfterReturn
void setInelasticStrainIncrementAfterReturn(const RankTwoTensor &stress_trial, Real gaE, const yieldAndFlow &smoothed_q, const RankFourTensor &elasticity_tensor, const RankTwoTensor &returned_stress, RankTwoTensor &inelastic_strain_increment) const override
Sets inelastic strain increment from the returned configuration This is called after the return-map p...
Definition: TwoParameterPlasticityStressUpdate.C:279
TwoParameterPlasticityStressUpdate::consistentTangentOperatorV
void consistentTangentOperatorV(const RankTwoTensor &stress_trial, const std::vector< Real > &trial_stress_params, const RankTwoTensor &stress, const std::vector< Real > &stress_params, Real gaE, const yieldAndFlow &smoothed_q, const RankFourTensor &Eijkl, bool compute_full_tangent_operator, const std::vector< std::vector< Real >> &dvar_dtrial, RankFourTensor &cto) override
Calculates the consistent tangent operator.
Definition: TwoParameterPlasticityStressUpdate.C:167
RankTwoTensorTempl< Real >
TwoParameterPlasticityStressUpdate::_Eqq
Real _Eqq
elasticity tensor in q direction
Definition: TwoParameterPlasticityStressUpdate.h:52
MultiParameterPlasticityStressUpdate
MultiParameterPlasticityStressUpdate performs the return-map algorithm and associated stress updates ...
Definition: MultiParameterPlasticityStressUpdate.h:99
RankFourTensor
RankFourTensorTempl< Real > RankFourTensor
Definition: ACGrGrElasticDrivingForce.h:20
TwoParameterPlasticityStressUpdate::computeAllQ
virtual void computeAllQ(Real p, Real q, const std::vector< Real > &intnl, std::vector< yieldAndFlow > &all_q) const =0
Completely fills all_q with correct values.
MultiParameterPlasticityStressUpdate::_tensor_dimensionality
constexpr static unsigned _tensor_dimensionality
Internal dimensionality of tensors (currently this is 3 throughout tensor_mechanics)
Definition: MultiParameterPlasticityStressUpdate.h:129
TwoParameterPlasticityStressUpdate::consistentTangentOperator
virtual void consistentTangentOperator(const RankTwoTensor &stress_trial, Real p_trial, Real q_trial, const RankTwoTensor &stress, Real p, Real q, Real gaE, const yieldAndFlow &smoothed_q, const RankFourTensor &Eijkl, bool compute_full_tangent_operator, RankFourTensor &cto) const
Calculates the consistent tangent operator.
Definition: TwoParameterPlasticityStressUpdate.C:203
TwoParameterPlasticityStressUpdate::dqdstress
virtual RankTwoTensor dqdstress(const RankTwoTensor &stress) const =0
d(q)/d(stress) Derived classes must override this
TwoParameterPlasticityStressUpdate.h