https://mooseframework.inl.gov
Public Member Functions | Static Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | List of all members
HillElastoPlasticityStressUpdateTempl< is_ad > Class Template Reference

This class uses the stress update material in an anisotropic return mapping. More...

#include <HillElastoPlasticityStressUpdate.h>

Inheritance diagram for HillElastoPlasticityStressUpdateTempl< is_ad >:
[legend]

Public Member Functions

 HillElastoPlasticityStressUpdateTempl (const InputParameters &parameters)
 
virtual Real computeStrainEnergyRateDensity (const GenericMaterialProperty< RankTwoTensor, is_ad > &stress, const GenericMaterialProperty< RankTwoTensor, is_ad > &strain_rate) override
 

Static Public Member Functions

static InputParameters validParams ()
 

Protected Types

enum  Axis { Axis::X, Axis::Y, Axis::Z }
 

Protected Member Functions

virtual void initQpStatefulProperties () override
 
virtual void computeStressInitialize (const GenericDenseVector< is_ad > &stress_dev, const GenericDenseVector< is_ad > &stress, const GenericRankFourTensor< is_ad > &elasticity_tensor) override
 
void computeElasticityTensorEigenDecomposition ()
 Compute eigendecomposition of element-wise elasticity tensor. More...
 
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 GenericReal< is_ad > computeDerivative (const GenericDenseVector< is_ad > &effective_trial_stress, const GenericDenseVector< is_ad > &stress_new, const GenericReal< is_ad > &scalar) 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
 
virtual void propagateQpStatefulProperties () override
 
bool requiresIsotropicTensor () override
 Does the model require the elasticity tensor to be isotropic? No, this class does anisotropic elasto-plasticity More...
 
Real computeHardeningDerivative ()
 
GenericReal< is_ad > computeHardeningValue (const GenericReal< is_ad > &scalar, const GenericReal< is_ad > &omega)
 
void computeHillTensorEigenDecomposition (const DenseMatrix< GenericReal< is_ad >> &hill_tensor)
 Compute eigendecomposition of Hill's tensor for anisotropic plasticity. More...
 
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. More...
 
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. More...
 
GenericReal< is_ad > computeOmega (const GenericReal< is_ad > &delta_gamma, const GenericDenseVector< is_ad > &stress_trial)
 
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)
 
void paramError (const std::string &param, Args... args) const
 
virtual GenericReal< is_ad > computeStressDerivative (const Real, const Real) override
 Calculate the derivative of the strain increment with respect to the updated stress. More...
 

Protected Attributes

 usingTransientInterfaceMembers
 
GenericReal< is_ad > _qsigma
 Square of the q function for orthotropy. More...
 
GenericDenseVector< is_ad > _eigenvalues_hill
 
GenericDenseMatrix< is_ad > _eigenvectors_hill
 
GenericDenseMatrix< is_ad > _anisotropic_elastic_tensor
 
const std::string _elasticity_tensor_name
 Name of the elasticity tensor material property. More...
 
const GenericMaterialProperty< RankFourTensor, is_ad > & _elasticity_tensor
 Anisotropic elasticity tensor material property. More...
 
const Real _hardening_constant
 
const Real _hardening_exponent
 
GenericMaterialProperty< Real, is_ad > & _hardening_variable
 
const MaterialProperty< Real > & _hardening_variable_old
 
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _elasticity_eigenvectors
 
GenericMaterialProperty< DenseVector< Real >, is_ad > & _elasticity_eigenvalues
 
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _b_eigenvectors
 
GenericMaterialProperty< DenseVector< Real >, is_ad > & _b_eigenvalues
 
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _alpha_matrix
 
GenericMaterialProperty< DenseVector< Real >, is_ad > & _sigma_tilde
 
GenericMaterialProperty< DenseVector< Real >, is_ad > & _sigma_tilde_rotated
 
GenericReal< is_ad > _hardening_derivative
 
GenericReal< is_ad > _yield_condition
 
GenericReal< is_ad > _yield_stress
 
const MaterialProperty< DenseMatrix< Real > > & _hill_tensor
 Hill tensor, when global axes do not (somehow) align with those of the material Example: Large rotation due to rigid body and/or large deformation kinematics. More...
 
MaterialProperty< RankTwoTensor > & _rotation_matrix
 
MaterialProperty< RankTwoTensor > & _rotation_matrix_transpose
 
const MaterialProperty< RankTwoTensor > & _rotation_matrix_old
 
const MaterialProperty< RankTwoTensor > & _rotation_matrix_transpose_old
 
const bool _local_cylindrical_csys
 
enum HillElastoPlasticityStressUpdateTempl::Axis _axis
 
const Elem *const & _current_elem
 
const MooseArray< Point > & _q_point
 
unsigned int _qp
 
GenericMaterialProperty< RankTwoTensor, is_ad > & _plasticity_strain
 Plasticity strain tensor material property. More...
 
const MaterialProperty< RankTwoTensor > & _plasticity_strain_old
 

Detailed Description

template<bool is_ad>
class HillElastoPlasticityStressUpdateTempl< is_ad >

This class uses the stress update material in an anisotropic return mapping.

This class is one of the generalized radial return constitutive models based on Hill's criterion; it assumes and anisotropic elasticity tensor and an anisotropic plastic yield surface. Constitutive models that combine creep and plasticity can be used.

Definition at line 22 of file HillElastoPlasticityStressUpdate.h.

Member Enumeration Documentation

◆ Axis

template<bool is_ad>
enum HillElastoPlasticityStressUpdateTempl::Axis
strongprotected
Enumerator

Definition at line 166 of file HillElastoPlasticityStressUpdate.h.

166 { X, Y, Z } _axis;
enum HillElastoPlasticityStressUpdateTempl::Axis _axis
static const std::string Z
Definition: NS.h:169

Constructor & Destructor Documentation

◆ HillElastoPlasticityStressUpdateTempl()

Definition at line 47 of file HillElastoPlasticityStressUpdate.C.

