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

FiniteStrainCrystalPlasticity 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 <FiniteStrainCrystalPlasticity.h>

Inheritance diagram for FiniteStrainCrystalPlasticity:
[legend]

Public Member Functions

 FiniteStrainCrystalPlasticity (const InputParameters &parameters)
 

Protected Member Functions

virtual void computeQpStress ()
 This function updates the stress at a quadrature point. More...
 
virtual void computeQpElasticityTensor ()
 This function updates the elasticity tensor at a quadrature point. More...
 
virtual void initQpStatefulProperties ()
 This function initializes the stateful properties such as stress, plastic deformation gradient, slip system resistances, etc. More...
 
virtual void calc_resid_jacob (RankTwoTensor &, RankFourTensor &)
 This function calls the residual and jacobian functions used in the stress update algorithm. More...
 
virtual void getSlipIncrements ()
 This function updates the slip increments. More...
 
virtual void update_slip_system_resistance ()
 This function updates the slip system resistances. More...
 
virtual void updateGss ()
 This function updates the slip system resistances. More...
 
virtual void getSlipSystems ()
 This function reads slip system from file - see test. More...
 
virtual void assignSlipSysRes ()
 This function assign initial values of slip system resistances/internal variables read from getSlipSystems(). More...
 
virtual void readFileInitSlipSysRes ()
 This function read slip system resistances from file - see test. More...
 
virtual void getInitSlipSysRes ()
 This function assign slip system resistances - see test. More...
 
virtual void readFileFlowRateParams ()
 This function read flow rate parameters from file - see test. More...
 
virtual void getFlowRateParams ()
 This function assign flow rate parameters - see test. More...
 
virtual void readFileHardnessParams ()
 This function read hardness parameters from file. More...
 
virtual void getHardnessParams ()
 This function assign flow rate parameters from .i file - see test. More...
 
virtual void initSlipSysProps ()
 This function initializes slip system resistances. More...
 
virtual void initAdditionalProps ()
 This function initializes additional parameters. More...
 
virtual void preSolveQp ()
 This function set variables for stress and internal variable solve. More...
 
virtual void solveQp ()
 This function solves stress and internal variables. More...
 
virtual void postSolveQp ()
 This function update stress and internal variable after solve. More...
 
virtual void preSolveStatevar ()
 This function set variables for internal variable solve. More...
 
virtual void solveStatevar ()
 This function solves internal variables. More...
 
virtual void postSolveStatevar ()
 This function update internal variable after solve. More...
 
virtual void preSolveStress ()
 This function set variables for stress solve. More...
 
virtual void solveStress ()
 This function solves for stress, updates plastic deformation gradient. More...
 
virtual void postSolveStress ()
 This function update stress and plastic deformation gradient after solve. More...
 
virtual void calcResidual (RankTwoTensor &)
 This function calculate stress residual. More...
 
virtual void calcJacobian (RankFourTensor &)
 This function calculate jacobian. More...
 
virtual RankFourTensor calcTangentModuli ()
 This function calculate the tangent moduli for preconditioner. More...
 
virtual RankFourTensor elasticTangentModuli ()
 This function calculate the elastic tangent moduli for preconditioner. More...
 
virtual RankFourTensor elastoPlasticTangentModuli ()
 This function calculate the exact tangent moduli for preconditioner. More...
 
RankTwoTensor get_current_rotation (const RankTwoTensor &a)
 This function perform RU decomposition to obtain the rotation tensor. More...
 
RankTwoTensor getMatRot (const RankTwoTensor &a)
 This function perform RU decomposition to obtain the rotation tensor. More...
 
void calc_schmid_tensor ()
 This function calculate the Schmid tensor. More...
 
bool line_search_update (const Real rnorm_prev, const RankTwoTensor)
 This function performs the line search update. More...
 
void internalVariableUpdateNRiteration ()
 This function updates internal variables after each NewTon Raphson iteration (_fp_inv) More...
 
virtual void computeQpProperties () override
 

Protected Attributes

const unsigned int _nss
 Number of slip system resistance. More...
 
std::vector< Real > _gprops
 
std::vector< Real > _hprops
 
std::vector< Real > _flowprops
 
std::string _slip_sys_file_name
 File should contain slip plane normal and direction. See test. More...
 
std::string _slip_sys_res_prop_file_name
 File should contain initial values of the slip system resistances. More...
 
std::string _slip_sys_flow_prop_file_name
 File should contain values of the flow rate equation parameters. More...
 
std::string _slip_sys_hard_prop_file_name
 The hardening parameters in this class are read from .i file. The user can override to read from file. More...
 
Real _rtol
 Stress residual equation relative tolerance. More...
 
Real _abs_tol
 Stress residual equation absolute tolerance. More...
 
Real _gtol
 Internal variable update equation tolerance. More...
 
Real _slip_incr_tol
 Slip increment tolerance. More...
 
unsigned int _maxiter
 Maximum number of iterations for stress update. More...
 
unsigned int _maxiterg
 Maximum number of iterations for internal variable update. More...
 
unsigned int _num_slip_sys_flowrate_props
 Number of slip system flow rate parameters. More...
 
MooseEnum _tan_mod_type
 Type of tangent moduli calculation. More...
 
MooseEnum _intvar_read_type
 Read from options for initial values of internal variables. More...
 
unsigned int _num_slip_sys_props
 Number of slip system specific properties provided in the file containing slip system normals and directions. More...
 
bool _gen_rndm_stress_flag
 
bool _input_rndm_scale_var
 Input option for scaling variable to generate random stress when convergence fails. More...
 
Real _rndm_scale_var
 Scaling value. More...
 
unsigned int _rndm_seed
 Seed value. 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
 
MaterialProperty< RankTwoTensor > & _fp
 
const MaterialProperty< RankTwoTensor > & _fp_old
 
MaterialProperty< RankTwoTensor > & _pk2
 
const MaterialProperty< RankTwoTensor > & _pk2_old
 
MaterialProperty< RankTwoTensor > & _lag_e
 
const MaterialProperty< RankTwoTensor > & _lag_e_old
 
MaterialProperty< std::vector< Real > > & _gss
 
const MaterialProperty< std::vector< Real > > & _gss_old
 
MaterialProperty< Real > & _acc_slip
 
const MaterialProperty< Real > & _acc_slip_old
 
MaterialProperty< RankTwoTensor > & _update_rot
 
const MaterialProperty< RankTwoTensor > & _deformation_gradient
 
const MaterialProperty< RankTwoTensor > & _deformation_gradient_old
 
const MaterialProperty< RankFourTensor > & _elasticity_tensor
 
const MaterialProperty< RankTwoTensor > & _crysrot
 
DenseVector< Real > _mo
 
DenseVector< Real > _no
 
DenseVector< Real > _a0
 
DenseVector< Real > _xm
 
Real _h0
 
Real _tau_sat
 
Real _tau_init
 
Real _r
 
RankTwoTensor _dfgrd_tmp
 
RankTwoTensor _fe
 
RankTwoTensor _fp_old_inv
 
RankTwoTensor _fp_inv
 
RankTwoTensor _fp_prev_inv
 
DenseVector< Real > _slip_incr
 
DenseVector< Real > _tau
 
DenseVector< Real > _dslipdtau
 
std::vector< RankTwoTensor > _s0
 
RankTwoTensor _pk2_tmp
 
RankTwoTensor _pk2_tmp_old
 
Real _accslip_tmp
 
Real _accslip_tmp_old
 
std::vector< Real > _gss_tmp
 
std::vector< Real > _gss_tmp_old
 
DenseVector< Real > _slip_sys_props
 
DenseMatrix< Real > _dgss_dsliprate
 
bool _read_from_slip_sys_file
 
bool _err_tol
 
RankTwoTensor _delta_dfgrd
 Flag to check whether convergence is achieved. More...
 
RankTwoTensor _dfgrd_tmp_old
 
Real _dfgrd_scale_factor
 Scales the substepping increment to obtain deformation gradient at a substep iteration. More...
 
bool _first_step_iter
 Flags to reset variables and reinitialize variables. More...
 
bool _last_step_iter
 
bool _first_substep
 
const std::string _base_name
 
const std::string _elasticity_tensor_name
 
const MaterialProperty< RankTwoTensor > & _mechanical_strain
 
MaterialProperty< RankTwoTensor > & _stress
 
MaterialProperty< RankTwoTensor > & _elastic_strain
 
const MaterialProperty< RankTwoTensor > & _extra_stress
 Extra stress tensor. More...
 
std::vector< Function * > _initial_stress_fcn
 initial stress components More...
 
MaterialProperty< RankFourTensor > & _Jacobian_mult
 derivative of stress w.r.t. strain (_dstress_dstrain) More...
 

Detailed Description

FiniteStrainCrystalPlasticity 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.

Definition at line 27 of file FiniteStrainCrystalPlasticity.h.

Constructor & Destructor Documentation

◆ FiniteStrainCrystalPlasticity()

FiniteStrainCrystalPlasticity::FiniteStrainCrystalPlasticity ( const InputParameters &  parameters)

Definition at line 100 of file FiniteStrainCrystalPlasticity.C.

