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

RateDepSmearCrackModel is the base class for rate dependent continuum damage model. More...

#include <RateDepSmearCrackModel.h>

Inheritance diagram for RateDepSmearCrackModel:
[legend]

Public Member Functions

 RateDepSmearCrackModel (const InputParameters &parameters)
 
virtual ~RateDepSmearCrackModel ()
 
void setQp (unsigned int qp)
 Sets the value of the variable _qp for inheriting classes. More...
 
virtual bool modifyStrainIncrement (const Elem &, SymmTensor &strain_increment, SymmTensor &d_strain_dT)
 
virtual bool updateElasticityTensor (SymmElasticityTensor &)
 
virtual bool applyThermalStrain (SymmTensor &strain_increment, SymmTensor &d_strain_dT)
 

Protected Member Functions

virtual void computeStress (const Elem &current_elem, const SymmElasticityTensor &elasticity_tensor, const SymmTensor &stress_old, SymmTensor &strain_increment, SymmTensor &stress_new)
 
virtual void initQpStatefulProperties ()
 
virtual void initVariables ()
 
virtual void solve ()
 This function solves the state variables. More...
 
virtual void postSolveVariables ()
 This function updates the internal variables after solve. More...
 
virtual void postSolveStress ()
 This function updates the stress after solve. More...
 
virtual void updateVariables ()
 This function updates variables during solve a(i+1) = a(i) + da(i+1) More...
 
bool getConvergeVar ()
 This function returns true if convergence is not achieved. More...
 
virtual void calcResidual ()
 This function calculates the residual as r = v - v_old - dv. More...
 
virtual void calcStateIncr ()
 This function calculated thw increment of the state variables (dv) used to form the residual. More...
 
virtual void calcJacobian ()
 This function calculates the Jacobian. More...
 
int matrixInversion (std::vector< Real > &A, int n) const
 

Protected Attributes

Real _ref_damage_rate
 
unsigned int _nstate
 reference damage rate More...
 
Real _exponent
 Number of state variables. More...
 
unsigned int _maxiter
 
Real _tol
 Maximum number of Newton Raphson iteration. More...
 
Real _zero_tol
 Relative tolerance factor for convergence of the Newton Raphson solve. More...
 
Real _intvar_incr_tol
 Tolerance for zero. More...
 
bool _input_rndm_scale_var
 Allowable relative increment size of state variables (dv) More...
 
Real _rndm_scale_var
 Flag to specify scaling parameter to generate random stress. More...
 
MaterialProperty< std::vector< Real > > & _intvar
 Variable value. More...
 
const MaterialProperty< std::vector< Real > > & _intvar_old
 
MaterialProperty< SymmTensor > & _stress_undamaged
 
const MaterialProperty< SymmTensor > & _stress_undamaged_old
 
std::vector< Real > _intvar_incr
 
std::vector< Real > _intvar_tmp
 
std::vector< Real > _intvar_old_tmp
 
std::vector< Real > _resid
 
std::vector< Real > _jac
 
std::vector< Real > _dvar
 
SymmElasticityTensor _elasticity
 
SymmTensor _stress_old
 
SymmTensor _dstrain
 
SymmTensor _stress_new
 
SymmTensor _stress0
 
SymmTensor _dstress0
 
bool _nconv
 
bool _err_tol
 Convergence flag. More...
 
const bool _has_temp
 
const VariableValue & _temperature
 
const VariableValue & _temperature_old
 
const Real _alpha
 
Function * _alpha_function
 
bool _has_stress_free_temp
 
Real _stress_free_temp
 
bool _mean_alpha_function
 
Real _ref_temp
 
bool & _step_zero_cm
 Restartable data to check for the zeroth and first time steps. More...
 
bool & _step_one_cm
 

Detailed Description

RateDepSmearCrackModel is the base class for rate dependent continuum damage model.

The model is local and hence mesh sensitive.

Definition at line 33 of file RateDepSmearCrackModel.h.

Constructor & Destructor Documentation

◆ RateDepSmearCrackModel()

RateDepSmearCrackModel::RateDepSmearCrackModel ( const InputParameters &  parameters)

Definition at line 41 of file RateDepSmearCrackModel.C.

