www.mooseframework.org
FiniteStrainPlasticMaterial.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 #include "libmesh/utility.h"
12 
14 
16 
17 InputParameters
19 {
20  InputParameters params = ComputeStressBase::validParams();
21 
22  params.addRequiredParam<std::vector<Real>>(
23  "yield_stress",
24  "Input data as pairs of equivalent plastic strain and yield stress: Should "
25  "start with equivalent plastic strain 0");
26  params.addParam<Real>("rtol", 1e-8, "Plastic strain NR tolerance");
27  params.addParam<Real>("ftol", 1e-4, "Consistency condition NR tolerance");
28  params.addParam<Real>("eptol", 1e-7, "Equivalent plastic strain NR tolerance");
29  params.addClassDescription("Associative J2 plasticity with isotropic hardening.");
30 
31  return params;
32 }
33 
35  : ComputeStressBase(parameters),
36  _yield_stress_vector(getParam<std::vector<Real>>("yield_stress")), // Read from input file
37  _plastic_strain(declareProperty<RankTwoTensor>(_base_name + "plastic_strain")),
38  _plastic_strain_old(getMaterialPropertyOld<RankTwoTensor>(_base_name + "plastic_strain")),
39  _eqv_plastic_strain(declareProperty<Real>(_base_name + "eqv_plastic_strain")),
40  _eqv_plastic_strain_old(getMaterialPropertyOld<Real>(_base_name + "eqv_plastic_strain")),
41  _stress_old(getMaterialPropertyOld<RankTwoTensor>(_base_name + "stress")),
42  _strain_increment(getMaterialProperty<RankTwoTensor>(_base_name + "strain_increment")),
43  _rotation_increment(getMaterialProperty<RankTwoTensor>(_base_name + "rotation_increment")),
44  _elasticity_tensor_name(_base_name + "elasticity_tensor"),
45  _elasticity_tensor(getMaterialPropertyByName<RankFourTensor>(_elasticity_tensor_name)),
46  _rtol(getParam<Real>("rtol")),
47  _ftol(getParam<Real>("ftol")),
48  _eptol(getParam<Real>("eptol")),
49  _deltaOuter(RankTwoTensor::Identity().outerProduct(RankTwoTensor::Identity())),
50  _deltaMixed(RankTwoTensor::Identity().mixedProductIkJl(RankTwoTensor::Identity()))
51 {
52 }
53 
54 void
56 {
58  _plastic_strain[_qp].zero();
59  _eqv_plastic_strain[_qp] = 0.0;
60 }
61 
62 void
64 {
65 
66  // perform the return-mapping algorithm
70  _strain_increment[_qp],
71  _elasticity_tensor[_qp],
72  _stress[_qp],
74  _plastic_strain[_qp]);
75 
76  // Rotate the stress tensor to the current configuration
77  _stress[_qp] = _rotation_increment[_qp] * _stress[_qp] * _rotation_increment[_qp].transpose();
78 
79  // Rotate plastic strain tensor to the current configuration
80  _plastic_strain[_qp] =
81  _rotation_increment[_qp] * _plastic_strain[_qp] * _rotation_increment[_qp].transpose();
82 
83  // Calculate the elastic strain_increment
85 
87 }
88 
150 void
152  const Real eqvpstrain_old,
153  const RankTwoTensor & plastic_strain_old,
154  const RankTwoTensor & delta_d,
155  const RankFourTensor & E_ijkl,
156  RankTwoTensor & sig,
157  Real & eqvpstrain,
158  RankTwoTensor & plastic_strain)
159 {
160  // the yield function, must be non-positive
161  // Newton-Raphson sets this to zero if trial stress enters inadmissible region
162  Real f;
163 
164  // the consistency parameter, must be non-negative
165  // change in plastic strain in this timestep = flow_incr*flow_potential
166  Real flow_incr = 0.0;
167 
168  // direction of flow defined by the potential
169  RankTwoTensor flow_dirn;
170 
171  // Newton-Raphson sets this zero
172  // resid_ij = flow_incr*flow_dirn_ij - (plastic_strain - plastic_strain_old)
173  RankTwoTensor resid;
174 
175  // Newton-Raphson sets this zero
176  // rep = -flow_incr*internalPotential - (eqvpstrain - eqvpstrain_old)
177  Real rep;
178 
179  // change in the stress (sig) in a Newton-Raphson iteration
180  RankTwoTensor ddsig;
181 
182  // change in the consistency parameter in a Newton-Raphson iteration
183  Real dflow_incr = 0.0;
184 
185  // change in equivalent plastic strain in one Newton-Raphson iteration
186  Real deqvpstrain = 0.0;
187 
188  // convenience variable that holds the change in plastic strain incurred during the return
189  // delta_dp = plastic_strain - plastic_strain_old
190  // delta_dp = E^{-1}*(trial_stress - sig), where trial_stress = E*(strain - plastic_strain_old)
191  RankTwoTensor delta_dp;
192 
193  // d(yieldFunction)/d(stress)
194  RankTwoTensor df_dsig;
195 
196  // d(resid_ij)/d(sigma_kl)
197  RankFourTensor dr_dsig;
198 
199  // dr_dsig_inv_ijkl*dr_dsig_klmn = 0.5*(de_ij de_jn + de_ij + de_jm), where de_ij = 1 if i=j, but
200  // zero otherwise
201  RankFourTensor dr_dsig_inv;
202 
203  // d(yieldFunction)/d(eqvpstrain)
204  Real fq;
205 
206  // yield stress at the start of this "time step" (ie, evaluated with
207  // eqvpstrain_old). It is held fixed during the Newton-Raphson return,
208  // even if eqvpstrain != eqvpstrain_old.
209  Real yield_stress;
210 
211  // measures of whether the Newton-Raphson process has converged
212  Real err1, err2, err3;
213 
214  // number of Newton-Raphson iterations performed
215  unsigned int iter = 0;
216 
217  // maximum number of Newton-Raphson iterations allowed
218  unsigned int maxiter = 100;
219 
220  // plastic loading occurs if yieldFunction > toly
221  Real toly = 1.0e-8;
222 
223  // Assume this strain increment does not induce any plasticity
224  // This is the elastic-predictor
225  sig = sig_old + E_ijkl * delta_d; // the trial stress
226  eqvpstrain = eqvpstrain_old;
227  plastic_strain = plastic_strain_old;
228 
229  yield_stress = getYieldStress(eqvpstrain); // yield stress at this equivalent plastic strain
230  if (yieldFunction(sig, yield_stress) > toly)
231  {
232  // the sig just calculated is inadmissable. We must return to the yield surface.
233  // This is done iteratively, using a Newton-Raphson process.
234 
235  delta_dp.zero();
236 
237  sig = sig_old + E_ijkl * delta_d; // this is the elastic predictor
238 
239  flow_dirn = flowPotential(sig);
240 
241  resid = flow_dirn * flow_incr - delta_dp; // Residual 1 - refer Hughes Simo
242  f = yieldFunction(sig, yield_stress);
243  rep = -eqvpstrain + eqvpstrain_old - flow_incr * internalPotential(); // Residual 3 rep=0
244 
245  err1 = resid.L2norm();
246  err2 = std::abs(f);
247  err3 = std::abs(rep);
248 
249  while ((err1 > _rtol || err2 > _ftol || err3 > _eptol) &&
250  iter < maxiter) // Stress update iteration (hardness fixed)
251  {
252  iter++;
253 
254  df_dsig = dyieldFunction_dstress(sig);
255  getJac(sig, E_ijkl, flow_incr, dr_dsig); // gets dr_dsig = d(resid_ij)/d(sig_kl)
256  fq = dyieldFunction_dinternal(eqvpstrain); // d(f)/d(eqvpstrain)
257 
269  dr_dsig_inv = dr_dsig.invSymm();
270 
278  dflow_incr = (f - df_dsig.doubleContraction(dr_dsig_inv * resid) + fq * rep) /
279  (df_dsig.doubleContraction(dr_dsig_inv * flow_dirn) - fq);
280  ddsig =
281  dr_dsig_inv *
282  (-resid -
283  flow_dirn * dflow_incr); // from solving the top row of linear system, given dflow_incr
284  deqvpstrain =
285  rep + dflow_incr; // from solving the bottom row of linear system, given dflow_incr
286 
287  // update the variables
288  flow_incr += dflow_incr;
289  delta_dp -= E_ijkl.invSymm() * ddsig;
290  sig += ddsig;
291  eqvpstrain += deqvpstrain;
292 
293  // evaluate the RHS equations ready for next Newton-Raphson iteration
294  flow_dirn = flowPotential(sig);
295  resid = flow_dirn * flow_incr - delta_dp;
296  f = yieldFunction(sig, yield_stress);
297  rep = -eqvpstrain + eqvpstrain_old - flow_incr * internalPotential();
298 
299  err1 = resid.L2norm();
300  err2 = std::abs(f);
301  err3 = std::abs(rep);
302  }
303 
304  if (iter >= maxiter)
305  mooseError("Constitutive failure");
306 
307  plastic_strain += delta_dp;
308  }
309 }
310 
311 Real
312 FiniteStrainPlasticMaterial::yieldFunction(const RankTwoTensor & stress, const Real yield_stress)
313 {
314  return getSigEqv(stress) - yield_stress;
315 }
316 
319 {
320  RankTwoTensor deriv = sig.dsecondInvariant();
321  deriv *= std::sqrt(3.0 / sig.secondInvariant()) / 2.0;
322  return deriv;
323 }
324 
325 Real
326 FiniteStrainPlasticMaterial::dyieldFunction_dinternal(const Real equivalent_plastic_strain)
327 {
328  return -getdYieldStressdPlasticStrain(equivalent_plastic_strain);
329 }
330 
333 {
334  return dyieldFunction_dstress(sig); // this plasticity model assumes associative flow
335 }
336 
337 Real
339 {
340  return -1;
341 }
342 
343 Real
345 {
346  return std::sqrt(3 * stress.secondInvariant());
347 }
348 
349 // Jacobian for stress update algorithm
350 void
352  const RankFourTensor & E_ijkl,
353  Real flow_incr,
354  RankFourTensor & dresid_dsig)
355 {
356  RankTwoTensor sig_dev, df_dsig, flow_dirn;
357  RankTwoTensor dfi_dft, dfi_dsig;
358  RankFourTensor dft_dsig, dfd_dft, dfd_dsig;
359  Real sig_eqv;
360  Real f1, f2, f3;
361  RankFourTensor temp;
362 
363  sig_dev = sig.deviatoric();
364  sig_eqv = getSigEqv(sig);
365  df_dsig = dyieldFunction_dstress(sig);
366  flow_dirn = flowPotential(sig);
367 
368  f1 = 3.0 / (2.0 * sig_eqv);
369  f2 = f1 / 3.0;
370  f3 = 9.0 / (4.0 * Utility::pow<3>(sig_eqv));
371 
372  dft_dsig = f1 * _deltaMixed - f2 * _deltaOuter - f3 * sig_dev.outerProduct(sig_dev);
373 
374  dfd_dsig = dft_dsig;
375  dresid_dsig = E_ijkl.invSymm() + dfd_dsig * flow_incr;
376 }
377 
378 // Obtain yield stress for a given equivalent plastic strain (input)
379 Real
381 {
382  unsigned nsize;
383 
384  nsize = _yield_stress_vector.size();
385 
386  if (_yield_stress_vector[0] > 0.0 || nsize % 2 > 0) // Error check for input inconsitency
387  mooseError("Error in yield stress input: Should be a vector with eqv plastic strain and yield "
388  "stress pair values.\n");
389 
390  unsigned int ind = 0;
391  Real tol = 1e-8;
392 
393  while (ind < nsize)
394  {
395  if (std::abs(eqpe - _yield_stress_vector[ind]) < tol)
396  return _yield_stress_vector[ind + 1];
397 
398  if (ind + 2 < nsize)
399  {
400  if (eqpe > _yield_stress_vector[ind] && eqpe < _yield_stress_vector[ind + 2])
401  return _yield_stress_vector[ind + 1] +
402  (eqpe - _yield_stress_vector[ind]) /
403  (_yield_stress_vector[ind + 2] - _yield_stress_vector[ind]) *
404  (_yield_stress_vector[ind + 3] - _yield_stress_vector[ind + 1]);
405  }
406  else
407  return _yield_stress_vector[nsize - 1];
408 
409  ind += 2;
410  }
411 
412  return 0.0;
413 }
414 
415 Real
417 {
418  unsigned nsize;
419 
420  nsize = _yield_stress_vector.size();
421 
422  if (_yield_stress_vector[0] > 0.0 || nsize % 2 > 0) // Error check for input inconsitency
423  mooseError("Error in yield stress input: Should be a vector with eqv plastic strain and yield "
424  "stress pair values.\n");
425 
426  unsigned int ind = 0;
427 
428  while (ind < nsize)
429  {
430  if (ind + 2 < nsize)
431  {
432  if (eqpe >= _yield_stress_vector[ind] && eqpe < _yield_stress_vector[ind + 2])
433  return (_yield_stress_vector[ind + 3] - _yield_stress_vector[ind + 1]) /
434  (_yield_stress_vector[ind + 2] - _yield_stress_vector[ind]);
435  }
436  else
437  return 0.0;
438 
439  ind += 2;
440  }
441 
442  return 0.0;
443 }
FiniteStrainPlasticMaterial::_plastic_strain
MaterialProperty< RankTwoTensor > & _plastic_strain
Definition: FiniteStrainPlasticMaterial.h:40
FiniteStrainPlasticMaterial::_elasticity_tensor
const MaterialProperty< RankFourTensor > & _elasticity_tensor
Elasticity tensor material property.
Definition: FiniteStrainPlasticMaterial.h:50
registerMooseObject
registerMooseObject("TensorMechanicsApp", FiniteStrainPlasticMaterial)
ComputeStressBase::_stress
MaterialProperty< RankTwoTensor > & _stress
Stress material property.
Definition: ComputeStressBase.h:50
ComputeStressBase::_Jacobian_mult
MaterialProperty< RankFourTensor > & _Jacobian_mult
derivative of stress w.r.t. strain (_dstress_dstrain)
Definition: ComputeStressBase.h:61
FiniteStrainPlasticMaterial::_rotation_increment
const MaterialProperty< RankTwoTensor > & _rotation_increment
Definition: FiniteStrainPlasticMaterial.h:46
FiniteStrainPlasticMaterial::_deltaMixed
RankFourTensor _deltaMixed
Definition: FiniteStrainPlasticMaterial.h:56
FiniteStrainPlasticMaterial
FiniteStrainPlasticMaterial implements rate-independent associative J2 plasticity with isotropic hard...
Definition: FiniteStrainPlasticMaterial.h:29
FiniteStrainPlasticMaterial::_stress_old
const MaterialProperty< RankTwoTensor > & _stress_old
Definition: FiniteStrainPlasticMaterial.h:44
FiniteStrainPlasticMaterial::_yield_stress_vector
std::vector< Real > _yield_stress_vector
Definition: FiniteStrainPlasticMaterial.h:39
FiniteStrainPlasticMaterial::returnMap
virtual void returnMap(const RankTwoTensor &sig_old, const Real eqvpstrain_old, const RankTwoTensor &plastic_strain_old, const RankTwoTensor &delta_d, const RankFourTensor &E_ijkl, RankTwoTensor &sig, Real &eqvpstrain, RankTwoTensor &plastic_strain)
Implements the return map.
Definition: FiniteStrainPlasticMaterial.C:151
FiniteStrainPlasticMaterial::FiniteStrainPlasticMaterial
FiniteStrainPlasticMaterial(const InputParameters &parameters)
Definition: FiniteStrainPlasticMaterial.C:34
FiniteStrainPlasticMaterial::_eqv_plastic_strain_old
const MaterialProperty< Real > & _eqv_plastic_strain_old
Definition: FiniteStrainPlasticMaterial.h:43
FiniteStrainPlasticMaterial::dyieldFunction_dinternal
virtual Real dyieldFunction_dinternal(const Real equivalent_plastic_strain)
Derivative of yieldFunction with respect to the equivalent plastic strain.
Definition: FiniteStrainPlasticMaterial.C:326
FiniteStrainPlasticMaterial.h
defineLegacyParams
defineLegacyParams(FiniteStrainPlasticMaterial)
ComputeStressBase
ComputeStressBase is the base class for stress tensors.
Definition: ComputeStressBase.h:26
FiniteStrainPlasticMaterial::yieldFunction
virtual Real yieldFunction(const RankTwoTensor &stress, const Real yield_stress)
Calculates the yield function.
Definition: FiniteStrainPlasticMaterial.C:312
FiniteStrainPlasticMaterial::computeQpStress
virtual void computeQpStress()
Compute the stress and store it in the _stress material property for the current quadrature point.
Definition: FiniteStrainPlasticMaterial.C:63
FiniteStrainPlasticMaterial::getYieldStress
Real getYieldStress(const Real equivalent_plastic_strain)
yield stress as a function of equivalent plastic strain.
Definition: FiniteStrainPlasticMaterial.C:380
ComputeStressBase::validParams
static InputParameters validParams()
Definition: ComputeStressBase.C:17
FiniteStrainPlasticMaterial::getdYieldStressdPlasticStrain
Real getdYieldStressdPlasticStrain(const Real equivalent_plastic_strain)
d(yieldstress)/d(equivalent plastic strain)
Definition: FiniteStrainPlasticMaterial.C:416
FiniteStrainPlasticMaterial::getJac
virtual void getJac(const RankTwoTensor &sig, const RankFourTensor &E_ijkl, Real flow_incr, RankFourTensor &dresid_dsig)
Evaluates the derivative d(resid_ij)/d(sig_kl), where resid_ij = flow_incr*flowPotential_ij - (E^{-1}...
Definition: FiniteStrainPlasticMaterial.C:351
FiniteStrainPlasticMaterial::dyieldFunction_dstress
virtual RankTwoTensor dyieldFunction_dstress(const RankTwoTensor &stress)
Derivative of yieldFunction with respect to the stress.
Definition: FiniteStrainPlasticMaterial.C:318
FiniteStrainPlasticMaterial::_strain_increment
const MaterialProperty< RankTwoTensor > & _strain_increment
Definition: FiniteStrainPlasticMaterial.h:45
ComputeStressBase::initQpStatefulProperties
virtual void initQpStatefulProperties() override
Definition: ComputeStressBase.C:43
FiniteStrainPlasticMaterial::_rtol
Real _rtol
Definition: FiniteStrainPlasticMaterial.h:51
tol
const double tol
Definition: Setup.h:18
FiniteStrainPlasticMaterial::getSigEqv
Real getSigEqv(const RankTwoTensor &stress)
Equivalent stress.
Definition: FiniteStrainPlasticMaterial.C:344
FiniteStrainPlasticMaterial::_plastic_strain_old
const MaterialProperty< RankTwoTensor > & _plastic_strain_old
Definition: FiniteStrainPlasticMaterial.h:41
FiniteStrainPlasticMaterial::_ftol
Real _ftol
Definition: FiniteStrainPlasticMaterial.h:52
ComputeStressBase::_elastic_strain
MaterialProperty< RankTwoTensor > & _elastic_strain
Elastic strain material property.
Definition: ComputeStressBase.h:52
FiniteStrainPlasticMaterial::initQpStatefulProperties
virtual void initQpStatefulProperties()
Definition: FiniteStrainPlasticMaterial.C:55
FiniteStrainPlasticMaterial::internalPotential
virtual Real internalPotential()
The internal potential.
Definition: FiniteStrainPlasticMaterial.C:338
RankFourTensorTempl< Real >
FiniteStrainPlasticMaterial::flowPotential
virtual RankTwoTensor flowPotential(const RankTwoTensor &stress)
Flow potential, which in this case is just dyieldFunction_dstress because we are doing associative fl...
Definition: FiniteStrainPlasticMaterial.C:332
FiniteStrainPlasticMaterial::validParams
static InputParameters validParams()
Definition: FiniteStrainPlasticMaterial.C:18
FiniteStrainPlasticMaterial::_eqv_plastic_strain
MaterialProperty< Real > & _eqv_plastic_strain
Definition: FiniteStrainPlasticMaterial.h:42
ComputeStressBase::_mechanical_strain
const MaterialProperty< RankTwoTensor > & _mechanical_strain
Mechanical strain material property.
Definition: ComputeStressBase.h:48
RankTwoTensorTempl< Real >
FiniteStrainPlasticMaterial::_deltaOuter
RankFourTensor _deltaOuter
Definition: FiniteStrainPlasticMaterial.h:56
FiniteStrainPlasticMaterial::_eptol
Real _eptol
Definition: FiniteStrainPlasticMaterial.h:53