www.mooseframework.org
RadialReturnStressUpdate.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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"
15 
29 template <bool is_ad>
32 {
33 public:
35 
37 
39  using Material::_dt;
40  using Material::_q_point;
41  using Material::_qp;
42 
43  enum class SubsteppingType
44  {
45  NONE,
48  };
49 
64  virtual void updateState(
65  GenericRankTwoTensor<is_ad> & strain_increment,
66  GenericRankTwoTensor<is_ad> & inelastic_strain_increment,
67  const GenericRankTwoTensor<is_ad> & rotation_increment,
68  GenericRankTwoTensor<is_ad> & stress_new,
69  const RankTwoTensor & stress_old,
70  const GenericRankFourTensor<is_ad> & elasticity_tensor,
71  const RankTwoTensor & elastic_strain_old,
72  bool compute_full_tangent_operator = false,
74 
75  virtual void updateStateSubstepInternal(
76  GenericRankTwoTensor<is_ad> & /*strain_increment*/,
77  GenericRankTwoTensor<is_ad> & /*inelastic_strain_increment*/,
78  const GenericRankTwoTensor<is_ad> & /*rotation_increment*/,
79  GenericRankTwoTensor<is_ad> & /*stress_new*/,
80  const RankTwoTensor & /*stress_old*/,
81  const GenericRankFourTensor<is_ad> & /*elasticity_tensor*/,
82  const RankTwoTensor & /*elastic_strain_old*/,
83  unsigned int total_number_substeps,
84  bool compute_full_tangent_operator = false,
86 
90  virtual void updateStateSubstep(
91  GenericRankTwoTensor<is_ad> & /*strain_increment*/,
92  GenericRankTwoTensor<is_ad> & /*inelastic_strain_increment*/,
93  const GenericRankTwoTensor<is_ad> & /*rotation_increment*/,
94  GenericRankTwoTensor<is_ad> & /*stress_new*/,
95  const RankTwoTensor & /*stress_old*/,
96  const GenericRankFourTensor<is_ad> & /*elasticity_tensor*/,
97  const RankTwoTensor & /*elastic_strain_old*/,
98  bool compute_full_tangent_operator = false,
100 
101  virtual Real
102  computeReferenceResidual(const GenericReal<is_ad> & effective_trial_stress,
103  const GenericReal<is_ad> & scalar_effective_inelastic_strain) override;
104 
105  virtual GenericReal<is_ad>
106  minimumPermissibleValue(const GenericReal<is_ad> & /*effective_trial_stress*/) const override
107  {
108  return 0.0;
109  }
110 
111  virtual GenericReal<is_ad>
112  maximumPermissibleValue(const GenericReal<is_ad> & effective_trial_stress) const override;
113 
118  virtual Real computeTimeStepLimit() override;
119 
123  bool requiresIsotropicTensor() override { return true; }
124 
128  bool isIsotropic() override { return true; };
129 
137  virtual int
138  calculateNumberSubsteps(const GenericRankTwoTensor<is_ad> & strain_increment) override;
139 
145  virtual bool substeppingCapabilityRequested() override
146  {
148  }
149 
152  {
154  }
155 
157  {
159  }
160 
162  {
164  }
165 
169  void computeTangentOperator(Real effective_trial_stress,
170  const RankTwoTensor & stress_new,
171  RankFourTensor & tangent_operator);
172 
173 protected:
174  virtual void initQpStatefulProperties() override;
175 
185 
191  virtual void computeStressInitialize(const GenericReal<is_ad> & effective_trial_stress,
192  const GenericRankFourTensor<is_ad> & elasticity_tensor);
193 
199  virtual Real computeStressDerivative(const Real /*effective_trial_stress*/, const Real /*scalar*/)
200  {
201  return 0.0;
202  }
203 
208  virtual void
209  computeStressFinalize(const GenericRankTwoTensor<is_ad> & /*inelasticStrainIncrement*/)
210  {
211  }
212 
213  void outputIterationSummary(std::stringstream * iter_output,
214  const unsigned int total_it) override;
215 
218 
221 
224 
230 
237 
242 
247 
252 
254  const bool _apply_strain;
255 
258 
261 
263  const unsigned int _maximum_number_substeps;
264 
267 };
268 
271 
272 template <>
273 void
275  const RankTwoTensor & stress_new,
276  RankFourTensor & tangent_operator);
const MooseArray< Point > & _q_point
SubsteppingType _use_substepping
Whether user has requested the use of substepping technique to improve convergence [make const later]...
virtual void computeStressFinalize(const GenericRankTwoTensor< is_ad > &)
Perform any necessary steps to finalize state after return mapping iterations.
typename Moose::GenericType< RankFourTensor, is_ad > GenericRankFourTensor
void outputIterationSummary(std::stringstream *iter_output, const unsigned int total_it) override
Output summary information for the convergence history of the model.
constexpr auto increment(std::index_sequence< first, tail... >)
virtual void updateState(GenericRankTwoTensor< is_ad > &strain_increment, GenericRankTwoTensor< is_ad > &inelastic_strain_increment, const GenericRankTwoTensor< is_ad > &rotation_increment, GenericRankTwoTensor< is_ad > &stress_new, const RankTwoTensor &stress_old, const GenericRankFourTensor< is_ad > &elasticity_tensor, const RankTwoTensor &elastic_strain_old, bool compute_full_tangent_operator=false, RankFourTensor &tangent_operator=StressUpdateBaseTempl< is_ad >::_identityTensor) override
A radial return (J2) mapping method is performed with return mapping iterations.
static InputParameters validParams()
const RankFourTensor _deviatoric_projection_four
Rank four deviatoric projection tensor.
virtual GenericReal< is_ad > minimumPermissibleValue(const GenericReal< is_ad > &) const override
Compute the minimum permissible value of the scalar.
bool requiresIsotropicTensor() override
Does the model require the elasticity tensor to be isotropic?
Real _dt_original
original timestep (to be restored after substepping is completed)
GenericMaterialProperty< Real, is_ad > & _effective_inelastic_strain
virtual void computeStressInitialize(const GenericReal< is_ad > &effective_trial_stress, const GenericRankFourTensor< is_ad > &elasticity_tensor)
Perform any necessary initialization before return mapping iterations.
bool isIsotropic() override
Radial return mapped models should be isotropic by default!
RadialReturnStressUpdate computes the radial return stress increment for an isotropic elastic-viscopl...
virtual Real computeStressDerivative(const Real, const Real)
Calculate the derivative of the strain increment with respect to the updated stress.
virtual int calculateNumberSubsteps(const GenericRankTwoTensor< is_ad > &strain_increment) override
If substepping is enabled, calculate the number of substeps as a function of the elastic strain incre...
Real & _dt
GenericReal< is_ad > _three_shear_modulus
3 * shear modulus
void updateEffectiveInelasticStrain(const GenericReal< is_ad > &increment)
virtual void updateStateSubstepInternal(GenericRankTwoTensor< is_ad > &, GenericRankTwoTensor< is_ad > &, const GenericRankTwoTensor< is_ad > &, GenericRankTwoTensor< is_ad > &, const RankTwoTensor &, const GenericRankFourTensor< is_ad > &, const RankTwoTensor &, unsigned int total_number_substeps, bool compute_full_tangent_operator=false, RankFourTensor &tangent_operator=StressUpdateBaseTempl< is_ad >::_identityTensor)
unsigned int _qp
virtual void initQpStatefulProperties() override
const GenericReal< is_ad > & effectiveInelasticStrainIncrement() const
Current value of scalar inelastic strain.
typename GenericMaterialPropertyStruct< T, is_ad >::type GenericMaterialProperty
const MaterialProperty< Real > & _effective_inelastic_strain_old
Base class that provides capability for Newton return mapping iterations on a single variable...
virtual Real computeTimeStepLimit() override
Compute the limiting value of the time step for this material.
StressUpdateBase is a material that is not called by MOOSE because of the compute=false flag set in t...
const bool _apply_strain
Debugging option to enable specifying instead of calculating strain.
virtual bool substeppingCapabilityRequested() override
Has the user requested usage of (possibly) implemented substepping capability for inelastic models...
void computeTangentOperator(Real effective_trial_stress, const RankTwoTensor &stress_new, RankFourTensor &tangent_operator)
Calculate the tangent_operator.
Real _max_inelastic_increment
Maximum allowable scalar inelastic strain increment, used to control the timestep size in conjunction...
virtual GenericReal< is_ad > maximumPermissibleValue(const GenericReal< is_ad > &effective_trial_stress) const override
Compute the maximum permissible value of the scalar.
RadialReturnStressUpdateTempl(const InputParameters &parameters)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void propagateQpStatefulPropertiesRadialReturn()
Propagate the properties pertaining to this intermediate class.
GenericReal< is_ad > _effective_inelastic_strain_increment
Stores the scalar effective inelastic strain increment from Newton iteration.
const unsigned int _maximum_number_substeps
Maximum number of substeps. If the calculation results in a larger number, cut overall time step...
typename Moose::GenericType< RankTwoTensor, is_ad > GenericRankTwoTensor
const InputParameters & parameters() const
RadialReturnStressUpdateTempl< true > ADRadialReturnStressUpdate
typename Moose::GenericType< Real, is_ad > GenericReal
void updateEffectiveInelasticStrainIncrement(const GenericReal< is_ad > &eisi)
const RankFourTensor _identity_symmetric_four
Rank four symmetric identity tensor.
virtual void updateStateSubstep(GenericRankTwoTensor< is_ad > &, GenericRankTwoTensor< is_ad > &, const GenericRankTwoTensor< is_ad > &, GenericRankTwoTensor< is_ad > &, const RankTwoTensor &, const GenericRankFourTensor< is_ad > &, const RankTwoTensor &, bool compute_full_tangent_operator=false, RankFourTensor &tangent_operator=StressUpdateBaseTempl< is_ad >::_identityTensor) override
Similar to the updateState function, this method updates the strain and stress for one substep...
const bool _adaptive_substepping
Use adaptive substepping, cutting substep sizes until convergence is achieved.
RadialReturnStressUpdateTempl< false > RadialReturnStressUpdate
const Elem *const & _current_elem
const RankTwoTensor _identity_two
Rank two identity tensor.
virtual Real computeReferenceResidual(const GenericReal< is_ad > &effective_trial_stress, const GenericReal< is_ad > &scalar_effective_inelastic_strain) override
Compute a reference quantity to be used for checking relative convergence.
const Real _substep_tolerance
Used to calculate the number of substeps taken in the radial return algorithm, when substepping is en...