www.mooseframework.org
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
FiniteStrainPlasticMaterial Class Reference

FiniteStrainPlasticMaterial implements rate-independent associative J2 plasticity with isotropic hardening in the finite-strain framework. More...

#include <FiniteStrainPlasticMaterial.h>

Inheritance diagram for FiniteStrainPlasticMaterial:
[legend]

Public Member Functions

 FiniteStrainPlasticMaterial (const InputParameters &parameters)
 

Static Public Member Functions

static InputParameters validParams ()
 

Protected Member Functions

virtual void computeQpStress ()
 Compute the stress and store it in the _stress material property for the current quadrature point. More...
 
virtual void initQpStatefulProperties ()
 
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. More...
 
virtual Real yieldFunction (const RankTwoTensor &stress, const Real yield_stress)
 Calculates the yield function. More...
 
virtual RankTwoTensor dyieldFunction_dstress (const RankTwoTensor &stress)
 Derivative of yieldFunction with respect to the stress. More...
 
virtual Real dyieldFunction_dinternal (const Real equivalent_plastic_strain)
 Derivative of yieldFunction with respect to the equivalent plastic strain. More...
 
virtual RankTwoTensor flowPotential (const RankTwoTensor &stress)
 Flow potential, which in this case is just dyieldFunction_dstress because we are doing associative flow, and hence does not depend on the internal hardening parameter equivalent_plastic_strain. More...
 
virtual Real internalPotential ()
 The internal potential. More...
 
Real getSigEqv (const RankTwoTensor &stress)
 Equivalent stress. More...
 
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}(trial_stress - sig))_ij. More...
 
Real getYieldStress (const Real equivalent_plastic_strain)
 yield stress as a function of equivalent plastic strain. More...
 
Real getdYieldStressdPlasticStrain (const Real equivalent_plastic_strain)
 d(yieldstress)/d(equivalent plastic strain) More...
 
virtual void computeQpProperties () override
 

Protected Attributes

std::vector< Real > _yield_stress_vector
 
MaterialProperty< RankTwoTensor > & _plastic_strain
 
const MaterialProperty< RankTwoTensor > & _plastic_strain_old
 
MaterialProperty< Real > & _eqv_plastic_strain
 
const MaterialProperty< Real > & _eqv_plastic_strain_old
 
const MaterialProperty< RankTwoTensor > & _stress_old
 
const MaterialProperty< RankTwoTensor > & _strain_increment
 
const MaterialProperty< RankTwoTensor > & _rotation_increment
 
const std::string _elasticity_tensor_name
 Name of the elasticity tensor material property. More...
 
const MaterialProperty< RankFourTensor > & _elasticity_tensor
 Elasticity tensor material property. More...
 
Real _rtol
 
Real _ftol
 
Real _eptol
 
RankFourTensor _deltaOuter
 
RankFourTensor _deltaMixed
 
const std::string _base_name
 Base name prepended to all material property names to allow for multi-material systems. More...
 
const MaterialProperty< RankTwoTensor > & _mechanical_strain
 Mechanical strain material property. More...
 
MaterialProperty< RankTwoTensor > & _stress
 Stress material property. More...
 
MaterialProperty< RankTwoTensor > & _elastic_strain
 Elastic strain material property. More...
 
const MaterialProperty< RankTwoTensor > & _extra_stress
 Extra stress tensor. More...
 
std::vector< const Function * > _initial_stress_fcn
 initial stress components More...
 
MaterialProperty< RankFourTensor > & _Jacobian_mult
 derivative of stress w.r.t. strain (_dstress_dstrain) More...
 

Detailed Description

FiniteStrainPlasticMaterial implements rate-independent associative J2 plasticity with isotropic hardening in the finite-strain framework.

Yield function = sqrt(3*s_ij*s_ij/2) - K(equivalent plastic strain) where s_ij = stress_ij - delta_ij*trace(stress)/3 is the deviatoric stress and K is the yield stress, specified as a piecewise-linear function by the user. Integration is performed in an incremental manner using Newton Raphson.

Definition at line 29 of file FiniteStrainPlasticMaterial.h.

Constructor & Destructor Documentation

◆ FiniteStrainPlasticMaterial()

FiniteStrainPlasticMaterial::FiniteStrainPlasticMaterial ( const InputParameters &  parameters)

Definition at line 34 of file FiniteStrainPlasticMaterial.C.

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 }

Member Function Documentation

◆ computeQpProperties()

void ComputeStressBase::computeQpProperties ( )
overrideprotectedvirtualinherited

Definition at line 50 of file ComputeStressBase.C.