50  _qsigma(0.0),
52  _eigenvectors_hill(6, 6),
54  _elasticity_tensor_name(this->_base_name + "elasticity_tensor"),
56  this->template getGenericMaterialProperty<RankFourTensor, is_ad>(_elasticity_tensor_name)),
57  _hardening_constant(this->template getParam<Real>("hardening_constant")),
58  _hardening_exponent(this->template getParam<Real>("hardening_exponent")),
59  _hardening_variable(this->template declareGenericProperty<Real, is_ad>(this->_base_name +
60  "hardening_variable")),
62  this->template getMaterialPropertyOld<Real>(this->_base_name + "hardening_variable")),
63  _elasticity_eigenvectors(this->template declareGenericProperty<DenseMatrix<Real>, is_ad>(
64  this->_base_name + "elasticity_eigenvectors")),
65  _elasticity_eigenvalues(this->template declareGenericProperty<DenseVector<Real>, is_ad>(
66  this->_base_name + "elasticity_eigenvalues")),
67  _b_eigenvectors(this->template declareGenericProperty<DenseMatrix<Real>, is_ad>(
68  this->_base_name + "b_eigenvectors")),
69  _b_eigenvalues(this->template declareGenericProperty<DenseVector<Real>, is_ad>(
70  this->_base_name + "b_eigenvalues")),
71  _alpha_matrix(this->template declareGenericProperty<DenseMatrix<Real>, is_ad>(this->_base_name +
72  "alpha_matrix")),
73  _sigma_tilde(this->template declareGenericProperty<DenseVector<Real>, is_ad>(this->_base_name +
74  "sigma_tilde")),
75  _sigma_tilde_rotated(this->template declareGenericProperty<DenseVector<Real>, is_ad>(
76  this->_base_name + "sigma_tilde_rotated")),
78  _yield_condition(1.0),
79  _yield_stress(this->template getParam<Real>("yield_stress")),
80  _hill_tensor(this->template getMaterialPropertyByName<DenseMatrix<Real>>(this->_base_name +
81  "hill_tensor")),
83  this->template declareProperty<RankTwoTensor>(this->_base_name + "rotation_matrix")),
84  _rotation_matrix_transpose(this->template declareProperty<RankTwoTensor>(
85  this->_base_name + "rotation_matrix_transpose")),
87  this->template getMaterialPropertyOld<RankTwoTensor>(this->_base_name + "rotation_matrix")),
88  _rotation_matrix_transpose_old(this->template getMaterialPropertyOld<RankTwoTensor>(
89  this->_base_name + "rotation_matrix_transpose")),
90  _local_cylindrical_csys(this->template getParam<bool>("local_cylindrical_csys")),
91  _axis(this->template getParam<MooseEnum>("axis").template getEnum<Axis>())
92 {
93  if (_local_cylindrical_csys && !parameters.isParamSetByUser("axis"))
94  {
95  this->paramError(
96  "local_cylindrical_csys",
97  "\nIf parameter local_cylindrical_csys is provided then parameter axis should be also "
98  "provided.");
99  }
100 
101  if (!_local_cylindrical_csys && parameters.isParamSetByUser("axis"))
102  {
103  this->paramError("axis",
104  "\nIf parameter axis is provided then the parameter local_cylindrical_csys "
105  "should be set to true");
106  }
107 }
GenericReal< is_ad > _qsigma
Square of the q function for orthotropy.
GenericMaterialProperty< Real, is_ad > & _hardening_variable
const MaterialProperty< RankTwoTensor > & _rotation_matrix_transpose_old
const MaterialProperty< Real > & _hardening_variable_old
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _elasticity_eigenvectors
enum HillElastoPlasticityStressUpdateTempl::Axis _axis
MaterialProperty< RankTwoTensor > & _rotation_matrix_transpose
GenericMaterialProperty< DenseVector< Real >, is_ad > & _sigma_tilde_rotated
const std::string _elasticity_tensor_name
Name of the elasticity tensor material property.
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _b_eigenvectors
GenericMaterialProperty< DenseVector< Real >, is_ad > & _sigma_tilde
GenericMaterialProperty< DenseVector< Real >, is_ad > & _elasticity_eigenvalues
const MaterialProperty< DenseMatrix< Real > > & _hill_tensor
Hill tensor, when global axes do not (somehow) align with those of the material Example: Large rotati...
GenericMaterialProperty< DenseVector< Real >, is_ad > & _b_eigenvalues
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _alpha_matrix
bool isParamSetByUser(const std::string &name) const
MaterialProperty< RankTwoTensor > & _rotation_matrix
const MaterialProperty< RankTwoTensor > & _rotation_matrix_old
This class provides baseline functionality for anisotropic (Hill-like) plasticity models based on the...
void paramError(const std::string &param, Args... args) const
const GenericMaterialProperty< RankFourTensor, is_ad > & _elasticity_tensor
Anisotropic elasticity tensor material property.

Member Function Documentation

◆ computeDeltaDerivatives()

template<bool is_ad>
void HillElastoPlasticityStressUpdateTempl< is_ad >::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 
)
protected

Definition at line 461 of file HillElastoPlasticityStressUpdate.C.

468 {
469  omega_gamma = 0.0;
470  sy_gamma = 0.0;
471 
472  GenericDenseVector<is_ad> K_deltaGamma(6);
473 
475  omega = computeOmega(delta_gamma, _sigma_tilde_rotated[_qp]);
476  else
477  omega = computeOmega(delta_gamma, _sigma_tilde[_qp]);
478 
480  for (unsigned int i = 0; i < 6; i++)
481  K(i) = _b_eigenvalues[_qp](i) / Utility::pow<2>(1 + delta_gamma * _b_eigenvalues[_qp](i));
482 
483  for (unsigned int i = 0; i < 6; i++)
484  K_deltaGamma(i) =
485  -2.0 * _b_eigenvalues[_qp](i) * K(i) / (1 + _b_eigenvalues[_qp](i) * delta_gamma);
486 
487  for (unsigned int i = 0; i < 6; i++)
488  {
490  omega_gamma += K_deltaGamma(i) * _sigma_tilde_rotated[_qp](i) * _sigma_tilde_rotated[_qp](i);
491  else
492  omega_gamma += K_deltaGamma(i) * _sigma_tilde[_qp](i) * _sigma_tilde[_qp](i);
493  }
494 
495  omega_gamma /= 4.0 * omega;
496  sy_gamma = 2.0 * sy_alpha * (omega + delta_gamma * omega_gamma);
497 }
static const std::string K
Definition: NS.h:170
GenericMaterialProperty< DenseVector< Real >, is_ad > & _sigma_tilde_rotated
Moose::GenericType< DenseVector< Real >, is_ad > GenericDenseVector
GenericMaterialProperty< DenseVector< Real >, is_ad > & _sigma_tilde
GenericReal< is_ad > computeOmega(const GenericReal< is_ad > &delta_gamma, const GenericDenseVector< is_ad > &stress_trial)
GenericMaterialProperty< DenseVector< Real >, is_ad > & _b_eigenvalues

◆ computeDerivative()

template<bool is_ad>
GenericReal< is_ad > HillElastoPlasticityStressUpdateTempl< is_ad >::computeDerivative ( const GenericDenseVector< is_ad > &  effective_trial_stress,
const GenericDenseVector< is_ad > &  stress_new,
const GenericReal< is_ad > &  scalar 
)
overrideprotectedvirtual

Definition at line 432 of file HillElastoPlasticityStressUpdate.C.

