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 "AnisotropicReturnPlasticityStressUpdateBase.h" 13 : 14 : /** 15 : * This class uses the stress update material in an anisotropic return mapping. 16 : * This class is one of the generalized radial return constitutive models based on Hill's criterion; 17 : * it assumes and isotropic elasticity tensor and an anisotropic plastic yield surface. 18 : * Constitutive models that combine creep and plasticity can be used. 19 : */ 20 : 21 : template <bool is_ad> 22 : class HillPlasticityStressUpdateTempl 23 : : public AnisotropicReturnPlasticityStressUpdateBaseTempl<is_ad> 24 : { 25 : public: 26 : static InputParameters validParams(); 27 : 28 : HillPlasticityStressUpdateTempl(const InputParameters & parameters); 29 : 30 : protected: 31 : using AnisotropicReturnPlasticityStressUpdateBaseTempl<is_ad>::_effective_inelastic_strain; 32 : using AnisotropicReturnPlasticityStressUpdateBaseTempl<is_ad>::_effective_inelastic_strain_old; 33 : using AnisotropicReturnPlasticityStressUpdateBaseTempl<is_ad>::_plasticity_strain; 34 : using AnisotropicReturnPlasticityStressUpdateBaseTempl<is_ad>::_plasticity_strain_old; 35 : using AnisotropicReturnPlasticityStressUpdateBaseTempl<is_ad>::isBlockDiagonal; 36 : 37 : using Material::_qp; 38 : 39 : virtual void 40 : computeStressInitialize(const GenericDenseVector<is_ad> & stress_dev, 41 : const GenericDenseVector<is_ad> & stress, 42 : const GenericRankFourTensor<is_ad> & elasticity_tensor) override; 43 : virtual GenericReal<is_ad> 44 : computeResidual(const GenericDenseVector<is_ad> & effective_trial_stress, 45 : const GenericDenseVector<is_ad> & stress_new, 46 : const GenericReal<is_ad> & scalar) override; 47 : virtual GenericReal<is_ad> 48 : computeDerivative(const GenericDenseVector<is_ad> & effective_trial_stress, 49 : const GenericDenseVector<is_ad> & stress_new, 50 : const GenericReal<is_ad> & scalar) override; 51 : 52 : virtual Real 53 : computeReferenceResidual(const GenericDenseVector<is_ad> & effective_trial_stress, 54 : const GenericDenseVector<is_ad> & stress_new, 55 : const GenericReal<is_ad> & residual, 56 : const GenericReal<is_ad> & scalar_effective_inelastic_strain) override; 57 : virtual void propagateQpStatefulProperties() override; 58 : /** 59 : * Does the model require the elasticity tensor to be isotropic? Yes, this class only does 60 : * anisotropic *plasticity* 61 : */ 62 18 : bool requiresIsotropicTensor() override { return true; } 63 : 64 : Real computeHardeningDerivative(); 65 : GenericReal<is_ad> computeHardeningValue(const GenericReal<is_ad> & scalar, 66 : const GenericReal<is_ad> & omega); 67 : /** 68 : * Compute eigendecomposition of Hill's tensor for anisotropic plasticity 69 : * @param hill_tensor 6x6 matrix representing fourth order Hill's tensor describing anisotropy 70 : */ 71 : void computeHillTensorEigenDecomposition(const DenseMatrix<Real> & hill_tensor); 72 : 73 : /** 74 : * Perform any necessary steps to finalize strain increment after return mapping iterations 75 : * @param inelasticStrainIncrement Inelastic strain increment 76 : */ 77 : virtual void computeStrainFinalize(GenericRankTwoTensor<is_ad> & /*inelasticStrainIncrement*/, 78 : const GenericRankTwoTensor<is_ad> & /*stress*/, 79 : const GenericDenseVector<is_ad> & /*stress_dev*/, 80 : const GenericReal<is_ad> & /*delta_gamma*/) override; 81 : 82 : /** 83 : * Perform any necessary steps to finalize state after return mapping iterations 84 : * @param inelasticStrainIncrement Inelastic strain increment 85 : */ 86 : virtual void 87 : computeStressFinalize(const GenericRankTwoTensor<is_ad> & inelasticStrainIncrement, 88 : const GenericReal<is_ad> & delta_gamma, 89 : GenericRankTwoTensor<is_ad> & stress, 90 : const GenericDenseVector<is_ad> & /*stress_dev*/, 91 : const GenericRankTwoTensor<is_ad> & /*stress_old*/, 92 : const GenericRankFourTensor<is_ad> & /*elasticity_tensor*/) override; 93 : 94 : GenericReal<is_ad> computeOmega(const GenericReal<is_ad> & delta_gamma, 95 : const GenericDenseVector<is_ad> & stress_trial); 96 : 97 : void computeDeltaDerivatives(const GenericReal<is_ad> & delta_gamma, 98 : const GenericDenseVector<is_ad> & stress_trial, 99 : const GenericReal<is_ad> & sy_alpha, 100 : GenericReal<is_ad> & omega, 101 : GenericReal<is_ad> & omega_gamma, 102 : GenericReal<is_ad> & sy_gamma); 103 : 104 : /// Square of the q function for orthotropy 105 : GenericReal<is_ad> _qsigma; 106 : 107 : GenericDenseVector<is_ad> _eigenvalues_hill; 108 : GenericDenseMatrix<is_ad> _eigenvectors_hill; 109 : 110 : const Real _hardening_constant; 111 : const Real _hardening_exponent; 112 : 113 : GenericMaterialProperty<Real, is_ad> & _hardening_variable; 114 : const MaterialProperty<Real> & _hardening_variable_old; 115 : GenericReal<is_ad> _hardening_derivative; 116 : GenericReal<is_ad> _yield_condition; 117 : GenericReal<is_ad> _yield_stress; 118 : 119 : /// Hill tensor, when global axes do not (somehow) align with those of the material 120 : /// Example: Large rotation due to rigid body and/or large deformation kinematics 121 : const MaterialProperty<DenseMatrix<Real>> & _hill_tensor; 122 : 123 : GenericDenseVector<is_ad> _stress_np1; 124 : /// 2 * shear modulus 125 : GenericReal<is_ad> _two_shear_modulus; 126 : }; 127 : 128 : typedef HillPlasticityStressUpdateTempl<false> HillPlasticityStressUpdate; 129 : typedef HillPlasticityStressUpdateTempl<true> ADHillPlasticityStressUpdate;