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

FiniteStrainUObasedCP uses the multiplicative decomposition of deformation gradient and solves the PK2 stress residual equation at the intermediate configuration to evolve the material state. More...

#include <FiniteStrainUObasedCP.h>

Inheritance diagram for FiniteStrainUObasedCP:
[legend]

Public Member Functions

 FiniteStrainUObasedCP (const InputParameters &parameters)
 

Static Public Member Functions

static InputParameters validParams ()
 

Protected Member Functions

virtual void computeQpStress ()
 updates the stress at a quadrature point. More...
 
virtual void initQpStatefulProperties ()
 initializes the stateful properties such as stress, plastic deformation gradient, slip system resistances, etc. More...
 
virtual void calcResidJacob ()
 calls the residual and jacobian functions used in the stress update algorithm. More...
 
virtual void updateSlipSystemResistanceAndStateVariable ()
 updates the slip system resistances and state variables. More...
 
virtual void preSolveQp ()
 set variables for stress and internal variable solve. More...
 
virtual void solveQp ()
 solve stress and internal variables. More...
 
virtual void postSolveQp ()
 update stress and internal variable after solve. More...
 
virtual void preSolveStatevar ()
 set variables for internal variable solve. More...
 
virtual void solveStatevar ()
 solve internal variables. More...
 
virtual void postSolveStatevar ()
 update internal variable after solve. More...
 
virtual void preSolveStress ()
 set variables for stress solve. More...
 
virtual void solveStress ()
 solves for stress, updates plastic deformation gradient. More...
 
virtual void postSolveStress ()
 update stress and plastic deformation gradient after solve. More...
 
virtual void calcResidual ()
 calculate stress residual. More...
 
virtual void calcJacobian ()
 calculate jacobian. More...
 
virtual void getSlipRates ()
 updates the slip rates. More...
 
virtual void calcTangentModuli ()
 calculate the tangent moduli for preconditioner. More...
 
virtual void elasticTangentModuli ()
 calculate the elastic tangent moduli for preconditioner. More...
 
virtual void elastoPlasticTangentModuli ()
 calculate the exact tangent moduli for preconditioner. More...
 
bool lineSearchUpdate (const Real rnorm_prev, const RankTwoTensor)
 performs the line search update More...
 
virtual bool isStateVariablesConverged ()
 evaluates convergence of state variables. More...
 
virtual void computeQpProperties () override
 

Protected Attributes

std::vector< const CrystalPlasticitySlipRate * > _uo_slip_rates
 User objects that define the slip rate. More...
 
std::vector< const CrystalPlasticitySlipResistance * > _uo_slip_resistances
 User objects that define the slip resistance. More...
 
std::vector< const CrystalPlasticityStateVariable * > _uo_state_vars
 User objects that define the state variable. More...
 
std::vector< const CrystalPlasticityStateVarRateComponent * > _uo_state_var_evol_rate_comps
 User objects that define the state variable evolution rate component. More...
 
std::vector< MaterialProperty< std::vector< Real > > * > _mat_prop_slip_rates
 Slip rates material property. More...
 
std::vector< MaterialProperty< std::vector< Real > > * > _mat_prop_slip_resistances
 Slip resistance material property. More...
 
std::vector< MaterialProperty< std::vector< Real > > * > _mat_prop_state_vars
 State variable material property. More...
 
std::vector< const MaterialProperty< std::vector< Real > > * > _mat_prop_state_vars_old
 Old state variable material property. More...
 
std::vector< MaterialProperty< std::vector< Real > > * > _mat_prop_state_var_evol_rate_comps
 State variable evolution rate component material property. More...
 
unsigned int _num_uo_slip_rates
 Number of slip rate user objects. More...
 
unsigned int _num_uo_slip_resistances
 Number of slip resistance user objects. More...
 
unsigned int _num_uo_state_vars
 Number of state variable user objects. More...
 
unsigned int _num_uo_state_var_evol_rate_comps
 Number of state variable evolution rate component user objects. More...
 
std::vector< std::vector< Real > > _state_vars_old
 Local state variable. More...
 
std::vector< std::vector< Real > > _state_vars_old_stored
 Local stored state variable (for sub-stepping) More...
 
std::vector< std::vector< Real > > _state_vars_prev
 Local old state variable. More...
 
Real _rtol
 Stress residual equation relative tolerance. More...
 
Real _abs_tol
 Stress residual equation absolute tolerance. More...
 
Real _stol
 Internal variable update equation tolerance. More...
 
Real _zero_tol
 Residual tolerance when variable value is zero. Default 1e-12. More...
 
RankTwoTensor _resid
 Residual tensor. More...
 
RankFourTensor _jac
 Jacobian tensor. More...
 
unsigned int _maxiter
 Maximum number of iterations for stress update. More...
 
unsigned int _maxiterg
 Maximum number of iterations for internal variable update. More...
 
MooseEnum _tan_mod_type
 Type of tangent moduli calculation. More...
 
unsigned int _max_substep_iter
 Maximum number of substep iterations. More...
 
bool _use_line_search
 Flag to activate line serach. More...
 
Real _min_lsrch_step
 Minimum line search step size. More...
 
Real _lsrch_tol
 Line search bisection method tolerance. More...
 
unsigned int _lsrch_max_iter
 Line search bisection method maximum iteration number. More...
 
MooseEnum _lsrch_method
 Line search method. More...
 
MaterialProperty< RankTwoTensor > & _fp
 
const MaterialProperty< RankTwoTensor > & _fp_old
 
MaterialProperty< RankTwoTensor > & _pk2
 
const MaterialProperty< RankTwoTensor > & _pk2_old
 
MaterialProperty< RankTwoTensor > & _lag_e
 
MaterialProperty< RankTwoTensor > & _update_rot
 
const MaterialProperty< RankTwoTensor > & _update_rot_old
 
const std::string _elasticity_tensor_name
 Name of the elasticity tensor material property. More...
 
const MaterialProperty< RankFourTensor > & _elasticity_tensor
 Elasticity tensor material property. More...
 
const MaterialProperty< RankTwoTensor > & _deformation_gradient
 
const MaterialProperty< RankTwoTensor > & _deformation_gradient_old
 
const MaterialProperty< RankTwoTensor > & _crysrot
 Crystal rotation. More...
 
RankTwoTensor _dfgrd_tmp
 
RankTwoTensor _fe
 
RankTwoTensor _fp_old_inv
 
RankTwoTensor _fp_inv
 
DenseVector< Real > _tau
 
std::vector< MaterialProperty< std::vector< RankTwoTensor > > * > _flow_direction
 
bool _err_tol
 Flag to check whether convergence is achieved. More...
 
RankTwoTensor _delta_dfgrd
 Used for substepping; Uniformly divides the increment in deformation gradient. More...
 
RankTwoTensor _dfgrd_tmp_old
 
Real _dfgrd_scale_factor
 Scales the substepping increment to obtain deformation gradient at a substep iteration. More...
 
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

FiniteStrainUObasedCP uses the multiplicative decomposition of deformation gradient and solves the PK2 stress residual equation at the intermediate configuration to evolve the material state.

The internal variables are updated using an interative predictor-corrector algorithm. Backward Euler integration rule is used for the rate equations.

Involves 4 different types of user objects that calculates: State variables - update state variable (derive from CrystalPlasticityStateVariable) State variable evolution compoment - individual component of state variable incremental rate (derive from CrystalPlasticityStateVariableEvolutionRateComponent) Slip resistance - calcuate slip resistances (derive from CrystalPlasticitySlipResistances) Slip rates - calcuate flow direction and slip rates (derive from CrystalPlasticitySlipRates)

Definition at line 38 of file FiniteStrainUObasedCP.h.

Constructor & Destructor Documentation

◆ FiniteStrainUObasedCP()

FiniteStrainUObasedCP::FiniteStrainUObasedCP ( const InputParameters &  parameters)

Definition at line 66 of file FiniteStrainUObasedCP.C.