101  : ComputeStressBase(parameters),
102  _nss(getParam<int>("nss")),
103  _gprops(getParam<std::vector<Real>>("gprops")),
104  _hprops(getParam<std::vector<Real>>("hprops")),
105  _flowprops(getParam<std::vector<Real>>("flowprops")),
106  _slip_sys_file_name(getParam<FileName>("slip_sys_file_name")),
107  _slip_sys_res_prop_file_name(getParam<FileName>("slip_sys_res_prop_file_name")),
108  _slip_sys_flow_prop_file_name(getParam<FileName>("slip_sys_flow_prop_file_name")),
109  _slip_sys_hard_prop_file_name(getParam<FileName>("slip_sys_hard_prop_file_name")),
110  _rtol(getParam<Real>("rtol")),
111  _abs_tol(getParam<Real>("abs_tol")),
112  _gtol(getParam<Real>("gtol")),
113  _slip_incr_tol(getParam<Real>("slip_incr_tol")),
114  _maxiter(getParam<unsigned int>("maxiter")),
115  _maxiterg(getParam<unsigned int>("maxitergss")),
116  _num_slip_sys_flowrate_props(getParam<unsigned int>("num_slip_sys_flowrate_props")),
117  _tan_mod_type(getParam<MooseEnum>("tan_mod_type")),
118  _intvar_read_type(getParam<MooseEnum>("intvar_read_type")),
119  _num_slip_sys_props(getParam<unsigned int>("num_slip_sys_props")),
120  _gen_rndm_stress_flag(getParam<bool>("gen_random_stress_flag")),
121  _input_rndm_scale_var(getParam<bool>("input_random_scaling_var")),
122  _rndm_scale_var(getParam<Real>("random_scaling_var")),
123  _rndm_seed(getParam<unsigned int>("random_seed")),
124  _max_substep_iter(getParam<unsigned int>("maximum_substep_iteration")),
125  _use_line_search(getParam<bool>("use_line_search")),
126  _min_lsrch_step(getParam<Real>("min_line_search_step_size")),
127  _lsrch_tol(getParam<Real>("line_search_tol")),
128  _lsrch_max_iter(getParam<unsigned int>("line_search_maxiter")),
129  _lsrch_method(getParam<MooseEnum>("line_search_method")),
130  _fp(declareProperty<RankTwoTensor>("fp")), // Plastic deformation gradient
131  _fp_old(getMaterialPropertyOld<RankTwoTensor>(
132  "fp")), // Plastic deformation gradient of previous increment
133  _pk2(declareProperty<RankTwoTensor>("pk2")), // 2nd Piola Kirchoff Stress
134  _pk2_old(getMaterialPropertyOld<RankTwoTensor>(
135  "pk2")), // 2nd Piola Kirchoff Stress of previous increment
136  _lag_e(declareProperty<RankTwoTensor>("lage")), // Lagrangian strain
137  _lag_e_old(
138  getMaterialPropertyOld<RankTwoTensor>("lage")), // Lagrangian strain of previous increment
139  _gss(declareProperty<std::vector<Real>>("gss")), // Slip system resistances
140  _gss_old(getMaterialPropertyOld<std::vector<Real>>(
141  "gss")), // Slip system resistances of previous increment
142  _acc_slip(declareProperty<Real>("acc_slip")), // Accumulated slip
144  getMaterialPropertyOld<Real>("acc_slip")), // Accumulated alip of previous increment
145  _update_rot(declareProperty<RankTwoTensor>(
146  "update_rot")), // Rotation tensor considering material rotation and crystal orientation
147  _deformation_gradient(getMaterialProperty<RankTwoTensor>("deformation_gradient")),
148  _deformation_gradient_old(getMaterialPropertyOld<RankTwoTensor>("deformation_gradient")),
149  _elasticity_tensor(getMaterialProperty<RankFourTensor>("elasticity_tensor")),
150  _crysrot(getMaterialProperty<RankTwoTensor>("crysrot")),
151  _mo(_nss * LIBMESH_DIM),
152  _no(_nss * LIBMESH_DIM),
153  _slip_incr(_nss),
154  _tau(_nss),
155  _dslipdtau(_nss),
156  _s0(_nss),
157  _gss_tmp(_nss),
160 {
161  _err_tol = false;
162 
163  if (_num_slip_sys_props > 0)
165 
166  _pk2_tmp.zero();
167  _delta_dfgrd.zero();
168 
169  _first_step_iter = false;
170  _last_step_iter = false;
171  // Initialize variables in the first iteration of substepping
172  _first_substep = true;
173 
174  _read_from_slip_sys_file = false;
175  if (_intvar_read_type == "slip_sys_file")
177 
179  mooseError("Crystal Plasticity Error: Specify number of internal variable's initial values to "
180  "be read from slip system file");
181 
182  getSlipSystems();
183 
184  RankTwoTensor::initRandom(_rndm_seed);
185 }
const MaterialProperty< RankTwoTensor > & _pk2_old
MaterialProperty< RankTwoTensor > & _lag_e
const MaterialProperty< RankTwoTensor > & _deformation_gradient
bool _input_rndm_scale_var
Input option for scaling variable to generate random stress when convergence fails.
bool _use_line_search
Flag to activate line serach.
const MaterialProperty< std::vector< Real > > & _gss_old
const MaterialProperty< RankFourTensor > & _elasticity_tensor
MaterialProperty< RankTwoTensor > & _pk2
std::string _slip_sys_flow_prop_file_name
File should contain values of the flow rate equation parameters.
unsigned int _num_slip_sys_flowrate_props
Number of slip system flow rate parameters.
const MaterialProperty< RankTwoTensor > & _fp_old
Real _rtol
Stress residual equation relative tolerance.
std::string _slip_sys_hard_prop_file_name
The hardening parameters in this class are read from .i file. The user can override to read from file...
unsigned int _maxiter
Maximum number of iterations for stress update.
const unsigned int _nss
Number of slip system resistance.
RankTwoTensor _delta_dfgrd
Flag to check whether convergence is achieved.
MaterialProperty< std::vector< Real > > & _gss
const MaterialProperty< RankTwoTensor > & _deformation_gradient_old
const MaterialProperty< RankTwoTensor > & _crysrot
MaterialProperty< RankTwoTensor > & _fp
unsigned int _lsrch_max_iter
Line search bisection method maximum iteration number.
Real _gtol
Internal variable update equation tolerance.
std::string _slip_sys_res_prop_file_name
File should contain initial values of the slip system resistances.
const MaterialProperty< RankTwoTensor > & _lag_e_old
ComputeStressBase(const InputParameters &parameters)
MooseEnum _tan_mod_type
Type of tangent moduli calculation.
Real _lsrch_tol
Line search bisection method tolerance.
virtual void getSlipSystems()
This function reads slip system from file - see test.
unsigned int _num_slip_sys_props
Number of slip system specific properties provided in the file containing slip system normals and dir...
unsigned int _max_substep_iter
Maximum number of substep iterations.
bool _first_step_iter
Flags to reset variables and reinitialize variables.
MaterialProperty< RankTwoTensor > & _update_rot
Real _abs_tol
Stress residual equation absolute tolerance.
unsigned int _maxiterg
Maximum number of iterations for internal variable update.
const MaterialProperty< Real > & _acc_slip_old
MooseEnum _intvar_read_type
Read from options for initial values of internal variables.
Real _min_lsrch_step
Minimum line search step size.
std::string _slip_sys_file_name
File should contain slip plane normal and direction. See test.
Real _slip_incr_tol
Slip increment tolerance.

Member Function Documentation

◆ assignSlipSysRes()

void FiniteStrainCrystalPlasticity::assignSlipSysRes ( )
protectedvirtual

This function assign initial values of slip system resistances/internal variables read from getSlipSystems().

Definition at line 233 of file FiniteStrainCrystalPlasticity.C.

Referenced by initSlipSysProps().

234 {
235  _gss[_qp].resize(_nss);
236 
237  for (unsigned int i = 0; i < _nss; ++i)
238  _gss[_qp][i] = _slip_sys_props(i);
239 }
const unsigned int _nss
Number of slip system resistance.
MaterialProperty< std::vector< Real > > & _gss

◆ calc_resid_jacob()

void FiniteStrainCrystalPlasticity::calc_resid_jacob ( RankTwoTensor &  resid,
RankFourTensor &  jac 
)
protectedvirtual

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

Definition at line 881 of file FiniteStrainCrystalPlasticity.C.

Referenced by solveStress().

882 {
883  calcResidual(resid);
884  if (_err_tol)
885  return;
886  calcJacobian(jac);
887 }
virtual void calcJacobian(RankFourTensor &)
This function calculate jacobian.
virtual void calcResidual(RankTwoTensor &)
This function calculate stress residual.

◆ calc_schmid_tensor()

void FiniteStrainCrystalPlasticity::calc_schmid_tensor ( )
protected

This function calculate the Schmid tensor.

Definition at line 1052 of file FiniteStrainCrystalPlasticity.C.

Referenced by preSolveQp().

1053 {
1054  DenseVector<Real> mo(LIBMESH_DIM * _nss), no(LIBMESH_DIM * _nss);
1055 
1056  // Update slip direction and normal with crystal orientation
1057  for (unsigned int i = 0; i < _nss; ++i)
1058  {
1059  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
1060  {
1061  mo(i * LIBMESH_DIM + j) = 0.0;
1062  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
1063  mo(i * LIBMESH_DIM + j) =
1064  mo(i * LIBMESH_DIM + j) + _crysrot[_qp](j, k) * _mo(i * LIBMESH_DIM + k);
1065  }
1066 
1067  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
1068  {
1069  no(i * LIBMESH_DIM + j) = 0.0;
1070  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
1071  no(i * LIBMESH_DIM + j) =
1072  no(i * LIBMESH_DIM + j) + _crysrot[_qp](j, k) * _no(i * LIBMESH_DIM + k);
1073  }
1074  }
1075 
1076  // Calculate Schmid tensor and resolved shear stresses
1077  for (unsigned int i = 0; i < _nss; ++i)
1078  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
1079  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
1080  _s0[i](j, k) = mo(i * LIBMESH_DIM + j) * no(i * LIBMESH_DIM + k);
1081 }
const unsigned int _nss
Number of slip system resistance.
const MaterialProperty< RankTwoTensor > & _crysrot

◆ calcJacobian()

void FiniteStrainCrystalPlasticity::calcJacobian ( RankFourTensor &  jac)
protectedvirtual

This function calculate jacobian.

Definition at line 930 of file FiniteStrainCrystalPlasticity.C.

Referenced by calc_resid_jacob().

931 {
932  RankFourTensor dfedfpinv, deedfe, dfpinvdpk2;
933 
934  std::vector<RankTwoTensor> dtaudpk2(_nss), dfpinvdslip(_nss);
935 
936  for (unsigned int i = 0; i < _nss; ++i)
937  {
938  dtaudpk2[i] = _s0[i];
939  dfpinvdslip[i] = -_fp_old_inv * _s0[i];
940  }
941 
942  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
943  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
944  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
945  dfedfpinv(i, j, k, j) = _dfgrd_tmp(i, k);
946 
947  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
948  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
949  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
950  {
951  deedfe(i, j, k, i) = deedfe(i, j, k, i) + _fe(k, j) * 0.5;
952  deedfe(i, j, k, j) = deedfe(i, j, k, j) + _fe(k, i) * 0.5;
953  }
954 
955  for (unsigned int i = 0; i < _nss; ++i)
956  dfpinvdpk2 += (dfpinvdslip[i] * _dslipdtau(i)).outerProduct(dtaudpk2[i]);
957 
958  jac =
959  RankFourTensor::IdentityFour() - (_elasticity_tensor[_qp] * deedfe * dfedfpinv * dfpinvdpk2);
960 }
const MaterialProperty< RankFourTensor > & _elasticity_tensor
const unsigned int _nss
Number of slip system resistance.