42  : ConstitutiveModel(parameters),
43  _ref_damage_rate(getParam<Real>("ref_damage_rate")),
44  _nstate(getParam<unsigned int>("nstate")),
45  _exponent(getParam<Real>("exponent")),
46  _maxiter(getParam<unsigned int>("maxiter")),
47  _tol(getParam<Real>("tol")),
48  _zero_tol(getParam<Real>("zero_tol")),
49  _intvar_incr_tol(getParam<Real>("intvar_incr_tol")),
50  _input_rndm_scale_var(getParam<bool>("input_random_scaling_var")),
51  _rndm_scale_var(getParam<Real>("random_scaling_var")),
52  _intvar(declareProperty<std::vector<Real>>("intvar")),
53  _intvar_old(getMaterialPropertyOld<std::vector<Real>>("intvar")),
54  _stress_undamaged(declareProperty<SymmTensor>("stress_undamaged")),
55  _stress_undamaged_old(getMaterialPropertyOld<SymmTensor>("stress_undamaged"))
56 {
57 
58  _intvar_incr.resize(_nstate, 0.0);
59  _intvar_tmp.resize(_nstate, 0.0);
60  _intvar_old_tmp.resize(_nstate, 0.0);
61  _resid.resize(_nstate, 0.0);
62  _jac.resize(_nstate * _nstate, 0.0);
63  _dvar.resize(_nstate, 0.0);
64 
66 }
static void initRandom()
Definition: SymmTensor.h:443
Real _tol
Maximum number of Newton Raphson iteration.
Real _rndm_scale_var
Flag to specify scaling parameter to generate random stress.
ConstitutiveModel(const InputParameters &parameters)
std::vector< Real > _intvar_tmp
std::vector< Real > _intvar_incr
Real _exponent
Number of state variables.
const MaterialProperty< SymmTensor > & _stress_undamaged_old
const MaterialProperty< std::vector< Real > > & _intvar_old
MaterialProperty< SymmTensor > & _stress_undamaged
unsigned int _nstate
reference damage rate
bool _input_rndm_scale_var
Allowable relative increment size of state variables (dv)
MaterialProperty< std::vector< Real > > & _intvar
Variable value.
std::vector< Real > _intvar_old_tmp
Real _intvar_incr_tol
Tolerance for zero.
Real _zero_tol
Relative tolerance factor for convergence of the Newton Raphson solve.

◆ ~RateDepSmearCrackModel()

RateDepSmearCrackModel::~RateDepSmearCrackModel ( )
virtual

Definition at line 246 of file RateDepSmearCrackModel.C.

246 {}

Member Function Documentation

◆ applyThermalStrain()

bool ConstitutiveModel::applyThermalStrain ( SymmTensor strain_increment,
SymmTensor d_strain_dT 
)
virtualinherited

Definition at line 106 of file ConstitutiveModel.C.

Referenced by ConstitutiveModel::modifyStrainIncrement().

107 {
108  if (_t_step >= 1)
109  _step_zero_cm = false;
110 
111  if (_t_step >= 2)
112  _step_one_cm = false;
113 
114  if (_has_temp && !_step_zero_cm)
115  {
116  Real inc_thermal_strain;
117  Real d_thermal_strain_d_temp;
118 
119  Real old_temp;
121  old_temp = _stress_free_temp;
122  else
123  old_temp = _temperature_old[_qp];
124 
125  Real current_temp = _temperature[_qp];
126 
127  Real delta_t = current_temp - old_temp;
128 
129  Real alpha = _alpha;
130 
131  if (_alpha_function)
132  {
133  Point p;
134  Real alpha_current_temp = _alpha_function->value(current_temp, p);
135  Real alpha_old_temp = _alpha_function->value(old_temp, p);
136  Real alpha_stress_free_temperature = _alpha_function->value(_stress_free_temp, p);
137 
139  {
140  Real small(1e-6);
141 
142  Real numerator = alpha_current_temp * (current_temp - _ref_temp) -
143  alpha_old_temp * (old_temp - _ref_temp);
144  Real denominator = 1.0 + alpha_stress_free_temperature * (_stress_free_temp - _ref_temp);
145  if (denominator < small)
146  mooseError("Denominator too small in thermal strain calculation");
147  inc_thermal_strain = numerator / denominator;
148  d_thermal_strain_d_temp = alpha_current_temp * (current_temp - _ref_temp);
149  }
150  else
151  {
152  inc_thermal_strain = delta_t * 0.5 * (alpha_current_temp + alpha_old_temp);
153  d_thermal_strain_d_temp = alpha_current_temp;
154  }
155  }
156  else
157  {
158  inc_thermal_strain = delta_t * alpha;
159  d_thermal_strain_d_temp = alpha;
160  }
161 
162  strain_increment.addDiag(-inc_thermal_strain);
163  d_strain_dT.addDiag(-d_thermal_strain_d_temp);
164  }
165 
166  bool modified = true;
167  return modified;
168 }
const VariableValue & _temperature
Function * _alpha_function
bool & _step_zero_cm
Restartable data to check for the zeroth and first time steps.
void addDiag(Real value)
Definition: SymmTensor.h:282
const VariableValue & _temperature_old

◆ calcJacobian()

void RateDepSmearCrackModel::calcJacobian ( )
protectedvirtual

This function calculates the Jacobian.

Reimplemented in RateDepSmearIsoCrackModel.

Definition at line 216 of file RateDepSmearCrackModel.C.

Referenced by solve().

217 {
218 }

◆ calcResidual()

void RateDepSmearCrackModel::calcResidual ( )
protectedvirtual

This function calculates the residual as r = v - v_old - dv.

Definition at line 199 of file RateDepSmearCrackModel.C.

Referenced by solve().

200 {
201  calcStateIncr();
202 
203  if (_err_tol)
204  return;
205 
206  for (unsigned int i = 0; i < _nstate; ++i)
208 }
std::vector< Real > _intvar_tmp
std::vector< Real > _intvar_incr
bool _err_tol
Convergence flag.
unsigned int _nstate
reference damage rate
virtual void calcStateIncr()
This function calculated thw increment of the state variables (dv) used to form the residual...
std::vector< Real > _intvar_old_tmp

◆ calcStateIncr()

void RateDepSmearCrackModel::calcStateIncr ( )
protectedvirtual

This function calculated thw increment of the state variables (dv) used to form the residual.

Reimplemented in RateDepSmearIsoCrackModel.

Definition at line 211 of file RateDepSmearCrackModel.C.

Referenced by calcResidual().

212 {
213 }

◆ computeStress()

void RateDepSmearCrackModel::computeStress ( const Elem &  current_elem,
const SymmElasticityTensor elasticity_tensor,
const SymmTensor stress_old,
SymmTensor strain_increment,
SymmTensor stress_new 
)
protectedvirtual

Reimplemented from ConstitutiveModel.

Definition at line 76 of file RateDepSmearCrackModel.C.

81 {
82  _elasticity = elasticityTensor;
83  _stress_old = stress_old;
84  _dstrain = strain_increment;
85 
86  initVariables();
87  solve();
88 
90  _rndm_scale_var = elasticityTensor.valueAtIndex(0);
91 
92  if (_nconv || _err_tol)
93  {
94  mooseWarning("RateDepSmearCrackModel: Constitutive cutback");
96  }
97  else
98  {
100  postSolveStress();
101 
102  stress_new = _stress_new;
104  }
105 }
Real _rndm_scale_var
Flag to specify scaling parameter to generate random stress.
static SymmTensor genRandomSymmTensor(Real scalefactor)
Definition: SymmTensor.h:450
virtual void solve()
This function solves the state variables.
bool _err_tol
Convergence flag.
MaterialProperty< SymmTensor > & _stress_undamaged
virtual void postSolveStress()
This function updates the stress after solve.
bool _input_rndm_scale_var
Allowable relative increment size of state variables (dv)
virtual void postSolveVariables()
This function updates the internal variables after solve.
SymmElasticityTensor _elasticity

◆ getConvergeVar()

bool RateDepSmearCrackModel::getConvergeVar ( )
protected

This function returns true if convergence is not achieved.

Definition at line 163 of file RateDepSmearCrackModel.C.

Referenced by solve().

164 {
165  Real vold, r;
166 
167  for (unsigned int i = 0; i < _nstate; i++)
168  {
169  vold = std::abs(_intvar_old_tmp[i]);
170  r = std::abs(_resid[i]);
171 
172  if (vold > _zero_tol)
173  {
174  if (r > _tol * vold)
175  return true;
176  }
177  else
178  {
179  if (r > _zero_tol)
180  return true;
181  }
182  }
183  return false;
184 }
Real _tol
Maximum number of Newton Raphson iteration.
unsigned int _nstate
reference damage rate
std::vector< Real > _intvar_old_tmp
Real _zero_tol
Relative tolerance factor for convergence of the Newton Raphson solve.

◆ initQpStatefulProperties()

void RateDepSmearCrackModel::initQpStatefulProperties ( )
protectedvirtual

Reimplemented in RateDepSmearIsoCrackModel.

Definition at line 69 of file RateDepSmearCrackModel.C.

Referenced by RateDepSmearIsoCrackModel::initQpStatefulProperties().

70 {
71  _intvar[_qp].resize(_nstate, 0.0);
72  const_cast<MaterialProperty<std::vector<Real>> &>(_intvar_old)[_qp].resize(_nstate, 0.0);
73 }
const MaterialProperty< std::vector< Real > > & _intvar_old
unsigned int _nstate
reference damage rate
MaterialProperty< std::vector< Real > > & _intvar
Variable value.

◆ initVariables()

void RateDepSmearCrackModel::initVariables ( )
protectedvirtual

Reimplemented in RateDepSmearIsoCrackModel.

Definition at line 108 of file RateDepSmearCrackModel.C.

Referenced by computeStress(), and RateDepSmearIsoCrackModel::initVariables().

109 {
110 
113 
114  for (unsigned int i = 0; i < _nstate; ++i)
115  {
116  _intvar_tmp[i] = _intvar_old[_qp][i];
117  _intvar_old_tmp[i] = _intvar_old[_qp][i];
118  }
119 }
std::vector< Real > _intvar_tmp
const MaterialProperty< SymmTensor > & _stress_undamaged_old
const MaterialProperty< std::vector< Real > > & _intvar_old
unsigned int _nstate
reference damage rate
std::vector< Real > _intvar_old_tmp
SymmElasticityTensor _elasticity

◆ matrixInversion()

int RateDepSmearCrackModel::matrixInversion ( std::vector< Real > &  A,
int  n 
) const
protected

Definition at line 221 of file RateDepSmearCrackModel.C.

Referenced by updateVariables().

222 {
223  int return_value, buffer_size = n * 64;
224  std::vector<PetscBLASInt> ipiv(n);
225  std::vector<PetscScalar> buffer(buffer_size);
226 
227  // Following does a LU decomposition of "square matrix A"
228  // upon return "A = P*L*U" if return_value == 0
229  // Here i use quotes because A is actually an array of length n^2, not a matrix of size n-by-n
230  LAPACKgetrf_(&n, &n, &A[0], &n, &ipiv[0], &return_value);
231 
232  if (return_value != 0)
233  // couldn't LU decompose because: illegal value in A; or, A singular
234  return return_value;
235 
236 // get the inverse of A
237 #if PETSC_VERSION_LESS_THAN(3, 5, 0)
238  FORTRAN_CALL(dgetri)(&n, &A[0], &n, &ipiv[0], &buffer[0], &buffer_size, &return_value);
239 #else
240  LAPACKgetri_(&n, &A[0], &n, &ipiv[0], &buffer[0], &buffer_size, &return_value);
241 #endif
242 
243  return return_value;
244 }
void FORTRAN_CALL() dgetri(...)

◆ modifyStrainIncrement()

virtual bool ConstitutiveModel::modifyStrainIncrement ( const Elem &  ,
SymmTensor strain_increment,
SymmTensor d_strain_dT 
)
inlinevirtualinherited

Reimplemented in CombinedCreepPlasticity.

Definition at line 39 of file ConstitutiveModel.h.

42  {
43  return applyThermalStrain(strain_increment, d_strain_dT);
44  }
virtual bool applyThermalStrain(SymmTensor &strain_increment, SymmTensor &d_strain_dT)

◆ postSolveStress()

void RateDepSmearCrackModel::postSolveStress ( )
protectedvirtual

This function updates the stress after solve.

In the base class it is defined as s = exp ( -d) * s0.

Reimplemented in RateDepSmearIsoCrackModel.

Definition at line 194 of file RateDepSmearCrackModel.C.

Referenced by computeStress().

195 {
196 }

◆ postSolveVariables()

void RateDepSmearCrackModel::postSolveVariables ( )
protectedvirtual

This function updates the internal variables after solve.

Definition at line 187 of file RateDepSmearCrackModel.C.

Referenced by computeStress().

188 {
189  for (unsigned int i = 0; i < _nstate; ++i)
190  _intvar[_qp][i] = _intvar_tmp[i];
191 }
std::vector< Real > _intvar_tmp
unsigned int _nstate
reference damage rate
MaterialProperty< std::vector< Real > > & _intvar
Variable value.

◆ setQp()

void ConstitutiveModel::setQp ( unsigned int  qp)
inherited

Sets the value of the variable _qp for inheriting classes.

Definition at line 89 of file ConstitutiveModel.C.

Referenced by CombinedCreepPlasticity::computeStress().

90 {
91  _qp = qp;
92 }

◆ solve()