67  : ComputeStressBase(parameters),
68  _num_uo_slip_rates(parameters.get<std::vector<UserObjectName>>("uo_slip_rates").size()),
70  parameters.get<std::vector<UserObjectName>>("uo_slip_resistances").size()),
71  _num_uo_state_vars(parameters.get<std::vector<UserObjectName>>("uo_state_vars").size()),
73  parameters.get<std::vector<UserObjectName>>("uo_state_var_evol_rate_comps").size()),
74  _rtol(getParam<Real>("rtol")),
75  _abs_tol(getParam<Real>("abs_tol")),
76  _stol(getParam<Real>("stol")),
77  _zero_tol(getParam<Real>("zero_tol")),
78  _maxiter(getParam<unsigned int>("maxiter")),
79  _maxiterg(getParam<unsigned int>("maxiter_state_variable")),
80  _tan_mod_type(getParam<MooseEnum>("tan_mod_type")),
81  _max_substep_iter(getParam<unsigned int>("maximum_substep_iteration")),
82  _use_line_search(getParam<bool>("use_line_search")),
83  _min_lsrch_step(getParam<Real>("min_line_search_step_size")),
84  _lsrch_tol(getParam<Real>("line_search_tol")),
85  _lsrch_max_iter(getParam<unsigned int>("line_search_maxiter")),
86  _lsrch_method(getParam<MooseEnum>("line_search_method")),
87  _fp(declareProperty<RankTwoTensor>("fp")), // Plastic deformation gradient
88  _fp_old(getMaterialPropertyOld<RankTwoTensor>(
89  "fp")), // Plastic deformation gradient of previous increment
90  _pk2(declareProperty<RankTwoTensor>("pk2")), // 2nd Piola-Kirchoff Stress
91  _pk2_old(getMaterialPropertyOld<RankTwoTensor>(
92  "pk2")), // 2nd Piola Kirchoff Stress of previous increment
93  _lag_e(declareProperty<RankTwoTensor>("lage")), // Lagrangian strain
94  _update_rot(declareProperty<RankTwoTensor>(
95  "update_rot")), // Rotation tensor considering material rotation and crystal orientation
96  _update_rot_old(getMaterialPropertyOld<RankTwoTensor>("update_rot")),
97  _elasticity_tensor_name(_base_name + "elasticity_tensor"),
98  _elasticity_tensor(getMaterialPropertyByName<RankFourTensor>(_elasticity_tensor_name)),
99  _deformation_gradient(getMaterialProperty<RankTwoTensor>("deformation_gradient")),
100  _deformation_gradient_old(getMaterialPropertyOld<RankTwoTensor>("deformation_gradient")),
101  _crysrot(getMaterialProperty<RankTwoTensor>("crysrot"))
102 {
103  _err_tol = false;
104 
105  _delta_dfgrd.zero();
106 
107  // resize the material properties for each userobject
113 
114  // resize the flow direction
116 
117  // resize local state variables
121 
122  // resize user objects
127 
128  // assign the user objects
129  for (unsigned int i = 0; i < _num_uo_slip_rates; ++i)
130  {
131  _uo_slip_rates[i] = &getUserObjectByName<CrystalPlasticitySlipRate>(
132  parameters.get<std::vector<UserObjectName>>("uo_slip_rates")[i]);
133  _mat_prop_slip_rates[i] = &declareProperty<std::vector<Real>>(
134  parameters.get<std::vector<UserObjectName>>("uo_slip_rates")[i]);
135  _flow_direction[i] = &declareProperty<std::vector<RankTwoTensor>>(
136  parameters.get<std::vector<UserObjectName>>("uo_slip_rates")[i] + "_flow_direction");
137  }
138 
139  for (unsigned int i = 0; i < _num_uo_slip_resistances; ++i)
140  {
141  _uo_slip_resistances[i] = &getUserObjectByName<CrystalPlasticitySlipResistance>(
142  parameters.get<std::vector<UserObjectName>>("uo_slip_resistances")[i]);
143  _mat_prop_slip_resistances[i] = &declareProperty<std::vector<Real>>(
144  parameters.get<std::vector<UserObjectName>>("uo_slip_resistances")[i]);
145  }
146 
147  for (unsigned int i = 0; i < _num_uo_state_vars; ++i)
148  {
149  _uo_state_vars[i] = &getUserObjectByName<CrystalPlasticityStateVariable>(
150  parameters.get<std::vector<UserObjectName>>("uo_state_vars")[i]);
151  _mat_prop_state_vars[i] = &declareProperty<std::vector<Real>>(
152  parameters.get<std::vector<UserObjectName>>("uo_state_vars")[i]);
153  _mat_prop_state_vars_old[i] = &getMaterialPropertyOld<std::vector<Real>>(
154  parameters.get<std::vector<UserObjectName>>("uo_state_vars")[i]);
155  }
156 
157  for (unsigned int i = 0; i < _num_uo_state_var_evol_rate_comps; ++i)
158  {
159  _uo_state_var_evol_rate_comps[i] = &getUserObjectByName<CrystalPlasticityStateVarRateComponent>(
160  parameters.get<std::vector<UserObjectName>>("uo_state_var_evol_rate_comps")[i]);
161  _mat_prop_state_var_evol_rate_comps[i] = &declareProperty<std::vector<Real>>(
162  parameters.get<std::vector<UserObjectName>>("uo_state_var_evol_rate_comps")[i]);
163  }
164 }

Member Function Documentation

◆ calcJacobian()

void FiniteStrainUObasedCP::calcJacobian ( )
protectedvirtual

calculate jacobian.

Definition at line 531 of file FiniteStrainUObasedCP.C.

532 {
533  RankFourTensor dfedfpinv, deedfe, dfpinvdpk2;
534 
535  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
536  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
537  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
538  dfedfpinv(i, j, k, j) = _dfgrd_tmp(i, k);
539 
540  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
541  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
542  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
543  {
544  deedfe(i, j, k, i) = deedfe(i, j, k, i) + _fe(k, j) * 0.5;
545  deedfe(i, j, k, j) = deedfe(i, j, k, j) + _fe(k, i) * 0.5;
546  }
547 
548  for (unsigned int i = 0; i < _num_uo_slip_rates; ++i)
549  {
550  unsigned int nss = _uo_slip_rates[i]->variableSize();
551  std::vector<RankTwoTensor> dtaudpk2(nss), dfpinvdslip(nss);
552  std::vector<Real> dslipdtau;
553  dslipdtau.resize(nss);
554  _uo_slip_rates[i]->calcSlipRateDerivative(_qp, _dt, dslipdtau);
555  for (unsigned int j = 0; j < nss; j++)
556  {
557  dtaudpk2[j] = (*_flow_direction[i])[_qp][j];
558  dfpinvdslip[j] = -_fp_old_inv * (*_flow_direction[i])[_qp][j];
559  }
560 
561  for (unsigned int j = 0; j < nss; j++)
562  dfpinvdpk2 += (dfpinvdslip[j] * dslipdtau[j] * _dt).outerProduct(dtaudpk2[j]);
563  }
564  _jac =
565  RankFourTensor::IdentityFour() - (_elasticity_tensor[_qp] * deedfe * dfedfpinv * dfpinvdpk2);
566 }

Referenced by calcResidJacob().

◆ calcResidJacob()

void FiniteStrainUObasedCP::calcResidJacob ( )
protectedvirtual

calls the residual and jacobian functions used in the stress update algorithm.

Definition at line 483 of file FiniteStrainUObasedCP.C.

484 {
485  calcResidual();
486  if (_err_tol)
487  return;
488  calcJacobian();
489 }

Referenced by solveStress().

◆ calcResidual()

void FiniteStrainUObasedCP::calcResidual ( )
protectedvirtual

calculate stress residual.

Definition at line 505 of file FiniteStrainUObasedCP.C.

506 {
507  RankTwoTensor iden(RankTwoTensor::initIdentity), ce, ee, ce_pk2, eqv_slip_incr, pk2_new;
508 
509  getSlipRates();
510  if (_err_tol)
511  return;
512 
513  for (unsigned int i = 0; i < _num_uo_slip_rates; ++i)
514  for (unsigned int j = 0; j < _uo_slip_rates[i]->variableSize(); ++j)
515  eqv_slip_incr += (*_flow_direction[i])[_qp][j] * (*_mat_prop_slip_rates[i])[_qp][j] * _dt;
516 
517  eqv_slip_incr = iden - eqv_slip_incr;
518  _fp_inv = _fp_old_inv * eqv_slip_incr;
519  _fe = _dfgrd_tmp * _fp_inv;
520 
521  ce = _fe.transpose() * _fe;
522  ee = ce - iden;
523  ee *= 0.5;
524 
525  pk2_new = _elasticity_tensor[_qp] * ee;
526 
527  _resid = _pk2[_qp] - pk2_new;
528 }

Referenced by calcResidJacob(), and lineSearchUpdate().

◆ calcTangentModuli()

void FiniteStrainUObasedCP::calcTangentModuli ( )
protectedvirtual