436 {
437  // If in elastic regime, return unit derivative
438  if (_yield_condition <= 0.0)
439  return 1.0;
440 
441  GenericReal<is_ad> omega = computeOmega(delta_gamma, stress_sigma);
443 
444  GenericReal<is_ad> hardeningVariable = computeHardeningValue(delta_gamma, omega);
445  GenericReal<is_ad> sy =
446  _hardening_constant * std::pow(hardeningVariable + 1.0e-30, _hardening_exponent) +
449 
450  GenericReal<is_ad> omega_gamma;
451  GenericReal<is_ad> sy_gamma;
452 
453  computeDeltaDerivatives(delta_gamma, stress_sigma, sy_alpha, omega, omega_gamma, sy_gamma);
454  GenericReal<is_ad> residual_derivative = 1 / omega * (sy_gamma - 1 / omega * omega_gamma * sy);
455 
456  return residual_derivative;
457 }
Moose::GenericType< Real, is_ad > GenericReal
GenericReal< is_ad > computeHardeningValue(const GenericReal< is_ad > &scalar, const GenericReal< is_ad > &omega)
GenericReal< is_ad > computeOmega(const GenericReal< is_ad > &delta_gamma, const GenericDenseVector< is_ad > &stress_trial)
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)
MooseUnits pow(const MooseUnits &, int)

◆ computeElasticityTensorEigenDecomposition()

template<bool is_ad>
void HillElastoPlasticityStressUpdateTempl< is_ad >::computeElasticityTensorEigenDecomposition ( )
protected

Compute eigendecomposition of element-wise elasticity tensor.

Definition at line 218 of file HillElastoPlasticityStressUpdate.C.

219 {
220  // Helper method to compute qp matrices required for the elasto-plasticity algorithm.
221 
223 
225  {
227  }
228 
229  ElasticityTensorTools::toVoigtNotation<is_ad>(_anisotropic_elastic_tensor, elasticity_tensor);
230 
231  const unsigned int dimension = _anisotropic_elastic_tensor.n();
232 
233  AnisotropyMatrixReal A = AnisotropyMatrixReal::Zero();
234  for (unsigned int index_i = 0; index_i < dimension; index_i++)
235  for (unsigned int index_j = 0; index_j < dimension; index_j++)
236  A(index_i, index_j) = MetaPhysicL::raw_value(_anisotropic_elastic_tensor(index_i, index_j));
237 
238  Eigen::SelfAdjointEigenSolver<AnisotropyMatrixReal> es(A);
239 
240  auto lambda = es.eigenvalues();
241  auto v = es.eigenvectors();
242 
243  for (unsigned int index_i = 0; index_i < dimension; index_i++)
244  _elasticity_eigenvalues[_qp](index_i) = lambda(index_i);
245 
246  for (unsigned int index_i = 0; index_i < dimension; index_i++)
247  for (unsigned int index_j = 0; index_j < dimension; index_j++)
248  _elasticity_eigenvectors[_qp](index_i, index_j) = v(index_i, index_j);
249 
250  // Compute sqrt(Delta_c) * QcT * A * Qc * sqrt(Delta_c)
251  GenericDenseMatrix<is_ad> sqrt_Delta(6, 6);
252  for (unsigned int i = 0; i < 6; i++)
253  {
254  sqrt_Delta(i, i) = std::sqrt(_elasticity_eigenvalues[_qp](i));
255  }
256 
257  GenericDenseMatrix<is_ad> eigenvectors_elasticity_transpose(6, 6);
258  _elasticity_eigenvectors[_qp].get_transpose(eigenvectors_elasticity_transpose);
259 
260  GenericDenseMatrix<is_ad> b_matrix(6, 6);
261  b_matrix = _hill_tensor[_qp];
262 
263  // Right multiply by matrix of eigenvectors transpose
264  b_matrix.right_multiply(_elasticity_eigenvectors[_qp]);
265  b_matrix.right_multiply(sqrt_Delta);
266  b_matrix.left_multiply(eigenvectors_elasticity_transpose);
267  b_matrix.left_multiply(sqrt_Delta);
268 
269  AnisotropyMatrixReal B_eigen = AnisotropyMatrixReal::Zero();
270  for (unsigned int index_i = 0; index_i < dimension; index_i++)
271  for (unsigned int index_j = 0; index_j < dimension; index_j++)
272  B_eigen(index_i, index_j) = MetaPhysicL::raw_value(b_matrix(index_i, index_j));
273 
274  if (isBlockDiagonal(B_eigen))
275  {
276  Eigen::SelfAdjointEigenSolver<AnisotropyMatrixRealBlock> es_b_left(B_eigen.block<3, 3>(0, 0));
277  Eigen::SelfAdjointEigenSolver<AnisotropyMatrixRealBlock> es_b_right(B_eigen.block<3, 3>(3, 3));
278 
279  auto lambda_b_left = es_b_left.eigenvalues();
280  auto v_b_left = es_b_left.eigenvectors();
281 
282  auto lambda_b_right = es_b_right.eigenvalues();
283  auto v_b_right = es_b_right.eigenvectors();
284 
285  _b_eigenvalues[_qp](0) = lambda_b_left(0);
286  _b_eigenvalues[_qp](1) = lambda_b_left(1);
287  _b_eigenvalues[_qp](2) = lambda_b_left(2);
288  _b_eigenvalues[_qp](3) = lambda_b_right(0);
289  _b_eigenvalues[_qp](4) = lambda_b_right(1);
290  _b_eigenvalues[_qp](5) = lambda_b_right(2);
291 
292  _b_eigenvectors[_qp](0, 0) = v_b_left(0, 0);
293  _b_eigenvectors[_qp](0, 1) = v_b_left(0, 1);
294  _b_eigenvectors[_qp](0, 2) = v_b_left(0, 2);
295  _b_eigenvectors[_qp](1, 0) = v_b_left(1, 0);
296  _b_eigenvectors[_qp](1, 1) = v_b_left(1, 1);
297  _b_eigenvectors[_qp](1, 2) = v_b_left(1, 2);
298  _b_eigenvectors[_qp](2, 0) = v_b_left(2, 0);
299  _b_eigenvectors[_qp](2, 1) = v_b_left(2, 1);
300  _b_eigenvectors[_qp](2, 2) = v_b_left(2, 2);
301  _b_eigenvectors[_qp](3, 3) = v_b_right(0, 0);
302  _b_eigenvectors[_qp](3, 4) = v_b_right(0, 1);
303  _b_eigenvectors[_qp](3, 5) = v_b_right(0, 2);
304  _b_eigenvectors[_qp](4, 3) = v_b_right(1, 0);
305  _b_eigenvectors[_qp](4, 4) = v_b_right(1, 1);
306  _b_eigenvectors[_qp](4, 5) = v_b_right(1, 2);
307  _b_eigenvectors[_qp](5, 3) = v_b_right(2, 0);
308  _b_eigenvectors[_qp](5, 4) = v_b_right(2, 1);
309  _b_eigenvectors[_qp](5, 5) = v_b_right(2, 2);
310  }
311  else
312  {
313  Eigen::SelfAdjointEigenSolver<AnisotropyMatrixReal> es_b(B_eigen);
314 
315  auto lambda_b = es_b.eigenvalues();
316  auto v_b = es_b.eigenvectors();
317  for (unsigned int index_i = 0; index_i < dimension; index_i++)
318  _b_eigenvalues[_qp](index_i) = lambda_b(index_i);
319 
320  for (unsigned int index_i = 0; index_i < dimension; index_i++)
321  for (unsigned int index_j = 0; index_j < dimension; index_j++)
322  _b_eigenvectors[_qp](index_i, index_j) = v_b(index_i, index_j);
323  }
324 
325  _alpha_matrix[_qp] = sqrt_Delta;
326  _alpha_matrix[_qp].right_multiply(_b_eigenvectors[_qp]);
327  _alpha_matrix[_qp].left_multiply(_elasticity_eigenvectors[_qp]);
328 }
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _elasticity_eigenvectors
Eigen::Matrix< Real, 6, 6 > AnisotropyMatrixReal
MaterialProperty< RankTwoTensor > & _rotation_matrix_transpose
auto raw_value(const Eigen::Map< T > &in)
Moose::GenericType< RankFourTensor, is_ad > GenericRankFourTensor
Real elasticity_tensor(unsigned int i, unsigned int j, unsigned int k, unsigned int l)
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _b_eigenvectors
GenericMaterialProperty< DenseVector< Real >, is_ad > & _elasticity_eigenvalues
const MaterialProperty< DenseMatrix< Real > > & _hill_tensor
Hill tensor, when global axes do not (somehow) align with those of the material Example: Large rotati...
GenericMaterialProperty< DenseVector< Real >, is_ad > & _b_eigenvalues
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _alpha_matrix
static const std::string v
Definition: NS.h:84
const GenericMaterialProperty< RankFourTensor, is_ad > & _elasticity_tensor
Anisotropic elasticity tensor material property.
Moose::GenericType< DenseMatrix< Real >, is_ad > GenericDenseMatrix