◆ calcResidual()

void FiniteStrainCrystalPlasticity::calcResidual ( RankTwoTensor &  resid)
protectedvirtual

This function calculate stress residual.

Definition at line 890 of file FiniteStrainCrystalPlasticity.C.

Referenced by calc_resid_jacob(), and line_search_update().

891 {
892  RankTwoTensor iden, ce, ee, ce_pk2, eqv_slip_incr, pk2_new;
893 
894  iden.zero();
895  iden.addIa(1.0);
896 
897  _fe = _dfgrd_tmp * _fp_prev_inv; // _fp_inv ==> _fp_prev_inv
898 
899  ce = _fe.transpose() * _fe;
900  ce_pk2 = ce * _pk2_tmp;
901  ce_pk2 = ce_pk2 / _fe.det();
902 
903  // Calculate Schmid tensor and resolved shear stresses
904  for (unsigned int i = 0; i < _nss; ++i)
905  _tau(i) = ce_pk2.doubleContraction(_s0[i]);
906 
907  getSlipIncrements(); // Calculate dslip,dslipdtau
908 
909  if (_err_tol)
910  return;
911 
912  eqv_slip_incr.zero();
913  for (unsigned int i = 0; i < _nss; ++i)
914  eqv_slip_incr += _s0[i] * _slip_incr(i);
915 
916  eqv_slip_incr = iden - eqv_slip_incr;
917  _fp_inv = _fp_old_inv * eqv_slip_incr;
918  _fe = _dfgrd_tmp * _fp_inv;
919 
920  ce = _fe.transpose() * _fe;
921  ee = ce - iden;
922  ee *= 0.5;
923 
924  pk2_new = _elasticity_tensor[_qp] * ee;
925 
926  resid = _pk2_tmp - pk2_new;
927 }
const MaterialProperty< RankFourTensor > & _elasticity_tensor
const unsigned int _nss
Number of slip system resistance.
virtual void getSlipIncrements()
This function updates the slip increments.

◆ calcTangentModuli()

RankFourTensor FiniteStrainCrystalPlasticity::calcTangentModuli ( )
protectedvirtual

This function 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 1035 of file FiniteStrainCrystalPlasticity.C.

Referenced by postSolveQp().

1036 {
1037  RankFourTensor tan_mod;
1038 
1039  switch (_tan_mod_type)
1040  {
1041  case 0:
1042  tan_mod = elastoPlasticTangentModuli();
1043  break;
1044  default:
1045  tan_mod = elasticTangentModuli();
1046  }
1047 
1048  return tan_mod;
1049 }
virtual RankFourTensor elastoPlasticTangentModuli()
This function calculate the exact tangent moduli for preconditioner.
virtual RankFourTensor elasticTangentModuli()
This function calculate the elastic tangent moduli for preconditioner.
MooseEnum _tan_mod_type
Type of tangent moduli calculation.

◆ computeQpElasticityTensor()

void FiniteStrainCrystalPlasticity::computeQpElasticityTensor ( )
protectedvirtual

This function updates the elasticity tensor at a quadrature point.

Presently void.

Definition at line 1030 of file FiniteStrainCrystalPlasticity.C.

1031 {
1032 }

◆ computeQpProperties()

void ComputeStressBase::computeQpProperties ( )
overrideprotectedvirtualinherited

Definition at line 51 of file ComputeStressBase.C.

52 {
54 
55  // Add in extra stress
56  _stress[_qp] += _extra_stress[_qp];
57 }
virtual void computeQpStress()=0
MaterialProperty< RankTwoTensor > & _stress
const MaterialProperty< RankTwoTensor > & _extra_stress
Extra stress tensor.

◆ computeQpStress()

void FiniteStrainCrystalPlasticity::computeQpStress ( )
protectedvirtual

This function updates the stress at a quadrature point.

Solves stress residual equation using NR.

Updates slip system resistances iteratively.

Implements ComputeStressBase.

Definition at line 492 of file FiniteStrainCrystalPlasticity.C.

493 {
494  unsigned int substep_iter = 1; // Depth of substepping; Limited to maximum substep iteration
495  unsigned int num_substep = 1; // Calculated from substep_iter as 2^substep_iter
496  Real dt_original = _dt; // Stores original _dt; Reset at the end of solve
497  _first_substep = true; // Initialize variables at substep_iter = 1
498 
499  if (_max_substep_iter > 1)
500  {
502  if (_dfgrd_tmp_old.det() == 0)
503  _dfgrd_tmp_old.addIa(1.0);
504 
506  _err_tol = true; // Indicator to continue substepping
507  }
508 
509  // Substepping loop
510  while (_err_tol && _max_substep_iter > 1)
511  {
512  _dt = dt_original / num_substep;
513 
514  for (unsigned int istep = 0; istep < num_substep; ++istep)
515  {
516  _first_step_iter = false;
517  if (istep == 0)
518  _first_step_iter = true;
519 
520  _last_step_iter = false;
521  if (istep == num_substep - 1)
522  _last_step_iter = true;
523 
524  _dfgrd_scale_factor = (static_cast<Real>(istep) + 1) / num_substep;
525 
526  preSolveQp();
527  solveQp();
528 
529  if (_err_tol)
530  {
531  substep_iter++;
532  num_substep *= 2;
533  break;
534  }
535  }
536 
537  _first_substep = false; // Prevents reinitialization
538  _dt = dt_original; // Resets dt
539 
540 #ifdef DEBUG
541  if (substep_iter > _max_substep_iter)
542  mooseWarning("FiniteStrainCrystalPlasticity: Failure with substepping");
543 #endif
544 
545  if (!_err_tol || substep_iter > _max_substep_iter)
546  postSolveQp(); // Evaluate variables after successful solve or indicate failure
547  }
548 
549  // No substepping
550  if (_max_substep_iter == 1)
551  {
552  preSolveQp();
553  solveQp();
554  postSolveQp();
555  }
556 }
const MaterialProperty< RankTwoTensor > & _deformation_gradient
Real _dfgrd_scale_factor
Scales the substepping increment to obtain deformation gradient at a substep iteration.
RankTwoTensor _delta_dfgrd
Flag to check whether convergence is achieved.
const MaterialProperty< RankTwoTensor > & _deformation_gradient_old
virtual void postSolveQp()
This function update stress and internal variable after solve.
unsigned int _max_substep_iter
Maximum number of substep iterations.
bool _first_step_iter
Flags to reset variables and reinitialize variables.
virtual void solveQp()
This function solves stress and internal variables.
virtual void preSolveQp()
This function set variables for stress and internal variable solve.

◆ elasticTangentModuli()

RankFourTensor FiniteStrainCrystalPlasticity::elasticTangentModuli ( )
protectedvirtual

This function calculate the elastic tangent moduli for preconditioner.

Definition at line 1123 of file FiniteStrainCrystalPlasticity.C.

Referenced by calcTangentModuli().

1124 {
1125  return _elasticity_tensor[_qp]; // update jacobian_mult
1126 }
const MaterialProperty< RankFourTensor > & _elasticity_tensor

◆ elastoPlasticTangentModuli()

RankFourTensor FiniteStrainCrystalPlasticity::elastoPlasticTangentModuli ( )
protectedvirtual

This function calculate the exact tangent moduli for preconditioner.

Definition at line 1084 of file FiniteStrainCrystalPlasticity.C.

Referenced by calcTangentModuli().

1085 {
1086  RankFourTensor tan_mod;
1087  RankTwoTensor pk2fet, fepk2;
1088  RankFourTensor deedfe, dsigdpk2dfe;
1089 
1090  // Fill in the matrix stiffness material property
1091 
1092  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
1093  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
1094  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
1095  {
1096  deedfe(i, j, k, i) = deedfe(i, j, k, i) + _fe(k, j) * 0.5;
1097  deedfe(i, j, k, j) = deedfe(i, j, k, j) + _fe(k, i) * 0.5;
1098  }
1099 
1100  dsigdpk2dfe = _fe.mixedProductIkJl(_fe) * _elasticity_tensor[_qp] * deedfe;
1101 
1102  pk2fet = _pk2_tmp * _fe.transpose();
1103  fepk2 = _fe * _pk2_tmp;
1104 
1105  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
1106  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
1107  for (unsigned int l = 0; l < LIBMESH_DIM; ++l)
1108  {
1109  tan_mod(i, j, i, l) = tan_mod(i, j, i, l) + pk2fet(l, j);
1110  tan_mod(i, j, j, l) = tan_mod(i, j, j, l) + fepk2(i, l);
1111  }
1112 
1113  tan_mod += dsigdpk2dfe;
1114 
1115  Real je = _fe.det();
1116  if (je > 0.0)
1117  tan_mod /= je;
1118 
1119  return tan_mod;
1120 }
const MaterialProperty< RankFourTensor > & _elasticity_tensor

◆ get_current_rotation()

RankTwoTensor FiniteStrainCrystalPlasticity::get_current_rotation ( const RankTwoTensor &  a)
protected

This function perform RU decomposition to obtain the rotation tensor.

Definition at line 988 of file FiniteStrainCrystalPlasticity.C.

Referenced by postSolveQp().

989 {
990  return getMatRot(a);
991 }
RankTwoTensor getMatRot(const RankTwoTensor &a)
This function perform RU decomposition to obtain the rotation tensor.

◆ getFlowRateParams()

void FiniteStrainCrystalPlasticity::getFlowRateParams ( )
protectedvirtual

This function assign flow rate parameters - see test.

.i input file format start_slip_sys_num, end_slip_sys_num, value1, value2

Definition at line 343 of file FiniteStrainCrystalPlasticity.C.

Referenced by initSlipSysProps().

