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 : #pragma once 11 : 12 : #include "RadialReturnStressUpdate.h" 13 : 14 : /** 15 : * This class uses the Discrete material in a radial return isotropic plasticity 16 : * model. This class is one of the basic radial return constitutive models; 17 : * more complex constitutive models combine creep and plasticity. 18 : * 19 : * This class inherits from RadialReturnStressUpdate and must be used 20 : * in conjunction with ComputeReturnMappingStress. This class calculates 21 : * an effective trial stress, an effective scalar plastic strain 22 : * increment, and the derivative of the scalar effective plastic strain increment; 23 : * these values are passed to the RadialReturnStressUpdate to compute 24 : * the radial return stress increment. This isotropic plasticity class also 25 : * computes the plastic strain as a stateful material property. 26 : * 27 : * This class is based on the implicit integration algorithm in F. Dunne and N. 28 : * Petrinic's Introduction to Computational Plasticity (2004) Oxford University 29 : * Press, pg. 146 - 149. 30 : */ 31 : 32 : template <bool is_ad> 33 : class IsotropicPlasticityStressUpdateTempl : public RadialReturnStressUpdateTempl<is_ad> 34 : { 35 : public: 36 : static InputParameters validParams(); 37 : 38 : IsotropicPlasticityStressUpdateTempl(const InputParameters & parameters); 39 : 40 : using Material::_qp; 41 : using RadialReturnStressUpdateTempl<is_ad>::_base_name; 42 : using RadialReturnStressUpdateTempl<is_ad>::_three_shear_modulus; 43 : 44 : virtual void 45 : computeStressInitialize(const GenericReal<is_ad> & effective_trial_stress, 46 : const GenericRankFourTensor<is_ad> & elasticity_tensor) override; 47 : virtual GenericReal<is_ad> computeResidual(const GenericReal<is_ad> & effective_trial_stress, 48 : const GenericReal<is_ad> & scalar) override; 49 : virtual GenericReal<is_ad> computeDerivative(const GenericReal<is_ad> & effective_trial_stress, 50 : const GenericReal<is_ad> & scalar) override; 51 : 52 : virtual void computeYieldStress(const GenericRankFourTensor<is_ad> & elasticity_tensor); 53 : 54 238688 : GenericReal<is_ad> yieldCondition() const { return _yield_condition; } 55 : 56 : virtual void 57 : computeStressFinalize(const GenericRankTwoTensor<is_ad> & plastic_strain_increment) override; 58 : 59 : protected: 60 : virtual void initQpStatefulProperties() override; 61 : virtual void propagateQpStatefulProperties() override; 62 : 63 : virtual void iterationFinalize(const GenericReal<is_ad> & scalar) override; 64 : 65 : virtual GenericReal<is_ad> computeHardeningValue(const GenericReal<is_ad> & scalar); 66 : virtual GenericReal<is_ad> computeHardeningDerivative(const GenericReal<is_ad> & scalar); 67 : 68 : /// a string to prepend to the plastic strain Material Property name 69 : const std::string _plastic_prepend; 70 : 71 : const Function * _yield_stress_function; 72 : GenericReal<is_ad> _yield_stress; 73 : const Real _hardening_constant; 74 : const Function * const _hardening_function; 75 : 76 : GenericReal<is_ad> _yield_condition; 77 : GenericReal<is_ad> _hardening_slope; 78 : 79 : /// plastic strain in this model 80 : GenericMaterialProperty<RankTwoTensor, is_ad> & _plastic_strain; 81 : 82 : /// old value of plastic strain 83 : const MaterialProperty<RankTwoTensor> & _plastic_strain_old; 84 : 85 : GenericMaterialProperty<Real, is_ad> & _hardening_variable; 86 : const MaterialProperty<Real> & _hardening_variable_old; 87 : const GenericVariableValue<is_ad> & _temperature; 88 : }; 89 : 90 : typedef IsotropicPlasticityStressUpdateTempl<false> IsotropicPlasticityStressUpdate; 91 : typedef IsotropicPlasticityStressUpdateTempl<true> ADIsotropicPlasticityStressUpdate;