51 {
53 
54  // Add in extra stress
55  _stress[_qp] += _extra_stress[_qp];
56 }

◆ computeQpStress()

void FiniteStrainPlasticMaterial::computeQpStress ( )
protectedvirtual

Compute the stress and store it in the _stress material property for the current quadrature point.

Implements ComputeStressBase.

Definition at line 63 of file FiniteStrainPlasticMaterial.C.

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 }

◆ dyieldFunction_dinternal()

Real FiniteStrainPlasticMaterial::dyieldFunction_dinternal ( const Real  equivalent_plastic_strain)
protectedvirtual

Derivative of yieldFunction with respect to the equivalent plastic strain.

Definition at line 326 of file FiniteStrainPlasticMaterial.C.

327 {
328  return -getdYieldStressdPlasticStrain(equivalent_plastic_strain);
329 }

Referenced by returnMap().

◆ dyieldFunction_dstress()

RankTwoTensor FiniteStrainPlasticMaterial::dyieldFunction_dstress ( const RankTwoTensor stress)
protectedvirtual

Derivative of yieldFunction with respect to the stress.

Definition at line 318 of file FiniteStrainPlasticMaterial.C.

319 {
320  RankTwoTensor deriv = sig.dsecondInvariant();
321  deriv *= std::sqrt(3.0 / sig.secondInvariant()) / 2.0;
322  return deriv;
323 }

Referenced by flowPotential(), getJac(), and returnMap().

◆ flowPotential()

RankTwoTensor FiniteStrainPlasticMaterial::flowPotential ( const RankTwoTensor stress)
protectedvirtual

Flow potential, which in this case is just dyieldFunction_dstress because we are doing associative flow, and hence does not depend on the internal hardening parameter equivalent_plastic_strain.

Definition at line 332 of file FiniteStrainPlasticMaterial.C.

333 {
334  return dyieldFunction_dstress(sig); // this plasticity model assumes associative flow
335 }

Referenced by getJac(), and returnMap().

◆ getdYieldStressdPlasticStrain()

Real FiniteStrainPlasticMaterial::getdYieldStressdPlasticStrain ( const Real  equivalent_plastic_strain)
protected

d(yieldstress)/d(equivalent plastic strain)

Definition at line 416 of file FiniteStrainPlasticMaterial.C.

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 }

Referenced by dyieldFunction_dinternal().

◆ getJac()

void FiniteStrainPlasticMaterial::getJac ( const RankTwoTensor sig,
const RankFourTensor E_ijkl,
Real  flow_incr,
RankFourTensor dresid_dsig 
)
protectedvirtual

Evaluates the derivative d(resid_ij)/d(sig_kl), where resid_ij = flow_incr*flowPotential_ij - (E^{-1}(trial_stress - sig))_ij.

Parameters
sigstress
E_ijklelasticity tensor (sig = E*(strain - plastic_strain))
flow_incrconsistency parameter
dresid_dsigthe required derivative (this is an output variable)

Definition at line 351 of file FiniteStrainPlasticMaterial.C.

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 }

Referenced by returnMap().

◆ getSigEqv()

Real FiniteStrainPlasticMaterial::getSigEqv ( const RankTwoTensor stress)
protected

Equivalent stress.

Definition at line 344 of file FiniteStrainPlasticMaterial.C.

345 {
346  return std::sqrt(3 * stress.secondInvariant());
347 }

Referenced by getJac(), and yieldFunction().

◆ getYieldStress()

Real FiniteStrainPlasticMaterial::getYieldStress ( const Real  equivalent_plastic_strain)
protected

yield stress as a function of equivalent plastic strain.

This is a piecewise linear function entered by the user in the yield_stress vector

Definition at line 380 of file FiniteStrainPlasticMaterial.C.

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 }

Referenced by returnMap().

◆ initQpStatefulProperties()

void FiniteStrainPlasticMaterial::initQpStatefulProperties ( )
protectedvirtual

Reimplemented from ComputeStressBase.

Definition at line 55 of file FiniteStrainPlasticMaterial.C.

56 {
58  _plastic_strain[_qp].zero();
59  _eqv_plastic_strain[_qp] = 0.0;
60 }

◆ internalPotential()

Real FiniteStrainPlasticMaterial::internalPotential ( )
protectedvirtual

The internal potential.

For associative J2 plasticity this is just -1

Definition at line 338 of file FiniteStrainPlasticMaterial.C.

339 {
340  return -1;
341 }

Referenced by returnMap().

◆ returnMap()

void FiniteStrainPlasticMaterial::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 
)
protectedvirtual

Implements the return map.