◆ computeHardeningDerivative()

template<bool is_ad>
Real HillElastoPlasticityStressUpdateTempl< is_ad >::computeHardeningDerivative ( )
protected

Definition at line 509 of file HillElastoPlasticityStressUpdate.C.

510 {
513  std::pow(_hardening_variable[_qp] + 1.0e-30, _hardening_exponent - 1));
514 }
GenericMaterialProperty< Real, is_ad > & _hardening_variable
auto raw_value(const Eigen::Map< T > &in)
MooseUnits pow(const MooseUnits &, int)

◆ computeHardeningValue()

template<bool is_ad>
GenericReal< is_ad > HillElastoPlasticityStressUpdateTempl< is_ad >::computeHardeningValue ( const GenericReal< is_ad > &  scalar,
const GenericReal< is_ad > &  omega 
)
protected

Definition at line 501 of file HillElastoPlasticityStressUpdate.C.

503 {
504  return _hardening_variable_old[_qp] + 2.0 * delta_gamma * omega;
505 }
const MaterialProperty< Real > & _hardening_variable_old

◆ computeHillTensorEigenDecomposition()

template<bool is_ad>
void HillElastoPlasticityStressUpdateTempl< is_ad >::computeHillTensorEigenDecomposition ( const DenseMatrix< GenericReal< is_ad >> &  hill_tensor)
protected

Compute eigendecomposition of Hill's tensor for anisotropic plasticity.

Parameters
hill_tensor6x6 matrix representing fourth order Hill's tensor describing anisotropy

◆ computeOmega()

template<bool is_ad>
GenericReal< is_ad > HillElastoPlasticityStressUpdateTempl< is_ad >::computeOmega ( const GenericReal< is_ad > &  delta_gamma,
const GenericDenseVector< is_ad > &  stress_trial 
)
protected

Definition at line 332 of file HillElastoPlasticityStressUpdate.C.

334 {
336  GenericReal<is_ad> omega = 0.0;
337 
338  for (unsigned int i = 0; i < 6; i++)
339  {
340  K(i) = _b_eigenvalues[_qp](i) / Utility::pow<2>(1 + delta_gamma * _b_eigenvalues[_qp](i));
341 
343  omega += K(i) * _sigma_tilde_rotated[_qp](i) * _sigma_tilde_rotated[_qp](i);
344  else
345  omega += K(i) * _sigma_tilde[_qp](i) * _sigma_tilde[_qp](i);
346  }
347  omega *= 0.5;
348 
349  if (omega == 0.0)
350  omega = 1.0e-36;
351 
352  return std::sqrt(omega);
353 }
Moose::GenericType< Real, is_ad > GenericReal
static const std::string K
Definition: NS.h:170
GenericMaterialProperty< DenseVector< Real >, is_ad > & _sigma_tilde_rotated
Moose::GenericType< DenseVector< Real >, is_ad > GenericDenseVector
GenericMaterialProperty< DenseVector< Real >, is_ad > & _sigma_tilde
GenericMaterialProperty< DenseVector< Real >, is_ad > & _b_eigenvalues

◆ computeReferenceResidual()

template<bool is_ad>
Real HillElastoPlasticityStressUpdateTempl< is_ad >::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 
)
overrideprotectedvirtual

Definition at line 357 of file HillElastoPlasticityStressUpdate.C.

362 {
363  // Equation is already normalized.
364  return 1.0;
365 }

◆ computeResidual()

template<bool is_ad>
GenericReal< is_ad > HillElastoPlasticityStressUpdateTempl< is_ad >::computeResidual ( const GenericDenseVector< is_ad > &  effective_trial_stress,
const GenericDenseVector< is_ad > &  stress_new,
const GenericReal< is_ad > &  scalar 
)
overrideprotectedvirtual

Definition at line 369 of file HillElastoPlasticityStressUpdate.C.

