https://mooseframework.inl.gov
ScalarDamageBase.h
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 #pragma once
11 
12 #include "DamageBase.h"
13 
14 // Forward declaration
15 
19 template <bool is_ad>
21 {
22 public:
24 
26 
27  virtual void initQpStatefulProperties() override;
28 
29  virtual void updateDamage() override;
30 
31  virtual void updateStressForDamage(GenericRankTwoTensor<is_ad> & stress_new) override;
32 
33  virtual void updateJacobianMultForDamage(RankFourTensor & jacobian_mult) override;
34 
35  virtual void computeUndamagedOldStress(RankTwoTensor & stress_old) override;
36 
37  virtual Real computeTimeStepLimit() override;
38 
42  const GenericReal<is_ad> & getQpDamageIndex(unsigned int qp);
43 
47  const std::string getDamageIndexName() const { return _damage_index_name; }
48 
49 protected:
51  const MaterialPropertyName _damage_index_name;
52 
54  virtual void updateQpDamageIndex() = 0;
55 
61 
63  const bool _use_old_damage;
64 
67 
70 
73 
78 };
79 
Moose::GenericType< Real, is_ad > GenericReal
virtual void updateJacobianMultForDamage(RankFourTensor &jacobian_mult) override
Update the material constitutive matrix.
const bool _use_old_damage
If true, use the damage index from the old state (rather than the current state)
virtual void updateQpDamageIndex()=0
Update the damage index at the current qpoint.
virtual void updateStressForDamage(GenericRankTwoTensor< is_ad > &stress_new) override
Update the current stress tensor for effects of damage.
const InputParameters & parameters() const
static InputParameters validParams()
const Real & _maximum_damage_increment
Maximum damage increment allowed for the time step.
virtual Real computeTimeStepLimit() override
Compute the limiting value of the time step for this material.
GenericMaterialProperty< Real, is_ad > & _damage_index
Material property that provides the damage index.
ScalarDamageBaseTempl< false > ScalarDamageBase
const GenericReal< is_ad > & getQpDamageIndex(unsigned int qp)
Get the value of the damage index for the current quadrature point.
const Real & _residual_stiffness_fraction
Residual fraction of stiffness used for material that is fully damaged.
const MaterialPropertyName _damage_index_name
Name of the material property where the damage index is stored.
virtual void initQpStatefulProperties() override
typename GenericMaterialPropertyStruct< T, is_ad >::type GenericMaterialProperty
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)
ScalarDamageBaseTempl< true > ADScalarDamageBase
const Real & _maximum_damage
Maximum allowed value for the damage index.
const MaterialProperty< Real > & _damage_index_older
virtual void computeUndamagedOldStress(RankTwoTensor &stress_old) override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const MaterialProperty< Real > & _damage_index_old
const std::string getDamageIndexName() const
Get the name of the material property containing the damage index.
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