344 {
345  if (_flowprops.size() <= 0)
346  mooseError("FiniteStrainCrystalPLasticity: Error in reading flow rate properties: Specify "
347  "input in .i file or a slip_sys_flow_prop_file_name");
348 
349  _a0.resize(_nss);
350  _xm.resize(_nss);
351 
352  unsigned int num_data_grp = 2 + _num_slip_sys_flowrate_props; // Number of data per group e.g.
353  // start_slip_sys, end_slip_sys,
354  // value1, value2, ..
355 
356  for (unsigned int i = 0; i < _flowprops.size() / num_data_grp; ++i)
357  {
358  Real vs, ve;
359  unsigned int is, ie;
360 
361  vs = _flowprops[i * num_data_grp];
362  ve = _flowprops[i * num_data_grp + 1];
363 
364  if (vs <= 0 || ve <= 0)
365  mooseError("FiniteStrainCrystalPLasticity: Indices in flow rate parameter read must be "
366  "positive integers: is = ",
367  vs,
368  " ie = ",
369  ve);
370 
371  if (vs != floor(vs) || ve != floor(ve))
372  mooseError("FiniteStrainCrystalPLasticity: Error in reading flow props: Values specifying "
373  "start and end number of slip system groups should be integer");
374 
375  is = static_cast<unsigned int>(vs);
376  ie = static_cast<unsigned int>(ve);
377 
378  if (is > ie)
379  mooseError("FiniteStrainCrystalPLasticity: Start index is = ",
380  is,
381  " should be greater than end index ie = ",
382  ie,
383  " in flow rate parameter read");
384 
385  for (unsigned int j = is; j <= ie; ++j)
386  {
387  _a0(j - 1) = _flowprops[i * num_data_grp + 2];
388  _xm(j - 1) = _flowprops[i * num_data_grp + 3];
389  }
390  }
391 
392  for (unsigned int i = 0; i < _nss; ++i)
393  {
394  if (!(_a0(i) > 0.0 && _xm(i) > 0.0))
395  {
396  mooseWarning(
397  "FiniteStrainCrystalPlasticity: Non-positive flow rate parameters ", _a0(i), ",", _xm(i));
398  break;
399  }
400  }
401 }
unsigned int _num_slip_sys_flowrate_props
Number of slip system flow rate parameters.
const unsigned int _nss
Number of slip system resistance.

◆ getHardnessParams()

void FiniteStrainCrystalPlasticity::getHardnessParams ( )
protectedvirtual

This function assign flow rate parameters from .i file - see test.

Definition at line 411 of file FiniteStrainCrystalPlasticity.C.

Referenced by initSlipSysProps().

412 {
413  if (_hprops.size() <= 0)
414  mooseError("FiniteStrainCrystalPLasticity: Error in reading hardness properties: Specify input "
415  "in .i file or a slip_sys_hard_prop_file_name");
416 
417  _r = _hprops[0];
418  _h0 = _hprops[1];
419  _tau_init = _hprops[2];
420  _tau_sat = _hprops[3];
421 }

◆ getInitSlipSysRes()

void FiniteStrainCrystalPlasticity::getInitSlipSysRes ( )
protectedvirtual

This function assign slip system resistances - see test.

.i input file format start_slip_sys_num, end_slip_sys_num, value.

Definition at line 261 of file FiniteStrainCrystalPlasticity.C.

Referenced by initSlipSysProps().

262 {
263  if (_gprops.size() <= 0)
264  mooseError("FiniteStrainCrystalPLasticity: Error in reading slip system resistance properties: "
265  "Specify input in .i file or in slip_sys_res_prop_file or in slip_sys_file");
266 
267  _gss[_qp].resize(_nss, 0.0);
268 
269  unsigned int num_data_grp = 3; // Number of data per group e.g. start_slip_sys, end_slip_sys,
270  // value
271 
272  for (unsigned int i = 0; i < _gprops.size() / num_data_grp; ++i)
273  {
274  Real vs, ve;
275  unsigned int is, ie;
276 
277  vs = _gprops[i * num_data_grp];
278  ve = _gprops[i * num_data_grp + 1];
279 
280  if (vs <= 0 || ve <= 0)
281  mooseError("FiniteStrainCrystalPLasticity: Indices in gss property read must be positive "
282  "integers: is = ",
283  vs,
284  " ie = ",
285  ve);
286 
287  if (vs != floor(vs) || ve != floor(ve))
288  mooseError("FiniteStrainCrystalPLasticity: Error in reading slip system resistances: Values "
289  "specifying start and end number of slip system groups should be integer");
290 
291  is = static_cast<unsigned int>(vs);
292  ie = static_cast<unsigned int>(ve);
293 
294  if (is > ie)
295  mooseError("FiniteStrainCrystalPLasticity: Start index is = ",
296  is,
297  " should be greater than end index ie = ",
298  ie,
299  " in slip system resistance property read");
300 
301  for (unsigned int j = is; j <= ie; ++j)
302  _gss[_qp][j - 1] = _gprops[i * num_data_grp + 2];
303  }
304 
305  for (unsigned int i = 0; i < _nss; ++i)
306  if (_gss[_qp][i] <= 0.0)
307  mooseError("FiniteStrainCrystalPLasticity: Value of resistance for slip system ",
308  i + 1,
309  " non positive");
310 }
const unsigned int _nss
Number of slip system resistance.
MaterialProperty< std::vector< Real > > & _gss

◆ getMatRot()

RankTwoTensor FiniteStrainCrystalPlasticity::getMatRot ( const RankTwoTensor &  a)
protected

This function perform RU decomposition to obtain the rotation tensor.

Definition at line 995 of file FiniteStrainCrystalPlasticity.C.

Referenced by get_current_rotation().

996 {
997  RankTwoTensor rot;
998  RankTwoTensor c, diag, evec;
999  PetscScalar cmat[LIBMESH_DIM][LIBMESH_DIM], work[10];
1000  PetscReal w[LIBMESH_DIM];
1001  PetscBLASInt nd = LIBMESH_DIM, lwork = 10, info;
1002 
1003  c = a.transpose() * a;
1004 
1005  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
1006  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
1007  cmat[i][j] = c(i, j);
1008 
1009  LAPACKsyev_("V", "U", &nd, &cmat[0][0], &nd, w, work, &lwork, &info);
1010 
1011  if (info != 0)
1012  mooseError("FiniteStrainCrystalPLasticity: DSYEV function call in getMatRot function failed");
1013 
1014  diag.zero();
1015 
1016  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
1017  diag(i, i) = std::sqrt(w[i]);
1018 
1019  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
1020  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
1021  evec(i, j) = cmat[i][j];
1022 
1023  rot = a * ((evec.transpose() * diag * evec).inverse());
1024 
1025  return rot;
1026 }

◆ getSlipIncrements()

void FiniteStrainCrystalPlasticity::getSlipIncrements ( )
protectedvirtual

This function updates the slip increments.

And derivative of slip w.r.t. resolved shear stress.

Reimplemented in FiniteStrainCPSlipRateRes.

Definition at line 964 of file FiniteStrainCrystalPlasticity.C.

Referenced by calcResidual(), and FiniteStrainCPSlipRateRes::getSlipIncrements().

965 {
966  for (unsigned int i = 0; i < _nss; ++i)
967  {
968  _slip_incr(i) = _a0(i) * std::pow(std::abs(_tau(i) / _gss_tmp[i]), 1.0 / _xm(i)) *
969  std::copysign(1.0, _tau(i)) * _dt;
970  if (std::abs(_slip_incr(i)) > _slip_incr_tol)
971  {
972  _err_tol = true;
973 #ifdef DEBUG
974  mooseWarning("Maximum allowable slip increment exceeded ", std::abs(_slip_incr(i)));
975 #endif
976  return;
977  }
978  }
979 
980  for (unsigned int i = 0; i < _nss; ++i)
981  _dslipdtau(i) = _a0(i) / _xm(i) *
982  std::pow(std::abs(_tau(i) / _gss_tmp[i]), 1.0 / _xm(i) - 1.0) / _gss_tmp[i] *
983  _dt;
984 }
const unsigned int _nss
Number of slip system resistance.
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
Real _slip_incr_tol
Slip increment tolerance.

◆ getSlipSystems()

void FiniteStrainCrystalPlasticity::getSlipSystems ( )
protectedvirtual

This function reads slip system from file - see test.

Definition at line 425 of file FiniteStrainCrystalPlasticity.C.

Referenced by FiniteStrainCrystalPlasticity().

426 {
427  Real vec[LIBMESH_DIM];
428  std::ifstream fileslipsys;
429 
430  MooseUtils::checkFileReadable(_slip_sys_file_name);
431 
432  fileslipsys.open(_slip_sys_file_name.c_str());
433 
434  for (unsigned int i = 0; i < _nss; ++i)
435  {
436  // Read the slip normal
437  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
438  if (!(fileslipsys >> vec[j]))
439  mooseError("Crystal Plasticity Error: Premature end of file reading slip system file \n");
440 
441  // Normalize the vectors
442  Real mag;
443  mag = Utility::pow<2>(vec[0]) + Utility::pow<2>(vec[1]) + Utility::pow<2>(vec[2]);
444  mag = std::sqrt(mag);
445 
446  for (unsigned j = 0; j < LIBMESH_DIM; ++j)
447  _no(i * LIBMESH_DIM + j) = vec[j] / mag;
448 
449  // Read the slip direction
450  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
451  if (!(fileslipsys >> vec[j]))
452  mooseError("Crystal Plasticity Error: Premature end of file reading slip system file \n");
453 
454  // Normalize the vectors
455  mag = Utility::pow<2>(vec[0]) + Utility::pow<2>(vec[1]) + Utility::pow<2>(vec[2]);
456  mag = std::sqrt(mag);
457 
458  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
459  _mo(i * LIBMESH_DIM + j) = vec[j] / mag;
460 
461  mag = 0.0;
462  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
463  mag += _mo(i * LIBMESH_DIM + j) * _no(i * LIBMESH_DIM + j);
464 
465  if (std::abs(mag) > 1e-8)
466  mooseError(
467  "Crystal Plasicity Error: Slip direction and normal not orthonormal, System number = ",
468  i,
469  "\n");
470 
472  for (unsigned int j = 0; j < _num_slip_sys_props; ++j)
473  if (!(fileslipsys >> _slip_sys_props(i * _num_slip_sys_props + j)))
474  mooseError("Crystal Plasticity Error: Premature end of file reading slip system file - "
475  "check in slip system file read input options/values\n");
476  }
477 
478  fileslipsys.close();
479 }
const unsigned int _nss
Number of slip system resistance.
unsigned int _num_slip_sys_props
Number of slip system specific properties provided in the file containing slip system normals and dir...
std::string _slip_sys_file_name
File should contain slip plane normal and direction. See test.

◆ initAdditionalProps()

void FiniteStrainCrystalPlasticity::initAdditionalProps ( )
protectedvirtual

This function initializes additional parameters.

Definition at line 483 of file FiniteStrainCrystalPlasticity.C.

Referenced by initQpStatefulProperties().

484 {
485 }

◆ initQpStatefulProperties()

void FiniteStrainCrystalPlasticity::initQpStatefulProperties ( )
protectedvirtual

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

Reimplemented from ComputeStressBase.

Definition at line 188 of file FiniteStrainCrystalPlasticity.C.