373 {
374 
375  // If in elastic regime, just return
376  if (_yield_condition <= 0.0)
377  return 0.0;
378 
379  GenericReal<is_ad> omega;
380 
382  {
383  // Get stress_tilde_rotated
385  RankTwoScalarTools::VoigtVectorToRankTwoTensor<is_ad>(stress_sigma, stress);
386 
387  stress.rotate(_rotation_matrix_transpose[_qp]);
388 
389  GenericDenseVector<is_ad> stress_sigma_rotated(6);
390  RankTwoScalarTools::RankTwoTensorToVoigtVector<is_ad>(stress, stress_sigma_rotated);
391 
392  GenericDenseVector<is_ad> stress_tilde_rotated(6);
393  GenericDenseMatrix<is_ad> alpha_temp_rotated(_alpha_matrix[_qp]);
394  alpha_temp_rotated.lu_solve(stress_sigma_rotated, stress_tilde_rotated);
395 
396  // Material property used in computeStressFinalize
397  _sigma_tilde_rotated[_qp] = stress_tilde_rotated;
398  omega = computeOmega(delta_gamma, stress_tilde_rotated);
399  }
400 
401  else
402  {
403  // Get stress_tilde
404  GenericDenseVector<is_ad> stress_tilde(6);
405  GenericDenseMatrix<is_ad> alpha_temp(_alpha_matrix[_qp]);
406  alpha_temp.lu_solve(stress_sigma, stress_tilde);
407 
408  // Material property used in computeStressFinalize
409  _sigma_tilde[_qp] = stress_tilde;
410  omega = computeOmega(delta_gamma, stress_tilde);
411  }
412 
413  // Hardening variable is \alpha isotropic hardening for now.
414  _hardening_variable[_qp] = computeHardeningValue(delta_gamma, omega);
415 
416  // A small value of 1.0e-30 is added to the hardening variable to avoid numerical issues
417  // related to hardening variable becoming negative early in the iteration leading to non-
418  // convergence. Note that std::pow(x,n) requires x to be positive if n is less than 1.
419 
420  GenericReal<is_ad> s_y =
423 
424  GenericReal<is_ad> residual = 0.0;
425  residual = s_y / omega - 1.0;
426 
427  return residual;
428 }
Moose::GenericType< Real, is_ad > GenericReal
GenericMaterialProperty< Real, is_ad > & _hardening_variable
MaterialProperty< RankTwoTensor > & _rotation_matrix_transpose
GenericMaterialProperty< DenseVector< Real >, is_ad > & _sigma_tilde_rotated
Moose::GenericType< DenseVector< Real >, is_ad > GenericDenseVector
GenericReal< is_ad > computeHardeningValue(const GenericReal< is_ad > &scalar, const GenericReal< is_ad > &omega)
GenericMaterialProperty< DenseVector< Real >, is_ad > & _sigma_tilde
GenericReal< is_ad > computeOmega(const GenericReal< is_ad > &delta_gamma, const GenericDenseVector< is_ad > &stress_trial)
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _alpha_matrix
MooseUnits pow(const MooseUnits &, int)
Moose::GenericType< DenseMatrix< Real >, is_ad > GenericDenseMatrix
Moose::GenericType< RankTwoTensor, is_ad > GenericRankTwoTensor

◆ computeStrainEnergyRateDensity()

template<bool is_ad>
Real HillElastoPlasticityStressUpdateTempl< is_ad >::computeStrainEnergyRateDensity ( const GenericMaterialProperty< RankTwoTensor, is_ad > &  stress,
const GenericMaterialProperty< RankTwoTensor, is_ad > &  strain_rate 
)
overridevirtual

Definition at line 610 of file HillElastoPlasticityStressUpdate.C.

613 {
614  mooseError("computeStrainEnergyRateDensity not implemented for anisotropic plasticity.");
615  return 0.0;
616 }
void mooseError(Args &&... args)

◆ computeStrainFinalize()

template<bool is_ad>
void HillElastoPlasticityStressUpdateTempl< is_ad >::computeStrainFinalize ( GenericRankTwoTensor< is_ad > &  inelasticStrainIncrement,
const GenericRankTwoTensor< is_ad > &  stress,
const GenericDenseVector< is_ad > &  stress_dev,
const GenericReal< is_ad > &  delta_gamma 
)
overrideprotectedvirtual

Perform any necessary steps to finalize strain increment after return mapping iterations.

Parameters
inelasticStrainIncrementInelastic strain increment
stressCauchy stresss tensor
stress_devDeviatoric partt of the Cauchy stresss tensor
delta_gammaGeneralized radial return's plastic multiplier

Reimplemented from AnisotropicReturnPlasticityStressUpdateBaseTempl< is_ad >.

Definition at line 518 of file HillElastoPlasticityStressUpdate.C.

523 {
524  // e^P = delta_gamma * hill_tensor * stress
525  GenericDenseVector<is_ad> inelasticStrainIncrement_vector(6);
526  GenericDenseVector<is_ad> hill_stress(6);
527  GenericDenseVector<is_ad> stress_vector(6);
528 
530  {
531  GenericRankTwoTensor<is_ad> stress_rotated(stress);
532  stress_rotated.rotate(_rotation_matrix_transpose[_qp]);
533 
534  RankTwoScalarTools::RankTwoTensorToVoigtVector<is_ad>(stress_rotated, stress_vector);
535  }
536  else
537  {
538  RankTwoScalarTools::RankTwoTensorToVoigtVector<is_ad>(stress, stress_vector);
539  }
540 
541  _hill_tensor[_qp].vector_mult(hill_stress, stress_vector);
542  hill_stress.scale(delta_gamma);
543  inelasticStrainIncrement_vector = hill_stress;
544 
545  inelasticStrainIncrement(0, 0) = inelasticStrainIncrement_vector(0);
546  inelasticStrainIncrement(1, 1) = inelasticStrainIncrement_vector(1);
547  inelasticStrainIncrement(2, 2) = inelasticStrainIncrement_vector(2);
548  inelasticStrainIncrement(0, 1) = inelasticStrainIncrement(1, 0) =
549  inelasticStrainIncrement_vector(3) / 2.0;
550  inelasticStrainIncrement(1, 2) = inelasticStrainIncrement(2, 1) =
551  inelasticStrainIncrement_vector(4) / 2.0;
552  inelasticStrainIncrement(0, 2) = inelasticStrainIncrement(2, 0) =
553  inelasticStrainIncrement_vector(5) / 2.0;
554 
556  inelasticStrainIncrement.rotate(_rotation_matrix[_qp]);
557 
558  // Calculate equivalent plastic strain
559  GenericDenseVector<is_ad> Mepsilon(6);
560  _hill_tensor[_qp].vector_mult(Mepsilon, inelasticStrainIncrement_vector);
561  GenericReal<is_ad> eq_plastic_strain_inc = Mepsilon.dot(inelasticStrainIncrement_vector);
562 
563  if (eq_plastic_strain_inc > 0.0)
564  eq_plastic_strain_inc = std::sqrt(eq_plastic_strain_inc);
565 
566  _effective_inelastic_strain[_qp] = _effective_inelastic_strain_old[_qp] + eq_plastic_strain_inc;
567 
569  inelasticStrainIncrement, stress, stress_dev, delta_gamma);
570 }
Moose::GenericType< Real, is_ad > GenericReal
MaterialProperty< RankTwoTensor > & _rotation_matrix_transpose
Moose::GenericType< DenseVector< Real >, is_ad > GenericDenseVector
const MaterialProperty< DenseMatrix< Real > > & _hill_tensor
Hill tensor, when global axes do not (somehow) align with those of the material Example: Large rotati...
virtual void computeStrainFinalize(GenericRankTwoTensor< is_ad > &, const GenericRankTwoTensor< is_ad > &, const GenericDenseVector< is_ad > &, const GenericReal< is_ad > &) override
Perform any necessary steps to finalize strain increment after return mapping iterations.
MaterialProperty< RankTwoTensor > & _rotation_matrix
Moose::GenericType< RankTwoTensor, is_ad > GenericRankTwoTensor

