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 "StressUpdateBase.h" 13 : #include "GeneralizedReturnMappingSolution.h" 14 : #include "MooseTypes.h" 15 : 16 : #include "libmesh/ignore_warnings.h" 17 : #include "Eigen/Dense" 18 : #include "Eigen/Eigenvalues" 19 : #include "libmesh/restore_warnings.h" 20 : 21 : #include <limits> 22 : 23 : namespace Eigen 24 : { 25 : namespace internal 26 : { 27 : template <> 28 : struct cast_impl<ADReal, int> 29 : { 30 : static inline int run(const ADReal & x) { return static_cast<int>(MetaPhysicL::raw_value(x)); } 31 : }; 32 : } // namespace internal 33 : } // namespace Eigen 34 : 35 : // Instantiation of ADReal will be enabled shortly, once 36 : // https://github.com/roystgnr/MetaPhysicL/issues/70 makes it to MOOSE 37 : // typedef Eigen::Matrix<ADReal, 6, 6> AnisotropyMatrix; 38 : 39 : typedef Eigen::Matrix<Real, 6, 6> AnisotropyMatrixReal; 40 : typedef Eigen::Matrix<Real, 3, 3> AnisotropyMatrixRealBlock; 41 : 42 : /** 43 : * ADGeneralizedRadialReturnStressUpdate computes the generalized radial return stress increment for 44 : * anisotropic (Hill-like) creep and plasticity. This generalized radial return mapping class 45 : * acts as a base class for the radial return of anisotropic creep and plasticity 46 : * classes/combinations. The stress increment computed by ADGeneralizedRadialReturnStressUpdate is 47 : * used by ComputeMultipleInelasticStress which computes the elastic stress for finite strains. This 48 : * class is based on the algorithm in Versino, D, Bennett, KC. Generalized radial return mapping 49 : * algorithm for anisotropic von Mises plasticity framed in material eigenspace. Int J Numer Methods 50 : * Eng. 2018. 116. 202 222 51 : */ 52 : template <bool is_ad> 53 : class GeneralizedRadialReturnStressUpdateTempl : public StressUpdateBaseTempl<is_ad>, 54 : public GeneralizedReturnMappingSolutionTempl<is_ad> 55 : { 56 : public: 57 : static InputParameters validParams(); 58 : 59 : GeneralizedRadialReturnStressUpdateTempl(const InputParameters & parameters); 60 : 61 : /** 62 : * A radial return (J2) mapping method is performed with return mapping 63 : * iterations. 64 : * @param strain_increment Sum of elastic and inelastic strain increments 65 : * @param inelastic_strain_increment Inelastic strain increment calculated by this class 66 : * @param rotation increment Not used by this class 67 : * @param stress_new New trial stress from pure elastic calculation 68 : * @param stress_old Old state of stress 69 : * @param elasticity_tensor Rank 4 C_{ijkl}, must be isotropic 70 : * @param elastic_strain_old Old state of total elastic strain 71 : */ 72 : virtual void 73 : updateState(GenericRankTwoTensor<is_ad> & strain_increment, 74 : GenericRankTwoTensor<is_ad> & inelastic_strain_increment, 75 : const GenericRankTwoTensor<is_ad> & rotation_increment, 76 : GenericRankTwoTensor<is_ad> & stress_new, 77 : const RankTwoTensor & stress_old, 78 : const GenericRankFourTensor<is_ad> & elasticity_tensor, 79 : const RankTwoTensor & elastic_strain_old, 80 : bool /*compute_full_tangent_operator = false*/, 81 : RankFourTensor & /*tangent_operator = StressUpdateBaseTempl<is_ad>::_identityTensor*/) 82 : override; 83 : 84 : virtual Real 85 : computeReferenceResidual(const GenericDenseVector<is_ad> & effective_trial_stress, 86 : const GenericDenseVector<is_ad> & stress_new, 87 : const GenericReal<is_ad> & residual, 88 : const GenericReal<is_ad> & scalar_effective_inelastic_strain) override; 89 : 90 1888576 : virtual GenericReal<is_ad> minimumPermissibleValue( 91 : const GenericDenseVector<is_ad> & /*effective_trial_stress*/) const override 92 : { 93 1888576 : return 0.0; 94 : } 95 : 96 : virtual GenericReal<is_ad> 97 : maximumPermissibleValue(const GenericDenseVector<is_ad> & effective_trial_stress) const override; 98 : 99 : /** 100 : * Compute the limiting value of the time step for this material 101 : * @return Limiting time step 102 : */ 103 : virtual Real computeTimeStepLimit() override; 104 : 105 : /** 106 : * Compute the limiting value of the time step for this material according to the numerical 107 : * integration error 108 : * @return Limiting time step 109 : */ 110 28880 : virtual Real computeIntegrationErrorTimeStep() { return std::numeric_limits<Real>::max(); } 111 : 112 : /** 113 : * Does the model require the elasticity tensor to be isotropic? 114 : */ 115 0 : bool requiresIsotropicTensor() override { return true; } 116 : 117 : /** 118 : * Check if an anisotropic matrix is block diagonal 119 : */ 120 : bool isBlockDiagonal(const AnisotropyMatrixReal & A); 121 : 122 : protected: 123 : using Material::_current_elem; 124 : using Material::_dt; 125 : using Material::_q_point; 126 : using Material::_qp; 127 : 128 : virtual void initQpStatefulProperties() override; 129 : 130 : /** 131 : * Propagate the properties pertaining to this intermediate class. This 132 : * is intended to be called from propagateQpStatefulProperties() in 133 : * classes that inherit from this one. 134 : * This is intentionally named uniquely because almost all models that derive 135 : * from this class have their own stateful properties, and this forces them 136 : * to define their own implementations of propagateQpStatefulProperties(). 137 : */ 138 : void propagateQpStatefulPropertiesRadialReturn(); 139 : 140 : /** 141 : * Perform any necessary initialization before return mapping iterations 142 : * @param effective_trial_stress Effective trial stress 143 : * @param elasticityTensor Elasticity tensor 144 : */ 145 : virtual void 146 : computeStressInitialize(const GenericDenseVector<is_ad> & /*stress_dev*/, 147 : const GenericDenseVector<is_ad> & /*stress*/, 148 : const GenericRankFourTensor<is_ad> & /*elasticity_tensor*/) = 0; 149 : 150 : /** 151 : * Calculate the derivative of the strain increment with respect to the updated stress. 152 : * @param effective_trial_stress Effective trial stress 153 : * @param scalar Inelastic strain increment magnitude being solved for 154 : */ 155 0 : virtual GenericReal<is_ad> computeStressDerivative(const Real /*effective_trial_stress*/, 156 : const Real /*scalar*/) 157 : { 158 0 : return 0.0; 159 : } 160 : 161 : /** 162 : * Perform any necessary steps to finalize state after return mapping iterations 163 : * @param inelasticStrainIncrement Inelastic strain increment 164 : * @param delta_gamma Plastic multiplier 165 : * @param stress Cauchy stress 166 : * @param stress_dev Deviatoric part of Cauchy stress 167 : */ 168 : virtual void computeStressFinalize(const GenericRankTwoTensor<is_ad> & inelasticStrainIncrement, 169 : const GenericReal<is_ad> & delta_gamma, 170 : GenericRankTwoTensor<is_ad> & stress, 171 : const GenericDenseVector<is_ad> & /*stress_dev*/, 172 : const GenericRankTwoTensor<is_ad> & stress_old, 173 : const GenericRankFourTensor<is_ad> & elasticity_tensor) = 0; 174 : /** 175 : * Perform any necessary steps to finalize strain increment after return mapping iterations 176 : * @param inelasticStrainIncrement Inelastic strain increment 177 : * @param stress Cauchy stress 178 : * @param stress_dev Deviatoric part of Cauchy stress 179 : * @param delta_gamma Plastic multiplier 180 : */ 181 : virtual void computeStrainFinalize(GenericRankTwoTensor<is_ad> & inelasticStrainIncrement, 182 : const GenericRankTwoTensor<is_ad> & stress, 183 : const GenericDenseVector<is_ad> & stress_dev, 184 : const GenericReal<is_ad> & delta_gamma) = 0; 185 : 186 : /* 187 : * Method that determines the tangent calculation method. For anisotropic creep only models, the 188 : * tangent calculation method is ELASTIC for now 189 : */ 190 96 : virtual TangentCalculationMethod getTangentCalculationMethod() override 191 : { 192 96 : return TangentCalculationMethod::ELASTIC; 193 : } 194 : 195 : void outputIterationSummary(std::stringstream * iter_output, 196 : const unsigned int total_it) override; 197 : 198 : /** 199 : * Calculate the tangent_operator. 200 : */ 201 : void computeTangentOperator(Real /*effective_trial_stress*/, 202 : RankTwoTensor & /*stress_new*/, 203 : bool /*compute_full_tangent_operator*/, 204 : RankFourTensor & /*tangent_operator*/); 205 : 206 : /// 3 * shear modulus 207 : GenericReal<is_ad> _three_shear_modulus; 208 : 209 : /// Equivalent creep/plastic strain 210 : GenericMaterialProperty<Real, is_ad> & _effective_inelastic_strain; 211 : const MaterialProperty<Real> & _effective_inelastic_strain_old; 212 : 213 : /// Equivalent creep/plastic strain rate: facilitates user's control of integration errors 214 : MaterialProperty<Real> & _inelastic_strain_rate; 215 : const MaterialProperty<Real> & _inelastic_strain_rate_old; 216 : 217 : /// Maximum inelastic strain increment for (next) time step prescription 218 : Real _max_inelastic_increment; 219 : 220 : /// Maximum integration error for creep 221 : Real _max_integration_error; 222 : 223 : /// Maximum integration error time step 224 : Real _max_integration_error_time_step; 225 : 226 : /// Whether to use fully transformed Hill's tensor due to rigid body or large deformation kinematic rotation 227 : const bool _use_transformation; 228 : }; 229 : 230 : typedef GeneralizedRadialReturnStressUpdateTempl<false> GeneralizedRadialReturnStressUpdate; 231 : typedef GeneralizedRadialReturnStressUpdateTempl<true> ADGeneralizedRadialReturnStressUpdate;