calculate the tangent moduli for preconditioner.

Default is the elastic stiffness matrix. Exact jacobian is currently implemented. tan_mod_type can be modified to exact in .i file to turn it on.

Definition at line 569 of file FiniteStrainUObasedCP.C.

570 {
571  switch (_tan_mod_type)
572  {
573  case 0:
575  break;
576  default:
578  }
579 }

Referenced by postSolveQp().

◆ 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 FiniteStrainUObasedCP::computeQpStress ( )
protectedvirtual

updates the stress at a quadrature point.

Solves stress residual equation using NR.

Updates slip system resistances iteratively.

Implements ComputeStressBase.

Definition at line 207 of file FiniteStrainUObasedCP.C.

208 {
209  // Userobject based crystal plasticity does not support face/boundary material property
210  // calculation.
211  if (isBoundaryMaterial())
212  return;
213  // Depth of substepping; Limited to maximum substep iteration
214  unsigned int substep_iter = 1;
215  // Calculated from substep_iter as 2^substep_iter
216  unsigned int num_substep = 1;
217  // Store original _dt; Reset at the end of solve
218  Real dt_original = _dt;
219 
221  if (_dfgrd_tmp_old.det() == 0)
222  _dfgrd_tmp_old.addIa(1.0);
223 
225 
226  // Saves the old stateful properties that is modified during sub stepping
227  for (unsigned int i = 0; i < _num_uo_state_vars; ++i)
229 
230  for (unsigned int i = 0; i < _num_uo_slip_rates; ++i)
231  _uo_slip_rates[i]->calcFlowDirection(_qp, (*_flow_direction[i])[_qp]);
232 
233  do
234  {
235  _err_tol = false;
236 
237  preSolveQp();
238 
239  _dt = dt_original / num_substep;
240 
241  for (unsigned int istep = 0; istep < num_substep; ++istep)
242  {
243  _dfgrd_tmp = (static_cast<Real>(istep) + 1) / num_substep * _delta_dfgrd + _dfgrd_tmp_old;
244 
245  solveQp();
246 
247  if (_err_tol)
248  {
249  substep_iter++;
250  num_substep *= 2;
251  break;
252  }
253  }
254  if (substep_iter > _max_substep_iter && _err_tol)
255  throw MooseException("FiniteStrainUObasedCP: Constitutive failure.");
256  } while (_err_tol);
257 
258  _dt = dt_original;
259 
260  postSolveQp();
261 }

◆ elasticTangentModuli()

void FiniteStrainUObasedCP::elasticTangentModuli ( )
protectedvirtual

calculate the elastic tangent moduli for preconditioner.

Definition at line 625 of file FiniteStrainUObasedCP.C.

626 {
627  // update jacobian_mult
629 }

Referenced by calcTangentModuli().

◆ elastoPlasticTangentModuli()

void FiniteStrainUObasedCP::elastoPlasticTangentModuli ( )
protectedvirtual

calculate the exact tangent moduli for preconditioner.

Definition at line 582 of file FiniteStrainUObasedCP.C.

583 {
584  RankFourTensor tan_mod;
585  RankTwoTensor pk2fet, fepk2;
586  RankFourTensor deedfe, dsigdpk2dfe, dfedf;
587 
588  // Fill in the matrix stiffness material property
589  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
590  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
591  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
592  {
593  deedfe(i, j, k, i) = deedfe(i, j, k, i) + _fe(k, j) * 0.5;
594  deedfe(i, j, k, j) = deedfe(i, j, k, j) + _fe(k, i) * 0.5;
595  }
596 
597  dsigdpk2dfe = _fe.mixedProductIkJl(_fe) * _elasticity_tensor[_qp] * deedfe;
598 
599  pk2fet = _pk2[_qp] * _fe.transpose();
600  fepk2 = _fe * _pk2[_qp];
601 
602  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
603  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
604  for (unsigned int l = 0; l < LIBMESH_DIM; ++l)
605  {
606  tan_mod(i, j, i, l) += pk2fet(l, j);
607  tan_mod(i, j, j, l) += fepk2(i, l);
608  }
609 
610  tan_mod += dsigdpk2dfe;
611 
612  Real je = _fe.det();
613  if (je > 0.0)
614  tan_mod /= je;
615 
616  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
617  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
618  for (unsigned int l = 0; l < LIBMESH_DIM; ++l)
619  dfedf(i, j, i, l) = _fp_inv(l, j);
620 
621  _Jacobian_mult[_qp] = tan_mod * dfedf;
622 }

Referenced by calcTangentModuli().

◆ getSlipRates()

void FiniteStrainUObasedCP::getSlipRates ( )
protectedvirtual

updates the slip rates.

Definition at line 492 of file FiniteStrainUObasedCP.C.

493 {
494  for (unsigned int i = 0; i < _num_uo_slip_rates; ++i)
495  {
496  if (!_uo_slip_rates[i]->calcSlipRate(_qp, _dt, (*_mat_prop_slip_rates[i])[_qp]))
497  {
498  _err_tol = true;
499  return;
500  }
501  }
502 }

Referenced by calcResidual().

◆ initQpStatefulProperties()

void FiniteStrainUObasedCP::initQpStatefulProperties ( )
protectedvirtual

initializes the stateful properties such as stress, plastic deformation gradient, slip system resistances, etc.

Reimplemented from ComputeStressBase.

Definition at line 167 of file FiniteStrainUObasedCP.C.

168 {
169  for (unsigned int i = 0; i < _num_uo_slip_rates; ++i)
170  {
171  (*_mat_prop_slip_rates[i])[_qp].resize(_uo_slip_rates[i]->variableSize());
172  (*_flow_direction[i])[_qp].resize(_uo_slip_rates[i]->variableSize());
173  }
174 
175  for (unsigned int i = 0; i < _num_uo_slip_resistances; ++i)
176  (*_mat_prop_slip_resistances[i])[_qp].resize(_uo_slip_resistances[i]->variableSize());
177 
178  for (unsigned int i = 0; i < _num_uo_state_vars; ++i)
179  {
180  (*_mat_prop_state_vars[i])[_qp].resize(_uo_state_vars[i]->variableSize());
181  _state_vars_old[i].resize(_uo_state_vars[i]->variableSize());
182  _state_vars_old_stored[i].resize(_uo_state_vars[i]->variableSize());
183  _state_vars_prev[i].resize(_uo_state_vars[i]->variableSize());
184  }
185 
186  for (unsigned int i = 0; i < _num_uo_state_var_evol_rate_comps; ++i)
187  (*_mat_prop_state_var_evol_rate_comps[i])[_qp].resize(
188  _uo_state_var_evol_rate_comps[i]->variableSize());
189 
190  _stress[_qp].zero();
191  _pk2[_qp].zero();
192  _lag_e[_qp].zero();
193 
194  _fp[_qp].setToIdentity();
195  _update_rot[_qp].setToIdentity();
196 
197  for (unsigned int i = 0; i < _num_uo_state_vars; ++i)
198  // Initializes slip system related properties
199  _uo_state_vars[i]->initSlipSysProps((*_mat_prop_state_vars[i])[_qp], _q_point[_qp]);
200 }

◆ isStateVariablesConverged()

bool FiniteStrainUObasedCP::isStateVariablesConverged ( )
protectedvirtual

evaluates convergence of state variables.

Definition at line 350 of file FiniteStrainUObasedCP.C.

351 {
352  Real diff;
353 
354  for (unsigned int i = 0; i < _num_uo_state_vars; ++i)
355  {
356  unsigned int n = (*_mat_prop_state_vars[i])[_qp].size();
357  for (unsigned j = 0; j < n; j++)
358  {
359  diff = std::abs((*_mat_prop_state_vars[i])[_qp][j] -
360  _state_vars_prev[i][j]); // Calculate increment size
361  if (std::abs(_state_vars_old_stored[i][j]) < _zero_tol && diff > _zero_tol)
362  return true;
363  if (std::abs(_state_vars_old_stored[i][j]) > _zero_tol &&
364  diff > _stol * std::abs(_state_vars_old_stored[i][j]))
365  return true;
366  }
367  }
368  return false;
369 }

Referenced by solveStatevar().

◆ lineSearchUpdate()

bool FiniteStrainUObasedCP::lineSearchUpdate ( const Real  rnorm_prev,
const RankTwoTensor  dpk2 
)
protected

performs the line search update

Definition at line 632 of file FiniteStrainUObasedCP.C.