◆ computeStressDerivative()

template<bool is_ad>
virtual GenericReal<is_ad> AnisotropicReturnPlasticityStressUpdateBaseTempl< is_ad >::computeStressDerivative ( const Real  ,
const Real   
)
inlineoverrideprotectedvirtualinherited

Calculate the derivative of the strain increment with respect to the updated stress.

Parameters
effective_trial_stressEffective trial stress
scalarInelastic strain increment magnitude being solved for

Definition at line 47 of file AnisotropicReturnPlasticityStressUpdateBase.h.

49  {
50  return 0.0;
51  }

◆ computeStressFinalize()

template<bool is_ad>
void HillElastoPlasticityStressUpdateTempl< is_ad >::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 
)
overrideprotectedvirtual

Perform any necessary steps to finalize state after return mapping iterations.

Parameters
inelasticStrainIncrementInelastic strain increment
delta_gammaGeneralized radial return's plastic multiplier
stressCauchy stresss tensor
stress_devDeviatoric partt of the Cauchy stresss tensor

Definition at line 574 of file HillElastoPlasticityStressUpdate.C.

581 {
582  // Need to compute this iteration's stress tensor based on the scalar variable
583  // For deviatoric
584  // sigma(n+1) = {Alpha [I + delta_gamma*Delta_b]^(-1) A^-1} sigma(trial)
585 
586  if (_yield_condition <= 0.0)
587  return;
588  GenericDenseMatrix<is_ad> inv_matrix(6, 6);
589  inv_matrix.zero();
590 
591  for (unsigned int i = 0; i < 6; i++)
592  inv_matrix(i, i) = 1 / (1 + delta_gamma * _b_eigenvalues[_qp](i));
593 
594  _alpha_matrix[_qp].right_multiply(inv_matrix);
595 
596  GenericDenseVector<is_ad> stress_output(6);
598  _alpha_matrix[_qp].vector_mult(stress_output, _sigma_tilde_rotated[_qp]);
599  else
600  _alpha_matrix[_qp].vector_mult(stress_output, _sigma_tilde[_qp]);
601 
602  RankTwoScalarTools::VoigtVectorToRankTwoTensor<is_ad>(stress_output, stress_new);
603 
605  stress_new.rotate(_rotation_matrix[_qp]);
606 }
GenericMaterialProperty< DenseVector< Real >, is_ad > & _sigma_tilde_rotated
Moose::GenericType< DenseVector< Real >, is_ad > GenericDenseVector
GenericMaterialProperty< DenseVector< Real >, is_ad > & _sigma_tilde
GenericMaterialProperty< DenseVector< Real >, is_ad > & _b_eigenvalues
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _alpha_matrix
MaterialProperty< RankTwoTensor > & _rotation_matrix
Moose::GenericType< DenseMatrix< Real >, is_ad > GenericDenseMatrix

◆ computeStressInitialize()

template<bool is_ad>
void HillElastoPlasticityStressUpdateTempl< is_ad >::computeStressInitialize ( const GenericDenseVector< is_ad > &  stress_dev,
const GenericDenseVector< is_ad > &  stress,
const GenericRankFourTensor< is_ad > &  elasticity_tensor 
)
overrideprotectedvirtual

Definition at line 190 of file HillElastoPlasticityStressUpdate.C.

194 {
197  _effective_inelastic_strain[_qp] = _effective_inelastic_strain_old[_qp];
198 
199  _elasticity_eigenvectors[_qp].resize(6, 6);
200  _elasticity_eigenvalues[_qp].resize(6);
201 
202  _b_eigenvectors[_qp].resize(6, 6);
203  _b_eigenvalues[_qp].resize(6);
204 
205  _alpha_matrix[_qp].resize(6, 6);
206  _sigma_tilde[_qp].resize(6);
207  _sigma_tilde_rotated[_qp].resize(6);
208 
209  // Algebra needed for the generalized return mapping of anisotropic elastoplasticity
211 
212  _yield_condition = 1.0; // Some positive value
213  _yield_condition = -computeResidual(stress_dev, stress, 0.0);
214 }
GenericMaterialProperty< Real, is_ad > & _hardening_variable
void computeElasticityTensorEigenDecomposition()
Compute eigendecomposition of element-wise elasticity tensor.
const MaterialProperty< Real > & _hardening_variable_old
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _elasticity_eigenvectors
GenericMaterialProperty< DenseVector< Real >, is_ad > & _sigma_tilde_rotated
virtual GenericReal< is_ad > computeResidual(const GenericDenseVector< is_ad > &effective_trial_stress, const GenericDenseVector< is_ad > &stress_new, const GenericReal< is_ad > &scalar) override
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _b_eigenvectors
GenericMaterialProperty< DenseVector< Real >, is_ad > & _sigma_tilde
GenericMaterialProperty< DenseVector< Real >, is_ad > & _elasticity_eigenvalues
GenericMaterialProperty< DenseVector< Real >, is_ad > & _b_eigenvalues
GenericMaterialProperty< DenseMatrix< Real >, is_ad > & _alpha_matrix
GenericMaterialProperty< RankTwoTensor, is_ad > & _plasticity_strain
Plasticity strain tensor material property.

◆ initQpStatefulProperties()

template<bool is_ad>
void HillElastoPlasticityStressUpdateTempl< is_ad >::initQpStatefulProperties ( )
overrideprotectedvirtual

Reimplemented from AnisotropicReturnPlasticityStressUpdateBaseTempl< is_ad >.

Definition at line 111 of file HillElastoPlasticityStressUpdate.C.

