www.mooseframework.org
ScalarDamageBase.h
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 #pragma once
11 
12 #include "DamageBase.h"
13 
14 // Forward declaration
15 class ScalarDamageBase;
16 
17 template <>
18 InputParameters validParams<ScalarDamageBase>();
19 
24 {
25 public:
26  static InputParameters validParams();
27 
28  ScalarDamageBase(const InputParameters & parameters);
29 
30  virtual void initQpStatefulProperties() override;
31 
32  virtual void updateDamage() override;
33 
34  virtual void updateStressForDamage(RankTwoTensor & stress_new) override;
35 
36  virtual void updateJacobianMultForDamage(RankFourTensor & jacobian_mult) override;
37 
38  virtual void computeUndamagedOldStress(RankTwoTensor & stress_old) override;
39 
40  virtual Real computeTimeStepLimit() override;
41 
42  const Real & getQpDamageIndex(unsigned int qp);
43 
44  const std::string getDamageIndexName() const { return _damage_index_name; }
45 
46 protected:
48  const std::string _damage_index_name;
49 
51  virtual void updateQpDamageIndex() = 0;
52 
54  MaterialProperty<Real> & _damage_index;
55  const MaterialProperty<Real> & _damage_index_old;
56  const MaterialProperty<Real> & _damage_index_older;
58 
60  const bool _use_old_damage;
61 
64 
67 };
ScalarDamageBase::updateJacobianMultForDamage
virtual void updateJacobianMultForDamage(RankFourTensor &jacobian_mult) override
Update the material constitutive matrix.
Definition: ScalarDamageBase.C:93
validParams< ScalarDamageBase >
InputParameters validParams< ScalarDamageBase >()
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
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
DamageBase.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::getDamageIndexName
const std::string getDamageIndexName() const
Definition: ScalarDamageBase.h:44
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 >
ScalarDamageBase::_damage_index_name
const std::string _damage_index_name
Name of the material property where the damage index is stored.
Definition: ScalarDamageBase.h:48