Implements the return-map algorithm via a Newton-Raphson process.

Parameters
sig_oldThe stress at the previous "time" step
eqvpstrain_oldThe equivalent plastic strain at the previous "time" step
plastic_strain_oldThe value of plastic strain at the previous "time" step
delta_dThe total strain increment for this "time" step
E_ijklThe elasticity tensor. If no plasiticity then sig_new = sig_old + E_ijkl*delta_d
sigThe stress after returning to the yield surface (this is an output variable)
eqvpstrainThe equivalent plastic strain after returning to the yield surface (this is an output variable)
plastic_strainThe value of plastic strain after returning to the yield surface (this is an output variable) Note that this algorithm doesn't do any rotations. In order to find the final stress and plastic_strain, sig and plastic_strain must be rotated using _rotation_increment.

This idea is fully explained in Simo and Hughes "Computational Inelasticity" Springer 1997, for instance, as well as many other books on plasticity. The basic idea is as follows. Given: sig_old - the stress at the start of the "time step" plastic_strain_old - the plastic strain at the start of the "time step" eqvpstrain_old - equivalent plastic strain at the start of the "time step" (In general we would be given some number of internal parameters that the yield function depends upon.) delta_d - the prescribed strain increment for this "time step" we want to determine the following parameters at the end of this "time step": sig - the stress plastic_strain - the plastic strain eqvpstrain - the equivalent plastic strain (again, in general, we would have an arbitrary number of internal parameters).

To determine these parameters, introduce the "yield function", f the "consistency parameter", flow_incr the "flow potential", flow_dirn_ij the "internal potential", internalPotential (in general there are as many internalPotential functions as there are internal parameters). All three of f, flow_dirn_ij, and internalPotential, are functions of sig and eqvpstrain. To find sig, plastic_strain and eqvpstrain, we need to solve the following resid_ij = 0 f = 0 rep = 0 This is done by using Newton-Raphson. There are 8 equations here: six from resid_ij=0 (more generally there are nine but in this case resid is symmetric); one from f=0; one from rep=0 (more generally, for N internal parameters there are N of these equations).

resid_ij = flow_incr*flow_dirn_ij - (plastic_strain - plastic_strain_old)_ij = flow_incr*flow_dirn_ij - (E^{-1}(trial_stress - sig))_ij Here trial_stress = E*(strain - plastic_strain_old) sig = E*(strain - plastic_strain) Note: flow_dirn_ij is evaluated at sig and eqvpstrain (not the old values).

f is the yield function, evaluated at sig and eqvpstrain

rep = -flow_incr*internalPotential - (eqvpstrain - eqvpstrain_old) Here internalPotential are evaluated at sig and eqvpstrain (not the old values).

The Newton-Raphson procedure has sig, flow_incr, and eqvpstrain as its variables. Therefore we need the derivatives of resid_ij, f, and rep with respect to these parameters

In this associative J2 with isotropic hardening, things are a little more specialised. (1) f = sqrt(3*s_ij*s_ij/2) - K(eqvpstrain) (this is called "isotropic hardening") (2) associativity means that flow_dirn_ij = df/d(sig_ij) = s_ij*sqrt(3/2/(s_ij*s_ij)), and this means that flow_dirn_ij*flow_dirn_ij = 3/2, so when resid_ij=0, we get (plastic_strain_dot)_ij*(plastic_strain_dot)_ij = (3/2)*flow_incr^2, where plastic_strain_dot = plastic_strain - plastic_strain_old (3) The definition of equivalent plastic strain is through eqvpstrain_dot = sqrt(2*plastic_strain_dot_ij*plastic_strain_dot_ij/3), so using (2), we obtain eqvpstrain_dot = flow_incr, and this yields internalPotential = -1 in the "rep" equation.

The linear system is ( dr_dsig flow_dirn 0 )( ddsig ) ( - resid ) ( df_dsig 0 fq )( dflow_incr ) = ( - f ) ( 0 1 -1 )( deqvpstrain ) ( - rep ) The zeroes are: d(resid_ij)/d(eqvpstrain) = flow_dirn*d(df/d(sig_ij))/d(eqvpstrain) = 0 and df/d(flow_dirn) = 0 (this is always true, even for general hardening and non-associative) and d(rep)/d(sig_ij) = -flow_incr*d(internalPotential)/d(sig_ij) = 0

Because of the zeroes and ones, the linear system is not impossible to solve by hand. NOTE: andy believes there was originally a sign-error in the next line. The next line is unchanged, however andy's definition of fq is negative of the original definition of fq. andy can't see any difference in any tests!

