www.mooseframework.org
ScalarDamageBase.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 
10 #include "ScalarDamageBase.h"
11 #include "MooseUtils.h"
12 
14 
15 InputParameters
17 {
18  InputParameters params = DamageBase::validParams();
19  params.addClassDescription("Base class for damage model based on a scalar damage parameter");
20  params.addParam<bool>(
21  "use_old_damage",
22  false,
23  "Whether to use the damage index from the previous step in the stress computation");
24  params.addRangeCheckedParam<Real>(
25  "residual_stiffness_fraction",
26  1.e-8,
27  "residual_stiffness_fraction>=0 & residual_stiffness_fraction<1",
28  "Minimum fraction of original material stiffness retained for fully "
29  "damaged material (when damage_index=1)");
30  params.addRangeCheckedParam<Real>(
31  "maximum_damage_increment",
32  0.1,
33  "maximum_damage_increment>0 & maximum_damage_increment<1",
34  "maximum damage increment allowed for simulations with adaptive time step");
35  params.addParam<std::string>("damage_index_name",
36  "damage_index",
37  "name of the material property where the damage index is stored");
38  return params;
39 }
40 
41 ScalarDamageBase::ScalarDamageBase(const InputParameters & parameters)
42  : DamageBase(parameters),
43  _damage_index_name(getParam<std::string>("damage_index_name")),
44  _damage_index(declareProperty<Real>(_base_name + _damage_index_name)),
45  _damage_index_old(getMaterialPropertyOld<Real>(_base_name + _damage_index_name)),
46  _damage_index_older(getMaterialPropertyOlder<Real>(_base_name + _damage_index_name)),
47  _use_old_damage(getParam<bool>("use_old_damage")),
48  _residual_stiffness_fraction(getParam<Real>("residual_stiffness_fraction")),
49  _maximum_damage_increment(getParam<Real>("maximum_damage_increment"))
50 {
51 }
52 
53 void
55 {
56  _damage_index[_qp] = 0.0;
57 }
58 
59 const Real &
61 {
62  setQp(qp);
64  return _damage_index[_qp];
65 }
66 
67 void
69 {
71 }
72 
73 void
75 {
76  // Avoid multiplying by a small negative number, which could occur if damage_index
77  // is slightly greater than 1.0
78  stress_new *=
79  std::max((1.0 - (_use_old_damage ? _damage_index_old[_qp] : _damage_index[_qp])), 0.0);
80 }
81 
82 void
84 {
85  Real damage_index_old =
86  std::max((1.0 - (_use_old_damage ? _damage_index_older[_qp] : _damage_index_old[_qp])), 0.0);
87 
88  if (damage_index_old > 0.0)
89  stress_old /= damage_index_old;
90 }
91 
92 void
94 {
95  jacobian_mult *= std::max((1.0 - (_use_old_damage ? _damage_index_old[_qp] : _damage_index[_qp])),
97 }
98 
99 Real
101 {
102  Real current_damage_increment = (_damage_index[_qp] - _damage_index_old[_qp]);
103  if (MooseUtils::absoluteFuzzyEqual(current_damage_increment, 0.0))
104  return std::numeric_limits<Real>::max();
105 
106  return _dt * _maximum_damage_increment / current_damage_increment;
107 }
ScalarDamageBase::updateJacobianMultForDamage
virtual void updateJacobianMultForDamage(RankFourTensor &jacobian_mult) override
Update the material constitutive matrix.
Definition: ScalarDamageBase.C:93
ScalarDamageBase::updateDamage
virtual void updateDamage() override
Update the internal variable(s) that evolve the damage.
Definition: ScalarDamageBase.C:68
ScalarDamageBase::_residual_stiffness_fraction
const Real & _residual_stiffness_fraction
Residual fraction of stiffness used for material that is fully damaged.
Definition: ScalarDamageBase.h:63
ScalarDamageBase::initQpStatefulProperties
virtual void initQpStatefulProperties() override
Definition: ScalarDamageBase.C:54
ScalarDamageBase::_damage_index_old
const MaterialProperty< Real > & _damage_index_old
Definition: ScalarDamageBase.h:55
ScalarDamageBase::validParams
static InputParameters validParams()
Definition: ScalarDamageBase.C:16
ScalarDamageBase
Base class for scalar damage models.
Definition: ScalarDamageBase.h:23
ScalarDamageBase::_use_old_damage
const bool _use_old_damage
If true, use the damage index from the old state (rather than the current state)
Definition: ScalarDamageBase.h:60
ScalarDamageBase::ScalarDamageBase
ScalarDamageBase(const InputParameters &parameters)
Definition: ScalarDamageBase.C:41
DamageBase::validParams
static InputParameters validParams()
Definition: DamageBase.C:15
defineLegacyParams
defineLegacyParams(ScalarDamageBase)
DamageBase
DamageBase is a base class for damage models, which modify the stress tensor computed by another mode...
Definition: DamageBase.h:28
ScalarDamageBase::_maximum_damage_increment
const Real & _maximum_damage_increment
Maximum damage increment allowed for the time step.
Definition: ScalarDamageBase.h:66
ScalarDamageBase.h
ScalarDamageBase::_damage_index_older
const MaterialProperty< Real > & _damage_index_older
Definition: ScalarDamageBase.h:56
ScalarDamageBase::updateStressForDamage
virtual void updateStressForDamage(RankTwoTensor &stress_new) override
Update the current stress tensor for effects of damage.
Definition: ScalarDamageBase.C:74
ScalarDamageBase::getQpDamageIndex
const Real & getQpDamageIndex(unsigned int qp)
Definition: ScalarDamageBase.C:60
ScalarDamageBase::_damage_index
MaterialProperty< Real > & _damage_index
Material property that provides the damage index.
Definition: ScalarDamageBase.h:54
ScalarDamageBase::computeTimeStepLimit
virtual Real computeTimeStepLimit() override
Compute the limiting value of the time step for this material.
Definition: ScalarDamageBase.C:100
RankFourTensorTempl< Real >
ScalarDamageBase::computeUndamagedOldStress
virtual void computeUndamagedOldStress(RankTwoTensor &stress_old) override
Definition: ScalarDamageBase.C:83
ScalarDamageBase::updateQpDamageIndex
virtual void updateQpDamageIndex()=0
Update the damage index at the current qpoint.
RankTwoTensorTempl< Real >
DamageBase::setQp
void setQp(unsigned int qp)
Sets the value of the member variable _qp for use in inheriting classes.
Definition: DamageBase.C:40