https://mooseframework.inl.gov
HillCreepStressUpdate.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 
13 
29 template <bool is_ad>
31 {
32 public:
34 
35  HillCreepStressUpdateTempl(const InputParameters & parameters);
36 
37 protected:
40  using Material::_q_point;
41  using Material::_qp;
42 
43  virtual void
45  const GenericDenseVector<is_ad> & stress,
46  const GenericRankFourTensor<is_ad> & elasticity_tensor) override;
47  virtual GenericReal<is_ad>
48  computeResidual(const GenericDenseVector<is_ad> & effective_trial_stress,
49  const GenericDenseVector<is_ad> & stress_new,
50  const GenericReal<is_ad> & scalar) override;
51 
52  virtual GenericReal<is_ad>
53  computeDerivative(const GenericDenseVector<is_ad> & effective_trial_stress,
54  const GenericDenseVector<is_ad> & stress_new,
55  const GenericReal<is_ad> & scalar) override;
56 
57  virtual Real
58  computeReferenceResidual(const GenericDenseVector<is_ad> & effective_trial_stress,
59  const GenericDenseVector<is_ad> & stress_new,
60  const GenericReal<is_ad> & residual,
61  const GenericReal<is_ad> & scalar_effective_inelastic_strain) override;
62 
63  void computeQsigmaChanged(GenericReal<is_ad> & qsigma_changed,
64  const GenericDenseVector<is_ad> & stress_new,
65  const GenericReal<is_ad> & delta_gamma,
66  GenericRankTwoTensor<is_ad> & stress_changed);
67 
75  virtual void computeStrainFinalize(GenericRankTwoTensor<is_ad> & inelasticStrainIncrement,
76  const GenericRankTwoTensor<is_ad> & stress,
77  const GenericDenseVector<is_ad> & stress_dev,
78  const GenericReal<is_ad> & delta_gamma) override;
79 
87  virtual void
88  computeStressFinalize(const GenericRankTwoTensor<is_ad> & inelasticStrainIncrement,
89  const GenericReal<is_ad> & delta_gamma,
91  const GenericDenseVector<is_ad> & stress_dev,
92  const GenericRankTwoTensor<is_ad> & stress_old,
93  const GenericRankFourTensor<is_ad> & elasticity_tensor) override;
94 
95  virtual GenericReal<is_ad>
96  initialGuess(const GenericDenseVector<is_ad> & /*stress_dev*/) override;
97 
102  bool requiresIsotropicTensor() override { return false; }
103 
110  {
111  return this->_max_integration_error_time_step;
112  }
113 
115  const bool _has_temp;
116 
119 
122 
125 
128 
131 
134 
137 
140 
143 
146 
150 
153 
156 
158  const std::string _elasticity_tensor_name;
159 
162 
165 
168 };
169 
const MooseArray< Point > & _q_point
virtual GenericReal< is_ad > computeDerivative(const GenericDenseVector< is_ad > &effective_trial_stress, const GenericDenseVector< is_ad > &stress_new, const GenericReal< is_ad > &scalar) override
This class uses the stress update material for an anisotropic creep model.
virtual void computeStressFinalize(const GenericRankTwoTensor< is_ad > &inelasticStrainIncrement, const GenericReal< is_ad > &delta_gamma, GenericRankTwoTensor< is_ad > &stress, const GenericDenseVector< is_ad > &stress_dev, const GenericRankTwoTensor< is_ad > &stress_old, const GenericRankFourTensor< is_ad > &elasticity_tensor) override
Perform any necessary steps to finalize state after return mapping iterations.
Moose::GenericType< Real, is_ad > GenericReal
const Real _start_time
Simulation start time.
virtual Real computeIntegrationErrorTimeStep() override
Compute the limiting value of the time step for this material according to the numerical integration ...
const Real _m_exponent
Exponent on time.
bool requiresIsotropicTensor() override
Does the model require the elasticity tensor to be isotropic? Not in principle.
virtual GenericReal< is_ad > initialGuess(const GenericDenseVector< is_ad > &) override
bool _anisotropic_elasticity
Materials&#39;s elasticity tensor is anisotropic or not.
Moose::GenericType< DenseVector< Real >, is_ad > GenericDenseVector
Moose::GenericType< RankFourTensor, is_ad > GenericRankFourTensor
const bool _has_temp
Flag to determine if temperature is supplied by the user.
Real _exp_time
Exponential calculated from current time.
unsigned int _qp
const Real _n_exponent
Exponent on the effective stress.
const MaterialProperty< DenseMatrix< Real > > * _hill_tensor
Hill tensor, when global axes do not (somehow) align with those of the material Example: Large rotati...
typename GenericMaterialPropertyStruct< T, is_ad >::type GenericMaterialProperty
const VariableValue & _temperature
Temperature variable value.
const MaterialProperty< std::vector< Real > > & _hill_constants
Hill constant material.
static InputParameters validParams()
const Function *const _prefactor_function
Prefactor for scaling creep.
virtual void computeStressInitialize(const GenericDenseVector< is_ad > &stress_dev, const GenericDenseVector< is_ad > &stress, const GenericRankFourTensor< is_ad > &elasticity_tensor) override
virtual Real computeReferenceResidual(const GenericDenseVector< is_ad > &effective_trial_stress, const GenericDenseVector< is_ad > &stress_new, const GenericReal< is_ad > &residual, const GenericReal< is_ad > &scalar_effective_inelastic_strain) override
const std::string _elasticity_tensor_name
Name of the elasticity tensor material property.
GenericReal< is_ad > _two_shear_modulus
2 * shear modulus
OutputTools< Real >::VariableValue VariableValue
void computeQsigmaChanged(GenericReal< is_ad > &qsigma_changed, const GenericDenseVector< is_ad > &stress_new, const GenericReal< is_ad > &delta_gamma, GenericRankTwoTensor< is_ad > &stress_changed)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
HillCreepStressUpdateTempl< true > ADHillCreepStressUpdate
virtual GenericReal< is_ad > computeResidual(const GenericDenseVector< is_ad > &effective_trial_stress, const GenericDenseVector< is_ad > &stress_new, const GenericReal< is_ad > &scalar) override
virtual void computeStrainFinalize(GenericRankTwoTensor< is_ad > &inelasticStrainIncrement, const GenericRankTwoTensor< is_ad > &stress, const GenericDenseVector< is_ad > &stress_dev, const GenericReal< is_ad > &delta_gamma) override
Perform any necessary steps to finalize strain increment after return mapping iterations.
HillCreepStressUpdateTempl(const InputParameters &parameters)
const Real _gas_constant
Gas constant for exp term.
Real _exponential
Exponential calculated from activation, gas constant, and temperature.
GenericDenseMatrix< is_ad > _C
Matrix form of the elasticity tensor.
const Real _activation_energy
Activation energy for exp term.
HillCreepStressUpdateTempl< false > HillCreepStressUpdate
This class provides baseline functionality for anisotropic (Hill-like) plasticity and creep models ba...
Moose::GenericType< DenseMatrix< Real >, is_ad > GenericDenseMatrix
const GenericMaterialProperty< RankFourTensor, is_ad > & _elasticity_tensor
Anisotropic elasticity tensor material property.
const Elem *const & _current_elem
const Real _coefficient
Leading coefficient.
Moose::GenericType< RankTwoTensor, is_ad > GenericRankTwoTensor