https://mooseframework.inl.gov
ADViscoplasticityStressUpdate.C
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 
11 
12 #include "libmesh/utility.h"
13 
15 
18 {
21  params.addClassDescription(
22  "This material computes the non-linear homogenized gauge stress in order to compute the "
23  "viscoplastic responce due to creep in porous materials. This material must be used in "
24  "conjunction with ADComputeMultiplePorousInelasticStress");
25  MooseEnum viscoplasticity_model("LPS GTN", "LPS");
26  params.addParam<MooseEnum>(
27  "viscoplasticity_model", viscoplasticity_model, "Which viscoplastic model to use");
28  MooseEnum pore_shape_model("spherical cylindrical", "spherical");
29  params.addParam<MooseEnum>("pore_shape_model", pore_shape_model, "Which pore shape model to use");
30  params.addRequiredParam<MaterialPropertyName>(
31  "coefficient", "Material property name for the leading coefficient for Norton power law");
33  "power", "power>=1.0", "Stress exponent for Norton power law");
34  params.addParam<Real>(
35  "maximum_gauge_ratio",
36  1.0e6,
37  "Maximum ratio between the gauge stress and the equivalent stress. This "
38  "should be a high number. Note that this does not set an upper bound on the value, but "
39  "rather will help with convergence of the inner Newton loop");
40  params.addParam<Real>(
41  "minimum_equivalent_stress",
42  1.0e-3,
43  "Minimum value of equivalent stress below which viscoplasticiy is not calculated.");
44  params.addParam<Real>("maximum_equivalent_stress",
45  1.0e12,
46  "Maximum value of equivalent stress above which an exception is thrown "
47  "instead of calculating the properties in this material.");
48 
49  params.addParamNamesToGroup("verbose maximum_gauge_ratio maximum_equivalent_stress", "Advanced");
50  return params;
51 }
52 
56  _model(parameters.get<MooseEnum>("viscoplasticity_model").getEnum<ViscoplasticityModel>()),
57  _pore_shape(parameters.get<MooseEnum>("pore_shape_model").getEnum<PoreShapeModel>()),
58  _pore_shape_factor(_pore_shape == PoreShapeModel::SPHERICAL ? 1.5 : std::sqrt(3.0)),
59  _power(getParam<Real>("power")),
60  _power_factor(_model == ViscoplasticityModel::LPS ? (_power - 1.0) / (_power + 1.0) : 1.0),
61  _coefficient(getADMaterialProperty<Real>("coefficient")),
62  _gauge_stress(declareADProperty<Real>(_base_name + "gauge_stress")),
63  _maximum_gauge_ratio(getParam<Real>("maximum_gauge_ratio")),
64  _minimum_equivalent_stress(getParam<Real>("minimum_equivalent_stress")),
65  _maximum_equivalent_stress(getParam<Real>("maximum_equivalent_stress")),
66  _hydro_stress(0.0),
67  _identity_two(RankTwoTensor::initIdentity),
68  _dhydro_stress_dsigma(_identity_two / 3.0),
69  _derivative(0.0)
70 {
71  _check_range = true;
72 }
73 
74 void
76  ADRankTwoTensor & inelastic_strain_increment,
77  const ADRankTwoTensor & /*rotation_increment*/,
78  ADRankTwoTensor & stress,
79  const RankTwoTensor & /*stress_old*/,
80  const ADRankFourTensor & elasticity_tensor,
81  const RankTwoTensor & elastic_strain_old,
82  bool /*compute_full_tangent_operator = false*/,
83  RankFourTensor & /*tangent_operator = _identityTensor*/)
84 {
85  using std::sqrt;
86 
87  // Compute initial hydrostatic stress and porosity
89  _hydro_stress = (stress(0, 0) + stress(1, 1)) / 2.0;
90  else
91  _hydro_stress = stress.trace() / 3.0;
92 
93  updateIntermediatePorosity(elastic_strain_increment);
94 
95  // Compute intermediate equivalent stress
96  const ADRankTwoTensor dev_stress = stress.deviatoric();
97  const ADReal dev_stress_squared = dev_stress.doubleContraction(dev_stress);
98  const ADReal equiv_stress = dev_stress_squared == 0.0 ? 0.0 : sqrt(1.5 * dev_stress_squared);
99 
101 
102  // Prepare values
105  inelastic_strain_increment.zero();
106 
107  // Protect against extremely high values of stresses calculated by other viscoplastic materials
108  if (equiv_stress > _maximum_equivalent_stress)
109  mooseException("In ",
110  _name,
111  ": equivalent stress (",
112  equiv_stress,
113  ") is higher than maximum_equivalent_stress (",
115  ").\nCutting time step.");
116 
117  // If equivalent stress is present, calculate creep strain increment
118  if (equiv_stress > _minimum_equivalent_stress)
119  {
120  // Initalize stress potential
121  ADReal dpsi_dgauge(0);
122 
124  dpsi_dgauge,
125  inelastic_strain_increment,
126  equiv_stress,
127  dev_stress,
128  stress);
129 
130  // Update elastic strain increment due to inelastic strain calculated here
131  elastic_strain_increment -= inelastic_strain_increment;
132  // Update stress due to new strain
133  stress = elasticity_tensor * (elastic_strain_old + elastic_strain_increment);
134 
135  // Compute effective strain from the stress potential. Note that this is approximate and to be
136  // used qualitatively
137  _effective_inelastic_strain[_qp] += dpsi_dgauge * _dt;
138  // Update creep strain due to currently computed inelastic strain
139  _inelastic_strain[_qp] += inelastic_strain_increment;
140  }
141 
142  const ADRankTwoTensor new_dev_stress = stress.deviatoric();
143  const ADReal new_dev_stress_squared = new_dev_stress.doubleContraction(new_dev_stress);
144  const ADReal new_equiv_stress =
145  new_dev_stress_squared == 0.0 ? 0.0 : sqrt(1.5 * new_dev_stress_squared);
146 
147  if (MooseUtils::relativeFuzzyGreaterThan(new_equiv_stress, equiv_stress))
148  mooseException("In ",
149  _name,
150  ": updated equivalent stress (",
151  MetaPhysicL::raw_value(new_equiv_stress),
152  ") is greater than initial equivalent stress (",
153  MetaPhysicL::raw_value(equiv_stress),
154  "). Try decreasing max_inelastic_increment to avoid this exception.");
155 
156  computeStressFinalize(inelastic_strain_increment);
157 }
158 
159 ADReal
161 {
162  return effective_trial_stress;
163 }
164 
165 ADReal
167 {
168  return effective_trial_stress * _maximum_gauge_ratio;
169 }
170 
171 ADReal
173 {
174  return effective_trial_stress;
175 }
176 
177 ADReal
179  const ADReal & trial_gauge)
180 {
181  using std::abs, std::cosh, std::sinh;
182 
183  const ADReal M = abs(_hydro_stress) / trial_gauge;
184  const ADReal dM_dtrial_gauge = -M / trial_gauge;
185 
186  const ADReal residual_left = Utility::pow<2>(equiv_stress / trial_gauge);
187  const ADReal dresidual_left_dtrial_gauge = -2.0 * residual_left / trial_gauge;
188 
189  ADReal residual = residual_left;
190  _derivative = dresidual_left_dtrial_gauge;
191 
193  {
194  residual *= 1.0 + _intermediate_porosity / 1.5;
195  _derivative *= 1.0 + _intermediate_porosity / 1.5;
196  }
197 
199  {
200  residual += 2.0 * _intermediate_porosity * cosh(_pore_shape_factor * M) - 1.0 -
201  Utility::pow<2>(_intermediate_porosity);
203  _pore_shape_factor * dM_dtrial_gauge;
204  }
205  else
206  {
207  const ADReal h = computeH(_power, M);
208  const ADReal dh_dM = computeH(_power, M, true);
209 
210  residual += _intermediate_porosity * (h + _power_factor / h) - 1.0 -
211  _power_factor * Utility::pow<2>(_intermediate_porosity);
212  const ADReal dresidual_dh = _intermediate_porosity * (1.0 - _power_factor / Utility::pow<2>(h));
213  _derivative += dresidual_dh * dh_dM * dM_dtrial_gauge;
214  }
215 
216  if (_verbose)
217  {
218  Moose::out << "in computeResidual:\n"
219  << " position: " << _q_point[_qp] << " hydro_stress: " << _hydro_stress
220  << " equiv_stress: " << equiv_stress << " trial_grage: " << trial_gauge
221  << " M: " << M << std::endl;
222  Moose::out << " residual: " << residual << " derivative: " << _derivative << std::endl;
223  }
224 
225  return residual;
226 }
227 
228 ADReal
229 ADViscoplasticityStressUpdate::computeH(const Real n, const ADReal & M, const bool derivative)
230 {
231  using std::pow;
232 
233  const ADReal mod = pow(M * _pore_shape_factor, (n + 1.0) / n);
234 
235  // Calculate derivative with respect to M
236  if (derivative)
237  {
238  const ADReal dmod_dM = (n + 1.0) / n * mod / M;
239  return dmod_dM * pow(1.0 + mod / n, n - 1.0);
240  }
241 
242  return pow(1.0 + mod / n, n);
243 }
244 
247  const ADReal & equiv_stress,
248  const ADRankTwoTensor & dev_stress,
249  const ADRankTwoTensor & stress)
250 {
251  using std::abs, std::sinh;
252 
253  // Compute the derivative of the gauge stress with respect to the equivalent and hydrostatic
254  // stress components
255  const ADReal M = abs(_hydro_stress) / gauge_stress;
256  const ADReal h = computeH(_power, M);
257 
258  // Compute the derviative of the residual with respect to the hydrostatic stress
259  ADReal dresidual_dhydro_stress = 0.0;
260  if (_hydro_stress != 0.0)
261  {
262  const ADReal dM_dhydro_stress = M / _hydro_stress;
264  {
265  dresidual_dhydro_stress = 2.0 * _intermediate_porosity * sinh(_pore_shape_factor * M) *
266  _pore_shape_factor * dM_dhydro_stress;
267  }
268  else
269  {
270  const ADReal dresidual_dh =
271  _intermediate_porosity * (1.0 - _power_factor / Utility::pow<2>(h));
272  const ADReal dh_dM = computeH(_power, M, true);
273  dresidual_dhydro_stress = dresidual_dh * dh_dM * dM_dhydro_stress;
274  }
275  }
276 
277  // Compute the derivative of the residual with respect to the equivalent stress
278  ADReal dresidual_dequiv_stress = 2.0 * equiv_stress / Utility::pow<2>(gauge_stress);
280  dresidual_dequiv_stress *= 1.0 + _intermediate_porosity / 1.5;
281 
282  // Compute the derivative of the equilvalent stress to the deviatoric stress
283  const ADRankTwoTensor dequiv_stress_dsigma = 1.5 * dev_stress / equiv_stress;
284 
285  // Compute the derivative of the residual with the stress
286  const ADRankTwoTensor dresidual_dsigma = dresidual_dhydro_stress * _dhydro_stress_dsigma +
287  dresidual_dequiv_stress * dequiv_stress_dsigma;
288 
289  // Compute the deritative of the gauge stress with respect to the stress
290  const ADRankTwoTensor dgauge_dsigma =
291  dresidual_dsigma * (gauge_stress / dresidual_dsigma.doubleContraction(stress));
292 
293  return dgauge_dsigma;
294 }
295 
296 void
298  ADReal & gauge_stress,
299  ADReal & dpsi_dgauge,
300  ADRankTwoTensor & inelastic_strain_increment,
301  const ADReal & equiv_stress,
302  const ADRankTwoTensor & dev_stress,
303  const ADRankTwoTensor & stress)
304 {
305  using std::sqrt, std::pow;
306 
307  // If hydrostatic stress and porosity present, compute non-linear gauge stress
308  if (_intermediate_porosity == 0.0)
309  gauge_stress = equiv_stress;
310  else if (_hydro_stress == 0.0)
311  gauge_stress = equiv_stress / sqrt(1.0 - (1.0 + _power_factor) * _intermediate_porosity +
312  _power_factor * Utility::pow<2>(_intermediate_porosity));
313  else
314  returnMappingSolve(equiv_stress, gauge_stress, _console);
315 
316  mooseAssert(gauge_stress >= equiv_stress,
317  "Gauge stress calculated in inner Newton solve is less than the equivalent stress.");
318 
319  // Compute stress potential
320  dpsi_dgauge = _coefficient[_qp] * pow(gauge_stress, _power);
321 
322  // Compute strain increment from stress potential and the gauge stress derivative with respect
323  // to the stress stress. The current form is explicit, and should eventually be changed
324  inelastic_strain_increment =
325  _dt * dpsi_dgauge * computeDGaugeDSigma(gauge_stress, equiv_stress, dev_stress, stress);
326 }
327 
328 void
330  const unsigned int total_it)
331 {
332  if (iter_output)
333  {
334  *iter_output << "At element " << _current_elem->id() << " _qp=" << _qp << " Coordinates "
335  << _q_point[_qp] << " block=" << _current_elem->subdomain_id() << '\n';
336  }
338 }
339 
340 Real
342  const ADReal & gauge_stress)
343 {
344  // Use gauge stress for relative tolerance criteria, defined as:
345  // abs(residual / gauge_stress) <= _relative_tolerance
346  return MetaPhysicL::raw_value(gauge_stress);
347 }
virtual ADReal maximumPermissibleValue(const ADReal &effective_trial_stress) const override
MetaPhysicL::DualNumber< V, D, asd > abs(const MetaPhysicL::DualNumber< V, D, asd > &a)
void addRequiredRangeCheckedParam(const std::string &name, const std::string &parsed_function, const std::string &doc_string)
virtual void computeStressFinalize(const GenericRankTwoTensor< is_ad > &)
Perform any necessary steps to finalize state after return mapping iterations.
const std::string & _name
GenericMaterialProperty< Real, is_ad > & _effective_inelastic_strain
Effective inelastic strain material property.
const Real _minimum_equivalent_stress
Minimum value of equivalent stress below which viscoplasticiy is not calculated.
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
virtual ADReal computeResidual(const ADReal &effective_trial_stress, const ADReal &scalar) override
Perform any necessary steps to finalize state after return mapping iterations.
virtual ADReal minimumPermissibleValue(const ADReal &effective_trial_stress) const override
auto raw_value(const Eigen::Map< T > &in)
ADReal _hydro_stress
Container for hydrostatic stress.
ViscoplasticityModel
Enum to choose which viscoplastic model to use.
const Real _power
Exponent on the effective stress.
void updateIntermediatePorosity(const GenericRankTwoTensor< is_ad > &elastic_strain_increment)
DualNumber< Real, DNDerivativeType, true > ADReal
void addRequiredParam(const std::string &name, const std::string &doc_string)
GenericReal< is_ad > _intermediate_porosity
Container for the porosity calculated from all other intelastic models except the current model...
ADMaterialProperty< Real > & _gauge_stress
Gauge stress.
virtual ADReal initialGuess(const ADReal &effective_trial_stress) override
Compute an initial guess for the value of the scalar.
Real elasticity_tensor(unsigned int i, unsigned int j, unsigned int k, unsigned int l)
RankTwoTensorTempl< T > deviatoric() const
const MaterialProperty< RankTwoTensor > & _inelastic_strain_old
void returnMappingSolve(const GenericReal< is_ad > &effective_trial_stress, GenericReal< is_ad > &scalar, const ConsoleStream &console)
Perform the return mapping iterations.
void outputIterationSummary(std::stringstream *iter_output, const unsigned int total_it) override
Output summary information for the convergence history of the model.
virtual Real computeReferenceResidual(const ADReal &effective_trial_stress, const ADReal &scalar_effective_inelastic_strain) override
ADViscoplasticityStressUpdate(const InputParameters &parameters)
enum ADViscoplasticityStressUpdate::ViscoplasticityModel _model
virtual void outputIterationSummary(std::stringstream *iter_output, const unsigned int total_it)
Output summary information for the convergence history of the model.
Base class that provides capability for Newton return mapping iterations on a single variable...
T doubleContraction(const RankTwoTensorTempl< T > &a) const
ADReal _derivative
Container for _derivative.
const Real _maximum_gauge_ratio
Maximum ratio between the gauge stress and the equilvalent stress.
SPHERICAL
const Real _maximum_equivalent_stress
Maximum value of equivalent stress above which an exception is thrown.
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template sinh(_arg) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(erf
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
const bool _verbose
Flag to enable verbose output.
const Real _power_factor
Power factor used for LPS model.
const Elem *const & _current_elem
GenericMaterialProperty< RankTwoTensor, is_ad > & _inelastic_strain
Creep strain material property.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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, bool compute_full_tangent_operator=false, RankFourTensor &tangent_operator=_identityTensor) override
void computeInelasticStrainIncrement(ADReal &gauge_stress, ADReal &dpsi_dgauge, ADRankTwoTensor &creep_strain_increment, const ADReal &equiv_stress, const ADRankTwoTensor &dev_stress, const ADRankTwoTensor &stress)
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sqrt(_arg)) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tanh
const MooseArray< Point > & _q_point
enum ADViscoplasticityStressUpdate::PoreShapeModel _pore_shape
virtual void computeStressInitialize(const GenericReal< is_ad > &, const GenericRankFourTensor< is_ad > &)
Perform any necessary initialization before return mapping iterations.
ADRankTwoTensor computeDGaugeDSigma(const ADReal &gauge_stress, const ADReal &equiv_stress, const ADRankTwoTensor &dev_stress, const ADRankTwoTensor &stress)
void addClassDescription(const std::string &doc_string)
const Real _pore_shape_factor
Pore shape factor depending on pore shape model.
ADReal computeH(const Real n, const ADReal &gauge_stress, const bool derivative=false)
registerMooseObject("SolidMechanicsApp", ADViscoplasticityStressUpdate)
const ConsoleStream _console
const ADMaterialProperty< Real > & _coefficient
Leading coefficient.
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template cosh(_arg) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(cosh
MooseUnits pow(const MooseUnits &, int)
bool _check_range
Whether to check to see whether iterative solution is within admissible range, and set within that ra...
PoreShapeModel
Enum to choose which pore shape model to use.
const Elem & get(const ElemType type_in)
const MaterialProperty< Real > & _effective_inelastic_strain_old
const RankTwoTensor _dhydro_stress_dsigma
Derivative of hydrostatic stress with respect to the stress tensor.
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)