void RateDepSmearCrackModel::solve ( )
protectedvirtual

This function solves the state variables.

In the present formulation the damaged stress (s) is related to the undamaged stress (s0) as s = exp(-d) * s0 where d is a state variable describing damage. d can be scalar or vector depending on the model A Newton-Raphson is used.

Definition at line 122 of file RateDepSmearCrackModel.C.

Referenced by computeStress().

123 {
124 
125  unsigned int iter = 0;
126  _err_tol = false;
127 
128  calcResidual();
129 
131 
132  while (_nconv && iter < _maxiter && !_err_tol)
133  {
134  calcJacobian();
135  updateVariables();
136  calcResidual();
137 
139 
140  iter++;
141  }
142 }
virtual void updateVariables()
This function updates variables during solve a(i+1) = a(i) + da(i+1)
bool _err_tol
Convergence flag.
virtual void calcResidual()
This function calculates the residual as r = v - v_old - dv.
virtual void calcJacobian()
This function calculates the Jacobian.
bool getConvergeVar()
This function returns true if convergence is not achieved.

◆ updateElasticityTensor()

virtual bool ConstitutiveModel::updateElasticityTensor ( SymmElasticityTensor )
inlinevirtualinherited

Definition at line 45 of file ConstitutiveModel.h.

45 { return false; }

◆ updateVariables()

void RateDepSmearCrackModel::updateVariables ( )
protectedvirtual

This function updates variables during solve a(i+1) = a(i) + da(i+1)

Definition at line 145 of file RateDepSmearCrackModel.C.

Referenced by solve().

146 {
147  int error = matrixInversion(_jac, _nstate);
148  if (error != 0)
149  mooseError("Error in Matrix Inversion in RankFourTensor");
150 
151  for (unsigned int i = 0; i < _nstate; i++)
152  {
153  _dvar[i] = 0.0;
154  for (unsigned int j = 0; j < _nstate; j++)
155  _dvar[i] += _jac[i * _nstate + j] * _resid[j];
156  }
157 
158  for (unsigned int i = 0; i < _nstate; i++)
159  _intvar_tmp[i] -= _dvar[i];
160 }
std::vector< Real > _intvar_tmp
unsigned int _nstate
reference damage rate
int matrixInversion(std::vector< Real > &A, int n) const

Member Data Documentation

◆ _alpha

const Real ConstitutiveModel::_alpha
protectedinherited

Definition at line 53 of file ConstitutiveModel.h.

Referenced by ConstitutiveModel::applyThermalStrain().

◆ _alpha_function

Function* ConstitutiveModel::_alpha_function
protectedinherited

◆ _dstrain

SymmTensor RateDepSmearCrackModel::_dstrain
protected

◆ _dstress0

SymmTensor RateDepSmearCrackModel::_dstress0
protected

Definition at line 123 of file RateDepSmearCrackModel.h.

Referenced by initVariables().

◆ _dvar

std::vector<Real> RateDepSmearCrackModel::_dvar
protected

Definition at line 119 of file RateDepSmearCrackModel.h.

Referenced by RateDepSmearCrackModel(), and updateVariables().

◆ _elasticity

SymmElasticityTensor RateDepSmearCrackModel::_elasticity
protected

Definition at line 121 of file RateDepSmearCrackModel.h.

Referenced by computeStress(), and initVariables().

◆ _err_tol

bool RateDepSmearCrackModel::_err_tol
protected

Convergence flag.

Definition at line 125 of file RateDepSmearCrackModel.h.

Referenced by calcResidual(), RateDepSmearIsoCrackModel::calcStateIncr(), computeStress(), and solve().

◆ _exponent

Real RateDepSmearCrackModel::_exponent
protected

Number of state variables.

Definition at line 101 of file RateDepSmearCrackModel.h.

Referenced by RateDepSmearIsoCrackModel::damageRate().

◆ _has_stress_free_temp

bool ConstitutiveModel::_has_stress_free_temp
protectedinherited

Definition at line 55 of file ConstitutiveModel.h.

Referenced by ConstitutiveModel::applyThermalStrain().

◆ _has_temp

const bool ConstitutiveModel::_has_temp
protectedinherited

◆ _input_rndm_scale_var

bool RateDepSmearCrackModel::_input_rndm_scale_var
protected

Allowable relative increment size of state variables (dv)

Definition at line 106 of file RateDepSmearCrackModel.h.

Referenced by computeStress().

◆ _intvar

MaterialProperty<std::vector<Real> >& RateDepSmearCrackModel::_intvar
protected