189 {
190  _stress[_qp].zero();
191 
192  _fp[_qp].zero();
193  _fp[_qp].addIa(1.0);
194 
195  _pk2[_qp].zero();
196  _acc_slip[_qp] = 0.0;
197  _lag_e[_qp].zero();
198 
199  _update_rot[_qp].zero();
200  _update_rot[_qp].addIa(1.0);
201 
202  initSlipSysProps(); // Initializes slip system related properties
204 }
MaterialProperty< RankTwoTensor > & _lag_e
MaterialProperty< RankTwoTensor > & _pk2
MaterialProperty< RankTwoTensor > & _stress
MaterialProperty< RankTwoTensor > & _fp
virtual void initSlipSysProps()
This function initializes slip system resistances.
virtual void initAdditionalProps()
This function initializes additional parameters.
MaterialProperty< RankTwoTensor > & _update_rot

◆ initSlipSysProps()

void FiniteStrainCrystalPlasticity::initSlipSysProps ( )
protectedvirtual

This function initializes slip system resistances.

Definition at line 207 of file FiniteStrainCrystalPlasticity.C.

Referenced by initQpStatefulProperties().

208 {
209  switch (_intvar_read_type)
210  {
211  case 0:
213  break;
214  case 1:
216  break;
217  default:
219  }
220 
221  if (_slip_sys_flow_prop_file_name.length() != 0)
223  else
225 
226  if (_slip_sys_hard_prop_file_name.length() != 0)
228  else
230 }
virtual void getHardnessParams()
This function assign flow rate parameters from .i file - see test.
std::string _slip_sys_flow_prop_file_name
File should contain values of the flow rate equation parameters.
virtual void readFileHardnessParams()
This function read hardness parameters from file.
std::string _slip_sys_hard_prop_file_name
The hardening parameters in this class are read from .i file. The user can override to read from file...
virtual void assignSlipSysRes()
This function assign initial values of slip system resistances/internal variables read from getSlipSy...
virtual void getFlowRateParams()
This function assign flow rate parameters - see test.
virtual void readFileFlowRateParams()
This function read flow rate parameters from file - see test.
virtual void getInitSlipSysRes()
This function assign slip system resistances - see test.
virtual void readFileInitSlipSysRes()
This function read slip system resistances from file - see test.
MooseEnum _intvar_read_type
Read from options for initial values of internal variables.

◆ internalVariableUpdateNRiteration()

void FiniteStrainCrystalPlasticity::internalVariableUpdateNRiteration ( )
protected

This function updates internal variables after each NewTon Raphson iteration (_fp_inv)

Definition at line 1212 of file FiniteStrainCrystalPlasticity.C.

Referenced by solveStress().

1213 {
1214  _fp_prev_inv = _fp_inv; // update _fp_prev_inv
1215 }

◆ line_search_update()

bool FiniteStrainCrystalPlasticity::line_search_update ( const Real  rnorm_prev,
const RankTwoTensor  dpk2 
)
protected

This function performs the line search update.

Definition at line 1129 of file FiniteStrainCrystalPlasticity.C.

Referenced by solveStress().

1130 {
1131  if (_lsrch_method == "CUT_HALF")
1132  {
1133  Real rnorm;
1134  RankTwoTensor resid;
1135  Real step = 1.0;
1136 
1137  do
1138  {
1139  _pk2_tmp = _pk2_tmp - step * dpk2;
1140  step /= 2.0;
1141  _pk2_tmp = _pk2_tmp + step * dpk2;
1142 
1143  calcResidual(resid);
1144  rnorm = resid.L2norm();
1145  } while (rnorm > rnorm_prev && step > _min_lsrch_step);
1146 
1147  if (rnorm > rnorm_prev && step <= _min_lsrch_step)
1148  return false;
1149 
1150  return true;
1151  }
1152  else if (_lsrch_method == "BISECTION")
1153  {
1154  unsigned int count = 0;
1155  Real step_a = 0.0;
1156  Real step_b = 1.0;
1157  Real step = 1.0;
1158  Real s_m = 1000.0;
1159  Real rnorm = 1000.0;
1160 
1161  RankTwoTensor resid;
1162  calcResidual(resid);
1163  Real s_b = resid.doubleContraction(dpk2);
1164  Real rnorm1 = resid.L2norm();
1165  _pk2_tmp = _pk2_tmp - dpk2;
1166  calcResidual(resid);
1167  Real s_a = resid.doubleContraction(dpk2);
1168  Real rnorm0 = resid.L2norm();
1169  _pk2_tmp = _pk2_tmp + dpk2;
1170 
1171  if ((rnorm1 / rnorm0) < _lsrch_tol || s_a * s_b > 0)
1172  {
1173  calcResidual(resid);
1174  return true;
1175  }
1176 
1177  while ((rnorm / rnorm0) > _lsrch_tol && count < _lsrch_max_iter)
1178  {
1179  _pk2_tmp = _pk2_tmp - step * dpk2;
1180  step = 0.5 * (step_b + step_a);
1181  _pk2_tmp = _pk2_tmp + step * dpk2;
1182  calcResidual(resid);
1183  s_m = resid.doubleContraction(dpk2);
1184  rnorm = resid.L2norm();
1185 
1186  if (s_m * s_a < 0.0)
1187  {
1188  step_b = step;
1189  s_b = s_m;
1190  }
1191  if (s_m * s_b < 0.0)
1192  {
1193  step_a = step;
1194  s_a = s_m;
1195  }
1196  count++;
1197  }
1198 
1199  if ((rnorm / rnorm0) < _lsrch_tol && count < _lsrch_max_iter)
1200  return true;
1201 
1202  return false;
1203  }
1204  else
1205  {
1206  mooseError("Line search meothod is not provided.");
1207  return false;
1208  }
1209 }
virtual void calcResidual(RankTwoTensor &)
This function calculate stress residual.
unsigned int _lsrch_max_iter
Line search bisection method maximum iteration number.
Real _lsrch_tol
Line search bisection method tolerance.
Real _min_lsrch_step
Minimum line search step size.

◆ postSolveQp()

void FiniteStrainCrystalPlasticity::postSolveQp ( )
protectedvirtual

This function update stress and internal variable after solve.

Definition at line 587 of file FiniteStrainCrystalPlasticity.C.

Referenced by computeQpStress().

588 {
589  if (_err_tol)
590  {
591  _err_tol = false;
593  {
595  _rndm_scale_var = _elasticity_tensor[_qp](0, 0, 0, 0);
596 
597  _stress[_qp] = RankTwoTensor::genRandomSymmTensor(_rndm_scale_var, 1.0);
598  }
599  else
600  mooseError("FiniteStrainCrystalPlasticity: Constitutive failure");
601  }
602  else
603  {
604  _stress[_qp] = _fe * _pk2[_qp] * _fe.transpose() / _fe.det();
605 
606  _Jacobian_mult[_qp] += calcTangentModuli(); // Calculate jacobian for preconditioner
607 
608  RankTwoTensor iden;
609  iden.addIa(1.0);
610 
611  _lag_e[_qp] = _deformation_gradient[_qp].transpose() * _deformation_gradient[_qp] - iden;
612  _lag_e[_qp] = _lag_e[_qp] * 0.5;
613 
614  RankTwoTensor rot;
615  rot = get_current_rotation(_deformation_gradient[_qp]); // Calculate material rotation
616  _update_rot[_qp] = rot * _crysrot[_qp];
617  }
618 }
MaterialProperty< RankFourTensor > & _Jacobian_mult
derivative of stress w.r.t. strain (_dstress_dstrain)
MaterialProperty< RankTwoTensor > & _lag_e
const MaterialProperty< RankTwoTensor > & _deformation_gradient
bool _input_rndm_scale_var
Input option for scaling variable to generate random stress when convergence fails.
const MaterialProperty< RankFourTensor > & _elasticity_tensor
MaterialProperty< RankTwoTensor > & _pk2
MaterialProperty< RankTwoTensor > & _stress
virtual RankFourTensor calcTangentModuli()
This function calculate the tangent moduli for preconditioner.
const MaterialProperty< RankTwoTensor > & _crysrot
MaterialProperty< RankTwoTensor > & _update_rot
RankTwoTensor get_current_rotation(const RankTwoTensor &a)
This function perform RU decomposition to obtain the rotation tensor.

◆ postSolveStatevar()

void FiniteStrainCrystalPlasticity::postSolveStatevar ( )
protectedvirtual

This function update internal variable after solve.

Definition at line 683 of file FiniteStrainCrystalPlasticity.C.

Referenced by solveQp().

684 {
685  if (_max_substep_iter == 1) // No substepping
686  {
687  _gss[_qp] = _gss_tmp;
688  _acc_slip[_qp] = _accslip_tmp;
689  }
690  else
691  {
692  if (_last_step_iter)
693  {
694  _gss[_qp] = _gss_tmp;
695  _acc_slip[_qp] = _accslip_tmp;
696  }
697  else
698  {
701  }
702  }
703 }
MaterialProperty< std::vector< Real > > & _gss
unsigned int _max_substep_iter
Maximum number of substep iterations.

◆ postSolveStress()

void FiniteStrainCrystalPlasticity::postSolveStress ( )
protectedvirtual

This function update stress and plastic deformation gradient after solve.

Definition at line 802 of file FiniteStrainCrystalPlasticity.C.

Referenced by FiniteStrainCPSlipRateRes::solveStatevar(), and solveStatevar().

803 {
804  if (_max_substep_iter == 1) // No substepping
805  {
806  _fp[_qp] = _fp_inv.inverse();
807  _pk2[_qp] = _pk2_tmp;
808  }
809  else
810  {
811  if (_last_step_iter)
812  {
813  _fp[_qp] = _fp_inv.inverse();
814  _pk2[_qp] = _pk2_tmp;
815  }
816  else
817  {
820  }
821  }
822 }
MaterialProperty< RankTwoTensor > & _pk2
MaterialProperty< RankTwoTensor > & _fp
unsigned int _max_substep_iter
Maximum number of substep iterations.

◆ preSolveQp()

void FiniteStrainCrystalPlasticity::preSolveQp ( )
protectedvirtual

This function set variables for stress and internal variable solve.

Definition at line 559 of file FiniteStrainCrystalPlasticity.C.

Referenced by computeQpStress().