Definition at line 151 of file FiniteStrainPlasticMaterial.C.

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 }

Referenced by computeQpStress().

◆ validParams()

InputParameters FiniteStrainPlasticMaterial::validParams ( )
static

Definition at line 18 of file FiniteStrainPlasticMaterial.C.

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 }

◆ yieldFunction()

Real FiniteStrainPlasticMaterial::yieldFunction ( const RankTwoTensor stress,
const Real  yield_stress 
)
protectedvirtual

Calculates the yield function.

Parameters
stressthe stress at which to calculate the yield function
yield_stressthe current value of the yield stress
Returns
equivalentstress - yield_stress

Definition at line 312 of file FiniteStrainPlasticMaterial.C.

313 {
314  return getSigEqv(stress) - yield_stress;
315 }

Referenced by returnMap().

Member Data Documentation

◆ _base_name

const std::string ComputeStressBase::_base_name
protectedinherited

Base name prepended to all material property names to allow for multi-material systems.

Definition at line 45 of file ComputeStressBase.h.

Referenced by ComputeLinearElasticStress::initialSetup(), and ComputeCosseratLinearElasticStress::initialSetup().

◆ _deltaMixed

RankFourTensor FiniteStrainPlasticMaterial::_deltaMixed
protected

Definition at line 56 of file FiniteStrainPlasticMaterial.h.

Referenced by getJac().

◆ _deltaOuter

RankFourTensor FiniteStrainPlasticMaterial::_deltaOuter
protected

Definition at line 56 of file FiniteStrainPlasticMaterial.h.

Referenced by getJac().

◆ _elastic_strain

MaterialProperty<RankTwoTensor>& ComputeStressBase::_elastic_strain
protectedinherited

◆ _elasticity_tensor

const MaterialProperty<RankFourTensor>& FiniteStrainPlasticMaterial::_elasticity_tensor
protected

Elasticity tensor material property.

Definition at line 50 of file FiniteStrainPlasticMaterial.h.

Referenced by computeQpStress().

◆ _elasticity_tensor_name

const std::string FiniteStrainPlasticMaterial::_elasticity_tensor_name
protected

Name of the elasticity tensor material property.

Definition at line 48 of file FiniteStrainPlasticMaterial.h.

◆ _eptol

Real FiniteStrainPlasticMaterial::_eptol
protected

Definition at line 53 of file FiniteStrainPlasticMaterial.h.

Referenced by returnMap().

◆ _eqv_plastic_strain

MaterialProperty<Real>& FiniteStrainPlasticMaterial::_eqv_plastic_strain
protected

Definition at line 42 of file FiniteStrainPlasticMaterial.h.

Referenced by computeQpStress(), and initQpStatefulProperties().

◆ _eqv_plastic_strain_old

const MaterialProperty<Real>& FiniteStrainPlasticMaterial::_eqv_plastic_strain_old
protected

Definition at line 43 of file FiniteStrainPlasticMaterial.h.

Referenced by computeQpStress().

◆ _extra_stress

const MaterialProperty<RankTwoTensor>& ComputeStressBase::_extra_stress
protectedinherited

Extra stress tensor.

Definition at line 55 of file ComputeStressBase.h.

Referenced by ComputeStressBase::computeQpProperties().

◆ _ftol

Real FiniteStrainPlasticMaterial::_ftol
protected

Definition at line 52 of file FiniteStrainPlasticMaterial.h.

Referenced by returnMap().

◆ _initial_stress_fcn

std::vector<const Function *> ComputeStressBase::_initial_stress_fcn
protectedinherited

initial stress components

Definition at line 58 of file ComputeStressBase.h.

◆ _Jacobian_mult

MaterialProperty<RankFourTensor>& ComputeStressBase::_Jacobian_mult
protectedinherited

◆ _mechanical_strain

const MaterialProperty<RankTwoTensor>& ComputeStressBase::_mechanical_strain
protectedinherited

◆ _plastic_strain

MaterialProperty<RankTwoTensor>& FiniteStrainPlasticMaterial::_plastic_strain
protected

Definition at line 40 of file FiniteStrainPlasticMaterial.h.

Referenced by computeQpStress(), and initQpStatefulProperties().

◆ _plastic_strain_old

const MaterialProperty<RankTwoTensor>& FiniteStrainPlasticMaterial::_plastic_strain_old
protected

Definition at line 41 of file FiniteStrainPlasticMaterial.h.

Referenced by computeQpStress().

◆ _rotation_increment

const MaterialProperty<RankTwoTensor>& FiniteStrainPlasticMaterial::_rotation_increment
protected

Definition at line 46 of file FiniteStrainPlasticMaterial.h.