633 {
634  switch (_lsrch_method)
635  {
636  case 0: // CUT_HALF
637  {
638  Real rnorm;
639  Real step = 1.0;
640 
641  do
642  {
643  _pk2[_qp] = _pk2[_qp] - step * dpk2;
644  step /= 2.0;
645  _pk2[_qp] = _pk2[_qp] + step * dpk2;
646 
647  calcResidual();
648  rnorm = _resid.L2norm();
649  } while (rnorm > rnorm_prev && step > _min_lsrch_step);
650 
651  // has norm improved or is the step still above minumum search step size?
652  return (rnorm <= rnorm_prev || step > _min_lsrch_step);
653  }
654 
655  case 1: // BISECTION
656  {
657  unsigned int count = 0;
658  Real step_a = 0.0;
659  Real step_b = 1.0;
660  Real step = 1.0;
661  Real s_m = 1000.0;
662  Real rnorm = 1000.0;
663 
664  calcResidual();
665  Real s_b = _resid.doubleContraction(dpk2);
666  Real rnorm1 = _resid.L2norm();
667  _pk2[_qp] = _pk2[_qp] - dpk2;
668  calcResidual();
669  Real s_a = _resid.doubleContraction(dpk2);
670  Real rnorm0 = _resid.L2norm();
671  _pk2[_qp] = _pk2[_qp] + dpk2;
672 
673  if ((rnorm1 / rnorm0) < _lsrch_tol || s_a * s_b > 0)
674  {
675  calcResidual();
676  return true;
677  }
678 
679  while ((rnorm / rnorm0) > _lsrch_tol && count < _lsrch_max_iter)
680  {
681  _pk2[_qp] = _pk2[_qp] - step * dpk2;
682  step = 0.5 * (step_b + step_a);
683  _pk2[_qp] = _pk2[_qp] + step * dpk2;
684  calcResidual();
685  s_m = _resid.doubleContraction(dpk2);
686  rnorm = _resid.L2norm();
687 
688  if (s_m * s_a < 0.0)
689  {
690  step_b = step;
691  s_b = s_m;
692  }
693  if (s_m * s_b < 0.0)
694  {
695  step_a = step;
696  s_a = s_m;
697  }
698  count++;
699  }
700 
701  // below tolerance and max iterations?
702  return ((rnorm / rnorm0) < _lsrch_tol && count < _lsrch_max_iter);
703  }
704 
705  default:
706  mooseError("Line search method is not provided.");
707  }
708 }

Referenced by solveStress().

◆ postSolveQp()

void FiniteStrainUObasedCP::postSolveQp ( )
protectedvirtual

update stress and internal variable after solve.

Definition at line 284 of file FiniteStrainUObasedCP.C.

285 {
286  _stress[_qp] = _fe * _pk2[_qp] * _fe.transpose() / _fe.det();
287 
288  // Calculate jacobian for preconditioner
290 
291  RankTwoTensor iden(RankTwoTensor::initIdentity);
292 
293  _lag_e[_qp] = _deformation_gradient[_qp].transpose() * _deformation_gradient[_qp] - iden;
294  _lag_e[_qp] = _lag_e[_qp] * 0.5;
295 
296  RankTwoTensor rot;
297  // Calculate material rotation
298  _deformation_gradient[_qp].getRUDecompositionRotation(rot);
299  _update_rot[_qp] = rot * _crysrot[_qp];
300 }

Referenced by computeQpStress().

◆ postSolveStatevar()

void FiniteStrainUObasedCP::postSolveStatevar ( )
protectedvirtual

update internal variable after solve.

Definition at line 372 of file FiniteStrainUObasedCP.C.

373 {
374  for (unsigned int i = 0; i < _num_uo_state_vars; ++i)
376 
378 }

Referenced by solveQp().

◆ postSolveStress()

void FiniteStrainUObasedCP::postSolveStress ( )
protectedvirtual

update stress and plastic deformation gradient after solve.

Definition at line 455 of file FiniteStrainUObasedCP.C.

456 {
457  _fp[_qp] = _fp_inv.inverse();
458 }

Referenced by solveStatevar().

◆ preSolveQp()

void FiniteStrainUObasedCP::preSolveQp ( )
protectedvirtual

set variables for stress and internal variable solve.

Definition at line 264 of file FiniteStrainUObasedCP.C.

265 {
266  for (unsigned int i = 0; i < _num_uo_state_vars; ++i)
268 
269  _pk2[_qp] = _pk2_old[_qp];
270  _fp_old_inv = _fp_old[_qp].inverse();
271 }

Referenced by computeQpStress().

◆ preSolveStatevar()

void FiniteStrainUObasedCP::preSolveStatevar ( )
protectedvirtual

set variables for internal variable solve.

Definition at line 303 of file FiniteStrainUObasedCP.C.

304 {
305  for (unsigned int i = 0; i < _num_uo_state_vars; ++i)
307 
308  for (unsigned int i = 0; i < _num_uo_slip_resistances; ++i)
309  _uo_slip_resistances[i]->calcSlipResistance(_qp, (*_mat_prop_slip_resistances[i])[_qp]);
310 
312 }

Referenced by solveQp().

◆ preSolveStress()

void FiniteStrainUObasedCP::preSolveStress ( )
protectedvirtual

set variables for stress solve.

Definition at line 381 of file FiniteStrainUObasedCP.C.

382 {
383 }

Referenced by solveStatevar().

◆ solveQp()

void FiniteStrainUObasedCP::solveQp ( )
protectedvirtual

solve stress and internal variables.

Definition at line 274 of file FiniteStrainUObasedCP.C.

275 {
277  solveStatevar();
278  if (_err_tol)
279  return;
281 }

Referenced by computeQpStress().

◆ solveStatevar()

void FiniteStrainUObasedCP::solveStatevar ( )
protectedvirtual

solve internal variables.

Definition at line 315 of file FiniteStrainUObasedCP.C.

316 {
317  unsigned int iterg;
318  bool iter_flag = true;
319 
320  iterg = 0;
321  // Check for slip system resistance update tolerance
322  while (iter_flag && iterg < _maxiterg)
323  {
324  preSolveStress();
325  solveStress();
326  if (_err_tol)
327  return;
328  postSolveStress();
329 
330  // Update slip system resistance and state variable
332 
333  if (_err_tol)
334  return;
335 
336  iter_flag = isStateVariablesConverged();
337  iterg++;
338  }
339 
340  if (iterg == _maxiterg)
341  {
342 #ifdef DEBUG
343  mooseWarning("FiniteStrainUObasedCP: Hardness Integration error\n");
344 #endif
345  _err_tol = true;
346  }
347 }

Referenced by solveQp().

◆ solveStress()

void FiniteStrainUObasedCP::solveStress ( )
protectedvirtual

solves for stress, updates plastic deformation gradient.

Definition at line 386 of file FiniteStrainUObasedCP.C.

387 {
388  unsigned int iter = 0;
389  RankTwoTensor dpk2;
390  Real rnorm, rnorm0, rnorm_prev;
391 
392  // Calculate stress residual
393  calcResidJacob();
394  if (_err_tol)
395  {
396 #ifdef DEBUG
397  mooseWarning("FiniteStrainUObasedCP: Slip increment exceeds tolerance - Element number ",
398  _current_elem->id(),
399  " Gauss point = ",
400  _qp);
401 #endif
402  return;
403  }
404 
405  rnorm = _resid.L2norm();
406  rnorm0 = rnorm;
407 
408  // Check for stress residual tolerance
409  while (rnorm > _rtol * rnorm0 && rnorm0 > _abs_tol && iter < _maxiter)
410  {
411  // Calculate stress increment
412  dpk2 = -_jac.invSymm() * _resid;
413  _pk2[_qp] = _pk2[_qp] + dpk2;
414  calcResidJacob();
415 
416  if (_err_tol)
417  {
418 #ifdef DEBUG
419  mooseWarning("FiniteStrainUObasedCP: Slip increment exceeds tolerance - Element number ",
420  _current_elem->id(),
421  " Gauss point = ",
422  _qp);
423 #endif
424  return;
425  }
426 
427  rnorm_prev = rnorm;
428  rnorm = _resid.L2norm();
429 
430  if (_use_line_search && rnorm > rnorm_prev && !lineSearchUpdate(rnorm_prev, dpk2))
431  {
432 #ifdef DEBUG
433  mooseWarning("FiniteStrainUObasedCP: Failed with line search");
434 #endif
435  _err_tol = true;
436  return;
437  }
438 
439  if (_use_line_search)
440  rnorm = _resid.L2norm();
441 
442  iter++;
443  }
444 
445  if (iter >= _maxiter)
446  {
447 #ifdef DEBUG
448  mooseWarning("FiniteStrainUObasedCP: Stress Integration error rmax = ", rnorm);
449 #endif
450  _err_tol = true;
451  }
452 }

Referenced by solveStatevar().

◆ updateSlipSystemResistanceAndStateVariable()

void FiniteStrainUObasedCP::updateSlipSystemResistanceAndStateVariable ( )
protectedvirtual

updates the slip system resistances and state variables.

override to modify slip system resistance and state variable evolution.

Definition at line 461 of file FiniteStrainUObasedCP.C.

462 {
463  for (unsigned int i = 0; i < _num_uo_state_vars; ++i)
464  _state_vars_prev[i] = (*_mat_prop_state_vars[i])[_qp];
465 
466  for (unsigned int i = 0; i < _num_uo_state_var_evol_rate_comps; ++i)
467  _uo_state_var_evol_rate_comps[i]->calcStateVariableEvolutionRateComponent(
468  _qp, (*_mat_prop_state_var_evol_rate_comps[i])[_qp]);
469 
470  for (unsigned int i = 0; i < _num_uo_state_vars; ++i)
471  {
472  if (!_uo_state_vars[i]->updateStateVariable(
473  _qp, _dt, (*_mat_prop_state_vars[i])[_qp], _state_vars_old_stored[i]))
474  _err_tol = true;
475  }
476 
477  for (unsigned int i = 0; i < _num_uo_slip_resistances; ++i)
478  _uo_slip_resistances[i]->calcSlipResistance(_qp, (*_mat_prop_slip_resistances[i])[_qp]);
479 }

Referenced by solveStatevar().

◆ validParams()

InputParameters FiniteStrainUObasedCP::validParams ( )
static

Definition at line 23 of file FiniteStrainUObasedCP.C.

24 {
25  InputParameters params = ComputeStressBase::validParams();
26  params.addClassDescription("UserObject based Crystal Plasticity system.");
27  params.addParam<Real>("rtol", 1e-6, "Constitutive stress residue relative tolerance");
28  params.addParam<Real>("abs_tol", 1e-6, "Constitutive stress residue absolute tolerance");
29  params.addParam<Real>(
30  "stol", 1e-2, "Constitutive slip system resistance relative residual tolerance");
31  params.addParam<Real>(
32  "zero_tol", 1e-12, "Tolerance for residual check when variable value is zero");
33  params.addParam<unsigned int>("maxiter", 100, "Maximum number of iterations for stress update");
34  params.addParam<unsigned int>(
35  "maxiter_state_variable", 100, "Maximum number of iterations for state variable update");
36  MooseEnum tan_mod_options("exact none", "none"); // Type of read
37  params.addParam<MooseEnum>("tan_mod_type",
38  tan_mod_options,
39  "Type of tangent moduli for preconditioner: default elastic");
40  params.addParam<unsigned int>(
41  "maximum_substep_iteration", 1, "Maximum number of substep iteration");
42  params.addParam<bool>("use_line_search", false, "Use line search in constitutive update");
43  params.addParam<Real>("min_line_search_step_size", 0.01, "Minimum line search step size");
44  params.addParam<Real>("line_search_tol", 0.5, "Line search bisection method tolerance");
45  params.addParam<unsigned int>(
46  "line_search_maxiter", 20, "Line search bisection method maximum number of iteration");
47  MooseEnum line_search_method("CUT_HALF BISECTION", "CUT_HALF");
48  params.addParam<MooseEnum>(
49  "line_search_method", line_search_method, "The method used in line search");
50  params.addRequiredParam<std::vector<UserObjectName>>(
51  "uo_slip_rates",
52  "List of names of user objects that define the slip rates for this material.");
53  params.addRequiredParam<std::vector<UserObjectName>>(
54  "uo_slip_resistances",
55  "List of names of user objects that define the slip resistances for this material.");
56  params.addRequiredParam<std::vector<UserObjectName>>(
57  "uo_state_vars",
58  "List of names of user objects that define the state variable for this material.");
59  params.addRequiredParam<std::vector<UserObjectName>>(
60  "uo_state_var_evol_rate_comps",
61  "List of names of user objects that define the state "
62  "variable evolution rate components for this material.");
63  return params;
64 }

Member Data Documentation

◆ _abs_tol

Real FiniteStrainUObasedCP::_abs_tol
protected

Stress residual equation absolute tolerance.

Definition at line 208 of file FiniteStrainUObasedCP.h.

Referenced by solveStress().

◆ _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().

◆ _crysrot

const MaterialProperty<RankTwoTensor>& FiniteStrainUObasedCP::_crysrot
protected

Crystal rotation.

Definition at line 262 of file FiniteStrainUObasedCP.h.

Referenced by postSolveQp().

◆ _deformation_gradient

const MaterialProperty<RankTwoTensor>& FiniteStrainUObasedCP::_deformation_gradient
protected

Definition at line 258 of file FiniteStrainUObasedCP.h.

Referenced by computeQpStress(), and postSolveQp().

◆ _deformation_gradient_old

const MaterialProperty<RankTwoTensor>& FiniteStrainUObasedCP::_deformation_gradient_old
protected

Definition at line 259 of file FiniteStrainUObasedCP.h.

Referenced by computeQpStress().

◆ _delta_dfgrd

RankTwoTensor FiniteStrainUObasedCP::_delta_dfgrd
protected

Used for substepping; Uniformly divides the increment in deformation gradient.

Definition at line 273 of file FiniteStrainUObasedCP.h.

Referenced by computeQpStress(), and FiniteStrainUObasedCP().

◆ _dfgrd_scale_factor

Real FiniteStrainUObasedCP::_dfgrd_scale_factor
protected

Scales the substepping increment to obtain deformation gradient at a substep iteration.

Definition at line 275 of file FiniteStrainUObasedCP.h.

◆ _dfgrd_tmp

RankTwoTensor FiniteStrainUObasedCP::_dfgrd_tmp
protected

Definition at line 264 of file FiniteStrainUObasedCP.h.

Referenced by calcJacobian(), calcResidual(), and computeQpStress().

◆ _dfgrd_tmp_old

RankTwoTensor FiniteStrainUObasedCP::_dfgrd_tmp_old
protected

Definition at line 273 of file FiniteStrainUObasedCP.h.

Referenced by computeQpStress().

◆ _elastic_strain

MaterialProperty<RankTwoTensor>& ComputeStressBase::_elastic_strain
protectedinherited

◆ _elasticity_tensor

const MaterialProperty<RankFourTensor>& FiniteStrainUObasedCP::_elasticity_tensor
protected

Elasticity tensor material property.

Definition at line 257 of file FiniteStrainUObasedCP.h.

Referenced by calcJacobian(), calcResidual(), elasticTangentModuli(), and elastoPlasticTangentModuli().

◆ _elasticity_tensor_name

const std::string FiniteStrainUObasedCP::_elasticity_tensor_name
protected

Name of the elasticity tensor material property.

Definition at line 255 of file FiniteStrainUObasedCP.h.

◆ _err_tol

bool FiniteStrainUObasedCP::_err_tol
protected

◆ _extra_stress

const MaterialProperty<RankTwoTensor>& ComputeStressBase::_extra_stress
protectedinherited

Extra stress tensor.

Definition at line 55 of file ComputeStressBase.h.

Referenced by ComputeStressBase::computeQpProperties().

◆ _fe

RankTwoTensor FiniteStrainUObasedCP::_fe
protected

◆ _flow_direction

std::vector<MaterialProperty<std::vector<RankTwoTensor> > *> FiniteStrainUObasedCP::_flow_direction
protected

◆ _fp

MaterialProperty<RankTwoTensor>& FiniteStrainUObasedCP::_fp
protected

Definition at line 246 of file FiniteStrainUObasedCP.h.

Referenced by initQpStatefulProperties(), and postSolveStress().

◆ _fp_inv

RankTwoTensor FiniteStrainUObasedCP::_fp_inv
protected

◆ _fp_old

const MaterialProperty<RankTwoTensor>& FiniteStrainUObasedCP::_fp_old
protected

Definition at line 247 of file FiniteStrainUObasedCP.h.

Referenced by preSolveQp().

◆ _fp_old_inv

RankTwoTensor FiniteStrainUObasedCP::_fp_old_inv
protected

◆ _initial_stress_fcn

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

initial stress components

Definition at line 58 of file ComputeStressBase.h.

