12 #include "libmesh/utility.h"
18 template <ComputeStage compute_stage>
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",
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",
45 "Minimum value of equivalent stress below which viscoplasticiy is not calculated.");
46 params.addParam<Real>(
"maximum_equivalent_stress",
48 "Maximum value of equivalent stress above which an exception is thrown "
49 "instead of calculating the properties in this material.");
51 params.addParamNamesToGroup(
"verbose maximum_gauge_ratio maximum_equivalent_stress",
"Advanced");
55 template <ComputeStage compute_stage>
57 const InputParameters & parameters)
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")),
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")),
71 _dhydro_stress_dsigma(_identity_two / 3.0)
76 template <ComputeStage compute_stage>
79 ADRankTwoTensor & elastic_strain_increment,
80 ADRankTwoTensor & inelastic_strain_increment,
81 const ADRankTwoTensor & ,
82 ADRankTwoTensor & stress,
84 const ADRankFourTensor & elasticity_tensor,
88 if (_pore_shape == PoreShapeModel::CYLINDRICAL)
89 _hydro_stress = (stress(0, 0) + stress(1, 1)) / 2.0;
91 _hydro_stress = stress.trace() / 3.0;
93 updateIntermediatePorosity(elastic_strain_increment);
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);
100 computeStressInitialize(equiv_stress, elasticity_tensor);
103 _effective_inelastic_strain[_qp] = _effective_inelastic_strain_old[_qp];
104 _inelastic_strain[_qp] = _inelastic_strain_old[_qp];
105 inelastic_strain_increment.zero();
108 if (equiv_stress > _maximum_equivalent_stress)
109 mooseException(
"In ",
111 ": equivalent stress (",
113 ") is higher than maximum_equivalent_stress (",
114 _maximum_equivalent_stress,
115 ").\nCutting time step.");
118 if (equiv_stress > _minimum_equivalent_stress)
121 ADReal dpsi_dgauge(0);
123 computeInelasticStrainIncrement(_gauge_stress[_qp],
125 inelastic_strain_increment,
131 elastic_strain_increment -= inelastic_strain_increment;
133 stress = elasticity_tensor * (elastic_strain_old + elastic_strain_increment);
137 _effective_inelastic_strain[_qp] += dpsi_dgauge * _dt;
139 _inelastic_strain[_qp] += inelastic_strain_increment;
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);
147 if (MooseUtils::relativeFuzzyGreaterThan(new_equiv_stress, equiv_stress))
148 mooseException(
"In ",
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.");
156 computeStressFinalize(inelastic_strain_increment);
159 template <ComputeStage compute_stage>
163 return effective_trial_stress;
166 template <ComputeStage compute_stage>
169 const ADReal & effective_trial_stress)
const
171 return effective_trial_stress * _maximum_gauge_ratio;
174 template <ComputeStage compute_stage>
177 const ADReal & effective_trial_stress)
const
179 return effective_trial_stress;
182 template <ComputeStage compute_stage>
185 const ADReal & trial_gauge)
187 const ADReal M = std::abs(_hydro_stress) / trial_gauge;
188 const ADReal dM_dtrial_gauge = -M / trial_gauge;
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;
193 ADReal residual = residual_left;
194 _derivative = dresidual_left_dtrial_gauge;
196 if (_pore_shape == PoreShapeModel::SPHERICAL)
198 residual *= 1.0 + _intermediate_porosity / 1.5;
199 _derivative *= 1.0 + _intermediate_porosity / 1.5;
202 if (_model == ViscoplasticityModel::GTN)
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;
211 const ADReal h = computeH(_power, M);
212 const ADReal dh_dM = computeH(_power, M,
true);
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;
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;
232 template <ComputeStage compute_stage>
236 const bool derivative)
238 const ADReal mod =
std::pow(M * _pore_shape_factor, (n + 1.0) / n);
243 const ADReal dmod_dM = (n + 1.0) / n * mod / M;
244 return dmod_dM *
std::pow(1.0 + mod / n, n - 1.0);
250 template <ComputeStage compute_stage>
253 const ADReal & gauge_stress,
254 const ADReal & equiv_stress,
255 const ADRankTwoTensor & dev_stress,
256 const ADRankTwoTensor & stress)
260 const ADReal M = std::abs(_hydro_stress) / gauge_stress;
261 const ADReal h = computeH(_power, M);
264 ADReal dresidual_dhydro_stress = 0.0;
265 if (_hydro_stress != 0.0)
267 const ADReal dM_dhydro_stress = M / _hydro_stress;
268 if (_model == ViscoplasticityModel::GTN)
270 dresidual_dhydro_stress = 2.0 * _intermediate_porosity * std::sinh(_pore_shape_factor * M) *
271 _pore_shape_factor * dM_dhydro_stress;
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;
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;
288 const ADRankTwoTensor dequiv_stress_dsigma = 1.5 * dev_stress / equiv_stress;
291 const ADRankTwoTensor dresidual_dsigma = dresidual_dhydro_stress * _dhydro_stress_dsigma +
292 dresidual_dequiv_stress * dequiv_stress_dsigma;
295 const ADRankTwoTensor dgauge_dsigma =
296 dresidual_dsigma * (gauge_stress / dresidual_dsigma.doubleContraction(stress));
298 return dgauge_dsigma;
301 template <ComputeStage compute_stage>
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)
312 if (_intermediate_porosity == 0.0)
313 gauge_stress = equiv_stress;
314 else if (_hydro_stress == 0.0)
316 equiv_stress / std::sqrt(1.0 - (1.0 + _power_factor) * _intermediate_porosity +
317 _power_factor * Utility::pow<2>(_intermediate_porosity));
319 returnMappingSolve(equiv_stress, gauge_stress, _console);
321 mooseAssert(gauge_stress >= equiv_stress,
322 "Gauge stress calculated in inner Newton solve is less than the equivalent stress.");
325 dpsi_dgauge = _coefficient[_qp] *
std::pow(gauge_stress, _power);
329 inelastic_strain_increment =
330 _dt * dpsi_dgauge * computeDGaugeDSigma(gauge_stress, equiv_stress, dev_stress, stress);