Referenced by computeQpStress().

◆ _rtol

Real FiniteStrainPlasticMaterial::_rtol
protected

Definition at line 51 of file FiniteStrainPlasticMaterial.h.

Referenced by returnMap().

◆ _strain_increment

const MaterialProperty<RankTwoTensor>& FiniteStrainPlasticMaterial::_strain_increment
protected

Definition at line 45 of file FiniteStrainPlasticMaterial.h.

Referenced by computeQpStress().

◆ _stress

MaterialProperty<RankTwoTensor>& ComputeStressBase::_stress
protectedinherited

Stress material property.

Definition at line 50 of file ComputeStressBase.h.

Referenced by ComputeMultipleInelasticCosseratStress::computeAdmissibleState(), ComputeMultipleInelasticStress::computeAdmissibleState(), ComputeStressBase::computeQpProperties(), ComputeStrainIncrementBasedStress::computeQpStress(), ComputeLinearElasticStress::computeQpStress(), ComputeDamageStress::computeQpStress(), ComputeFiniteStrainElasticStress::computeQpStress(), ComputeCosseratLinearElasticStress::computeQpStress(), ComputeSmearedCrackingStress::computeQpStress(), ComputeLinearElasticPFFractureStress::computeQpStress(), computeQpStress(), ComputeMultiPlasticityStress::computeQpStress(), ComputeLinearViscoelasticStress::computeQpStress(), ComputeMultipleInelasticStress::computeQpStress(), ComputeMultipleInelasticStress::computeQpStressIntermediateConfiguration(), ComputeLinearElasticPFFractureStress::computeStrainSpectral(), ComputeLinearElasticPFFractureStress::computeStrainVolDev(), ComputeLinearElasticPFFractureStress::computeStressSpectral(), ComputeMultipleInelasticStress::finiteStrainRotation(), ComputeStressBase::initQpStatefulProperties(), FiniteStrainCrystalPlasticity::initQpStatefulProperties(), FiniteStrainUObasedCP::initQpStatefulProperties(), FiniteStrainHyperElasticViscoPlastic::initQpStatefulProperties(), ComputeMultiPlasticityStress::postReturnMap(), FiniteStrainUObasedCP::postSolveQp(), FiniteStrainHyperElasticViscoPlastic::postSolveQp(), FiniteStrainCrystalPlasticity::postSolveQp(), ComputeSmearedCrackingStress::updateCrackingStateAndStress(), ComputeMultipleInelasticStress::updateQpState(), and ComputeMultipleInelasticStress::updateQpStateSingleModel().

◆ _stress_old

const MaterialProperty<RankTwoTensor>& FiniteStrainPlasticMaterial::_stress_old
protected

Definition at line 44 of file FiniteStrainPlasticMaterial.h.

Referenced by computeQpStress().

◆ _yield_stress_vector

std::vector<Real> FiniteStrainPlasticMaterial::_yield_stress_vector
protected

Definition at line 39 of file FiniteStrainPlasticMaterial.h.

Referenced by getdYieldStressdPlasticStrain(), and getYieldStress().


The documentation for this class was generated from the following files:
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
ComputeStressBase::_stress
MaterialProperty< RankTwoTensor > & _stress
Stress material property.
Definition: ComputeStressBase.h:50
ComputeStressBase::_extra_stress
const MaterialProperty< RankTwoTensor > & _extra_stress
Extra stress tensor.
Definition: ComputeStressBase.h:55
FiniteStrainPlasticMaterial::_elasticity_tensor_name
const std::string _elasticity_tensor_name
Name of the elasticity tensor material property.
Definition: FiniteStrainPlasticMaterial.h:48
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::_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::_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
ComputeStressBase::computeQpStress
virtual void computeQpStress()=0
Compute the stress and store it in the _stress material property for the current quadrature point.
FiniteStrainPlasticMaterial::yieldFunction
virtual Real yieldFunction(const RankTwoTensor &stress, const Real yield_stress)
Calculates the yield function.
Definition: FiniteStrainPlasticMaterial.C:312
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
ComputeStressBase::_base_name
const std::string _base_name
Base name prepended to all material property names to allow for multi-material systems.
Definition: ComputeStressBase.h:45
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::internalPotential
virtual Real internalPotential()
The internal potential.
Definition: FiniteStrainPlasticMaterial.C:338
RankFourTensorTempl
Definition: ACGrGrElasticDrivingForce.h:20
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::_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
ComputeStressBase::ComputeStressBase
ComputeStressBase(const InputParameters &parameters)
Definition: ComputeStressBase.C:28