◆ _intvar_incr

std::vector<Real> RateDepSmearCrackModel::_intvar_incr
protected

◆ _intvar_incr_tol

Real RateDepSmearCrackModel::_intvar_incr_tol
protected

Tolerance for zero.

Definition at line 105 of file RateDepSmearCrackModel.h.

Referenced by RateDepSmearIsoCrackModel::calcStateIncr().

◆ _intvar_old

const MaterialProperty<std::vector<Real> >& RateDepSmearCrackModel::_intvar_old
protected

◆ _intvar_old_tmp

std::vector<Real> RateDepSmearCrackModel::_intvar_old_tmp
protected

◆ _intvar_tmp

std::vector<Real> RateDepSmearCrackModel::_intvar_tmp
protected

◆ _jac

std::vector<Real> RateDepSmearCrackModel::_jac
protected

◆ _maxiter

unsigned int RateDepSmearCrackModel::_maxiter
protected

Definition at line 102 of file RateDepSmearCrackModel.h.

Referenced by solve().

◆ _mean_alpha_function

bool ConstitutiveModel::_mean_alpha_function
protectedinherited

◆ _nconv

bool RateDepSmearCrackModel::_nconv
protected

Definition at line 124 of file RateDepSmearCrackModel.h.

Referenced by computeStress(), and solve().

◆ _nstate

unsigned int RateDepSmearCrackModel::_nstate
protected

◆ _ref_damage_rate

Real RateDepSmearCrackModel::_ref_damage_rate
protected

Definition at line 99 of file RateDepSmearCrackModel.h.

Referenced by RateDepSmearIsoCrackModel::damageRate().

◆ _ref_temp

Real ConstitutiveModel::_ref_temp
protectedinherited

◆ _resid

std::vector<Real> RateDepSmearCrackModel::_resid
protected

◆ _rndm_scale_var

Real RateDepSmearCrackModel::_rndm_scale_var
protected

Flag to specify scaling parameter to generate random stress.

Definition at line 107 of file RateDepSmearCrackModel.h.

Referenced by computeStress().

◆ _step_one_cm

bool& ConstitutiveModel::_step_one_cm
protectedinherited

Definition at line 62 of file ConstitutiveModel.h.

Referenced by ConstitutiveModel::applyThermalStrain().

◆ _step_zero_cm

bool& ConstitutiveModel::_step_zero_cm
protectedinherited

Restartable data to check for the zeroth and first time steps.

Definition at line 61 of file ConstitutiveModel.h.

Referenced by ConstitutiveModel::applyThermalStrain().

◆ _stress0

SymmTensor RateDepSmearCrackModel::_stress0
protected

◆ _stress_free_temp

Real ConstitutiveModel::_stress_free_temp
protectedinherited

Definition at line 56 of file ConstitutiveModel.h.

Referenced by ConstitutiveModel::applyThermalStrain().

◆ _stress_new

SymmTensor RateDepSmearCrackModel::_stress_new
protected

◆ _stress_old

SymmTensor RateDepSmearCrackModel::_stress_old
protected

Definition at line 122 of file RateDepSmearCrackModel.h.

Referenced by computeStress().

◆ _stress_undamaged

MaterialProperty<SymmTensor>& RateDepSmearCrackModel::_stress_undamaged
protected

Definition at line 112 of file RateDepSmearCrackModel.h.

Referenced by computeStress().

◆ _stress_undamaged_old

const MaterialProperty<SymmTensor>& RateDepSmearCrackModel::_stress_undamaged_old
protected

Definition at line 113 of file RateDepSmearCrackModel.h.

Referenced by initVariables().

◆ _temperature

const VariableValue& ConstitutiveModel::_temperature
protectedinherited

◆ _temperature_old

const VariableValue& ConstitutiveModel::_temperature_old
protectedinherited

Definition at line 52 of file ConstitutiveModel.h.

Referenced by ConstitutiveModel::applyThermalStrain().

◆ _tol

Real RateDepSmearCrackModel::_tol
protected

Maximum number of Newton Raphson iteration.

Definition at line 103 of file RateDepSmearCrackModel.h.

Referenced by getConvergeVar().

◆ _zero_tol

Real RateDepSmearCrackModel::_zero_tol
protected

Relative tolerance factor for convergence of the Newton Raphson solve.

Definition at line 104 of file RateDepSmearCrackModel.h.

Referenced by RateDepSmearIsoCrackModel::calcStateIncr(), and getConvergeVar().


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