112 {
113  RealVectorValue normal_vector(0, 0, 0);
114  RealVectorValue axial_vector(0, 0, 0);
115 
117  {
118  if (_axis == Axis::X)
119  {
120  normal_vector(0) = 0.0;
121  normal_vector(1) = _q_point[_qp](1);
122  normal_vector(2) = _q_point[_qp](2);
123  axial_vector(0) = 1.0;
124  }
125  else if (_axis == Axis::Y)
126  {
127  normal_vector(0) = _q_point[_qp](0);
128  normal_vector(1) = 0.0;
129  normal_vector(2) = _q_point[_qp](2);
130  axial_vector(1) = 1.0;
131  }
132  else if (_axis == Axis::Z)
133  {
134  normal_vector(0) = _q_point[_qp](0);
135  normal_vector(1) = _q_point[_qp](1);
136  normal_vector(2) = 0.0;
137  axial_vector(2) = 1.0;
138  }
139  else
140  mooseError("\nInvalid value for axis parameter provided in HillElastoPlasticityStressUpdate");
141  }
142 
144  {
145  Real nv_norm = normal_vector.norm();
146  Real av_norm = axial_vector.norm();
147 
148  if (!(MooseUtils::absoluteFuzzyEqual(nv_norm, 0.0)))
149  normal_vector /= nv_norm;
150  else
151  {
152  mooseError("The normal vector cannot be a zero vector in "
153  "HillElastoPlasticityStressUpdate");
154  }
155 
156  if (!(MooseUtils::absoluteFuzzyEqual(av_norm, 0.0)))
157  axial_vector /= av_norm;
158  else
159  {
160  mooseError("The axial vector cannot be a zero vector in "
161  "HillElastoPlasticityStressUpdate");
162  }
163 
165  normal_vector, axial_vector, _rotation_matrix[_qp], false);
167  }
168 
170 }
bool absoluteFuzzyEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
void mooseError(Args &&... args)
enum HillElastoPlasticityStressUpdateTempl::Axis _axis
MaterialProperty< RankTwoTensor > & _rotation_matrix_transpose
void setRotationMatrix(const RealVectorValue &outwardnormal, const RealVectorValue &axialVector, RankTwoTensor &rotationMatrix, const bool transpose)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
MaterialProperty< RankTwoTensor > & _rotation_matrix
const MooseArray< Point > & _q_point

◆ paramError()

template<bool is_ad>
void Material::paramError
protected

◆ propagateQpStatefulProperties()

template<bool is_ad>
void HillElastoPlasticityStressUpdateTempl< is_ad >::propagateQpStatefulProperties ( )
overrideprotectedvirtual

Reimplemented from AnisotropicReturnPlasticityStressUpdateBaseTempl< is_ad >.

Definition at line 174 of file HillElastoPlasticityStressUpdate.C.

175 {
178 
180  {
183  }
184 
186 }
GenericMaterialProperty< Real, is_ad > & _hardening_variable
const MaterialProperty< RankTwoTensor > & _rotation_matrix_transpose_old
const MaterialProperty< Real > & _hardening_variable_old
MaterialProperty< RankTwoTensor > & _rotation_matrix_transpose
MaterialProperty< RankTwoTensor > & _rotation_matrix
GenericMaterialProperty< RankTwoTensor, is_ad > & _plasticity_strain
Plasticity strain tensor material property.
const MaterialProperty< RankTwoTensor > & _rotation_matrix_old

◆ requiresIsotropicTensor()

template<bool is_ad>
bool HillElastoPlasticityStressUpdateTempl< is_ad >::requiresIsotropicTensor ( )
inlineoverrideprotected

Does the model require the elasticity tensor to be isotropic? No, this class does anisotropic elasto-plasticity

Definition at line 77 of file HillElastoPlasticityStressUpdate.h.

77 { return false; }

◆ validParams()

template<bool is_ad>
InputParameters HillElastoPlasticityStressUpdateTempl< is_ad >::validParams ( )
static

Definition at line 19 of file HillElastoPlasticityStressUpdate.C.

20 {
22  params.addClassDescription(
23  "This class uses the generalized radial return for anisotropic elasto-plasticity model."
24  "This class can be used in conjunction with other creep and plasticity materials for "
25  "more complex simulations.");
26 
27  params.addRequiredParam<Real>("hardening_constant",
28  "Hardening constant (H) for anisotropic plasticity");
29  params.addParam<Real>(
30  "hardening_exponent", 1.0, "Hardening exponent (n) for anisotropic plasticity");
31  params.addRequiredParam<Real>("yield_stress",
32  "Yield stress (constant value) for anisotropic plasticity");
33  params.addParam<bool>(
34  "local_cylindrical_csys",
35  false,
36  "Compute inelastic strain increment in local cylindrical coordinate system");
37 
38  params.addParam<MooseEnum>("axis",
39  MooseEnum("x y z", "y"),
40  "The axis of cylindrical component about which to rotate when "
41  "computing inelastic strain increment in local coordinate system");
42 
43  return params;
44 }
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
void addRequiredParam(const std::string &name, const std::string &doc_string)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void addClassDescription(const std::string &doc_string)

Member Data Documentation

◆ _alpha_matrix

template<bool is_ad>
GenericMaterialProperty<DenseMatrix<Real>, is_ad>& HillElastoPlasticityStressUpdateTempl< is_ad >::_alpha_matrix
protected

Definition at line 149 of file HillElastoPlasticityStressUpdate.h.

◆ _anisotropic_elastic_tensor

template<bool is_ad>
GenericDenseMatrix<is_ad> HillElastoPlasticityStressUpdateTempl< is_ad >::_anisotropic_elastic_tensor
protected

Definition at line 130 of file HillElastoPlasticityStressUpdate.h.

◆ _axis

template<bool is_ad>
enum HillElastoPlasticityStressUpdateTempl::Axis HillElastoPlasticityStressUpdateTempl< is_ad >::_axis
protected

◆ _b_eigenvalues

template<bool is_ad>
GenericMaterialProperty<DenseVector<Real>, is_ad>& HillElastoPlasticityStressUpdateTempl< is_ad >::_b_eigenvalues
protected

Definition at line 147 of file HillElastoPlasticityStressUpdate.h.

◆ _b_eigenvectors

template<bool is_ad>
GenericMaterialProperty<DenseMatrix<Real>, is_ad>& HillElastoPlasticityStressUpdateTempl< is_ad >::_b_eigenvectors
protected

Definition at line 146 of file HillElastoPlasticityStressUpdate.h.

◆ _current_elem

template<bool is_ad>
const Elem *const & Material::_current_elem
protected

◆ _eigenvalues_hill

template<bool is_ad>
GenericDenseVector<is_ad> HillElastoPlasticityStressUpdateTempl< is_ad >::_eigenvalues_hill
protected

Definition at line 128 of file HillElastoPlasticityStressUpdate.h.

◆ _eigenvectors_hill

template<bool is_ad>
GenericDenseMatrix<is_ad> HillElastoPlasticityStressUpdateTempl< is_ad >::_eigenvectors_hill
protected

