FiniteStrainPlasticMaterial implements rate-independent associative J2 plasticity with isotropic hardening in the finite-strain framework. More...
#include <FiniteStrainPlasticMaterial.h>
Public Member Functions | |
FiniteStrainPlasticMaterial (const InputParameters ¶meters) | |
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... | |
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.
FiniteStrainPlasticMaterial::FiniteStrainPlasticMaterial | ( | const InputParameters & | parameters | ) |
Definition at line 34 of file FiniteStrainPlasticMaterial.C.
|
overrideprotectedvirtualinherited |
Definition at line 50 of file ComputeStressBase.C.
|
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.
|
protectedvirtual |
Derivative of yieldFunction with respect to the equivalent plastic strain.
Definition at line 326 of file FiniteStrainPlasticMaterial.C.
Referenced by returnMap().
|
protectedvirtual |
Derivative of yieldFunction with respect to the stress.
Definition at line 318 of file FiniteStrainPlasticMaterial.C.
Referenced by flowPotential(), getJac(), and returnMap().
|
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.
Referenced by getJac(), and returnMap().
|
protected |
d(yieldstress)/d(equivalent plastic strain)
Definition at line 416 of file FiniteStrainPlasticMaterial.C.
Referenced by dyieldFunction_dinternal().
|
protectedvirtual |
Evaluates the derivative d(resid_ij)/d(sig_kl), where resid_ij = flow_incr*flowPotential_ij - (E^{-1}(trial_stress - sig))_ij.
sig | stress |
E_ijkl | elasticity tensor (sig = E*(strain - plastic_strain)) |
flow_incr | consistency parameter |
dresid_dsig | the required derivative (this is an output variable) |
Definition at line 351 of file FiniteStrainPlasticMaterial.C.
Referenced by returnMap().
|
protected |
Equivalent stress.
Definition at line 344 of file FiniteStrainPlasticMaterial.C.
Referenced by getJac(), and yieldFunction().
|
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.
Referenced by returnMap().
|
protectedvirtual |
|
protectedvirtual |
The internal potential.
For associative J2 plasticity this is just -1
Definition at line 338 of file FiniteStrainPlasticMaterial.C.
Referenced by returnMap().
|
protectedvirtual |
Implements the return map.
Implements the return-map algorithm via a Newton-Raphson process.
sig_old | The stress at the previous "time" step |
eqvpstrain_old | The equivalent plastic strain at the previous "time" step |
plastic_strain_old | The value of plastic strain at the previous "time" step |
delta_d | The total strain increment for this "time" step |
E_ijkl | The elasticity tensor. If no plasiticity then sig_new = sig_old + E_ijkl*delta_d |
sig | The stress after returning to the yield surface (this is an output variable) |
eqvpstrain | The equivalent plastic strain after returning to the yield surface (this is an output variable) |
plastic_strain | The 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.
Referenced by computeQpStress().
|
static |
Definition at line 18 of file FiniteStrainPlasticMaterial.C.
|
protectedvirtual |
Calculates the yield function.
stress | the stress at which to calculate the yield function |
yield_stress | the current value of the yield stress |
Definition at line 312 of file FiniteStrainPlasticMaterial.C.
Referenced by returnMap().
|
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().
|
protected |
Definition at line 56 of file FiniteStrainPlasticMaterial.h.
Referenced by getJac().
|
protected |
Definition at line 56 of file FiniteStrainPlasticMaterial.h.
Referenced by getJac().
|
protectedinherited |
Elastic strain material property.
Definition at line 52 of file ComputeStressBase.h.
Referenced by ComputeSmearedCrackingStress::computeCrackStrainAndOrientation(), ComputeLinearElasticStress::computeQpStress(), ComputeFiniteStrainElasticStress::computeQpStress(), ComputeCosseratLinearElasticStress::computeQpStress(), ComputeSmearedCrackingStress::computeQpStress(), computeQpStress(), ComputeMultiPlasticityStress::computeQpStress(), ComputeLinearViscoelasticStress::computeQpStress(), ComputeMultipleInelasticStress::computeQpStressIntermediateConfiguration(), ComputeMultipleInelasticStress::finiteStrainRotation(), and ComputeStressBase::initQpStatefulProperties().
|
protected |
Elasticity tensor material property.
Definition at line 50 of file FiniteStrainPlasticMaterial.h.
Referenced by computeQpStress().
|
protected |
Name of the elasticity tensor material property.
Definition at line 48 of file FiniteStrainPlasticMaterial.h.
|
protected |
Definition at line 53 of file FiniteStrainPlasticMaterial.h.
Referenced by returnMap().
|
protected |
Definition at line 42 of file FiniteStrainPlasticMaterial.h.
Referenced by computeQpStress(), and initQpStatefulProperties().
|
protected |
Definition at line 43 of file FiniteStrainPlasticMaterial.h.
Referenced by computeQpStress().
|
protectedinherited |
Extra stress tensor.
Definition at line 55 of file ComputeStressBase.h.
Referenced by ComputeStressBase::computeQpProperties().
|
protected |
Definition at line 52 of file FiniteStrainPlasticMaterial.h.
Referenced by returnMap().
|
protectedinherited |
initial stress components
Definition at line 58 of file ComputeStressBase.h.
|
protectedinherited |
derivative of stress w.r.t. strain (_dstress_dstrain)
Definition at line 61 of file ComputeStressBase.h.
Referenced by ComputeStrainIncrementBasedStress::computeQpJacobian(), FiniteStrainHyperElasticViscoPlastic::computeQpJacobian(), ComputeMultipleInelasticCosseratStress::computeQpJacobianMult(), ComputeMultipleInelasticStress::computeQpJacobianMult(), ComputeLinearElasticStress::computeQpStress(), ComputeDamageStress::computeQpStress(), ComputeFiniteStrainElasticStress::computeQpStress(), ComputeCosseratLinearElasticStress::computeQpStress(), ComputeSmearedCrackingStress::computeQpStress(), computeQpStress(), ComputeLinearElasticPFFractureStress::computeQpStress(), ComputeMultiPlasticityStress::computeQpStress(), ComputeLinearViscoelasticStress::computeQpStress(), ComputeMultipleInelasticStress::computeQpStress(), ComputeMultipleInelasticStress::computeQpStressIntermediateConfiguration(), ComputeLinearElasticPFFractureStress::computeStrainSpectral(), ComputeLinearElasticPFFractureStress::computeStrainVolDev(), ComputeLinearElasticPFFractureStress::computeStressSpectral(), FiniteStrainUObasedCP::elasticTangentModuli(), FiniteStrainUObasedCP::elastoPlasticTangentModuli(), ComputeMultipleInelasticStress::finiteStrainRotation(), ComputeMultiPlasticityStress::postReturnMap(), FiniteStrainCrystalPlasticity::postSolveQp(), FiniteStrainCrystalPlasticity::preSolveQp(), and ComputeMultipleInelasticStress::updateQpStateSingleModel().
|
protectedinherited |
Mechanical strain material property.
Definition at line 48 of file ComputeStressBase.h.
Referenced by ComputeLinearElasticStress::computeQpStress(), ComputeFiniteStrainElasticStress::computeQpStress(), ComputeCosseratLinearElasticStress::computeQpStress(), computeQpStress(), ComputeLinearElasticPFFractureStress::computeQpStress(), ComputeLinearViscoelasticStress::computeQpStress(), ComputeLinearElasticPFFractureStress::computeStrainSpectral(), ComputeLinearElasticPFFractureStress::computeStrainVolDev(), and ComputeLinearElasticPFFractureStress::computeStressSpectral().
|
protected |
Definition at line 40 of file FiniteStrainPlasticMaterial.h.
Referenced by computeQpStress(), and initQpStatefulProperties().
|
protected |
Definition at line 41 of file FiniteStrainPlasticMaterial.h.
Referenced by computeQpStress().
|
protected |
Definition at line 46 of file FiniteStrainPlasticMaterial.h.
Referenced by computeQpStress().
|
protected |
Definition at line 51 of file FiniteStrainPlasticMaterial.h.
Referenced by returnMap().
|
protected |
Definition at line 45 of file FiniteStrainPlasticMaterial.h.
Referenced by computeQpStress().
|
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().
|
protected |
Definition at line 44 of file FiniteStrainPlasticMaterial.h.
Referenced by computeQpStress().
|
protected |
Definition at line 39 of file FiniteStrainPlasticMaterial.h.
Referenced by getdYieldStressdPlasticStrain(), and getYieldStress().