https://mooseframework.inl.gov
HillElastoPlasticityStressUpdate.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 
21 template <bool is_ad>
24 {
25 public:
27 
29 
32  const GenericMaterialProperty<RankTwoTensor, is_ad> & strain_rate) override;
33 
34 protected:
42  using Material::_q_point;
43  using Material::_qp;
45 
46  virtual void initQpStatefulProperties() override;
47 
48  virtual void
50  const GenericDenseVector<is_ad> & stress,
51  const GenericRankFourTensor<is_ad> & elasticity_tensor) override;
52 
57 
58  virtual GenericReal<is_ad>
59  computeResidual(const GenericDenseVector<is_ad> & effective_trial_stress,
60  const GenericDenseVector<is_ad> & stress_new,
61  const GenericReal<is_ad> & scalar) override;
62  virtual GenericReal<is_ad>
63  computeDerivative(const GenericDenseVector<is_ad> & effective_trial_stress,
64  const GenericDenseVector<is_ad> & stress_new,
65  const GenericReal<is_ad> & scalar) override;
66 
67  virtual Real
68  computeReferenceResidual(const GenericDenseVector<is_ad> & effective_trial_stress,
69  const GenericDenseVector<is_ad> & stress_new,
70  const GenericReal<is_ad> & residual,
71  const GenericReal<is_ad> & scalar_effective_inelastic_strain) override;
72  virtual void propagateQpStatefulProperties() override;
77  bool requiresIsotropicTensor() override { return false; }
78 
81  const GenericReal<is_ad> & omega);
86  void computeHillTensorEigenDecomposition(const DenseMatrix<GenericReal<is_ad>> & hill_tensor);
87 
95  virtual void computeStrainFinalize(GenericRankTwoTensor<is_ad> & inelasticStrainIncrement,
96  const GenericRankTwoTensor<is_ad> & stress,
97  const GenericDenseVector<is_ad> & stress_dev,
98  const GenericReal<is_ad> & delta_gamma) override;
99 
107  virtual void
108  computeStressFinalize(const GenericRankTwoTensor<is_ad> & inelasticStrainIncrement,
109  const GenericReal<is_ad> & delta_gamma,
111  const GenericDenseVector<is_ad> & stress_dev,
112  const GenericRankTwoTensor<is_ad> & stress_old,
113  const GenericRankFourTensor<is_ad> & elasticity_tensor) override;
114 
116  const GenericDenseVector<is_ad> & stress_trial);
117 
118  void computeDeltaDerivatives(const GenericReal<is_ad> & delta_gamma,
119  const GenericDenseVector<is_ad> & stress_trial,
120  const GenericReal<is_ad> & sy_alpha,
121  GenericReal<is_ad> & omega,
122  GenericReal<is_ad> & omega_gamma,
123  GenericReal<is_ad> & sy_gamma);
124 
127 
131 
133  const std::string _elasticity_tensor_name;
136 
139 
142 
145 
148 
152 
156 
160 
166  const enum class Axis { X, Y, Z } _axis;
167 };
168 
const MooseArray< Point > & _q_point
GenericReal< is_ad > _qsigma
Square of the q function for orthotropy.
Moose::GenericType< Real, is_ad > GenericReal
GenericMaterialProperty< Real, is_ad > & _hardening_variable
const MaterialProperty< RankTwoTensor > & _rotation_matrix_transpose_old
void computeElasticityTensorEigenDecomposition()
Compute eigendecomposition of element-wise elasticity tensor.
const MaterialProperty< Real > & _hardening_variable_old
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _elasticity_eigenvectors
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.
enum HillElastoPlasticityStressUpdateTempl::Axis _axis
MaterialProperty< RankTwoTensor > & _rotation_matrix_transpose
HillElastoPlasticityStressUpdateTempl< true > ADHillElastoPlasticityStressUpdate
This class uses the stress update material in an anisotropic return mapping.
void computeHillTensorEigenDecomposition(const DenseMatrix< GenericReal< is_ad >> &hill_tensor)
Compute eigendecomposition of Hill&#39;s tensor for anisotropic plasticity.
GenericMaterialProperty< DenseVector< Real >, is_ad > & _sigma_tilde_rotated
Moose::GenericType< DenseVector< Real >, is_ad > GenericDenseVector
Moose::GenericType< RankFourTensor, is_ad > GenericRankFourTensor
unsigned int _qp
const std::string _elasticity_tensor_name
Name of the elasticity tensor material property.
virtual GenericReal< is_ad > computeResidual(const GenericDenseVector< is_ad > &effective_trial_stress, const GenericDenseVector< is_ad > &stress_new, const GenericReal< is_ad > &scalar) override
GenericReal< is_ad > computeHardeningValue(const GenericReal< is_ad > &scalar, const GenericReal< is_ad > &omega)
HillElastoPlasticityStressUpdateTempl< false > HillElastoPlasticityStressUpdate
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
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _b_eigenvectors
HillElastoPlasticityStressUpdateTempl(const InputParameters &parameters)
GenericMaterialProperty< DenseVector< Real >, is_ad > & _sigma_tilde
GenericMaterialProperty< DenseVector< Real >, is_ad > & _elasticity_eigenvalues
typename GenericMaterialPropertyStruct< T, is_ad >::type GenericMaterialProperty
virtual GenericReal< is_ad > computeDerivative(const GenericDenseVector< is_ad > &effective_trial_stress, const GenericDenseVector< is_ad > &stress_new, const GenericReal< is_ad > &scalar) override
void paramError(const std::string &param, Args... args) const
const MaterialProperty< DenseMatrix< Real > > & _hill_tensor
Hill tensor, when global axes do not (somehow) align with those of the material Example: Large rotati...
GenericReal< is_ad > computeOmega(const GenericReal< is_ad > &delta_gamma, const GenericDenseVector< is_ad > &stress_trial)
GenericMaterialProperty< DenseVector< Real >, is_ad > & _b_eigenvalues
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _alpha_matrix
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.
virtual void computeStressInitialize(const GenericDenseVector< is_ad > &stress_dev, const GenericDenseVector< is_ad > &stress, const GenericRankFourTensor< is_ad > &elasticity_tensor) override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
MaterialProperty< RankTwoTensor > & _rotation_matrix
const MaterialProperty< RankTwoTensor > & _rotation_matrix_old
bool requiresIsotropicTensor() override
Does the model require the elasticity tensor to be isotropic? No, this class does anisotropic elasto-...
This class provides baseline functionality for anisotropic (Hill-like) plasticity models based on the...
const GenericMaterialProperty< RankFourTensor, is_ad > & _elasticity_tensor
Anisotropic elasticity tensor material property.
void computeDeltaDerivatives(const GenericReal< is_ad > &delta_gamma, const GenericDenseVector< is_ad > &stress_trial, const GenericReal< is_ad > &sy_alpha, GenericReal< is_ad > &omega, GenericReal< is_ad > &omega_gamma, GenericReal< is_ad > &sy_gamma)
Moose::GenericType< DenseMatrix< Real >, is_ad > GenericDenseMatrix
const Elem *const & _current_elem
Moose::GenericType< RankTwoTensor, is_ad > GenericRankTwoTensor
virtual Real computeStrainEnergyRateDensity(const GenericMaterialProperty< RankTwoTensor, is_ad > &stress, const GenericMaterialProperty< RankTwoTensor, is_ad > &strain_rate) override