www.mooseframework.org
ADViscoplasticityStressUpdate.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 
13 
14 template <ComputeStage>
16 
18 
19 template <ComputeStage compute_stage>
21 {
22 public:
23  static InputParameters validParams();
24 
25  ADViscoplasticityStressUpdate(const InputParameters & parameters);
26 
27  virtual void updateState(ADRankTwoTensor & strain_increment,
28  ADRankTwoTensor & inelastic_strain_increment,
29  const ADRankTwoTensor & rotation_increment,
30  ADRankTwoTensor & stress_new,
31  const RankTwoTensor & stress_old,
32  const ADRankFourTensor & elasticity_tensor,
33  const RankTwoTensor & elastic_strain_old) override;
34 
35  virtual ADReal minimumPermissibleValue(const ADReal & effective_trial_stress) const override;
36 
37  virtual ADReal maximumPermissibleValue(const ADReal & effective_trial_stress) const override;
38 
39 protected:
47  virtual ADReal initialGuess(const ADReal & effective_trial_stress) override;
48 
53  virtual ADReal computeResidual(const ADReal & effective_trial_stress,
54  const ADReal & scalar) override;
55 
56  ADReal computeH(const Real n, const ADReal & gauge_stress, const bool derivative = false);
57 
58  ADRankTwoTensor computeDGaugeDSigma(const ADReal & gauge_stress,
59  const ADReal & equiv_stress,
60  const ADRankTwoTensor & dev_stress,
61  const ADRankTwoTensor & stress);
62 
63  void computeInelasticStrainIncrement(ADReal & gauge_stress,
64  ADReal & dpsi_dgauge,
65  ADRankTwoTensor & creep_strain_increment,
66  const ADReal & equiv_stress,
67  const ADRankTwoTensor & dev_stress,
68  const ADRankTwoTensor & stress);
69 
71  const enum class ViscoplasticityModel { LPS, GTN } _model;
72 
75 
77  const Real _pore_shape_factor;
78 
80  const Real _power;
81 
83  const Real _power_factor;
84 
86  const ADMaterialProperty(Real) & _coefficient;
87 
89  ADMaterialProperty(Real) & _gauge_stress;
90 
93 
96 
99 
102 
105 
108 
110 };
ADViscoplasticityStressUpdate::_pore_shape_factor
const Real _pore_shape_factor
Pore shape factor depending on pore shape model.
Definition: ADViscoplasticityStressUpdate.h:77
ADViscoplasticityStressUpdate::computeH
ADReal computeH(const Real n, const ADReal &gauge_stress, const bool derivative=false)
Definition: ADViscoplasticityStressUpdate.C:234
ADViscoplasticityStressUpdate::ViscoplasticityModel::LPS
ADViscoplasticityStressUpdate::PoreShapeModel
PoreShapeModel
Enum to choose which pore shape model to use.
Definition: ADViscoplasticityStressUpdate.h:74
ADViscoplasticityStressUpdate::_model
enum ADViscoplasticityStressUpdate::ViscoplasticityModel _model
ADViscoplasticityStressUpdate::minimumPermissibleValue
virtual ADReal minimumPermissibleValue(const ADReal &effective_trial_stress) const override
Compute the minimum permissible value of the scalar.
Definition: ADViscoplasticityStressUpdate.C:176
ADViscoplasticityStressUpdateBase.h
declareADValidParams
declareADValidParams(ADViscoplasticityStressUpdate)
ADViscoplasticityStressUpdate::PoreShapeModel::CYLINDRICAL
ADViscoplasticityStressUpdate::_identity_two
const RankTwoTensor _identity_two
Rank two identity tensor.
Definition: ADViscoplasticityStressUpdate.h:104
ADViscoplasticityStressUpdate::computeInelasticStrainIncrement
void computeInelasticStrainIncrement(ADReal &gauge_stress, ADReal &dpsi_dgauge, ADRankTwoTensor &creep_strain_increment, const ADReal &equiv_stress, const ADRankTwoTensor &dev_stress, const ADRankTwoTensor &stress)
Definition: ADViscoplasticityStressUpdate.C:303
ADViscoplasticityStressUpdate::_minimum_equivalent_stress
const Real _minimum_equivalent_stress
Minimum value of equivalent stress below which viscoplasticiy is not calculated.
Definition: ADViscoplasticityStressUpdate.h:95
ADViscoplasticityStressUpdate::initialGuess
virtual ADReal initialGuess(const ADReal &effective_trial_stress) override
Compute an initial guess for the value of the scalar.
Definition: ADViscoplasticityStressUpdate.C:161
ADViscoplasticityStressUpdate::_pore_shape
enum ADViscoplasticityStressUpdate::PoreShapeModel _pore_shape
ADViscoplasticityStressUpdate::_hydro_stress
ADReal _hydro_stress
Container for hydrostatic stress.
Definition: ADViscoplasticityStressUpdate.h:101
ADViscoplasticityStressUpdate::ADMaterialProperty
const ADMaterialProperty(Real) &_coefficient
Leading coefficient.
ADViscoplasticityStressUpdate::PoreShapeModel::SPHERICAL
ADViscoplasticityStressUpdate::ViscoplasticityModel
ViscoplasticityModel
Enum to choose which viscoplastic model to use.
Definition: ADViscoplasticityStressUpdate.h:71
ADViscoplasticityStressUpdate::_dhydro_stress_dsigma
const RankTwoTensor _dhydro_stress_dsigma
Derivative of hydrostatic stress with respect to the stress tensor.
Definition: ADViscoplasticityStressUpdate.h:107
ADViscoplasticityStressUpdate::_maximum_equivalent_stress
const Real _maximum_equivalent_stress
Maximum value of equivalent stress above which an exception is thrown.
Definition: ADViscoplasticityStressUpdate.h:98
ADViscoplasticityStressUpdate::_power_factor
const Real _power_factor
Power factor used for LPS model.
Definition: ADViscoplasticityStressUpdate.h:83
ADViscoplasticityStressUpdate::validParams
static InputParameters validParams()
Definition: ADViscoplasticityStressUpdate.C:20
ADViscoplasticityStressUpdate::_maximum_gauge_ratio
const Real _maximum_gauge_ratio
Maximum ratio between the gauge stress and the equilvalent stress.
Definition: ADViscoplasticityStressUpdate.h:92
ADViscoplasticityStressUpdate::_power
const Real _power
Exponent on the effective stress.
Definition: ADViscoplasticityStressUpdate.h:80
ADViscoplasticityStressUpdate::updateState
virtual void updateState(ADRankTwoTensor &strain_increment, ADRankTwoTensor &inelastic_strain_increment, const ADRankTwoTensor &rotation_increment, ADRankTwoTensor &stress_new, const RankTwoTensor &stress_old, const ADRankFourTensor &elasticity_tensor, const RankTwoTensor &elastic_strain_old) override
Given a strain increment that results in a trial stress, perform some procedure (such as an iterative...
Definition: ADViscoplasticityStressUpdate.C:78
ADViscoplasticityStressUpdate::computeResidual
virtual ADReal computeResidual(const ADReal &effective_trial_stress, const ADReal &scalar) override
Perform any necessary steps to finalize state after return mapping iterations.
Definition: ADViscoplasticityStressUpdate.C:184
ADViscoplasticityStressUpdate
Definition: ADViscoplasticityStressUpdate.h:15
ADViscoplasticityStressUpdateBase
Definition: ADViscoplasticityStressUpdateBase.h:32
ADViscoplasticityStressUpdate::usingViscoplasticityStressUpdateBaseMembers
usingViscoplasticityStressUpdateBaseMembers
Definition: ADViscoplasticityStressUpdate.h:109
ADViscoplasticityStressUpdate::ADViscoplasticityStressUpdate
ADViscoplasticityStressUpdate(const InputParameters &parameters)
Definition: ADViscoplasticityStressUpdate.C:56
RankTwoTensorTempl< Real >
ADViscoplasticityStressUpdate::ViscoplasticityModel::GTN
ADViscoplasticityStressUpdate::computeDGaugeDSigma
ADRankTwoTensor computeDGaugeDSigma(const ADReal &gauge_stress, const ADReal &equiv_stress, const ADRankTwoTensor &dev_stress, const ADRankTwoTensor &stress)
Definition: ADViscoplasticityStressUpdate.C:252
ADViscoplasticityStressUpdate::maximumPermissibleValue
virtual ADReal maximumPermissibleValue(const ADReal &effective_trial_stress) const override
Compute the maximum permissible value of the scalar.
Definition: ADViscoplasticityStressUpdate.C:168