560 {
561  // Initialize variable
562  if (_first_substep)
563  {
564  _Jacobian_mult[_qp].zero(); // Initializes jacobian for preconditioner
566  }
567 
568  if (_max_substep_iter == 1)
569  _dfgrd_tmp = _deformation_gradient[_qp]; // Without substepping
570  else
572 
573  _err_tol = false;
574 }
MaterialProperty< RankFourTensor > & _Jacobian_mult
derivative of stress w.r.t. strain (_dstress_dstrain)
const MaterialProperty< RankTwoTensor > & _deformation_gradient
Real _dfgrd_scale_factor
Scales the substepping increment to obtain deformation gradient at a substep iteration.
RankTwoTensor _delta_dfgrd
Flag to check whether convergence is achieved.
void calc_schmid_tensor()
This function calculate the Schmid tensor.
unsigned int _max_substep_iter
Maximum number of substep iterations.

◆ preSolveStatevar()

void FiniteStrainCrystalPlasticity::preSolveStatevar ( )
protectedvirtual

This function set variables for internal variable solve.

Definition at line 621 of file FiniteStrainCrystalPlasticity.C.

Referenced by solveQp().

622 {
623  if (_max_substep_iter == 1) // No substepping
624  {
625  _gss_tmp = _gss_old[_qp];
627  }
628  else
629  {
630  if (_first_step_iter)
631  {
632  _gss_tmp = _gss_tmp_old = _gss_old[_qp];
634  }
635  else
637  }
638 }
const MaterialProperty< std::vector< Real > > & _gss_old
unsigned int _max_substep_iter
Maximum number of substep iterations.
bool _first_step_iter
Flags to reset variables and reinitialize variables.
const MaterialProperty< Real > & _acc_slip_old

◆ preSolveStress()

void FiniteStrainCrystalPlasticity::preSolveStress ( )
protectedvirtual

This function set variables for stress solve.

Reimplemented in FiniteStrainCPSlipRateRes.

Definition at line 706 of file FiniteStrainCrystalPlasticity.C.

Referenced by FiniteStrainCPSlipRateRes::preSolveStress(), and solveStatevar().

707 {
708  if (_max_substep_iter == 1) // No substepping
709  {
710  _pk2_tmp = _pk2_old[_qp];
711  _fp_old_inv = _fp_old[_qp].inverse();
714  }
715  else
716  {
717  if (_first_step_iter)
718  {
719  _pk2_tmp = _pk2_tmp_old = _pk2_old[_qp];
720  _fp_old_inv = _fp_old[_qp].inverse();
721  }
722  else
724 
727  }
728 }
const MaterialProperty< RankTwoTensor > & _pk2_old
const MaterialProperty< RankTwoTensor > & _fp_old
unsigned int _max_substep_iter
Maximum number of substep iterations.
bool _first_step_iter
Flags to reset variables and reinitialize variables.

◆ readFileFlowRateParams()

void FiniteStrainCrystalPlasticity::readFileFlowRateParams ( )
protectedvirtual

This function read flow rate parameters from file - see test.

Definition at line 314 of file FiniteStrainCrystalPlasticity.C.

Referenced by initSlipSysProps().

315 {
316  _a0.resize(_nss);
317  _xm.resize(_nss);
318 
319  MooseUtils::checkFileReadable(_slip_sys_flow_prop_file_name);
320 
321  std::ifstream file;
322  file.open(_slip_sys_flow_prop_file_name.c_str());
323 
324  std::vector<Real> vec;
325  vec.resize(_num_slip_sys_flowrate_props);
326 
327  for (unsigned int i = 0; i < _nss; ++i)
328  {
329  for (unsigned int j = 0; j < _num_slip_sys_flowrate_props; ++j)
330  if (!(file >> vec[j]))
331  mooseError(
332  "Error FiniteStrainCrystalPlasticity: Premature end of slip_sys_flow_rate_param file");
333 
334  _a0(i) = vec[0];
335  _xm(i) = vec[1];
336  }
337 
338  file.close();
339 }
std::string _slip_sys_flow_prop_file_name
File should contain values of the flow rate equation parameters.
unsigned int _num_slip_sys_flowrate_props
Number of slip system flow rate parameters.
const unsigned int _nss
Number of slip system resistance.

◆ readFileHardnessParams()

void FiniteStrainCrystalPlasticity::readFileHardnessParams ( )
protectedvirtual

This function read hardness parameters from file.

Definition at line 405 of file FiniteStrainCrystalPlasticity.C.

Referenced by initSlipSysProps().

406 {
407 }

◆ readFileInitSlipSysRes()

void FiniteStrainCrystalPlasticity::readFileInitSlipSysRes ( )
protectedvirtual

This function read slip system resistances from file - see test.

Definition at line 243 of file FiniteStrainCrystalPlasticity.C.

Referenced by initSlipSysProps().

244 {
245  _gss[_qp].resize(_nss);
246 
247  MooseUtils::checkFileReadable(_slip_sys_res_prop_file_name);
248 
249  std::ifstream file;
250  file.open(_slip_sys_res_prop_file_name.c_str());
251 
252  for (unsigned int i = 0; i < _nss; ++i)
253  if (!(file >> _gss[_qp][i]))
254  mooseError("Error FiniteStrainCrystalPlasticity: Premature end of slip_sys_res_prop file");
255 
256  file.close();
257 }
const unsigned int _nss
Number of slip system resistance.
MaterialProperty< std::vector< Real > > & _gss
std::string _slip_sys_res_prop_file_name
File should contain initial values of the slip system resistances.

◆ solveQp()

void FiniteStrainCrystalPlasticity::solveQp ( )
protectedvirtual

This function solves stress and internal variables.

Definition at line 577 of file FiniteStrainCrystalPlasticity.C.

Referenced by computeQpStress().

578 {
580  solveStatevar();
581  if (_err_tol)
582  return;
584 }
virtual void solveStatevar()
This function solves internal variables.
virtual void postSolveStatevar()
This function update internal variable after solve.
virtual void preSolveStatevar()
This function set variables for internal variable solve.

◆ solveStatevar()

void FiniteStrainCrystalPlasticity::solveStatevar ( )
protectedvirtual

This function solves internal variables.

Reimplemented in FiniteStrainCPSlipRateRes.

Definition at line 641 of file FiniteStrainCrystalPlasticity.C.

Referenced by solveQp().

642 {
643  Real gmax, gdiff;
644  unsigned int iterg;
645  std::vector<Real> gss_prev(_nss);
646 
647  gmax = 1.1 * _gtol;
648  iterg = 0;
649 
650  while (gmax > _gtol && iterg < _maxiterg) // Check for slip system resistance update tolerance
651  {
652  preSolveStress();
653  solveStress();
654  if (_err_tol)
655  return;
656  postSolveStress();
657 
658  gss_prev = _gss_tmp;
659 
660  update_slip_system_resistance(); // Update slip system resistance
661 
662  gmax = 0.0;
663  for (unsigned i = 0; i < _nss; ++i)
664  {
665  gdiff = std::abs(gss_prev[i] - _gss_tmp[i]); // Calculate increment size
666 
667  if (gdiff > gmax)
668  gmax = gdiff;
669  }
670  iterg++;
671  }
672 
673  if (iterg == _maxiterg)
674  {
675 #ifdef DEBUG
676  mooseWarning("FiniteStrainCrystalPLasticity: Hardness Integration error gmax", gmax, "\n");
677 #endif
678  _err_tol = true;
679  }
680 }
virtual void solveStress()
This function solves for stress, updates plastic deformation gradient.
virtual void update_slip_system_resistance()
This function updates the slip system resistances.
const unsigned int _nss
Number of slip system resistance.
Real _gtol
Internal variable update equation tolerance.
virtual void preSolveStress()
This function set variables for stress solve.
virtual void postSolveStress()
This function update stress and plastic deformation gradient after solve.
unsigned int _maxiterg
Maximum number of iterations for internal variable update.

◆ solveStress()

void FiniteStrainCrystalPlasticity::solveStress ( )
protectedvirtual

This function solves for stress, updates plastic deformation gradient.

Reimplemented in FiniteStrainCPSlipRateRes.

Definition at line 731 of file FiniteStrainCrystalPlasticity.C.

Referenced by solveStatevar().

732 {
733  unsigned int iter = 0;
734  RankTwoTensor resid, dpk2;
735  RankFourTensor jac;
736  Real rnorm, rnorm0, rnorm_prev;
737 
738  calc_resid_jacob(resid, jac); // Calculate stress residual
739  if (_err_tol)
740  {
741 #ifdef DEBUG
742  mooseWarning(
743  "FiniteStrainCrystalPLasticity: Slip increment exceeds tolerance - Element number ",
744  _current_elem->id(),
745  " Gauss point = ",
746  _qp);
747 #endif
748  return;
749  }
750 
751  rnorm = resid.L2norm();
752  rnorm0 = rnorm;
753 
754  while (rnorm > _rtol * rnorm0 && rnorm0 > _abs_tol &&
755  iter < _maxiter) // Check for stress residual tolerance
756  {
757  dpk2 = -jac.invSymm() * resid; // Calculate stress increment
758  _pk2_tmp = _pk2_tmp + dpk2; // Update stress
759  calc_resid_jacob(resid, jac);
760  internalVariableUpdateNRiteration(); // update _fp_prev_inv
761 
762  if (_err_tol)
763  {
764 #ifdef DEBUG
765  mooseWarning(
766  "FiniteStrainCrystalPLasticity: Slip increment exceeds tolerance - Element number ",
767  _current_elem->id(),
768  " Gauss point = ",
769  _qp);
770 #endif
771  return;
772  }
773 
774  rnorm_prev = rnorm;
775  rnorm = resid.L2norm();
776 
777  if (_use_line_search && rnorm > rnorm_prev && !line_search_update(rnorm_prev, dpk2))
778  {
779 #ifdef DEBUG
780  mooseWarning("FiniteStrainCrystalPLasticity: Failed with line search");
781 #endif
782  _err_tol = true;
783  return;
784  }
785 
786  if (_use_line_search)
787  rnorm = resid.L2norm();
788 
789  iter++;
790  }
791 
792  if (iter >= _maxiter)
793  {
794 #ifdef DEBUG
795  mooseWarning("FiniteStrainCrystalPLasticity: Stress Integration error rmax = ", rnorm);
796 #endif
797  _err_tol = true;
798  }
799 }
void internalVariableUpdateNRiteration()
This function updates internal variables after each NewTon Raphson iteration (_fp_inv) ...
bool _use_line_search
Flag to activate line serach.
Real _rtol
Stress residual equation relative tolerance.
unsigned int _maxiter
Maximum number of iterations for stress update.
bool line_search_update(const Real rnorm_prev, const RankTwoTensor)
This function performs the line search update.
virtual void calc_resid_jacob(RankTwoTensor &, RankFourTensor &)
This function calls the residual and jacobian functions used in the stress update algorithm...
Real _abs_tol
Stress residual equation absolute tolerance.

◆ update_slip_system_resistance()

void FiniteStrainCrystalPlasticity::update_slip_system_resistance ( )
protectedvirtual

This function updates the slip system resistances.

Definition at line 826 of file FiniteStrainCrystalPlasticity.C.

Referenced by FiniteStrainCPSlipRateRes::calcResidualSlipRate(), and solveStatevar().

827 {
828  updateGss();
829 }
virtual void updateGss()
This function updates the slip system resistances.

◆ updateGss()

void FiniteStrainCrystalPlasticity::updateGss ( )
protectedvirtual

This function updates the slip system resistances.

Old function to update slip system resistances.

Kept to avoid code break at computeQpstress

Definition at line 836 of file FiniteStrainCrystalPlasticity.C.

Referenced by update_slip_system_resistance().

837 {
838  DenseVector<Real> hb(_nss);
839  Real qab;
840 
841  Real a = _hprops[4]; // Kalidindi
842 
844  for (unsigned int i = 0; i < _nss; ++i)
845  _accslip_tmp += std::abs(_slip_incr(i));
846 
847  // Real val = std::cosh(_h0 * _accslip_tmp / (_tau_sat - _tau_init)); // Karthik
848  // val = _h0 * std::pow(1.0/val,2.0); // Kalidindi
849 
850  for (unsigned int i = 0; i < _nss; ++i)
851  // hb(i)=val;
852  hb(i) = _h0 * std::pow(std::abs(1.0 - _gss_tmp[i] / _tau_sat), a) *
853  std::copysign(1.0, 1.0 - _gss_tmp[i] / _tau_sat);
854 
855  for (unsigned int i = 0; i < _nss; ++i)
856  {
857  if (_max_substep_iter == 1) // No substepping
858  _gss_tmp[i] = _gss_old[_qp][i];
859  else
860  _gss_tmp[i] = _gss_tmp_old[i];
861 
862  for (unsigned int j = 0; j < _nss; ++j)
863  {
864  unsigned int iplane, jplane;
865  iplane = i / 3;
866  jplane = j / 3;
867 
868  if (iplane == jplane) // Kalidindi
869  qab = 1.0;
870  else
871  qab = _r;
872 
873  _gss_tmp[i] += qab * hb(j) * std::abs(_slip_incr(j));
874  _dgss_dsliprate(i, j) = qab * hb(j) * std::copysign(1.0, _slip_incr(j)) * _dt;
875  }
876  }
877 }
const MaterialProperty< std::vector< Real > > & _gss_old
const unsigned int _nss
Number of slip system resistance.
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
unsigned int _max_substep_iter
Maximum number of substep iterations.

Member Data Documentation

◆ _a0

DenseVector<Real> FiniteStrainCrystalPlasticity::_a0
protected

◆ _abs_tol

Real FiniteStrainCrystalPlasticity::_abs_tol
protected

Stress residual equation absolute tolerance.

Definition at line 251 of file FiniteStrainCrystalPlasticity.h.

Referenced by FiniteStrainCPSlipRateRes::solveStress(), and solveStress().

◆ _acc_slip

MaterialProperty<Real>& FiniteStrainCrystalPlasticity::_acc_slip
protected

Definition at line 311 of file FiniteStrainCrystalPlasticity.h.

Referenced by initQpStatefulProperties(), and postSolveStatevar().

◆ _acc_slip_old

const MaterialProperty<Real>& FiniteStrainCrystalPlasticity::_acc_slip_old
protected

Definition at line 312 of file FiniteStrainCrystalPlasticity.h.

Referenced by preSolveStatevar().

◆ _accslip_tmp

Real FiniteStrainCrystalPlasticity::_accslip_tmp
protected

Definition at line 337 of file FiniteStrainCrystalPlasticity.h.

Referenced by postSolveStatevar(), and updateGss().

◆ _accslip_tmp_old

Real FiniteStrainCrystalPlasticity::_accslip_tmp_old
protected

Definition at line 337 of file FiniteStrainCrystalPlasticity.h.

Referenced by postSolveStatevar(), preSolveStatevar(), and updateGss().

◆ _base_name

const std::string ComputeStressBase::_base_name
protectedinherited

◆ _crysrot

const MaterialProperty<RankTwoTensor>& FiniteStrainCrystalPlasticity::_crysrot
protected

Definition at line 318 of file FiniteStrainCrystalPlasticity.h.

Referenced by calc_schmid_tensor(), and postSolveQp().

◆ _deformation_gradient

const MaterialProperty<RankTwoTensor>& FiniteStrainCrystalPlasticity::_deformation_gradient
protected

Definition at line 315 of file FiniteStrainCrystalPlasticity.h.

Referenced by computeQpStress(), postSolveQp(), and preSolveQp().

◆ _deformation_gradient_old

const MaterialProperty<RankTwoTensor>& FiniteStrainCrystalPlasticity::_deformation_gradient_old
protected

Definition at line 316 of file FiniteStrainCrystalPlasticity.h.

Referenced by computeQpStress().

◆ _delta_dfgrd

RankTwoTensor FiniteStrainCrystalPlasticity::_delta_dfgrd
protected

Flag to check whether convergence is achieved.

Used for substepping; Uniformly divides the increment in deformation gradient

Definition at line 350 of file FiniteStrainCrystalPlasticity.h.

Referenced by computeQpStress(), FiniteStrainCrystalPlasticity(), and preSolveQp().

◆ _dfgrd_scale_factor

Real FiniteStrainCrystalPlasticity::_dfgrd_scale_factor
protected

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

Definition at line 352 of file FiniteStrainCrystalPlasticity.h.

Referenced by computeQpStress(), and preSolveQp().

◆ _dfgrd_tmp

RankTwoTensor FiniteStrainCrystalPlasticity::_dfgrd_tmp
protected

◆ _dfgrd_tmp_old

RankTwoTensor FiniteStrainCrystalPlasticity::_dfgrd_tmp_old
protected

Definition at line 350 of file FiniteStrainCrystalPlasticity.h.

Referenced by computeQpStress(), and preSolveQp().

◆ _dgss_dsliprate

DenseMatrix<Real> FiniteStrainCrystalPlasticity::_dgss_dsliprate
protected

◆ _dslipdtau

DenseVector<Real> FiniteStrainCrystalPlasticity::_dslipdtau
protected

◆ _elastic_strain

MaterialProperty<RankTwoTensor>& ComputeStressBase::_elastic_strain
protectedinherited

◆ _elasticity_tensor

const MaterialProperty<RankFourTensor>& FiniteStrainCrystalPlasticity::_elasticity_tensor
protected

◆ _elasticity_tensor_name

const std::string ComputeStressBase::_elasticity_tensor_name
protectedinherited

◆ _err_tol

bool FiniteStrainCrystalPlasticity::_err_tol
protected

◆ _extra_stress

const MaterialProperty<RankTwoTensor>& ComputeStressBase::_extra_stress
protectedinherited

Extra stress tensor.

Definition at line 47 of file ComputeStressBase.h.

Referenced by ComputeStressBase::computeQpProperties().

◆ _fe

RankTwoTensor FiniteStrainCrystalPlasticity::_fe
protected

◆ _first_step_iter

bool FiniteStrainCrystalPlasticity::_first_step_iter
protected

Flags to reset variables and reinitialize variables.

Definition at line 354 of file FiniteStrainCrystalPlasticity.h.

Referenced by computeQpStress(), FiniteStrainCrystalPlasticity(), preSolveStatevar(), and preSolveStress().

◆ _first_substep

bool FiniteStrainCrystalPlasticity::_first_substep
protected

◆ _flowprops

std::vector<Real> FiniteStrainCrystalPlasticity::_flowprops
protected

Definition at line 231 of file FiniteStrainCrystalPlasticity.h.

Referenced by getFlowRateParams().

◆ _fp

MaterialProperty<RankTwoTensor>& FiniteStrainCrystalPlasticity::_fp
protected

Definition at line 303 of file FiniteStrainCrystalPlasticity.h.

Referenced by initQpStatefulProperties(), and postSolveStress().

◆ _fp_inv

RankTwoTensor FiniteStrainCrystalPlasticity::_fp_inv
protected

◆ _fp_old

const MaterialProperty<RankTwoTensor>& FiniteStrainCrystalPlasticity::_fp_old
protected

Definition at line 304 of file FiniteStrainCrystalPlasticity.h.

Referenced by preSolveStress().

◆ _fp_old_inv

RankTwoTensor FiniteStrainCrystalPlasticity::_fp_old_inv
protected

◆ _fp_prev_inv

RankTwoTensor FiniteStrainCrystalPlasticity::_fp_prev_inv
protected

◆ _gen_rndm_stress_flag

bool FiniteStrainCrystalPlasticity::_gen_rndm_stress_flag
protected

Definition at line 274 of file FiniteStrainCrystalPlasticity.h.

Referenced by postSolveQp().

◆ _gprops

std::vector<Real> FiniteStrainCrystalPlasticity::_gprops
protected

Definition at line 229 of file FiniteStrainCrystalPlasticity.h.

Referenced by getInitSlipSysRes().

◆ _gss

MaterialProperty<std::vector<Real> >& FiniteStrainCrystalPlasticity::_gss
protected

◆ _gss_old

const MaterialProperty<std::vector<Real> >& FiniteStrainCrystalPlasticity::_gss_old
protected

Definition at line 310 of file FiniteStrainCrystalPlasticity.h.

Referenced by preSolveStatevar(), and updateGss().

◆ _gss_tmp

std::vector<Real> FiniteStrainCrystalPlasticity::_gss_tmp
protected

◆ _gss_tmp_old

std::vector<Real> FiniteStrainCrystalPlasticity::_gss_tmp_old
protected

Definition at line 339 of file FiniteStrainCrystalPlasticity.h.

Referenced by postSolveStatevar(), preSolveStatevar(), and updateGss().

◆ _gtol

Real FiniteStrainCrystalPlasticity::_gtol
protected

Internal variable update equation tolerance.

Definition at line 253 of file FiniteStrainCrystalPlasticity.h.

Referenced by solveStatevar().

◆ _h0

Real FiniteStrainCrystalPlasticity::_h0
protected

Definition at line 326 of file FiniteStrainCrystalPlasticity.h.

Referenced by getHardnessParams(), and updateGss().

