www.mooseframework.org
RateDepSmearIsoCrackModel.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 
13 
14 template <>
15 InputParameters
17 {
18 
19  InputParameters params = validParams<RateDepSmearCrackModel>();
20 
21  params.addRequiredParam<Real>("critical_energy", "Critical Energy");
22  params.addParam<Real>("k_fail", 1e-6, "Post failure stiffness");
23  params.addParam<Real>("upper_limit_damage",
24  5.0,
25  "Upper limit of damage beyond which constitutive check is not performed");
26 
27  return params;
28 }
29 
30 RateDepSmearIsoCrackModel::RateDepSmearIsoCrackModel(const InputParameters & parameters)
31  : RateDepSmearCrackModel(parameters),
32  _crit_energy(getParam<Real>("critical_energy")),
33  _kfail(getParam<Real>("k_fail")),
34  _upper_lim_damage(getParam<Real>("upper_limit_damage")),
35  _energy(declareProperty<Real>("energy")),
36  _energy_old(getMaterialPropertyOld<Real>("energy"))
37 {
38 
39  if (_nstate != 2)
40  mooseError(" RateDpeSmearIsoCrackModel Error: Requires 2 state variables - nstate = 2");
41 }
42 
43 void
45 {
47 
48  _intvar[_qp][1] = const_cast<MaterialProperty<std::vector<Real>> &>(_intvar_old)[_qp][1] =
50 
51  _energy[_qp] = 0.0;
52 }
53 
54 void
56 {
57 
59 
61 
62  for (unsigned int i = 0; i < LIBMESH_DIM; i++)
63  {
64  _s0_diag_pos(i, i) = (std::abs(_s0_diag(i, 0)) + _s0_diag(i, 0)) / 2.0;
65  _s0_diag_neg(i, i) = (std::abs(_s0_diag(i, 0)) - _s0_diag(i, 0)) / 2.0;
66  }
67 
69 
70  for (unsigned int i = 0; i < LIBMESH_DIM; i++)
71  {
72  _dstrain_diag_pos(i, i) = (std::abs(_dstrain_diag(i, 0)) + _dstrain_diag(i, 0)) / 2.0;
73  _dstrain_diag_neg(i, i) = (std::abs(_dstrain_diag(i, 0)) - _dstrain_diag(i, 0)) / 2.0;
74  }
75 }
76 
77 Real
79 {
80 
81  Real den = 0.0;
82 
83  for (unsigned int i = 0; i < LIBMESH_DIM; i++)
84  den += _s0_diag_pos(i, i) * _dstrain_diag_pos(i, i);
85 
86  _energy[_qp] = _energy_old[_qp] + den;
87 
88  Real e = _energy[_qp];
89 
90  Real fac = e - _intvar_tmp[1];
91  Real val = 0.0;
92 
93  _ddamagerate_drs = 0.0;
94 
95  if (fac > 0.0)
96  {
97  val = std::pow(fac, _exponent);
98  val *= _ref_damage_rate * _dt;
99 
100  Real dval = -_exponent * std::pow(fac, _exponent - 1.0);
101 
102  _ddamagerate_drs = dval * _ref_damage_rate * _dt;
103  }
104 
105  return val;
106 }
107 
108 void
110 {
111 
112  Real val = damageRate();
113 
114  if (std::abs(_intvar_old_tmp[0]) < _zero_tol && std::abs(val) > _intvar_incr_tol)
115  {
116  _err_tol = true;
117  return;
118  }
119 
120  if (std::abs(_intvar_old_tmp[0]) > _zero_tol &&
121  std::abs(_intvar_old_tmp[0]) < _upper_lim_damage &&
122  std::abs(val) > _intvar_incr_tol * std::abs(_intvar_old_tmp[0]))
123  {
124  _err_tol = true;
125  return;
126  }
127 
128  for (unsigned int i = 0; i < _nstate; i++)
129  _intvar_incr[i] = val;
130 }
131 
132 void
134 {
135 
136  _jac[0] = 1.0;
137  _jac[1] = -_ddamagerate_drs;
138  _jac[2] = 0.0;
139  _jac[3] = 1 - _ddamagerate_drs;
140 }
141 
142 void
144 {
145  ColumnMajorMatrix stress_cmm = _s0_diag_pos * (std::exp(-_intvar_tmp[0]) + _kfail) - _s0_diag_neg;
146  _stress_new = _s0_evec * stress_cmm * _s0_evec.transpose();
147 }
148 
RateDepSmearIsoCrackModel::_ddamagerate_drs
Real _ddamagerate_drs
Definition: RateDepSmearIsoCrackModel.h:49
RateDepSmearCrackModel::_intvar_incr_tol
Real _intvar_incr_tol
Tolerance for zero.
Definition: RateDepSmearCrackModel.h:104
RateDepSmearIsoCrackModel::_s0_evec
ColumnMajorMatrix _s0_evec
Definition: RateDepSmearIsoCrackModel.h:51
registerMooseObject
registerMooseObject("SolidMechanicsApp", RateDepSmearIsoCrackModel)
RateDepSmearIsoCrackModel::initQpStatefulProperties
virtual void initQpStatefulProperties()
Definition: RateDepSmearIsoCrackModel.C:44
RateDepSmearIsoCrackModel::_dstrain_evec
ColumnMajorMatrix _dstrain_evec
Definition: RateDepSmearIsoCrackModel.h:52
RateDepSmearIsoCrackModel::_energy
MaterialProperty< Real > & _energy
Definition: RateDepSmearIsoCrackModel.h:45
RateDepSmearCrackModel::_stress_new
SymmTensor _stress_new
Definition: RateDepSmearCrackModel.h:121
RateDepSmearIsoCrackModel::_dstrain_diag_neg
ColumnMajorMatrix _dstrain_diag_neg
Definition: RateDepSmearIsoCrackModel.h:52
pow
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
Definition: ExpressionBuilder.h:673
RateDepSmearCrackModel::initVariables
virtual void initVariables()
Definition: RateDepSmearCrackModel.C:108
RateDepSmearCrackModel
RateDepSmearCrackModel is the base class for rate dependent continuum damage model.
Definition: RateDepSmearCrackModel.h:32
RateDepSmearIsoCrackModel::calcJacobian
virtual void calcJacobian()
This function calculates the Jacobian.
Definition: RateDepSmearIsoCrackModel.C:133
RateDepSmearCrackModel::_intvar_tmp
std::vector< Real > _intvar_tmp
Definition: RateDepSmearCrackModel.h:115
RateDepSmearIsoCrackModel::_s0_diag_pos
ColumnMajorMatrix _s0_diag_pos
Definition: RateDepSmearIsoCrackModel.h:51
RateDepSmearIsoCrackModel::damageRate
virtual Real damageRate()
This function calculates rate of damage based on energy.
Definition: RateDepSmearIsoCrackModel.C:78
RateDepSmearIsoCrackModel.h
RateDepSmearCrackModel::_zero_tol
Real _zero_tol
Relative tolerance factor for convergence of the Newton Raphson solve.
Definition: RateDepSmearCrackModel.h:103
RateDepSmearCrackModel::_dstrain
SymmTensor _dstrain
Definition: RateDepSmearCrackModel.h:121
RateDepSmearIsoCrackModel::postSolveStress
virtual void postSolveStress()
This function updates the stress after solve.
Definition: RateDepSmearIsoCrackModel.C:143
SymmTensor::columnMajorMatrix
ColumnMajorMatrix columnMajorMatrix() const
Definition: SymmTensor.h:425
RateDepSmearCrackModel::_err_tol
bool _err_tol
Convergence flag.
Definition: RateDepSmearCrackModel.h:124
RateDepSmearCrackModel::_intvar_incr
std::vector< Real > _intvar_incr
Definition: RateDepSmearCrackModel.h:114
RateDepSmearIsoCrackModel::_dstrain_diag
ColumnMajorMatrix _dstrain_diag
Definition: RateDepSmearIsoCrackModel.h:52
RateDepSmearCrackModel::_jac
std::vector< Real > _jac
Definition: RateDepSmearCrackModel.h:117
RateDepSmearCrackModel::_ref_damage_rate
Real _ref_damage_rate
Definition: RateDepSmearCrackModel.h:98
RateDepSmearIsoCrackModel::~RateDepSmearIsoCrackModel
virtual ~RateDepSmearIsoCrackModel()
Definition: RateDepSmearIsoCrackModel.C:149
RateDepSmearIsoCrackModel::_crit_energy
Real _crit_energy
Definition: RateDepSmearIsoCrackModel.h:41
validParams< RateDepSmearIsoCrackModel >
InputParameters validParams< RateDepSmearIsoCrackModel >()
Definition: RateDepSmearIsoCrackModel.C:16
validParams< RateDepSmearCrackModel >
InputParameters validParams< RateDepSmearCrackModel >()
Definition: RateDepSmearCrackModel.C:16
RateDepSmearIsoCrackModel::RateDepSmearIsoCrackModel
RateDepSmearIsoCrackModel(const InputParameters &parameters)
Definition: RateDepSmearIsoCrackModel.C:30
RateDepSmearCrackModel::_nstate
unsigned int _nstate
reference damage rate
Definition: RateDepSmearCrackModel.h:99
RateDepSmearIsoCrackModel::_s0_diag_neg
ColumnMajorMatrix _s0_diag_neg
Definition: RateDepSmearIsoCrackModel.h:51
RateDepSmearIsoCrackModel::_s0_diag
ColumnMajorMatrix _s0_diag
Definition: RateDepSmearIsoCrackModel.h:51
RateDepSmearIsoCrackModel::_kfail
Real _kfail
Required parameter.
Definition: RateDepSmearIsoCrackModel.h:42
RateDepSmearCrackModel::_stress0
SymmTensor _stress0
Definition: RateDepSmearCrackModel.h:122
RateDepSmearCrackModel::initQpStatefulProperties
virtual void initQpStatefulProperties()
Definition: RateDepSmearCrackModel.C:69
RateDepSmearCrackModel::_intvar_old_tmp
std::vector< Real > _intvar_old_tmp
Definition: RateDepSmearCrackModel.h:115
RateDepSmearCrackModel::_intvar
MaterialProperty< std::vector< Real > > & _intvar
Variable value.
Definition: RateDepSmearCrackModel.h:108
RateDepSmearCrackModel::_exponent
Real _exponent
Number of state variables.
Definition: RateDepSmearCrackModel.h:100
RateDepSmearIsoCrackModel
In this class a rate dependent isotropic damage model is implemented.
Definition: RateDepSmearIsoCrackModel.h:23
RateDepSmearIsoCrackModel::initVariables
virtual void initVariables()
Definition: RateDepSmearIsoCrackModel.C:55
RateDepSmearIsoCrackModel::_energy_old
const MaterialProperty< Real > & _energy_old
Definition: RateDepSmearIsoCrackModel.h:46
RateDepSmearIsoCrackModel::_dstrain_diag_pos
ColumnMajorMatrix _dstrain_diag_pos
Definition: RateDepSmearIsoCrackModel.h:52
RateDepSmearIsoCrackModel::calcStateIncr
virtual void calcStateIncr()
This function calculated thw increment of the state variables (dv) used to form the residual.
Definition: RateDepSmearIsoCrackModel.C:109
RateDepSmearCrackModel::_intvar_old
const MaterialProperty< std::vector< Real > > & _intvar_old
Definition: RateDepSmearCrackModel.h:109
RateDepSmearIsoCrackModel::_upper_lim_damage
Real _upper_lim_damage
Used to avoid non-positive definiteness.
Definition: RateDepSmearIsoCrackModel.h:43