◆ _jac

RankFourTensor FiniteStrainUObasedCP::_jac
protected

Jacobian tensor.

Definition at line 218 of file FiniteStrainUObasedCP.h.

Referenced by calcJacobian(), and solveStress().

◆ _Jacobian_mult

MaterialProperty<RankFourTensor>& ComputeStressBase::_Jacobian_mult
protectedinherited

◆ _lag_e

MaterialProperty<RankTwoTensor>& FiniteStrainUObasedCP::_lag_e
protected

Definition at line 250 of file FiniteStrainUObasedCP.h.

Referenced by initQpStatefulProperties(), and postSolveQp().

◆ _lsrch_max_iter

unsigned int FiniteStrainUObasedCP::_lsrch_max_iter
protected

Line search bisection method maximum iteration number.

Definition at line 241 of file FiniteStrainUObasedCP.h.

Referenced by lineSearchUpdate().

◆ _lsrch_method

MooseEnum FiniteStrainUObasedCP::_lsrch_method
protected

Line search method.

Definition at line 244 of file FiniteStrainUObasedCP.h.

Referenced by lineSearchUpdate().

◆ _lsrch_tol

Real FiniteStrainUObasedCP::_lsrch_tol
protected

Line search bisection method tolerance.

Definition at line 238 of file FiniteStrainUObasedCP.h.

Referenced by lineSearchUpdate().

◆ _mat_prop_slip_rates

std::vector<MaterialProperty<std::vector<Real> > *> FiniteStrainUObasedCP::_mat_prop_slip_rates
protected

Slip rates material property.

Definition at line 170 of file FiniteStrainUObasedCP.h.

Referenced by calcResidual(), FiniteStrainUObasedCP(), getSlipRates(), and initQpStatefulProperties().

◆ _mat_prop_slip_resistances

std::vector<MaterialProperty<std::vector<Real> > *> FiniteStrainUObasedCP::_mat_prop_slip_resistances
protected

Slip resistance material property.

Definition at line 173 of file FiniteStrainUObasedCP.h.

Referenced by FiniteStrainUObasedCP(), initQpStatefulProperties(), preSolveStatevar(), and updateSlipSystemResistanceAndStateVariable().

◆ _mat_prop_state_var_evol_rate_comps

std::vector<MaterialProperty<std::vector<Real> > *> FiniteStrainUObasedCP::_mat_prop_state_var_evol_rate_comps
protected

State variable evolution rate component material property.

Definition at line 182 of file FiniteStrainUObasedCP.h.

Referenced by FiniteStrainUObasedCP(), initQpStatefulProperties(), and updateSlipSystemResistanceAndStateVariable().

◆ _mat_prop_state_vars

std::vector<MaterialProperty<std::vector<Real> > *> FiniteStrainUObasedCP::_mat_prop_state_vars
protected

◆ _mat_prop_state_vars_old

std::vector<const MaterialProperty<std::vector<Real> > *> FiniteStrainUObasedCP::_mat_prop_state_vars_old
protected

Old state variable material property.

Definition at line 179 of file FiniteStrainUObasedCP.h.

Referenced by computeQpStress(), and FiniteStrainUObasedCP().

◆ _max_substep_iter

unsigned int FiniteStrainUObasedCP::_max_substep_iter
protected

Maximum number of substep iterations.

Definition at line 229 of file FiniteStrainUObasedCP.h.

Referenced by computeQpStress().

◆ _maxiter

unsigned int FiniteStrainUObasedCP::_maxiter
protected

Maximum number of iterations for stress update.

Definition at line 221 of file FiniteStrainUObasedCP.h.

Referenced by solveStress().

◆ _maxiterg

unsigned int FiniteStrainUObasedCP::_maxiterg
protected

Maximum number of iterations for internal variable update.

Definition at line 223 of file FiniteStrainUObasedCP.h.

Referenced by solveStatevar().

◆ _mechanical_strain

const MaterialProperty<RankTwoTensor>& ComputeStressBase::_mechanical_strain
protectedinherited

◆ _min_lsrch_step

Real FiniteStrainUObasedCP::_min_lsrch_step
protected

Minimum line search step size.

Definition at line 235 of file FiniteStrainUObasedCP.h.

Referenced by lineSearchUpdate().

◆ _num_uo_slip_rates

unsigned int FiniteStrainUObasedCP::_num_uo_slip_rates
protected

Number of slip rate user objects.

Definition at line 185 of file FiniteStrainUObasedCP.h.

Referenced by calcJacobian(), calcResidual(), computeQpStress(), FiniteStrainUObasedCP(), getSlipRates(), and initQpStatefulProperties().

◆ _num_uo_slip_resistances

unsigned int FiniteStrainUObasedCP::_num_uo_slip_resistances
protected

Number of slip resistance user objects.

Definition at line 188 of file FiniteStrainUObasedCP.h.

Referenced by FiniteStrainUObasedCP(), initQpStatefulProperties(), preSolveStatevar(), and updateSlipSystemResistanceAndStateVariable().

◆ _num_uo_state_var_evol_rate_comps

unsigned int FiniteStrainUObasedCP::_num_uo_state_var_evol_rate_comps
protected

Number of state variable evolution rate component user objects.

Definition at line 194 of file FiniteStrainUObasedCP.h.

Referenced by FiniteStrainUObasedCP(), initQpStatefulProperties(), and updateSlipSystemResistanceAndStateVariable().

◆ _num_uo_state_vars

unsigned int FiniteStrainUObasedCP::_num_uo_state_vars
protected

◆ _pk2

MaterialProperty<RankTwoTensor>& FiniteStrainUObasedCP::_pk2
protected

◆ _pk2_old

const MaterialProperty<RankTwoTensor>& FiniteStrainUObasedCP::_pk2_old
protected

Definition at line 249 of file FiniteStrainUObasedCP.h.

Referenced by preSolveQp().

◆ _resid

RankTwoTensor FiniteStrainUObasedCP::_resid
protected

Residual tensor.

Definition at line 215 of file FiniteStrainUObasedCP.h.

Referenced by calcResidual(), lineSearchUpdate(), and solveStress().

◆ _rtol

Real FiniteStrainUObasedCP::_rtol
protected

Stress residual equation relative tolerance.

Definition at line 206 of file FiniteStrainUObasedCP.h.

Referenced by solveStress().

◆ _state_vars_old

std::vector<std::vector<Real> > FiniteStrainUObasedCP::_state_vars_old
protected

Local state variable.

Definition at line 197 of file FiniteStrainUObasedCP.h.

Referenced by computeQpStress(), FiniteStrainUObasedCP(), initQpStatefulProperties(), and preSolveQp().

◆ _state_vars_old_stored

std::vector<std::vector<Real> > FiniteStrainUObasedCP::_state_vars_old_stored
protected

◆ _state_vars_prev

std::vector<std::vector<Real> > FiniteStrainUObasedCP::_state_vars_prev
protected

◆ _stol

Real FiniteStrainUObasedCP::_stol
protected

Internal variable update equation tolerance.

Definition at line 210 of file FiniteStrainUObasedCP.h.

Referenced by isStateVariablesConverged().

◆ _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(), FiniteStrainPlasticMaterial::computeQpStress(), ComputeMultiPlasticityStress::computeQpStress(), ComputeLinearViscoelasticStress::computeQpStress(), ComputeMultipleInelasticStress::computeQpStress(), ComputeMultipleInelasticStress::computeQpStressIntermediateConfiguration(), ComputeLinearElasticPFFractureStress::computeStrainSpectral(), ComputeLinearElasticPFFractureStress::computeStrainVolDev(), ComputeLinearElasticPFFractureStress::computeStressSpectral(), ComputeMultipleInelasticStress::finiteStrainRotation(), ComputeStressBase::initQpStatefulProperties(), FiniteStrainCrystalPlasticity::initQpStatefulProperties(), initQpStatefulProperties(), FiniteStrainHyperElasticViscoPlastic::initQpStatefulProperties(), ComputeMultiPlasticityStress::postReturnMap(), postSolveQp(), FiniteStrainHyperElasticViscoPlastic::postSolveQp(), FiniteStrainCrystalPlasticity::postSolveQp(), ComputeSmearedCrackingStress::updateCrackingStateAndStress(), ComputeMultipleInelasticStress::updateQpState(), and ComputeMultipleInelasticStress::updateQpStateSingleModel().

◆ _tan_mod_type

MooseEnum FiniteStrainUObasedCP::_tan_mod_type
protected

Type of tangent moduli calculation.