◆ _hprops

std::vector<Real> FiniteStrainCrystalPlasticity::_hprops
protected

Definition at line 230 of file FiniteStrainCrystalPlasticity.h.

Referenced by getHardnessParams(), and updateGss().

◆ _initial_stress_fcn

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

initial stress components

Definition at line 50 of file ComputeStressBase.h.

◆ _input_rndm_scale_var

bool FiniteStrainCrystalPlasticity::_input_rndm_scale_var
protected

Input option for scaling variable to generate random stress when convergence fails.

Definition at line 277 of file FiniteStrainCrystalPlasticity.h.

Referenced by postSolveQp().

◆ _intvar_read_type

MooseEnum FiniteStrainCrystalPlasticity::_intvar_read_type
protected

Read from options for initial values of internal variables.

Definition at line 269 of file FiniteStrainCrystalPlasticity.h.

Referenced by FiniteStrainCrystalPlasticity(), and initSlipSysProps().

◆ _Jacobian_mult

MaterialProperty<RankFourTensor>& ComputeStressBase::_Jacobian_mult
protectedinherited

◆ _lag_e

MaterialProperty<RankTwoTensor>& FiniteStrainCrystalPlasticity::_lag_e
protected

Definition at line 307 of file FiniteStrainCrystalPlasticity.h.

Referenced by initQpStatefulProperties(), and postSolveQp().

◆ _lag_e_old

const MaterialProperty<RankTwoTensor>& FiniteStrainCrystalPlasticity::_lag_e_old
protected

Definition at line 308 of file FiniteStrainCrystalPlasticity.h.

◆ _last_step_iter

bool FiniteStrainCrystalPlasticity::_last_step_iter
protected

◆ _lsrch_max_iter

unsigned int FiniteStrainCrystalPlasticity::_lsrch_max_iter
protected

Line search bisection method maximum iteration number.

Definition at line 298 of file FiniteStrainCrystalPlasticity.h.

Referenced by line_search_update(), and FiniteStrainCPSlipRateRes::lineSearchUpdateSlipRate().

◆ _lsrch_method

MooseEnum FiniteStrainCrystalPlasticity::_lsrch_method
protected

◆ _lsrch_tol

Real FiniteStrainCrystalPlasticity::_lsrch_tol
protected

Line search bisection method tolerance.

Definition at line 295 of file FiniteStrainCrystalPlasticity.h.

Referenced by line_search_update(), and FiniteStrainCPSlipRateRes::lineSearchUpdateSlipRate().

◆ _max_substep_iter

unsigned int FiniteStrainCrystalPlasticity::_max_substep_iter
protected

Maximum number of substep iterations.

Definition at line 286 of file FiniteStrainCrystalPlasticity.h.

Referenced by computeQpStress(), postSolveStatevar(), postSolveStress(), preSolveQp(), preSolveStatevar(), preSolveStress(), and updateGss().

◆ _maxiter

unsigned int FiniteStrainCrystalPlasticity::_maxiter
protected

Maximum number of iterations for stress update.

Definition at line 258 of file FiniteStrainCrystalPlasticity.h.

Referenced by FiniteStrainCPSlipRateRes::solveStress(), and solveStress().

◆ _maxiterg

unsigned int FiniteStrainCrystalPlasticity::_maxiterg
protected

Maximum number of iterations for internal variable update.

Definition at line 260 of file FiniteStrainCrystalPlasticity.h.

Referenced by solveStatevar().

◆ _mechanical_strain

const MaterialProperty<RankTwoTensor>& ComputeStressBase::_mechanical_strain
protectedinherited

◆ _min_lsrch_step

Real FiniteStrainCrystalPlasticity::_min_lsrch_step
protected

Minimum line search step size.

Definition at line 292 of file FiniteStrainCrystalPlasticity.h.

Referenced by line_search_update(), and FiniteStrainCPSlipRateRes::lineSearchUpdateSlipRate().

◆ _mo

DenseVector<Real> FiniteStrainCrystalPlasticity::_mo
protected

Definition at line 320 of file FiniteStrainCrystalPlasticity.h.

Referenced by calc_schmid_tensor(), and getSlipSystems().

◆ _no

DenseVector<Real> FiniteStrainCrystalPlasticity::_no
protected

Definition at line 321 of file FiniteStrainCrystalPlasticity.h.

Referenced by calc_schmid_tensor(), and getSlipSystems().

◆ _nss

const unsigned int FiniteStrainCrystalPlasticity::_nss
protected

◆ _num_slip_sys_flowrate_props

unsigned int FiniteStrainCrystalPlasticity::_num_slip_sys_flowrate_props
protected

Number of slip system flow rate parameters.

Definition at line 263 of file FiniteStrainCrystalPlasticity.h.

Referenced by getFlowRateParams(), and readFileFlowRateParams().

◆ _num_slip_sys_props

unsigned int FiniteStrainCrystalPlasticity::_num_slip_sys_props
protected

Number of slip system specific properties provided in the file containing slip system normals and directions.

Definition at line 272 of file FiniteStrainCrystalPlasticity.h.

Referenced by FiniteStrainCrystalPlasticity(), and getSlipSystems().

◆ _pk2

MaterialProperty<RankTwoTensor>& FiniteStrainCrystalPlasticity::_pk2
protected

◆ _pk2_old

const MaterialProperty<RankTwoTensor>& FiniteStrainCrystalPlasticity::_pk2_old
protected

Definition at line 306 of file FiniteStrainCrystalPlasticity.h.

Referenced by preSolveStress().

◆ _pk2_tmp

RankTwoTensor FiniteStrainCrystalPlasticity::_pk2_tmp
protected

◆ _pk2_tmp_old

RankTwoTensor FiniteStrainCrystalPlasticity::_pk2_tmp_old
protected

Definition at line 336 of file FiniteStrainCrystalPlasticity.h.

Referenced by postSolveStress(), and preSolveStress().

◆ _r

Real FiniteStrainCrystalPlasticity::_r
protected

Definition at line 329 of file FiniteStrainCrystalPlasticity.h.

Referenced by getHardnessParams(), and updateGss().

◆ _read_from_slip_sys_file

bool FiniteStrainCrystalPlasticity::_read_from_slip_sys_file
protected

◆ _rndm_scale_var

Real FiniteStrainCrystalPlasticity::_rndm_scale_var
protected

Scaling value.

Definition at line 280 of file FiniteStrainCrystalPlasticity.h.

Referenced by postSolveQp().

◆ _rndm_seed

unsigned int FiniteStrainCrystalPlasticity::_rndm_seed
protected

Seed value.

Definition at line 283 of file FiniteStrainCrystalPlasticity.h.

Referenced by FiniteStrainCrystalPlasticity().

◆ _rtol

Real FiniteStrainCrystalPlasticity::_rtol
protected

Stress residual equation relative tolerance.

Definition at line 249 of file FiniteStrainCrystalPlasticity.h.

Referenced by FiniteStrainCPSlipRateRes::solveStress(), and solveStress().

◆ _s0

std::vector<RankTwoTensor> FiniteStrainCrystalPlasticity::_s0
protected

◆ _slip_incr

DenseVector<Real> FiniteStrainCrystalPlasticity::_slip_incr
protected

◆ _slip_incr_tol

Real FiniteStrainCrystalPlasticity::_slip_incr_tol
protected

Slip increment tolerance.

Definition at line 255 of file FiniteStrainCrystalPlasticity.h.

Referenced by getSlipIncrements().

◆ _slip_sys_file_name

std::string FiniteStrainCrystalPlasticity::_slip_sys_file_name
protected

File should contain slip plane normal and direction. See test.

Definition at line 234 of file FiniteStrainCrystalPlasticity.h.

Referenced by getSlipSystems().

◆ _slip_sys_flow_prop_file_name

std::string FiniteStrainCrystalPlasticity::_slip_sys_flow_prop_file_name
protected

File should contain values of the flow rate equation parameters.

Values for every slip system must be provided. Should have the same order of slip systens as in slip_sys_file. See test. The option of reading all the properties from .i is still present.

Definition at line 243 of file FiniteStrainCrystalPlasticity.h.

Referenced by initSlipSysProps(), and readFileFlowRateParams().

◆ _slip_sys_hard_prop_file_name

std::string FiniteStrainCrystalPlasticity::_slip_sys_hard_prop_file_name
protected

The hardening parameters in this class are read from .i file. The user can override to read from file.

Definition at line 246 of file FiniteStrainCrystalPlasticity.h.

Referenced by initSlipSysProps().

◆ _slip_sys_props

DenseVector<Real> FiniteStrainCrystalPlasticity::_slip_sys_props
protected

◆ _slip_sys_res_prop_file_name

std::string FiniteStrainCrystalPlasticity::_slip_sys_res_prop_file_name
protected

File should contain initial values of the slip system resistances.

Definition at line 237 of file FiniteStrainCrystalPlasticity.h.

Referenced by readFileInitSlipSysRes().

◆ _stress

MaterialProperty<RankTwoTensor>& ComputeStressBase::_stress
protectedinherited

◆ _tan_mod_type

MooseEnum FiniteStrainCrystalPlasticity::_tan_mod_type
protected

Type of tangent moduli calculation.

Definition at line 266 of file FiniteStrainCrystalPlasticity.h.

Referenced by calcTangentModuli().

◆ _tau

DenseVector<Real> FiniteStrainCrystalPlasticity::_tau
protected

◆ _tau_init

Real FiniteStrainCrystalPlasticity::_tau_init
protected

Definition at line 328 of file FiniteStrainCrystalPlasticity.h.

Referenced by getHardnessParams().

◆ _tau_sat

Real FiniteStrainCrystalPlasticity::_tau_sat
protected

Definition at line 327 of file FiniteStrainCrystalPlasticity.h.

Referenced by getHardnessParams(), and updateGss().

◆ _update_rot

MaterialProperty<RankTwoTensor>& FiniteStrainCrystalPlasticity::_update_rot
protected

Definition at line 313 of file FiniteStrainCrystalPlasticity.h.

Referenced by initQpStatefulProperties(), and postSolveQp().

◆ _use_line_search

bool FiniteStrainCrystalPlasticity::_use_line_search
protected

Flag to activate line serach.

Definition at line 289 of file FiniteStrainCrystalPlasticity.h.

Referenced by FiniteStrainCPSlipRateRes::solveStress(), and solveStress().

◆ _xm

DenseVector<Real> FiniteStrainCrystalPlasticity::_xm
protected

The documentation for this class was generated from the following files: