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

In this class a rate dependent isotropic damage model is implemented. More...

#include <RateDepSmearIsoCrackModel.h>

Inheritance diagram for RateDepSmearIsoCrackModel:
[legend]

Public Member Functions

 RateDepSmearIsoCrackModel (const InputParameters &parameters)
 
virtual ~RateDepSmearIsoCrackModel ()
 
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 initQpStatefulProperties ()
 
virtual void initVariables ()
 
virtual Real damageRate ()
 This function calculates rate of damage based on energy. 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...
 
virtual void postSolveStress ()
 This function updates the stress after solve. More...
 
virtual void computeStress (const Elem &current_elem, const SymmElasticityTensor &elasticity_tensor, const SymmTensor &stress_old, SymmTensor &strain_increment, SymmTensor &stress_new)
 
virtual void solve ()
 This function solves the state variables. More...
 
virtual void postSolveVariables ()
 This function updates the internal variables 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...
 
int matrixInversion (std::vector< Real > &A, int n) const
 

Protected Attributes

Real _crit_energy
 
Real _kfail
 Required parameter. More...
 
Real _upper_lim_damage
 Used to avoid non-positive definiteness. More...
 
MaterialProperty< Real > & _energy
 
const MaterialProperty< Real > & _energy_old
 
Real _ddamage
 
Real _ddamagerate_drs
 
ColumnMajorMatrix _s0_diag
 
ColumnMajorMatrix _s0_diag_pos
 
ColumnMajorMatrix _s0_diag_neg
 
ColumnMajorMatrix _s0_evec
 
ColumnMajorMatrix _dstrain_diag
 
ColumnMajorMatrix _dstrain_diag_pos
 
ColumnMajorMatrix _dstrain_diag_neg
 
ColumnMajorMatrix _dstrain_evec
 
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

In this class a rate dependent isotropic damage model is implemented.

Definition at line 21 of file RateDepSmearIsoCrackModel.h.

Constructor & Destructor Documentation

RateDepSmearIsoCrackModel::RateDepSmearIsoCrackModel ( const InputParameters &  parameters)

Definition at line 25 of file RateDepSmearIsoCrackModel.C.

26  : RateDepSmearCrackModel(parameters),
27  _crit_energy(getParam<Real>("critical_energy")),
28  _kfail(getParam<Real>("k_fail")),
29  _upper_lim_damage(getParam<Real>("upper_limit_damage")),
30  _energy(declareProperty<Real>("energy")),
31  _energy_old(getMaterialPropertyOld<Real>("energy"))
32 {
33 
34  if (_nstate != 2)
35  mooseError(" RateDpeSmearIsoCrackModel Error: Requires 2 state variables - nstate = 2");
36 }
MaterialProperty< Real > & _energy
RateDepSmearCrackModel(const InputParameters &parameters)
Real _upper_lim_damage
Used to avoid non-positive definiteness.
Real _kfail
Required parameter.
unsigned int _nstate
reference damage rate
const MaterialProperty< Real > & _energy_old
RateDepSmearIsoCrackModel::~RateDepSmearIsoCrackModel ( )
virtual

Definition at line 144 of file RateDepSmearIsoCrackModel.C.

144 {}

Member Function Documentation

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

Definition at line 103 of file ConstitutiveModel.C.

Referenced by ConstitutiveModel::modifyStrainIncrement(), and ConstitutiveModel::updateElasticityTensor().

104 {
105  if (_t_step >= 1)
106  _step_zero_cm = false;
107 
108  if (_t_step >= 2)
109  _step_one_cm = false;
110 
111  if (_has_temp && !_step_zero_cm)
112  {
113  Real inc_thermal_strain;
114  Real d_thermal_strain_d_temp;
115 
116  Real old_temp;
118  old_temp = _stress_free_temp;
119  else
120  old_temp = _temperature_old[_qp];
121 
122  Real current_temp = _temperature[_qp];
123 
124  Real delta_t = current_temp - old_temp;
125 
126  Real alpha = _alpha;
127 
128  if (_alpha_function)
129  {
130  Point p;
131  Real alpha_current_temp = _alpha_function->value(current_temp, p);
132  Real alpha_old_temp = _alpha_function->value(old_temp, p);
133  Real alpha_stress_free_temperature = _alpha_function->value(_stress_free_temp, p);
134 
136  {
137  Real small(1e-6);
138 
139  Real numerator = alpha_current_temp * (current_temp - _ref_temp) -
140  alpha_old_temp * (old_temp - _ref_temp);
141  Real denominator = 1.0 + alpha_stress_free_temperature * (_stress_free_temp - _ref_temp);
142  if (denominator < small)
143  mooseError("Denominator too small in thermal strain calculation");
144  inc_thermal_strain = numerator / denominator;
145  d_thermal_strain_d_temp = alpha_current_temp * (current_temp - _ref_temp);
146  }
147  else
148  {
149  inc_thermal_strain = delta_t * 0.5 * (alpha_current_temp + alpha_old_temp);
150  d_thermal_strain_d_temp = alpha_current_temp;
151  }
152  }
153  else
154  {
155  inc_thermal_strain = delta_t * alpha;
156  d_thermal_strain_d_temp = alpha;
157  }
158 
159  strain_increment.addDiag(-inc_thermal_strain);
160  d_strain_dT.addDiag(-d_thermal_strain_d_temp);
161  }
162 
163  bool modified = true;
164  return modified;
165 }
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:279
const VariableValue & _temperature_old
void RateDepSmearIsoCrackModel::calcJacobian ( )
protectedvirtual

This function calculates the Jacobian.

Reimplemented from RateDepSmearCrackModel.

Definition at line 128 of file RateDepSmearIsoCrackModel.C.

129 {
130 
131  _jac[0] = 1.0;
132  _jac[1] = -_ddamagerate_drs;
133  _jac[2] = 0.0;
134  _jac[3] = 1 - _ddamagerate_drs;
135 }
void RateDepSmearCrackModel::calcResidual ( )
protectedvirtualinherited

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

Definition at line 194 of file RateDepSmearCrackModel.C.

Referenced by RateDepSmearCrackModel::solve().

195 {
196  calcStateIncr();
197 
198  if (_err_tol)
199  return;
200 
201  for (unsigned int i = 0; i < _nstate; ++i)
203 }
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
void RateDepSmearIsoCrackModel::calcStateIncr ( )
protectedvirtual

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

Reimplemented from RateDepSmearCrackModel.

Definition at line 104 of file RateDepSmearIsoCrackModel.C.

105 {
106 
107  Real val = damageRate();
108 
109  if (std::abs(_intvar_old_tmp[0]) < _zero_tol && std::abs(val) > _intvar_incr_tol)
110  {
111  _err_tol = true;
112  return;
113  }
114 
115  if (std::abs(_intvar_old_tmp[0]) > _zero_tol &&
116  std::abs(_intvar_old_tmp[0]) < _upper_lim_damage &&
117  std::abs(val) > _intvar_incr_tol * std::abs(_intvar_old_tmp[0]))
118  {
119  _err_tol = true;
120  return;
121  }
122 
123  for (unsigned int i = 0; i < _nstate; i++)
124  _intvar_incr[i] = val;
125 }
virtual Real damageRate()
This function calculates rate of damage based on energy.
std::vector< Real > _intvar_incr
Real _upper_lim_damage
Used to avoid non-positive definiteness.
bool _err_tol
Convergence flag.
unsigned int _nstate
reference damage rate
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.
void RateDepSmearCrackModel::computeStress ( const Elem &  current_elem,
const SymmElasticityTensor elasticity_tensor,
const SymmTensor stress_old,
SymmTensor strain_increment,
SymmTensor stress_new 
)
protectedvirtualinherited

Reimplemented from ConstitutiveModel.

Definition at line 71 of file RateDepSmearCrackModel.C.

76 {
77  _elasticity = elasticityTensor;
78  _stress_old = stress_old;
79  _dstrain = strain_increment;
80 
81  initVariables();
82  solve();
83 
85  _rndm_scale_var = elasticityTensor.valueAtIndex(0);
86 
87  if (_nconv || _err_tol)
88  {
89  mooseWarning("RateDepSmearCrackModel: Constitutive cutback");
91  }
92  else
93  {
96 
97  stress_new = _stress_new;
99  }
100 }
Real _rndm_scale_var
Flag to specify scaling parameter to generate random stress.
static SymmTensor genRandomSymmTensor(Real scalefactor)
Definition: SymmTensor.h:447
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
Real RateDepSmearIsoCrackModel::damageRate ( )
protectedvirtual

This function calculates rate of damage based on energy.

Definition at line 73 of file RateDepSmearIsoCrackModel.C.

Referenced by calcStateIncr().

74 {
75 
76  Real den = 0.0;
77 
78  for (unsigned int i = 0; i < LIBMESH_DIM; i++)
79  den += _s0_diag_pos(i, i) * _dstrain_diag_pos(i, i);
80 
81  _energy[_qp] = _energy_old[_qp] + den;
82 
83  Real e = _energy[_qp];
84 
85  Real fac = e - _intvar_tmp[1];
86  Real val = 0.0;
87 
88  _ddamagerate_drs = 0.0;
89 
90  if (fac > 0.0)
91  {
92  val = std::pow(fac, _exponent);
93  val *= _ref_damage_rate * _dt;
94 
95  Real dval = -_exponent * std::pow(fac, _exponent - 1.0);
96 
97  _ddamagerate_drs = dval * _ref_damage_rate * _dt;
98  }
99 
100  return val;
101 }
MaterialProperty< Real > & _energy
std::vector< Real > _intvar_tmp
Real _exponent
Number of state variables.
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
const MaterialProperty< Real > & _energy_old
bool RateDepSmearCrackModel::getConvergeVar ( )
protectedinherited

This function returns true if convergence is not achieved.

Definition at line 158 of file RateDepSmearCrackModel.C.

Referenced by RateDepSmearCrackModel::solve().

159 {
160  Real vold, r;
161 
162  for (unsigned int i = 0; i < _nstate; i++)
163  {
164  vold = std::abs(_intvar_old_tmp[i]);
165  r = std::abs(_resid[i]);
166 
167  if (vold > _zero_tol)
168  {
169  if (r > _tol * vold)
170  return true;
171  }
172  else
173  {
174  if (r > _zero_tol)
175  return true;
176  }
177  }
178  return false;
179 }
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.
void RateDepSmearIsoCrackModel::initQpStatefulProperties ( )
protectedvirtual

Reimplemented from RateDepSmearCrackModel.

Definition at line 39 of file RateDepSmearIsoCrackModel.C.

40 {
42 
43  _intvar[_qp][1] = const_cast<MaterialProperty<std::vector<Real>> &>(_intvar_old)[_qp][1] =
45 
46  _energy[_qp] = 0.0;
47 }
MaterialProperty< Real > & _energy
virtual void initQpStatefulProperties()
const MaterialProperty< std::vector< Real > > & _intvar_old
MaterialProperty< std::vector< Real > > & _intvar
Variable value.
void RateDepSmearIsoCrackModel::initVariables ( )
protectedvirtual

Reimplemented from RateDepSmearCrackModel.

Definition at line 50 of file RateDepSmearIsoCrackModel.C.

51 {
52 
54 
56 
57  for (unsigned int i = 0; i < LIBMESH_DIM; i++)
58  {
59  _s0_diag_pos(i, i) = (std::abs(_s0_diag(i, 0)) + _s0_diag(i, 0)) / 2.0;
60  _s0_diag_neg(i, i) = (std::abs(_s0_diag(i, 0)) - _s0_diag(i, 0)) / 2.0;
61  }
62 
64 
65  for (unsigned int i = 0; i < LIBMESH_DIM; i++)
66  {
67  _dstrain_diag_pos(i, i) = (std::abs(_dstrain_diag(i, 0)) + _dstrain_diag(i, 0)) / 2.0;
68  _dstrain_diag_neg(i, i) = (std::abs(_dstrain_diag(i, 0)) - _dstrain_diag(i, 0)) / 2.0;
69  }
70 }
ColumnMajorMatrix columnMajorMatrix() const
Definition: SymmTensor.h:423
int RateDepSmearCrackModel::matrixInversion ( std::vector< Real > &  A,
int  n 
) const
protectedinherited

Definition at line 216 of file RateDepSmearCrackModel.C.

Referenced by RateDepSmearCrackModel::updateVariables().

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

Reimplemented in CombinedCreepPlasticity.

Definition at line 34 of file ConstitutiveModel.h.

37  {
38  return applyThermalStrain(strain_increment, d_strain_dT);
39  }
virtual bool applyThermalStrain(SymmTensor &strain_increment, SymmTensor &d_strain_dT)
void RateDepSmearIsoCrackModel::postSolveStress ( )
protectedvirtual

This function updates the stress after solve.

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

Reimplemented from RateDepSmearCrackModel.

Definition at line 138 of file RateDepSmearIsoCrackModel.C.

139 {
140  ColumnMajorMatrix stress_cmm = _s0_diag_pos * (std::exp(-_intvar_tmp[0]) + _kfail) - _s0_diag_neg;
141  _stress_new = _s0_evec * stress_cmm * _s0_evec.transpose();
142 }
std::vector< Real > _intvar_tmp
Real _kfail
Required parameter.
void RateDepSmearCrackModel::postSolveVariables ( )
protectedvirtualinherited

This function updates the internal variables after solve.

Definition at line 182 of file RateDepSmearCrackModel.C.

Referenced by RateDepSmearCrackModel::computeStress().

183 {
184  for (unsigned int i = 0; i < _nstate; ++i)
185  _intvar[_qp][i] = _intvar_tmp[i];
186 }
std::vector< Real > _intvar_tmp
unsigned int _nstate
reference damage rate
MaterialProperty< std::vector< Real > > & _intvar
Variable value.
void ConstitutiveModel::setQp ( unsigned int  qp)
inherited

Sets the value of the variable _qp for inheriting classes.

Definition at line 86 of file ConstitutiveModel.C.

Referenced by CombinedCreepPlasticity::computeStress(), and ConstitutiveModel::~ConstitutiveModel().

87 {
88  _qp = qp;
89 }
void RateDepSmearCrackModel::solve ( )
protectedvirtualinherited

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 117 of file RateDepSmearCrackModel.C.

Referenced by RateDepSmearCrackModel::computeStress().

118 {
119 
120  unsigned int iter = 0;
121  _err_tol = false;
122 
123  calcResidual();
124 
126 
127  while (_nconv && iter < _maxiter && !_err_tol)
128  {
129  calcJacobian();
130  updateVariables();
131  calcResidual();
132 
134 
135  iter++;
136  }
137 }
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.
virtual bool ConstitutiveModel::updateElasticityTensor ( SymmElasticityTensor )
inlinevirtualinherited

Definition at line 40 of file ConstitutiveModel.h.

40 { return false; }
void RateDepSmearCrackModel::updateVariables ( )
protectedvirtualinherited

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

Definition at line 140 of file RateDepSmearCrackModel.C.

Referenced by RateDepSmearCrackModel::solve().

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

Member Data Documentation

const Real ConstitutiveModel::_alpha
protectedinherited

Definition at line 48 of file ConstitutiveModel.h.

Referenced by ConstitutiveModel::applyThermalStrain().

Function* ConstitutiveModel::_alpha_function
protectedinherited
Real RateDepSmearIsoCrackModel::_crit_energy
protected

Definition at line 39 of file RateDepSmearIsoCrackModel.h.

Referenced by initQpStatefulProperties().

Real RateDepSmearIsoCrackModel::_ddamage
protected