Definition at line 226 of file FiniteStrainUObasedCP.h.

Referenced by calcTangentModuli().

◆ _tau

DenseVector<Real> FiniteStrainUObasedCP::_tau
protected

Definition at line 266 of file FiniteStrainUObasedCP.h.

◆ _uo_slip_rates

std::vector<const CrystalPlasticitySlipRate *> FiniteStrainUObasedCP::_uo_slip_rates
protected

User objects that define the slip rate.

Definition at line 158 of file FiniteStrainUObasedCP.h.

Referenced by calcJacobian(), calcResidual(), computeQpStress(), FiniteStrainUObasedCP(), getSlipRates(), and initQpStatefulProperties().

◆ _uo_slip_resistances

std::vector<const CrystalPlasticitySlipResistance *> FiniteStrainUObasedCP::_uo_slip_resistances
protected

User objects that define the slip resistance.

Definition at line 161 of file FiniteStrainUObasedCP.h.

Referenced by FiniteStrainUObasedCP(), initQpStatefulProperties(), preSolveStatevar(), and updateSlipSystemResistanceAndStateVariable().

◆ _uo_state_var_evol_rate_comps

std::vector<const CrystalPlasticityStateVarRateComponent *> FiniteStrainUObasedCP::_uo_state_var_evol_rate_comps
protected

User objects that define the state variable evolution rate component.

Definition at line 167 of file FiniteStrainUObasedCP.h.

Referenced by FiniteStrainUObasedCP(), initQpStatefulProperties(), and updateSlipSystemResistanceAndStateVariable().

◆ _uo_state_vars

std::vector<const CrystalPlasticityStateVariable *> FiniteStrainUObasedCP::_uo_state_vars
protected

User objects that define the state variable.

Definition at line 164 of file FiniteStrainUObasedCP.h.

Referenced by FiniteStrainUObasedCP(), initQpStatefulProperties(), and updateSlipSystemResistanceAndStateVariable().

◆ _update_rot

MaterialProperty<RankTwoTensor>& FiniteStrainUObasedCP::_update_rot
protected

Definition at line 251 of file FiniteStrainUObasedCP.h.

Referenced by initQpStatefulProperties(), and postSolveQp().

◆ _update_rot_old

const MaterialProperty<RankTwoTensor>& FiniteStrainUObasedCP::_update_rot_old
protected

Definition at line 252 of file FiniteStrainUObasedCP.h.

◆ _use_line_search

bool FiniteStrainUObasedCP::_use_line_search
protected

Flag to activate line serach.

Definition at line 232 of file FiniteStrainUObasedCP.h.

Referenced by solveStress().

◆ _zero_tol

Real FiniteStrainUObasedCP::_zero_tol
protected

Residual tolerance when variable value is zero. Default 1e-12.

Definition at line 212 of file FiniteStrainUObasedCP.h.

Referenced by isStateVariablesConverged().


