www.mooseframework.org
ADViscoplasticityStressUpdate.C
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 
11 
12 #include "libmesh/utility.h"
13 
15 
17 
18 template <ComputeStage compute_stage>
19 InputParameters
21 {
23  params.addClassDescription(
24  "This material computes the non-linear homogenized gauge stress in order to compute the "
25  "viscoplastic responce due to creep in porous materials. This material must be used in "
26  "conjunction with ADComputeMultiplePorousInelasticStress");
27  MooseEnum viscoplasticity_model("LPS GTN", "LPS");
28  params.addParam<MooseEnum>(
29  "viscoplasticity_model", viscoplasticity_model, "Which viscoplastic model to use");
30  MooseEnum pore_shape_model("spherical cylindrical", "spherical");
31  params.addParam<MooseEnum>("pore_shape_model", pore_shape_model, "Which pore shape model to use");
32  params.addRequiredParam<MaterialPropertyName>(
33  "coefficient", "Material property name for the leading coefficient for Norton power law");
34  params.addRequiredRangeCheckedParam<Real>(
35  "power", "power>=1.0", "Stress exponent for Norton power law");
36  params.addParam<Real>(
37  "maximum_gauge_ratio",
38  1.0e6,
39  "Maximum ratio between the gauge stress and the equivalent stress. This "
40  "should be a high number. Note that this does not set an upper bound on the value, but "
41  "rather will help with convergence of the inner Newton loop");
42  params.addParam<Real>(
43  "minimum_equivalent_stress",
44  1.0e-3,
45  "Minimum value of equivalent stress below which viscoplasticiy is not calculated.");
46  params.addParam<Real>("maximum_equivalent_stress",
47  1.0e12,
48  "Maximum value of equivalent stress above which an exception is thrown "
49  "instead of calculating the properties in this material.");
50 
51  params.addParamNamesToGroup("verbose maximum_gauge_ratio maximum_equivalent_stress", "Advanced");
52  return params;
53 }
54 
55 template <ComputeStage compute_stage>
57  const InputParameters & parameters)
58  : ADViscoplasticityStressUpdateBase<compute_stage>(parameters),
59  _model(parameters.get<MooseEnum>("viscoplasticity_model").getEnum<ViscoplasticityModel>()),
60  _pore_shape(parameters.get<MooseEnum>("pore_shape_model").getEnum<PoreShapeModel>()),
61  _pore_shape_factor(_pore_shape == PoreShapeModel::SPHERICAL ? 1.5 : std::sqrt(3.0)),
62  _power(getParam<Real>("power")),
63  _power_factor(_model == ViscoplasticityModel::LPS ? (_power - 1.0) / (_power + 1.0) : 1.0),
64  _coefficient(getADMaterialProperty<Real>("coefficient")),
65  _gauge_stress(declareADProperty<Real>(_base_name + "gauge_stress")),
66  _maximum_gauge_ratio(getParam<Real>("maximum_gauge_ratio")),
67  _minimum_equivalent_stress(getParam<Real>("minimum_equivalent_stress")),
68  _maximum_equivalent_stress(getParam<Real>("maximum_equivalent_stress")),
69  _hydro_stress(0.0),
70  _identity_two(RankTwoTensor::initIdentity),
71  _dhydro_stress_dsigma(_identity_two / 3.0)
72 {
73  _check_range = true;
74 }
75 
76 template <ComputeStage compute_stage>
77 void
79  ADRankTwoTensor & elastic_strain_increment,
80  ADRankTwoTensor & inelastic_strain_increment,
81  const ADRankTwoTensor & /*rotation_increment*/,
82  ADRankTwoTensor & stress,
83  const RankTwoTensor & /*stress_old*/,
84  const ADRankFourTensor & elasticity_tensor,
85  const RankTwoTensor & elastic_strain_old)
86 {
87  // Compute initial hydrostatic stress and porosity
88  if (_pore_shape == PoreShapeModel::CYLINDRICAL)
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 : std::sqrt(1.5 * dev_stress_squared);
99 
100  computeStressInitialize(equiv_stress, elasticity_tensor);
101 
102  // Prepare values
103  _effective_inelastic_strain[_qp] = _effective_inelastic_strain_old[_qp];
104  _inelastic_strain[_qp] = _inelastic_strain_old[_qp];
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 (",
114  _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 
123  computeInelasticStrainIncrement(_gauge_stress[_qp],
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 : std::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 template <ComputeStage compute_stage>
160 ADReal
162 {
163  return effective_trial_stress;
164 }
165 
166 template <ComputeStage compute_stage>
167 ADReal
169  const ADReal & effective_trial_stress) const
170 {
171  return effective_trial_stress * _maximum_gauge_ratio;
172 }
173 
174 template <ComputeStage compute_stage>
175 ADReal
177  const ADReal & effective_trial_stress) const
178 {
179  return effective_trial_stress;
180 }
181 
182 template <ComputeStage compute_stage>
183 ADReal
185  const ADReal & trial_gauge)
186 {
187  const ADReal M = std::abs(_hydro_stress) / trial_gauge;
188  const ADReal dM_dtrial_gauge = -M / trial_gauge;
189 
190  const ADReal residual_left = Utility::pow<2>(equiv_stress / trial_gauge);
191  const ADReal dresidual_left_dtrial_gauge = -2.0 * residual_left / trial_gauge;
192 
193  ADReal residual = residual_left;
194  _derivative = dresidual_left_dtrial_gauge;
195 
196  if (_pore_shape == PoreShapeModel::SPHERICAL)
197  {
198  residual *= 1.0 + _intermediate_porosity / 1.5;
199  _derivative *= 1.0 + _intermediate_porosity / 1.5;
200  }
201 
202  if (_model == ViscoplasticityModel::GTN)
203  {
204  residual += 2.0 * _intermediate_porosity * std::cosh(_pore_shape_factor * M) - 1.0 -
205  Utility::pow<2>(_intermediate_porosity);
206  _derivative += 2.0 * _intermediate_porosity * std::sinh(_pore_shape_factor * M) *
207  _pore_shape_factor * dM_dtrial_gauge;
208  }
209  else
210  {
211  const ADReal h = computeH(_power, M);
212  const ADReal dh_dM = computeH(_power, M, true);
213 
214  residual += _intermediate_porosity * (h + _power_factor / h) - 1.0 -
215  _power_factor * Utility::pow<2>(_intermediate_porosity);
216  const ADReal dresidual_dh = _intermediate_porosity * (1.0 - _power_factor / Utility::pow<2>(h));
217  _derivative += dresidual_dh * dh_dM * dM_dtrial_gauge;
218  }
219 
220  if (_verbose)
221  {
222  Moose::out << "in computeResidual:\n"
223  << " position: " << _q_point[_qp] << " hydro_stress: " << _hydro_stress
224  << " equiv_stress: " << equiv_stress << " trial_grage: " << trial_gauge
225  << " M: " << M << std::endl;
226  Moose::out << " residual: " << residual << " derivative: " << _derivative << std::endl;
227  }
228 
229  return residual;
230 }
231 
232 template <ComputeStage compute_stage>
233 ADReal
235  const ADReal & M,
236  const bool derivative)
237 {
238  const ADReal mod = std::pow(M * _pore_shape_factor, (n + 1.0) / n);
239 
240  // Calculate derivative with respect to M
241  if (derivative)
242  {
243  const ADReal dmod_dM = (n + 1.0) / n * mod / M;
244  return dmod_dM * std::pow(1.0 + mod / n, n - 1.0);
245  }
246 
247  return std::pow(1.0 + mod / n, n);
248 }
249 
250 template <ComputeStage compute_stage>
251 ADRankTwoTensor
253  const ADReal & gauge_stress,
254  const ADReal & equiv_stress,
255  const ADRankTwoTensor & dev_stress,
256  const ADRankTwoTensor & stress)
257 {
258  // Compute the derivative of the gauge stress with respect to the equilvalent and hydrostatic
259  // stress components
260  const ADReal M = std::abs(_hydro_stress) / gauge_stress;
261  const ADReal h = computeH(_power, M);
262 
263  // Compute the derviative of the residual with respect to the hydrostatic stress
264  ADReal dresidual_dhydro_stress = 0.0;
265  if (_hydro_stress != 0.0)
266  {
267  const ADReal dM_dhydro_stress = M / _hydro_stress;
268  if (_model == ViscoplasticityModel::GTN)
269  {
270  dresidual_dhydro_stress = 2.0 * _intermediate_porosity * std::sinh(_pore_shape_factor * M) *
271  _pore_shape_factor * dM_dhydro_stress;
272  }
273  else
274  {
275  const ADReal dresidual_dh =
276  _intermediate_porosity * (1.0 - _power_factor / Utility::pow<2>(h));
277  const ADReal dh_dM = computeH(_power, M, true);
278  dresidual_dhydro_stress = dresidual_dh * dh_dM * dM_dhydro_stress;
279  }
280  }
281 
282  // Compute the derivative of the residual with respect to the equilvalent stress
283  ADReal dresidual_dequiv_stress = 2.0 * equiv_stress / Utility::pow<2>(gauge_stress);
284  if (_pore_shape == PoreShapeModel::SPHERICAL)
285  dresidual_dequiv_stress *= 1.0 + _intermediate_porosity / 1.5;
286 
287  // Compute the derivative of the equilvalent stress to the deviatoric stress
288  const ADRankTwoTensor dequiv_stress_dsigma = 1.5 * dev_stress / equiv_stress;
289 
290  // Compute the derivative of the residual with the stress
291  const ADRankTwoTensor dresidual_dsigma = dresidual_dhydro_stress * _dhydro_stress_dsigma +
292  dresidual_dequiv_stress * dequiv_stress_dsigma;
293 
294  // Compute the deritative of the gauge stress with respect to the stress
295  const ADRankTwoTensor dgauge_dsigma =
296  dresidual_dsigma * (gauge_stress / dresidual_dsigma.doubleContraction(stress));
297 
298  return dgauge_dsigma;
299 }
300 
301 template <ComputeStage compute_stage>
302 void
304  ADReal & gauge_stress,
305  ADReal & dpsi_dgauge,
306  ADRankTwoTensor & inelastic_strain_increment,
307  const ADReal & equiv_stress,
308  const ADRankTwoTensor & dev_stress,
309  const ADRankTwoTensor & stress)
310 {
311  // If hydrostatic stress and porosity present, compute non-linear gauge stress
312  if (_intermediate_porosity == 0.0)
313  gauge_stress = equiv_stress;
314  else if (_hydro_stress == 0.0)
315  gauge_stress =
316  equiv_stress / std::sqrt(1.0 - (1.0 + _power_factor) * _intermediate_porosity +
317  _power_factor * Utility::pow<2>(_intermediate_porosity));
318  else
319  returnMappingSolve(equiv_stress, gauge_stress, _console);
320 
321  mooseAssert(gauge_stress >= equiv_stress,
322  "Gauge stress calculated in inner Newton solve is less than the equivalent stress.");
323 
324  // Compute stress potential
325  dpsi_dgauge = _coefficient[_qp] * std::pow(gauge_stress, _power);
326 
327  // Compute strain increment from stress potential and the gauge stress derivative with respect
328  // to the stress stress. The current form is explicit, and should eventually be changed
329  inelastic_strain_increment =
330  _dt * dpsi_dgauge * computeDGaugeDSigma(gauge_stress, equiv_stress, dev_stress, stress);
331 }
ADViscoplasticityStressUpdate::computeH
ADReal computeH(const Real n, const ADReal &gauge_stress, const bool derivative=false)
Definition: ADViscoplasticityStressUpdate.C:234
registerADMooseObject
registerADMooseObject("TensorMechanicsApp", ADViscoplasticityStressUpdate)
ADViscoplasticityStressUpdate::PoreShapeModel
PoreShapeModel
Enum to choose which pore shape model to use.
Definition: ADViscoplasticityStressUpdate.h:74
ADViscoplasticityStressUpdate::minimumPermissibleValue
virtual ADReal minimumPermissibleValue(const ADReal &effective_trial_stress) const override
Compute the minimum permissible value of the scalar.
Definition: ADViscoplasticityStressUpdate.C:176
pow
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
Definition: ExpressionBuilder.h:673
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
ADViscoplasticityStressUpdateBase::validParams
static InputParameters validParams()
Definition: ADViscoplasticityStressUpdateBase.C:16
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::ViscoplasticityModel
ViscoplasticityModel
Enum to choose which viscoplastic model to use.
Definition: ADViscoplasticityStressUpdate.h:71
ADViscoplasticityStressUpdate::validParams
static InputParameters validParams()
Definition: ADViscoplasticityStressUpdate.C:20
defineADLegacyParams
defineADLegacyParams(ADViscoplasticityStressUpdate)
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
ADViscoplasticityStressUpdate.h
ADSingleVariableReturnMappingSolution::_check_range
bool _check_range
Whether to check to see whether iterative solution is within admissible range, and set within that ra...
Definition: ADSingleVariableReturnMappingSolution.h:124
ADViscoplasticityStressUpdateBase
Definition: ADViscoplasticityStressUpdateBase.h:32
ADViscoplasticityStressUpdate::ADViscoplasticityStressUpdate
ADViscoplasticityStressUpdate(const InputParameters &parameters)
Definition: ADViscoplasticityStressUpdate.C:56
RankTwoTensorTempl< Real >
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