https://mooseframework.inl.gov
ScalarDamageBase.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
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 
10 #include "ScalarDamageBase.h"
11 #include "MooseUtils.h"
12 
13 template <bool is_ad>
16 {
18  params.addClassDescription("Base class for damage model based on a scalar damage parameter");
19  params.addParam<bool>(
20  "use_old_damage",
21  false,
22  "Whether to use the damage index from the previous step in the stress computation");
23  params.addRangeCheckedParam<Real>(
24  "residual_stiffness_fraction",
25  1.e-8,
26  "residual_stiffness_fraction>=0 & residual_stiffness_fraction<1",
27  "Minimum fraction of original material stiffness retained for fully "
28  "damaged material (when damage_index=1)");
29  params.addRangeCheckedParam<Real>(
30  "maximum_damage_increment",
31  0.1,
32  "maximum_damage_increment>0 & maximum_damage_increment<1",
33  "maximum damage increment allowed for simulations with adaptive time step");
34  params.addRangeCheckedParam<Real>("maximum_damage",
35  1.0,
36  "maximum_damage>=0 & maximum_damage<=1",
37  "Maximum value allowed for damage index");
38  params.addParam<MaterialPropertyName>(
39  "damage_index_name",
40  "damage_index",
41  "name of the material property where the damage index is stored");
42  return params;
43 }
44 
45 template <bool is_ad>
47  : DamageBaseTempl<is_ad>(parameters),
48  _damage_index_name(this->template getParam<MaterialPropertyName>("damage_index_name")),
49  _damage_index(
50  this->template declareGenericPropertyByName<Real, is_ad>(_base_name + _damage_index_name)),
51  _damage_index_old(this->template getMaterialPropertyOld<Real>(_base_name + _damage_index_name)),
52  _damage_index_older(
53  this->template getMaterialPropertyOlder<Real>(_base_name + _damage_index_name)),
54  _use_old_damage(this->template getParam<bool>("use_old_damage")),
55  _residual_stiffness_fraction(this->template getParam<Real>("residual_stiffness_fraction")),
56  _maximum_damage_increment(this->template getParam<Real>("maximum_damage_increment")),
57  _maximum_damage(this->template getParam<Real>("maximum_damage"))
58 {
59 }
60 
61 template <bool is_ad>
62 void
64 {
65  _damage_index[_qp] = 0.0;
66 }
67 
68 template <bool is_ad>
69 const GenericReal<is_ad> &
71 {
72  setQp(qp);
73  updateDamage();
74  return _damage_index[_qp];
75 }
76 
77 template <bool is_ad>
78 void
80 {
81  using std::min;
82  updateQpDamageIndex();
83  _damage_index[_qp] = min(_maximum_damage, _damage_index[_qp]);
84 }
85 
86 template <bool is_ad>
87 void
89 {
90  using std::max;
91  // Avoid multiplying by a small negative number, which could occur if damage_index
92  // is slightly greater than 1.0
93  stress_new *= max((1.0 - (_use_old_damage ? _damage_index_old[_qp] : _damage_index[_qp])), 0.0);
94 }
95 
96 template <bool is_ad>
97 void
99 {
100  Real damage_index_old =
101  std::max((1.0 - (_use_old_damage ? _damage_index_older[_qp] : _damage_index_old[_qp])), 0.0);
102 
103  if (damage_index_old > 0.0)
104  stress_old /= damage_index_old;
105 }
106 
107 template <bool is_ad>
108 void
110 {
111  jacobian_mult *= std::max((1.0 - (_use_old_damage ? _damage_index_old[_qp]
112  : MetaPhysicL::raw_value(_damage_index[_qp]))),
113  _residual_stiffness_fraction);
114 }
115 
116 template <bool is_ad>
117 Real
119 {
120  Real current_damage_increment =
121  (MetaPhysicL::raw_value(_damage_index[_qp]) - _damage_index_old[_qp]);
122  if (MooseUtils::absoluteFuzzyEqual(current_damage_increment, 0.0))
123  return std::numeric_limits<Real>::max();
124 
125  return _dt * _maximum_damage_increment / current_damage_increment;
126 }
127 
128 template class ScalarDamageBaseTempl<false>;
129 template class ScalarDamageBaseTempl<true>;
Moose::GenericType< Real, is_ad > GenericReal
virtual void updateJacobianMultForDamage(RankFourTensor &jacobian_mult) override
Update the material constitutive matrix.
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
virtual void updateStressForDamage(GenericRankTwoTensor< is_ad > &stress_new) override
Update the current stress tensor for effects of damage.
auto raw_value(const Eigen::Map< T > &in)
static InputParameters validParams()
auto max(const L &left, const R &right)
virtual Real computeTimeStepLimit() override
Compute the limiting value of the time step for this material.
static InputParameters validParams()
Definition: DamageBase.C:14
const GenericReal< is_ad > & getQpDamageIndex(unsigned int qp)
Get the value of the damage index for the current quadrature point.
virtual void initQpStatefulProperties() override
DamageBase is a base class for damage models, which modify the stress tensor computed by another mode...
Definition: DamageBase.h:25
ScalarDamageBaseTempl(const InputParameters &parameters)
virtual void computeUndamagedOldStress(RankTwoTensor &stress_old) override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void addClassDescription(const std::string &doc_string)
void addRangeCheckedParam(const std::string &name, const T &value, const std::string &parsed_function, const std::string &doc_string)
virtual void updateDamage() override
Update the internal variable(s) that evolve the damage.
auto min(const L &left, const R &right)
Base class for scalar damage models.
Moose::GenericType< RankTwoTensor, is_ad > GenericRankTwoTensor