The documentation for this class was generated from the following files:
FiniteStrainUObasedCP::_uo_slip_rates
std::vector< const CrystalPlasticitySlipRate * > _uo_slip_rates
User objects that define the slip rate.
Definition: FiniteStrainUObasedCP.h:158
FiniteStrainUObasedCP::_uo_slip_resistances
std::vector< const CrystalPlasticitySlipResistance * > _uo_slip_resistances
User objects that define the slip resistance.
Definition: FiniteStrainUObasedCP.h:161
FiniteStrainUObasedCP::_dfgrd_tmp
RankTwoTensor _dfgrd_tmp
Definition: FiniteStrainUObasedCP.h:264
FiniteStrainUObasedCP::_tan_mod_type
MooseEnum _tan_mod_type
Type of tangent moduli calculation.
Definition: FiniteStrainUObasedCP.h:226
FiniteStrainUObasedCP::_fp_old
const MaterialProperty< RankTwoTensor > & _fp_old
Definition: FiniteStrainUObasedCP.h:247
FiniteStrainUObasedCP::_update_rot
MaterialProperty< RankTwoTensor > & _update_rot
Definition: FiniteStrainUObasedCP.h:251
FiniteStrainUObasedCP::_num_uo_slip_resistances
unsigned int _num_uo_slip_resistances
Number of slip resistance user objects.
Definition: FiniteStrainUObasedCP.h:188
FiniteStrainUObasedCP::_abs_tol
Real _abs_tol
Stress residual equation absolute tolerance.
Definition: FiniteStrainUObasedCP.h:208
ComputeStressBase::_stress
MaterialProperty< RankTwoTensor > & _stress
Stress material property.
Definition: ComputeStressBase.h:50
FiniteStrainUObasedCP::_fp
MaterialProperty< RankTwoTensor > & _fp
Definition: FiniteStrainUObasedCP.h:246
ComputeStressBase::_extra_stress
const MaterialProperty< RankTwoTensor > & _extra_stress
Extra stress tensor.
Definition: ComputeStressBase.h:55
FiniteStrainUObasedCP::_lsrch_tol
Real _lsrch_tol
Line search bisection method tolerance.
Definition: FiniteStrainUObasedCP.h:238
ComputeStressBase::_Jacobian_mult
MaterialProperty< RankFourTensor > & _Jacobian_mult
derivative of stress w.r.t. strain (_dstress_dstrain)
Definition: ComputeStressBase.h:61
FiniteStrainUObasedCP::_mat_prop_slip_rates
std::vector< MaterialProperty< std::vector< Real > > * > _mat_prop_slip_rates
Slip rates material property.
Definition: FiniteStrainUObasedCP.h:170
FiniteStrainUObasedCP::_stol
Real _stol
Internal variable update equation tolerance.
Definition: FiniteStrainUObasedCP.h:210
FiniteStrainUObasedCP::elastoPlasticTangentModuli
virtual void elastoPlasticTangentModuli()
calculate the exact tangent moduli for preconditioner.
Definition: FiniteStrainUObasedCP.C:582
FiniteStrainUObasedCP::_crysrot
const MaterialProperty< RankTwoTensor > & _crysrot
Crystal rotation.
Definition: FiniteStrainUObasedCP.h:262
FiniteStrainUObasedCP::elasticTangentModuli
virtual void elasticTangentModuli()
calculate the elastic tangent moduli for preconditioner.
Definition: FiniteStrainUObasedCP.C:625
FiniteStrainUObasedCP::_mat_prop_state_vars_old
std::vector< const MaterialProperty< std::vector< Real > > * > _mat_prop_state_vars_old
Old state variable material property.
Definition: FiniteStrainUObasedCP.h:179
FiniteStrainUObasedCP::_flow_direction
std::vector< MaterialProperty< std::vector< RankTwoTensor > > * > _flow_direction
Definition: FiniteStrainUObasedCP.h:267
FiniteStrainUObasedCP::_state_vars_prev
std::vector< std::vector< Real > > _state_vars_prev
Local old state variable.
Definition: FiniteStrainUObasedCP.h:203
FiniteStrainUObasedCP::_uo_state_vars
std::vector< const CrystalPlasticityStateVariable * > _uo_state_vars
User objects that define the state variable.
Definition: FiniteStrainUObasedCP.h:164
FiniteStrainUObasedCP::_update_rot_old
const MaterialProperty< RankTwoTensor > & _update_rot_old
Definition: FiniteStrainUObasedCP.h:252
FiniteStrainUObasedCP::_num_uo_slip_rates
unsigned int _num_uo_slip_rates
Number of slip rate user objects.
Definition: FiniteStrainUObasedCP.h:185
FiniteStrainUObasedCP::_mat_prop_slip_resistances
std::vector< MaterialProperty< std::vector< Real > > * > _mat_prop_slip_resistances
Slip resistance material property.
Definition: FiniteStrainUObasedCP.h:173
FiniteStrainUObasedCP::_deformation_gradient
const MaterialProperty< RankTwoTensor > & _deformation_gradient
Definition: FiniteStrainUObasedCP.h:258
FiniteStrainUObasedCP::solveQp
virtual void solveQp()
solve stress and internal variables.
Definition: FiniteStrainUObasedCP.C:274
FiniteStrainUObasedCP::_resid
RankTwoTensor _resid
Residual tensor.
Definition: FiniteStrainUObasedCP.h:215
FiniteStrainUObasedCP::_lag_e
MaterialProperty< RankTwoTensor > & _lag_e
Definition: FiniteStrainUObasedCP.h:250
FiniteStrainUObasedCP::_num_uo_state_var_evol_rate_comps
unsigned int _num_uo_state_var_evol_rate_comps
Number of state variable evolution rate component user objects.
Definition: FiniteStrainUObasedCP.h:194
ComputeStressBase::computeQpStress
virtual void computeQpStress()=0
Compute the stress and store it in the _stress material property for the current quadrature point.
FiniteStrainUObasedCP::postSolveQp
virtual void postSolveQp()
update stress and internal variable after solve.
Definition: FiniteStrainUObasedCP.C:284
FiniteStrainUObasedCP::getSlipRates
virtual void getSlipRates()
updates the slip rates.
Definition: FiniteStrainUObasedCP.C:492
FiniteStrainUObasedCP::_min_lsrch_step
Real _min_lsrch_step
Minimum line search step size.
Definition: FiniteStrainUObasedCP.h:235
ComputeStressBase::validParams
static InputParameters validParams()
Definition: ComputeStressBase.C:17
FiniteStrainUObasedCP::_uo_state_var_evol_rate_comps
std::vector< const CrystalPlasticityStateVarRateComponent * > _uo_state_var_evol_rate_comps
User objects that define the state variable evolution rate component.
Definition: FiniteStrainUObasedCP.h:167
FiniteStrainUObasedCP::_elasticity_tensor_name
const std::string _elasticity_tensor_name
Name of the elasticity tensor material property.
Definition: FiniteStrainUObasedCP.h:255
FiniteStrainUObasedCP::_pk2
MaterialProperty< RankTwoTensor > & _pk2
Definition: FiniteStrainUObasedCP.h:248
FiniteStrainUObasedCP::isStateVariablesConverged
virtual bool isStateVariablesConverged()
evaluates convergence of state variables.
Definition: FiniteStrainUObasedCP.C:350
FiniteStrainUObasedCP::_pk2_old
const MaterialProperty< RankTwoTensor > & _pk2_old
Definition: FiniteStrainUObasedCP.h:249
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
FiniteStrainUObasedCP::_delta_dfgrd
RankTwoTensor _delta_dfgrd
Used for substepping; Uniformly divides the increment in deformation gradient.
Definition: FiniteStrainUObasedCP.h:273
FiniteStrainUObasedCP::_deformation_gradient_old
const MaterialProperty< RankTwoTensor > & _deformation_gradient_old
Definition: FiniteStrainUObasedCP.h:259
FiniteStrainUObasedCP::_err_tol
bool _err_tol
Flag to check whether convergence is achieved.
Definition: FiniteStrainUObasedCP.h:270
FiniteStrainUObasedCP::_mat_prop_state_vars
std::vector< MaterialProperty< std::vector< Real > > * > _mat_prop_state_vars
State variable material property.
Definition: FiniteStrainUObasedCP.h:176
FiniteStrainUObasedCP::preSolveStatevar
virtual void preSolveStatevar()
set variables for internal variable solve.
Definition: FiniteStrainUObasedCP.C:303
FiniteStrainUObasedCP::_zero_tol
Real _zero_tol
Residual tolerance when variable value is zero. Default 1e-12.
Definition: FiniteStrainUObasedCP.h:212
FiniteStrainUObasedCP::_dfgrd_tmp_old
RankTwoTensor _dfgrd_tmp_old
Definition: FiniteStrainUObasedCP.h:273
FiniteStrainUObasedCP::_fp_inv
RankTwoTensor _fp_inv
Definition: FiniteStrainUObasedCP.h:265
FiniteStrainUObasedCP::_use_line_search
bool _use_line_search
Flag to activate line serach.
Definition: FiniteStrainUObasedCP.h:232
FiniteStrainUObasedCP::_max_substep_iter
unsigned int _max_substep_iter
Maximum number of substep iterations.
Definition: FiniteStrainUObasedCP.h:229
RankFourTensorTempl
Definition: ACGrGrElasticDrivingForce.h:20
FiniteStrainUObasedCP::solveStatevar
virtual void solveStatevar()
solve internal variables.
Definition: FiniteStrainUObasedCP.C:315
FiniteStrainUObasedCP::_fp_old_inv
RankTwoTensor _fp_old_inv
Definition: FiniteStrainUObasedCP.h:265
FiniteStrainUObasedCP::preSolveStress
virtual void preSolveStress()
set variables for stress solve.
Definition: FiniteStrainUObasedCP.C:381
FiniteStrainUObasedCP::postSolveStress
virtual void postSolveStress()
update stress and plastic deformation gradient after solve.
Definition: FiniteStrainUObasedCP.C:455
FiniteStrainUObasedCP::_state_vars_old
std::vector< std::vector< Real > > _state_vars_old
Local state variable.
Definition: FiniteStrainUObasedCP.h:197
FiniteStrainUObasedCP::_fe
RankTwoTensor _fe
Definition: FiniteStrainUObasedCP.h:265
FiniteStrainUObasedCP::_maxiter
unsigned int _maxiter
Maximum number of iterations for stress update.
Definition: FiniteStrainUObasedCP.h:221
FiniteStrainUObasedCP::calcTangentModuli
virtual void calcTangentModuli()
calculate the tangent moduli for preconditioner.
Definition: FiniteStrainUObasedCP.C:569
FiniteStrainUObasedCP::_maxiterg
unsigned int _maxiterg
Maximum number of iterations for internal variable update.
Definition: FiniteStrainUObasedCP.h:223
FiniteStrainUObasedCP::_mat_prop_state_var_evol_rate_comps
std::vector< MaterialProperty< std::vector< Real > > * > _mat_prop_state_var_evol_rate_comps
State variable evolution rate component material property.
Definition: FiniteStrainUObasedCP.h:182
FiniteStrainUObasedCP::_lsrch_method
MooseEnum _lsrch_method
Line search method.
Definition: FiniteStrainUObasedCP.h:244
FiniteStrainUObasedCP::_jac
RankFourTensor _jac
Jacobian tensor.
Definition: FiniteStrainUObasedCP.h:218
FiniteStrainUObasedCP::calcResidJacob
virtual void calcResidJacob()
calls the residual and jacobian functions used in the stress update algorithm.
Definition: FiniteStrainUObasedCP.C:483
RankTwoTensorTempl< Real >
FiniteStrainUObasedCP::_num_uo_state_vars
unsigned int _num_uo_state_vars
Number of state variable user objects.
Definition: FiniteStrainUObasedCP.h:191
FiniteStrainUObasedCP::_rtol
Real _rtol
Stress residual equation relative tolerance.
Definition: FiniteStrainUObasedCP.h:206
FiniteStrainUObasedCP::updateSlipSystemResistanceAndStateVariable
virtual void updateSlipSystemResistanceAndStateVariable()
updates the slip system resistances and state variables.
Definition: FiniteStrainUObasedCP.C:461
FiniteStrainUObasedCP::_lsrch_max_iter
unsigned int _lsrch_max_iter
Line search bisection method maximum iteration number.
Definition: FiniteStrainUObasedCP.h:241
FiniteStrainUObasedCP::calcResidual
virtual void calcResidual()
calculate stress residual.
Definition: FiniteStrainUObasedCP.C:505
FiniteStrainUObasedCP::calcJacobian
virtual void calcJacobian()
calculate jacobian.
Definition: FiniteStrainUObasedCP.C:531
FiniteStrainUObasedCP::_state_vars_old_stored
std::vector< std::vector< Real > > _state_vars_old_stored
Local stored state variable (for sub-stepping)
Definition: FiniteStrainUObasedCP.h:200
FiniteStrainUObasedCP::lineSearchUpdate
bool lineSearchUpdate(const Real rnorm_prev, const RankTwoTensor)
performs the line search update
Definition: FiniteStrainUObasedCP.C:632
FiniteStrainUObasedCP::_elasticity_tensor
const MaterialProperty< RankFourTensor > & _elasticity_tensor
Elasticity tensor material property.
Definition: FiniteStrainUObasedCP.h:257
ComputeStressBase::ComputeStressBase
ComputeStressBase(const InputParameters &parameters)
Definition: ComputeStressBase.C:28
FiniteStrainUObasedCP::preSolveQp
virtual void preSolveQp()
set variables for stress and internal variable solve.
Definition: FiniteStrainUObasedCP.C:264
FiniteStrainUObasedCP::solveStress
virtual void solveStress()
solves for stress, updates plastic deformation gradient.
Definition: FiniteStrainUObasedCP.C:386
FiniteStrainUObasedCP::postSolveStatevar
virtual void postSolveStatevar()
update internal variable after solve.
Definition: FiniteStrainUObasedCP.C:372