Definition at line 129 of file HillElastoPlasticityStressUpdate.h.

◆ _elasticity_eigenvalues

template<bool is_ad>
GenericMaterialProperty<DenseVector<Real>, is_ad>& HillElastoPlasticityStressUpdateTempl< is_ad >::_elasticity_eigenvalues
protected

Definition at line 144 of file HillElastoPlasticityStressUpdate.h.

◆ _elasticity_eigenvectors

template<bool is_ad>
GenericMaterialProperty<DenseMatrix<Real>, is_ad>& HillElastoPlasticityStressUpdateTempl< is_ad >::_elasticity_eigenvectors
protected

Definition at line 143 of file HillElastoPlasticityStressUpdate.h.

◆ _elasticity_tensor

template<bool is_ad>
const GenericMaterialProperty<RankFourTensor, is_ad>& HillElastoPlasticityStressUpdateTempl< is_ad >::_elasticity_tensor
protected

Anisotropic elasticity tensor material property.

Definition at line 135 of file HillElastoPlasticityStressUpdate.h.

◆ _elasticity_tensor_name

template<bool is_ad>
const std::string HillElastoPlasticityStressUpdateTempl< is_ad >::_elasticity_tensor_name
protected

Name of the elasticity tensor material property.

Definition at line 133 of file HillElastoPlasticityStressUpdate.h.

◆ _hardening_constant

template<bool is_ad>
const Real HillElastoPlasticityStressUpdateTempl< is_ad >::_hardening_constant
protected

Definition at line 137 of file HillElastoPlasticityStressUpdate.h.

◆ _hardening_derivative

template<bool is_ad>
GenericReal<is_ad> HillElastoPlasticityStressUpdateTempl< is_ad >::_hardening_derivative
protected

Definition at line 153 of file HillElastoPlasticityStressUpdate.h.

◆ _hardening_exponent

template<bool is_ad>
const Real HillElastoPlasticityStressUpdateTempl< is_ad >::_hardening_exponent
protected

Definition at line 138 of file HillElastoPlasticityStressUpdate.h.

◆ _hardening_variable

template<bool is_ad>
GenericMaterialProperty<Real, is_ad>& HillElastoPlasticityStressUpdateTempl< is_ad >::_hardening_variable
protected

Definition at line 140 of file HillElastoPlasticityStressUpdate.h.

◆ _hardening_variable_old

template<bool is_ad>
const MaterialProperty<Real>& HillElastoPlasticityStressUpdateTempl< is_ad >::_hardening_variable_old
protected

Definition at line 141 of file HillElastoPlasticityStressUpdate.h.

◆ _hill_tensor

template<bool is_ad>
const MaterialProperty<DenseMatrix<Real> >& HillElastoPlasticityStressUpdateTempl< is_ad >::_hill_tensor
protected

Hill tensor, when global axes do not (somehow) align with those of the material Example: Large rotation due to rigid body and/or large deformation kinematics.

Definition at line 159 of file HillElastoPlasticityStressUpdate.h.

◆ _local_cylindrical_csys

template<bool is_ad>
const bool HillElastoPlasticityStressUpdateTempl< is_ad >::_local_cylindrical_csys
protected

◆ _plasticity_strain

template<bool is_ad>
GenericMaterialProperty<RankTwoTensor, is_ad>& AnisotropicReturnPlasticityStressUpdateBaseTempl< is_ad >::_plasticity_strain
protectedinherited

Plasticity strain tensor material property.

Definition at line 66 of file AnisotropicReturnPlasticityStressUpdateBase.h.

◆ _plasticity_strain_old

template<bool is_ad>
const MaterialProperty<RankTwoTensor>& AnisotropicReturnPlasticityStressUpdateBaseTempl< is_ad >::_plasticity_strain_old
protectedinherited

◆ _q_point

template<bool is_ad>
const MooseArray< Point > & Material::_q_point
protected

◆ _qp

template<bool is_ad>
unsigned int Material::_qp
protected

◆ _qsigma

template<bool is_ad>
GenericReal<is_ad> HillElastoPlasticityStressUpdateTempl< is_ad >::_qsigma
protected

Square of the q function for orthotropy.

Definition at line 126 of file HillElastoPlasticityStressUpdate.h.

◆ _rotation_matrix

template<bool is_ad>
MaterialProperty<RankTwoTensor>& HillElastoPlasticityStressUpdateTempl< is_ad >::_rotation_matrix
protected

Definition at line 161 of file HillElastoPlasticityStressUpdate.h.

◆ _rotation_matrix_old

template<bool is_ad>
const MaterialProperty<RankTwoTensor>& HillElastoPlasticityStressUpdateTempl< is_ad >::_rotation_matrix_old
protected

Definition at line 163 of file HillElastoPlasticityStressUpdate.h.

◆ _rotation_matrix_transpose

template<bool is_ad>
MaterialProperty<RankTwoTensor>& HillElastoPlasticityStressUpdateTempl< is_ad >::_rotation_matrix_transpose
protected

Definition at line 162 of file HillElastoPlasticityStressUpdate.h.

◆ _rotation_matrix_transpose_old

template<bool is_ad>
const MaterialProperty<RankTwoTensor>& HillElastoPlasticityStressUpdateTempl< is_ad >::_rotation_matrix_transpose_old
protected

Definition at line 164 of file HillElastoPlasticityStressUpdate.h.

◆ _sigma_tilde

template<bool is_ad>
GenericMaterialProperty<DenseVector<Real>, is_ad>& HillElastoPlasticityStressUpdateTempl< is_ad >::_sigma_tilde
protected

Definition at line 150 of file HillElastoPlasticityStressUpdate.h.

◆ _sigma_tilde_rotated

template<bool is_ad>
GenericMaterialProperty<DenseVector<Real>, is_ad>& HillElastoPlasticityStressUpdateTempl< is_ad >::_sigma_tilde_rotated
protected

Definition at line 151 of file HillElastoPlasticityStressUpdate.h.

◆ _yield_condition

template<bool is_ad>
GenericReal<is_ad> HillElastoPlasticityStressUpdateTempl< is_ad >::_yield_condition
protected

Definition at line 154 of file HillElastoPlasticityStressUpdate.h.

◆ _yield_stress

template<bool is_ad>
GenericReal<is_ad> HillElastoPlasticityStressUpdateTempl< is_ad >::_yield_stress
protected

Definition at line 155 of file HillElastoPlasticityStressUpdate.h.

◆ usingTransientInterfaceMembers

template<bool is_ad>
HillElastoPlasticityStressUpdateTempl< is_ad >::usingTransientInterfaceMembers
protected

Definition at line 35 of file HillElastoPlasticityStressUpdate.h.


The documentation for this class was generated from the following files: