18 params.
addClassDescription(
"Base class for damage model based on a scalar damage parameter");
22 "Whether to use the damage index from the previous step in the stress computation");
24 "residual_stiffness_fraction",
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)");
30 "maximum_damage_increment",
32 "maximum_damage_increment>0 & maximum_damage_increment<1",
33 "maximum damage increment allowed for simulations with adaptive time step");
36 "maximum_damage>=0 & maximum_damage<=1",
37 "Maximum value allowed for damage index");
38 params.
addParam<MaterialPropertyName>(
41 "name of the material property where the damage index is stored");
48 _damage_index_name(this->template getParam<MaterialPropertyName>(
"damage_index_name")),
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)),
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"))
65 _damage_index[_qp] = 0.0;
74 return _damage_index[_qp];
81 updateQpDamageIndex();
82 _damage_index[_qp] = std::min(_maximum_damage, _damage_index[_qp]);
92 std::max((1.0 - (_use_old_damage ? _damage_index_old[_qp] : _damage_index[_qp])), 0.0);
99 Real damage_index_old =
100 std::max((1.0 - (_use_old_damage ? _damage_index_older[_qp] : _damage_index_old[_qp])), 0.0);
102 if (damage_index_old > 0.0)
103 stress_old /= damage_index_old;
106 template <
bool is_ad>
110 jacobian_mult *= std::max((1.0 - (_use_old_damage ? _damage_index_old[_qp]
112 _residual_stiffness_fraction);
115 template <
bool is_ad>
119 Real current_damage_increment =
122 return std::numeric_limits<Real>::max();
124 return _dt * _maximum_damage_increment / current_damage_increment;
Moose::GenericType< Real, is_ad > GenericReal
virtual void updateJacobianMultForDamage(RankFourTensor &jacobian_mult) override
Update the material constitutive matrix.
bool absoluteFuzzyEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
virtual void updateStressForDamage(GenericRankTwoTensor< is_ad > &stress_new) override
Update the current stress tensor for effects of damage.
static InputParameters validParams()
virtual Real computeTimeStepLimit() override
Compute the limiting value of the time step for this material.
static InputParameters validParams()
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...
ScalarDamageBaseTempl(const InputParameters ¶meters)
virtual void computeUndamagedOldStress(RankTwoTensor &stress_old) override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void updateDamage() override
Update the internal variable(s) that evolve the damage.
Base class for scalar damage models.
Moose::GenericType< RankTwoTensor, is_ad > GenericRankTwoTensor