Line data Source code
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> 14 : InputParameters 15 1104 : ScalarDamageBaseTempl<is_ad>::validParams() 16 : { 17 1104 : InputParameters params = DamageBaseTempl<is_ad>::validParams(); 18 1104 : params.addClassDescription("Base class for damage model based on a scalar damage parameter"); 19 2208 : params.addParam<bool>( 20 : "use_old_damage", 21 2208 : false, 22 : "Whether to use the damage index from the previous step in the stress computation"); 23 3312 : params.addRangeCheckedParam<Real>( 24 : "residual_stiffness_fraction", 25 2208 : 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 3312 : params.addRangeCheckedParam<Real>( 30 : "maximum_damage_increment", 31 2208 : 0.1, 32 : "maximum_damage_increment>0 & maximum_damage_increment<1", 33 : "maximum damage increment allowed for simulations with adaptive time step"); 34 3312 : params.addRangeCheckedParam<Real>("maximum_damage", 35 2208 : 1.0, 36 : "maximum_damage>=0 & maximum_damage<=1", 37 : "Maximum value allowed for damage index"); 38 2208 : params.addParam<MaterialPropertyName>( 39 : "damage_index_name", 40 : "damage_index", 41 : "name of the material property where the damage index is stored"); 42 1104 : return params; 43 0 : } 44 : 45 : template <bool is_ad> 46 828 : ScalarDamageBaseTempl<is_ad>::ScalarDamageBaseTempl(const InputParameters & parameters) 47 : : DamageBaseTempl<is_ad>(parameters), 48 1656 : _damage_index_name(this->template getParam<MaterialPropertyName>("damage_index_name")), 49 828 : _damage_index( 50 828 : this->template declareGenericPropertyByName<Real, is_ad>(_base_name + _damage_index_name)), 51 1656 : _damage_index_old(this->template getMaterialPropertyOld<Real>(_base_name + _damage_index_name)), 52 828 : _damage_index_older( 53 828 : this->template getMaterialPropertyOlder<Real>(_base_name + _damage_index_name)), 54 1656 : _use_old_damage(this->template getParam<bool>("use_old_damage")), 55 1656 : _residual_stiffness_fraction(this->template getParam<Real>("residual_stiffness_fraction")), 56 1656 : _maximum_damage_increment(this->template getParam<Real>("maximum_damage_increment")), 57 2484 : _maximum_damage(this->template getParam<Real>("maximum_damage")) 58 : { 59 828 : } 60 : 61 : template <bool is_ad> 62 : void 63 25920 : ScalarDamageBaseTempl<is_ad>::initQpStatefulProperties() 64 : { 65 33920 : _damage_index[_qp] = 0.0; 66 25920 : } 67 : 68 : template <bool is_ad> 69 : const GenericReal<is_ad> & 70 2118000 : ScalarDamageBaseTempl<is_ad>::getQpDamageIndex(unsigned int qp) 71 : { 72 : setQp(qp); 73 2145904 : updateDamage(); 74 2145904 : return _damage_index[_qp]; 75 : } 76 : 77 : template <bool is_ad> 78 : void 79 4513116 : ScalarDamageBaseTempl<is_ad>::updateDamage() 80 : { 81 : using std::min; 82 4513116 : updateQpDamageIndex(); 83 4513112 : _damage_index[_qp] = min(_maximum_damage, _damage_index[_qp]); 84 4513112 : } 85 : 86 : template <bool is_ad> 87 : void 88 2367208 : ScalarDamageBaseTempl<is_ad>::updateStressForDamage(GenericRankTwoTensor<is_ad> & stress_new) 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 3328456 : stress_new *= max((1.0 - (_use_old_damage ? _damage_index_old[_qp] : _damage_index[_qp])), 0.0); 94 2367208 : } 95 : 96 : template <bool is_ad> 97 : void 98 12640 : ScalarDamageBaseTempl<is_ad>::computeUndamagedOldStress(RankTwoTensor & stress_old) 99 : { 100 12640 : Real damage_index_old = 101 12640 : std::max((1.0 - (_use_old_damage ? _damage_index_older[_qp] : _damage_index_old[_qp])), 0.0); 102 : 103 12640 : if (damage_index_old > 0.0) 104 12640 : stress_old /= damage_index_old; 105 12640 : } 106 : 107 : template <bool is_ad> 108 : void 109 1405960 : ScalarDamageBaseTempl<is_ad>::updateJacobianMultForDamage(RankFourTensor & jacobian_mult) 110 : { 111 1427940 : jacobian_mult *= std::max((1.0 - (_use_old_damage ? _damage_index_old[_qp] 112 1191576 : : MetaPhysicL::raw_value(_damage_index[_qp]))), 113 1405960 : _residual_stiffness_fraction); 114 1405960 : } 115 : 116 : template <bool is_ad> 117 : Real 118 2367208 : ScalarDamageBaseTempl<is_ad>::computeTimeStepLimit() 119 : { 120 : Real current_damage_increment = 121 2367208 : (MetaPhysicL::raw_value(_damage_index[_qp]) - _damage_index_old[_qp]); 122 2367208 : if (MooseUtils::absoluteFuzzyEqual(current_damage_increment, 0.0)) 123 : return std::numeric_limits<Real>::max(); 124 : 125 1791632 : return _dt * _maximum_damage_increment / current_damage_increment; 126 : } 127 : 128 : template class ScalarDamageBaseTempl<false>; 129 : template class ScalarDamageBaseTempl<true>;