Definition at line 46 of file RateDepSmearIsoCrackModel.h.

Real RateDepSmearIsoCrackModel::_ddamagerate_drs
protected

Definition at line 47 of file RateDepSmearIsoCrackModel.h.

Referenced by calcJacobian(), and damageRate().

SymmTensor RateDepSmearCrackModel::_dstrain
protectedinherited
ColumnMajorMatrix RateDepSmearIsoCrackModel::_dstrain_diag
protected

Definition at line 50 of file RateDepSmearIsoCrackModel.h.

Referenced by initVariables().

ColumnMajorMatrix RateDepSmearIsoCrackModel::_dstrain_diag_neg
protected

Definition at line 50 of file RateDepSmearIsoCrackModel.h.

Referenced by initVariables().

ColumnMajorMatrix RateDepSmearIsoCrackModel::_dstrain_diag_pos
protected

Definition at line 50 of file RateDepSmearIsoCrackModel.h.

Referenced by damageRate(), and initVariables().

ColumnMajorMatrix RateDepSmearIsoCrackModel::_dstrain_evec
protected

Definition at line 50 of file RateDepSmearIsoCrackModel.h.

Referenced by initVariables().

SymmTensor RateDepSmearCrackModel::_dstress0
protectedinherited

Definition at line 120 of file RateDepSmearCrackModel.h.

Referenced by RateDepSmearCrackModel::initVariables().

std::vector<Real> RateDepSmearCrackModel::_dvar
protectedinherited
SymmElasticityTensor RateDepSmearCrackModel::_elasticity
protectedinherited
MaterialProperty<Real>& RateDepSmearIsoCrackModel::_energy
protected

Definition at line 43 of file RateDepSmearIsoCrackModel.h.

Referenced by damageRate(), and initQpStatefulProperties().

const MaterialProperty<Real>& RateDepSmearIsoCrackModel::_energy_old
protected

Definition at line 44 of file RateDepSmearIsoCrackModel.h.

Referenced by damageRate().

bool RateDepSmearCrackModel::_err_tol
protectedinherited
Real RateDepSmearCrackModel::_exponent
protectedinherited

Number of state variables.

Definition at line 98 of file RateDepSmearCrackModel.h.

Referenced by damageRate().

bool ConstitutiveModel::_has_stress_free_temp
protectedinherited

Definition at line 50 of file ConstitutiveModel.h.

Referenced by ConstitutiveModel::applyThermalStrain().

const bool ConstitutiveModel::_has_temp
protectedinherited
bool RateDepSmearCrackModel::_input_rndm_scale_var
protectedinherited

Allowable relative increment size of state variables (dv)

Definition at line 103 of file RateDepSmearCrackModel.h.

Referenced by RateDepSmearCrackModel::computeStress().

MaterialProperty<std::vector<Real> >& RateDepSmearCrackModel::_intvar
protectedinherited
std::vector<Real> RateDepSmearCrackModel::_intvar_incr
protectedinherited
Real RateDepSmearCrackModel::_intvar_incr_tol
protectedinherited

Tolerance for zero.

Definition at line 102 of file RateDepSmearCrackModel.h.

Referenced by calcStateIncr().

const MaterialProperty<std::vector<Real> >& RateDepSmearCrackModel::_intvar_old
protectedinherited
std::vector<Real> RateDepSmearCrackModel::_intvar_old_tmp
protectedinherited
std::vector<Real> RateDepSmearCrackModel::_intvar_tmp
protectedinherited
std::vector<Real> RateDepSmearCrackModel::_jac
protectedinherited
Real RateDepSmearIsoCrackModel::_kfail
protected

Required parameter.

Definition at line 40 of file RateDepSmearIsoCrackModel.h.

Referenced by postSolveStress().

unsigned int RateDepSmearCrackModel::_maxiter
protectedinherited

Definition at line 99 of file RateDepSmearCrackModel.h.

Referenced by RateDepSmearCrackModel::solve().

bool ConstitutiveModel::_mean_alpha_function
protectedinherited
bool RateDepSmearCrackModel::_nconv
protectedinherited
unsigned int RateDepSmearCrackModel::_nstate
protectedinherited
Real RateDepSmearCrackModel::_ref_damage_rate
protectedinherited

Definition at line 96 of file RateDepSmearCrackModel.h.

Referenced by damageRate().

Real ConstitutiveModel::_ref_temp
protectedinherited
std::vector<Real> RateDepSmearCrackModel::_resid
protectedinherited
Real RateDepSmearCrackModel::_rndm_scale_var
protectedinherited

Flag to specify scaling parameter to generate random stress.

Definition at line 104 of file RateDepSmearCrackModel.h.

Referenced by RateDepSmearCrackModel::computeStress().

ColumnMajorMatrix RateDepSmearIsoCrackModel::_s0_diag
protected

Definition at line 49 of file RateDepSmearIsoCrackModel.h.

Referenced by initVariables().

ColumnMajorMatrix RateDepSmearIsoCrackModel::_s0_diag_neg
protected

Definition at line 49 of file RateDepSmearIsoCrackModel.h.

Referenced by initVariables(), and postSolveStress().

ColumnMajorMatrix RateDepSmearIsoCrackModel::_s0_diag_pos
protected

Definition at line 49 of file RateDepSmearIsoCrackModel.h.

Referenced by damageRate(), initVariables(), and postSolveStress().

ColumnMajorMatrix RateDepSmearIsoCrackModel::_s0_evec
protected

Definition at line 49 of file RateDepSmearIsoCrackModel.h.

Referenced by initVariables(), and postSolveStress().

bool& ConstitutiveModel::_step_one_cm
protectedinherited

Definition at line 57 of file ConstitutiveModel.h.

Referenced by ConstitutiveModel::applyThermalStrain().

bool& ConstitutiveModel::_step_zero_cm
protectedinherited

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

Definition at line 56 of file ConstitutiveModel.h.

Referenced by ConstitutiveModel::applyThermalStrain().

SymmTensor RateDepSmearCrackModel::_stress0
protectedinherited
Real ConstitutiveModel::_stress_free_temp
protectedinherited

Definition at line 51 of file ConstitutiveModel.h.

Referenced by ConstitutiveModel::applyThermalStrain().

SymmTensor RateDepSmearCrackModel::_stress_new
protectedinherited
SymmTensor RateDepSmearCrackModel::_stress_old
protectedinherited

Definition at line 119 of file RateDepSmearCrackModel.h.

Referenced by RateDepSmearCrackModel::computeStress().

MaterialProperty<SymmTensor>& RateDepSmearCrackModel::_stress_undamaged
protectedinherited

Definition at line 109 of file RateDepSmearCrackModel.h.

Referenced by RateDepSmearCrackModel::computeStress().

const MaterialProperty<SymmTensor>& RateDepSmearCrackModel::_stress_undamaged_old
protectedinherited

Definition at line 110 of file RateDepSmearCrackModel.h.

Referenced by RateDepSmearCrackModel::initVariables().

const VariableValue& ConstitutiveModel::_temperature
protectedinherited
const VariableValue& ConstitutiveModel::_temperature_old
protectedinherited

Definition at line 47 of file ConstitutiveModel.h.

Referenced by ConstitutiveModel::applyThermalStrain().

Real RateDepSmearCrackModel::_tol
protectedinherited

Maximum number of Newton Raphson iteration.

Definition at line 100 of file RateDepSmearCrackModel.h.

Referenced by RateDepSmearCrackModel::getConvergeVar().

Real RateDepSmearIsoCrackModel::_upper_lim_damage
protected

Used to avoid non-positive definiteness.

Definition at line 41 of file RateDepSmearIsoCrackModel.h.

Referenced by calcStateIncr().

Real RateDepSmearCrackModel::_zero_tol
protectedinherited

Relative tolerance factor for convergence of the Newton Raphson solve.

Definition at line 101 of file RateDepSmearCrackModel.h.

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


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