#include <ADViscoplasticityStressUpdate.h>
|
| ADViscoplasticityStressUpdate (const InputParameters ¶meters) |
|
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 return-mapping process) to produce an admissible stress, an elastic strain increment and an inelastic strain increment. More...
|
|
virtual ADReal | minimumPermissibleValue (const ADReal &effective_trial_stress) const override |
| Compute the minimum permissible value of the scalar. More...
|
|
virtual ADReal | maximumPermissibleValue (const ADReal &effective_trial_stress) const override |
| Compute the maximum permissible value of the scalar. More...
|
|
virtual Real | computeReferenceResidual (const ADReal &effective_trial_stress, const ADReal &scalar_effective_inelastic_strain) override |
| Compute a reference quantity to be used for checking relative convergence. More...
|
|
virtual Real | computeTimeStepLimit () override |
| Compute the limiting value of the time step for this material. More...
|
|
bool | requiresIsotropicTensor () override |
| Does the model require the elasticity tensor to be isotropic? More...
|
|
void | setQp (unsigned int qp) |
| Sets the value of the global variable _qp for inheriting classes. More...
|
|
|
void | resetQpProperties () final |
| Retained as empty methods to avoid a warning from Material.C in framework. These methods are unused in all inheriting classes and should not be overwritten. More...
|
|
void | resetProperties () final |
|
|
virtual ADReal | initialGuess (const ADReal &effective_trial_stress) override |
| Compute an initial guess for the value of the scalar. More...
|
|
virtual ADReal | computeResidual (const ADReal &effective_trial_stress, const ADReal &scalar) override |
| Perform any necessary steps to finalize state after return mapping iterations. More...
|
|
ADReal | computeH (const Real n, const ADReal &gauge_stress, const bool derivative=false) |
|
ADRankTwoTensor | computeDGaugeDSigma (const ADReal &gauge_stress, const ADReal &equiv_stress, const ADRankTwoTensor &dev_stress, const ADRankTwoTensor &stress) |
|
void | computeInelasticStrainIncrement (ADReal &gauge_stress, ADReal &dpsi_dgauge, ADRankTwoTensor &creep_strain_increment, const ADReal &equiv_stress, const ADRankTwoTensor &dev_stress, const ADRankTwoTensor &stress) |
|
const | ADMaterialProperty (Real) &_coefficient |
| Leading coefficient. More...
|
|
| ADMaterialProperty (Real) &_gauge_stress |
| Gauge stress. More...
|
|
virtual void | initQpStatefulProperties () override |
|
virtual void | propagateQpStatefulProperties () override |
| If updateState is not called during a timestep, this will be. More...
|
|
virtual void | computeStressInitialize (const ADReal &, const ADRankFourTensor &) |
| Perform any necessary initialization before return mapping iterations. More...
|
|
virtual void | computeStressFinalize (const ADRankTwoTensor &) |
| Perform any necessary steps to finalize state after return mapping iterations. More...
|
|
virtual ADReal | computeDerivative (const ADReal &, const ADReal &) override |
| Compute the derivative of the residual as a function of the scalar variable. More...
|
|
void | updateIntermediatePorosity (const ADRankTwoTensor &elastic_strain_increment) |
|
void | outputIterationSummary (std::stringstream *iter_output, const unsigned int total_it) override |
| Output summary information for the convergence history of the model. More...
|
|
const | ADMaterialProperty (RankTwoTensor) &_strain_increment |
| Material property for the total strain increment. More...
|
|
void | returnMappingSolve (const ADReal &effective_trial_stress, ADReal &scalar, const ConsoleStream &console) |
| Perform the return mapping iterations. More...
|
|
virtual void | iterationFinalize (ADReal) |
| Finalize internal state variables for a model for a given iteration. More...
|
|
virtual void | outputIterationStep (std::stringstream *iter_output, const ADReal &effective_trial_stress, const ADReal &scalar, const Real reference_residual) |
| Output information for a single iteration step to build the convergence history of the model. More...
|
|
bool | converged (const ADReal &residual, const Real reference) |
| Check to see whether the residual is within the convergence limits. More...
|
|
|
SolveState | internalSolve (const ADReal effective_trial_stress, ADReal &scalar, std::stringstream *iter_output=nullptr) |
| Method called from within this class to perform the actual return mappping iterations. More...
|
|
bool | convergedAcceptable (const unsigned int it, const Real reference) |
| Check to see whether the residual is within acceptable convergence limits. More...
|
|
void | checkPermissibleRange (ADReal &scalar, ADReal &scalar_increment, const ADReal &scalar_old, const ADReal min_permissible_scalar, const ADReal max_permissible_scalar, std::stringstream *iter_output) |
| Check to see whether solution is within admissible range, and set it within that range if it is not. More...
|
|
void | updateBounds (const ADReal &scalar, const ADReal &residual, const Real init_resid_sign, ADReal &scalar_upper_bound, ADReal &scalar_lower_bound, std::stringstream *iter_output) |
| Update the upper and lower bounds of the root for the effective inelastic strain. More...
|
|
template<ComputeStage compute_stage>
class ADViscoplasticityStressUpdate< compute_stage >
Definition at line 15 of file ADViscoplasticityStressUpdate.h.
◆ InternalSolveOutput
template<ComputeStage compute_stage>
◆ PoreShapeModel
template<ComputeStage compute_stage>
◆ SolveState
template<ComputeStage compute_stage>
◆ ViscoplasticityModel
template<ComputeStage compute_stage>
◆ ADViscoplasticityStressUpdate()
template<ComputeStage compute_stage>
◆ ADMaterialProperty() [1/3]
template<ComputeStage compute_stage>
Material property for the total strain increment.
◆ ADMaterialProperty() [2/3]
template<ComputeStage compute_stage>
◆ ADMaterialProperty() [3/3]
template<ComputeStage compute_stage>
◆ checkPermissibleRange()
template<ComputeStage compute_stage>
void ADSingleVariableReturnMappingSolution< compute_stage >::checkPermissibleRange |
( |
ADReal & |
scalar, |
|
|
ADReal & |
scalar_increment, |
|
|
const ADReal & |
scalar_old, |
|
|
const ADReal |
min_permissible_scalar, |
|
|
const ADReal |
max_permissible_scalar, |
|
|
std::stringstream * |
iter_output |
|
) |
| |
|
privateinherited |
Check to see whether solution is within admissible range, and set it within that range if it is not.
- Parameters
-
scalar | Current value of the inelastic strain increment |
scalar_increment | Incremental change in scalar from the previous iteration |
scalar_old | Previous value of scalar |
min_permissible_scalar | Minimum permissible value of scalar |
max_permissible_scalar | Maximum permissible value of scalar |
iter_output | Output stream |
Definition at line 334 of file ADSingleVariableReturnMappingSolution.C.
342 if (scalar > max_permissible_scalar)
344 scalar_increment = (max_permissible_scalar - scalar_old) / 2.0;
345 scalar = scalar_old + scalar_increment;
347 *iter_output <<
"Scalar greater than maximum ("
348 << MetaPhysicL::raw_value(max_permissible_scalar)
349 <<
") adjusted scalar=" << MetaPhysicL::raw_value(scalar)
350 <<
" scalar_increment=" << MetaPhysicL::raw_value(scalar_increment) << std::endl;
352 else if (scalar < min_permissible_scalar)
354 scalar_increment = (min_permissible_scalar - scalar_old) / 2.0;
355 scalar = scalar_old + scalar_increment;
357 *iter_output <<
"Scalar less than minimum (" << MetaPhysicL::raw_value(min_permissible_scalar)
358 <<
") adjusted scalar=" << MetaPhysicL::raw_value(scalar)
359 <<
" scalar_increment=" << MetaPhysicL::raw_value(scalar_increment) << std::endl;
◆ computeDerivative()
template<ComputeStage compute_stage>
|
inlineoverrideprotectedvirtualinherited |
◆ computeDGaugeDSigma()
template<ComputeStage compute_stage>
ADRankTwoTensor ADViscoplasticityStressUpdate< compute_stage >::computeDGaugeDSigma |
( |
const ADReal & |
gauge_stress, |
|
|
const ADReal & |
equiv_stress, |
|
|
const ADRankTwoTensor & |
dev_stress, |
|
|
const ADRankTwoTensor & |
stress |
|
) |
| |
|
protected |
Definition at line 252 of file ADViscoplasticityStressUpdate.C.
264 ADReal dresidual_dhydro_stress = 0.0;
275 const ADReal dresidual_dh =
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);
288 const ADRankTwoTensor dequiv_stress_dsigma = 1.5 * dev_stress / equiv_stress;
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;
◆ computeH()
template<ComputeStage compute_stage>
◆ computeInelasticStrainIncrement()
template<ComputeStage compute_stage>
void ADViscoplasticityStressUpdate< compute_stage >::computeInelasticStrainIncrement |
( |
ADReal & |
gauge_stress, |
|
|
ADReal & |
dpsi_dgauge, |
|
|
ADRankTwoTensor & |
creep_strain_increment, |
|
|
const ADReal & |
equiv_stress, |
|
|
const ADRankTwoTensor & |
dev_stress, |
|
|
const ADRankTwoTensor & |
stress |
|
) |
| |
|
protected |
Definition at line 303 of file ADViscoplasticityStressUpdate.C.
313 gauge_stress = equiv_stress;
321 mooseAssert(gauge_stress >= equiv_stress,
322 "Gauge stress calculated in inner Newton solve is less than the equivalent stress.");
329 inelastic_strain_increment =
◆ computeReferenceResidual()
template<ComputeStage compute_stage>
◆ computeResidual()
template<ComputeStage compute_stage>
Perform any necessary steps to finalize state after return mapping iterations.
- Parameters
-
inelasticStrainIncrement | Inelastic strain increment |
Implements ADSingleVariableReturnMappingSolution< compute_stage >.
Definition at line 184 of file ADViscoplasticityStressUpdate.C.
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;
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;
◆ computeStressFinalize()
template<ComputeStage compute_stage>
|
inlineprotectedvirtualinherited |
Perform any necessary steps to finalize state after return mapping iterations.
- Parameters
-
inelasticStrainIncrement | Inelastic strain increment |
Definition at line 79 of file ADViscoplasticityStressUpdateBase.h.
◆ computeStressInitialize()
template<ComputeStage compute_stage>
|
inlineprotectedvirtualinherited |
Perform any necessary initialization before return mapping iterations.
- Parameters
-
effective_trial_stress | Effective trial stress |
elasticityTensor | Elasticity tensor |
Definition at line 70 of file ADViscoplasticityStressUpdateBase.h.
◆ computeTimeStepLimit()
template<ComputeStage compute_stage>
Compute the limiting value of the time step for this material.
- Returns
- Limiting time step
Reimplemented from ADStressUpdateBase< compute_stage >.
Definition at line 92 of file ADViscoplasticityStressUpdateBase.C.
94 const Real scalar_inelastic_strain_incr =
95 MetaPhysicL::raw_value(_effective_inelastic_strain[_qp]) -
98 if (!scalar_inelastic_strain_incr)
99 return std::numeric_limits<Real>::max();
◆ converged()
template<ComputeStage compute_stage>
Check to see whether the residual is within the convergence limits.
- Parameters
-
residual | Current value of the residual |
reference | Current value of the reference quantity |
- Returns
- Whether the model converged
Definition at line 301 of file ADSingleVariableReturnMappingSolution.C.
304 const Real residual = MetaPhysicL::raw_value(ad_residual);
◆ convergedAcceptable()
template<ComputeStage compute_stage>
Check to see whether the residual is within acceptable convergence limits.
This will only return true if it has been determined that progress is no longer being made and that the residual is within the acceptable limits.
- Parameters
-
residual | Current iteration count |
residual | Current value of the residual |
reference | Current value of the reference quantity |
- Returns
- Whether the model converged
Definition at line 311 of file ADSingleVariableReturnMappingSolution.C.
322 const Real convergence_history_factor = 10.0;
323 if (std::abs(
_residual * convergence_history_factor) <
◆ initialGuess()
template<ComputeStage compute_stage>
Compute an initial guess for the value of the scalar.
For some cases, an intellegent starting point can provide enhanced robustness in the Newton iterations. This is also an opportunity for classes that derive from this to perform initialization tasks.
- Parameters
-
effective_trial_stress | Effective trial stress |
Reimplemented from ADSingleVariableReturnMappingSolution< compute_stage >.
Definition at line 161 of file ADViscoplasticityStressUpdate.C.
163 return effective_trial_stress;
◆ initQpStatefulProperties()
template<ComputeStage compute_stage>
|
overrideprotectedvirtualinherited |
◆ internalSolve()
template<ComputeStage compute_stage>
Method called from within this class to perform the actual return mappping iterations.
- Parameters
-
effective_trial_stress | Effective trial stress |
scalar | Inelastic strain increment magnitude being solved for |
iter_output | Output stream – if null, no output is produced |
- Returns
- Whether the solution was successful
Definition at line 162 of file ADSingleVariableReturnMappingSolution.C.
166 ADReal scalar_old = scalar;
167 ADReal scalar_increment = 0.0;
170 ADReal scalar_upper_bound = max_permissible_scalar;
171 ADReal scalar_lower_bound = min_permissible_scalar;
177 Real init_resid_sign = MathUtils::sign(MetaPhysicL::raw_value(
_residual));
194 scalar = scalar_old + scalar_increment;
200 min_permissible_scalar,
201 max_permissible_scalar,
210 scalar,
_residual, init_resid_sign, scalar_upper_bound, scalar_lower_bound, iter_output);
219 bool modified_increment =
false;
226 ADReal alpha = residual_old / (residual_old -
_residual);
227 alpha = MathUtils::clamp(alpha, 1.0e-2, 1.0);
231 modified_increment =
true;
232 scalar_increment *= alpha;
234 *iter_output <<
" Line search alpha = " << MetaPhysicL::raw_value(alpha)
235 <<
" increment = " << MetaPhysicL::raw_value(scalar_increment)
245 if (scalar_old + scalar_increment >= scalar_upper_bound ||
246 scalar_old + scalar_increment <= scalar_lower_bound)
248 if (scalar_upper_bound != max_permissible_scalar &&
249 scalar_lower_bound != min_permissible_scalar)
251 const Real frac = 0.5;
253 (1.0 - frac) * scalar_lower_bound + frac * scalar_upper_bound - scalar_old;
254 modified_increment =
true;
256 *iter_output <<
" Trial scalar_increment exceeded bounds. Setting between "
257 "lower/upper bounds. frac: "
258 << frac << std::endl;
265 if (modified_increment)
267 scalar = scalar_old + scalar_increment;
◆ iterationFinalize()
template<ComputeStage compute_stage>
|
inlineprotectedvirtualinherited |
◆ maximumPermissibleValue()
template<ComputeStage compute_stage>
◆ minimumPermissibleValue()
template<ComputeStage compute_stage>
◆ outputIterationStep()
template<ComputeStage compute_stage>
void ADSingleVariableReturnMappingSolution< compute_stage >::outputIterationStep |
( |
std::stringstream * |
iter_output, |
|
|
const ADReal & |
effective_trial_stress, |
|
|
const ADReal & |
scalar, |
|
|
const Real |
reference_residual |
|
) |
| |
|
protectedvirtualinherited |
Output information for a single iteration step to build the convergence history of the model.
- Parameters
-
iter_output | Output stream |
it | Current iteration count |
effective_trial_stress | Effective trial stress |
scalar | Inelastic strain increment magnitude being solved for |
residual | Current value of the residual |
reference | Current value of the reference quantity |
Definition at line 393 of file ADSingleVariableReturnMappingSolution.C.
402 const Real residual = MetaPhysicL::raw_value(
_residual);
404 *iter_output <<
" iteration=" << it
405 <<
" trial_stress=" << MetaPhysicL::raw_value(effective_trial_stress)
406 <<
" scalar=" << MetaPhysicL::raw_value(scalar) <<
" residual=" << residual
407 <<
" ref_res=" << reference_residual
408 <<
" rel_res=" << std::abs(residual) / reference_residual
◆ outputIterationSummary()
template<ComputeStage compute_stage>
|
overrideprotectedvirtualinherited |
◆ propagateQpStatefulProperties()
template<ComputeStage compute_stage>
|
overrideprotectedvirtualinherited |
◆ requiresIsotropicTensor()
template<ComputeStage compute_stage>
|
inlineoverridevirtualinherited |
◆ resetProperties()
template<ComputeStage compute_stage>
◆ resetQpProperties()
template<ComputeStage compute_stage>
Retained as empty methods to avoid a warning from Material.C in framework. These methods are unused in all inheriting classes and should not be overwritten.
Definition at line 103 of file ADStressUpdateBase.h.
◆ returnMappingSolve()
template<ComputeStage compute_stage>
Perform the return mapping iterations.
- Parameters
-
effective_trial_stress | Effective trial stress |
scalar | Inelastic strain increment magnitude being solved for |
console | Console output |
Definition at line 101 of file ADSingleVariableReturnMappingSolution.C.
105 std::unique_ptr<std::stringstream> iter_output =
107 ? libmesh_make_unique<std::stringstream>()
121 throw MooseException(
"");
125 iter_output = libmesh_make_unique<std::stringstream>();
131 *iter_output <<
"Encountered inf or nan in material return mapping iterations.\n";
135 *iter_output <<
"Exceeded maximum iterations in material return mapping iterations.\n";
139 mooseError(
"Unhandled solver state");
145 internalSolve(effective_trial_stress, scalar, iter_output.get());
149 throw MooseException(iter_output->str());
156 console << iter_output->str();
◆ setQp()
template<ComputeStage compute_stage>
Sets the value of the global variable _qp for inheriting classes.
Definition at line 47 of file ADStressUpdateBase.C.
◆ updateBounds()
template<ComputeStage compute_stage>
void ADSingleVariableReturnMappingSolution< compute_stage >::updateBounds |
( |
const ADReal & |
scalar, |
|
|
const ADReal & |
residual, |
|
|
const Real |
init_resid_sign, |
|
|
ADReal & |
scalar_upper_bound, |
|
|
ADReal & |
scalar_lower_bound, |
|
|
std::stringstream * |
iter_output |
|
) |
| |
|
privateinherited |
Update the upper and lower bounds of the root for the effective inelastic strain.
- Parameters
-
scalar | Current value of the inelastic strain increment |
residual | Current value of the residual |
init_resid_sign | Sign of the initial value of the residual |
scalar_upper_bound | Upper bound value of scalar |
scalar_lower_bound | Lower bound value of scalar |
iter_output | Output stream |
Definition at line 365 of file ADSingleVariableReturnMappingSolution.C.
373 if (residual * init_resid_sign < 0.0 && scalar < scalar_upper_bound)
375 scalar_upper_bound = scalar;
376 if (scalar_upper_bound < scalar_lower_bound)
378 scalar_upper_bound = scalar_lower_bound;
379 scalar_lower_bound = 0.0;
381 *iter_output <<
" Corrected for scalar_upper_bound < scalar_lower_bound" << std::endl;
386 else if (residual * init_resid_sign > 0.0 && scalar > scalar_lower_bound &&
387 scalar < scalar_upper_bound)
388 scalar_lower_bound = scalar;
◆ updateIntermediatePorosity()
template<ComputeStage compute_stage>
◆ updateState()
template<ComputeStage compute_stage>
void ADViscoplasticityStressUpdate< compute_stage >::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 |
|
) |
| |
|
overridevirtual |
Given a strain increment that results in a trial stress, perform some procedure (such as an iterative return-mapping process) to produce an admissible stress, an elastic strain increment and an inelastic strain increment.
If _fe_problem.currentlyComputingJacobian() = true, then updateState also computes d(stress)/d(strain) (or some approximation to it).
This method is called by ComputeMultipleInelasticStress. This method is pure virtual: all inheriting classes must overwrite this method.
- Parameters
-
strain_increment | Upon input: the strain increment. Upon output: the elastic strain increment |
inelastic_strain_increment | The inelastic_strain resulting from the iterative procedure |
rotation_increment | The finite-strain rotation increment |
stress_new | Upon input: the trial stress that results from applying strain_increment as an elastic strain. Upon output: the admissible stress |
stress_old | The old value of stress |
elasticity_tensor | The elasticity tensor |
Implements ADStressUpdateBase< compute_stage >.
Definition at line 78 of file ADViscoplasticityStressUpdate.C.
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);
105 inelastic_strain_increment.zero();
109 mooseException(
"In ",
111 ": equivalent stress (",
113 ") is higher than maximum_equivalent_stress (",
115 ").\nCutting time step.");
121 ADReal dpsi_dgauge(0);
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.");
◆ validParams()
template<ComputeStage compute_stage>
Definition at line 20 of file ADViscoplasticityStressUpdate.C.
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");
◆ _absolute_tolerance
template<ComputeStage compute_stage>
◆ _acceptable_multiplier
template<ComputeStage compute_stage>
◆ _base_name
template<ComputeStage compute_stage>
Name used as a prefix for all material properties related to the stress update model.
Definition at line 109 of file ADStressUpdateBase.h.
◆ _bracket_solution
template<ComputeStage compute_stage>
◆ _check_range
template<ComputeStage compute_stage>
◆ _derivative
template<ComputeStage compute_stage>
◆ _dhydro_stress_dsigma
template<ComputeStage compute_stage>
◆ _effective_inelastic_strain_old
template<ComputeStage compute_stage>
◆ _hydro_stress
template<ComputeStage compute_stage>
◆ _identity_two
template<ComputeStage compute_stage>
◆ _inelastic_strain_old
template<ComputeStage compute_stage>
◆ _initial_residual
template<ComputeStage compute_stage>
◆ _intermediate_porosity
template<ComputeStage compute_stage>
◆ _internal_solve_full_iteration_history
template<ComputeStage compute_stage>
◆ _internal_solve_output_on
template<ComputeStage compute_stage>
◆ _iteration
template<ComputeStage compute_stage>
◆ _line_search
template<ComputeStage compute_stage>
◆ _max_inelastic_increment
template<ComputeStage compute_stage>
◆ _max_its
template<ComputeStage compute_stage>
Maximum number of return mapping iterations.
This exists only to avoid an infinite loop, and is is intended to be a large number that is not settable by the user.
Definition at line 172 of file ADSingleVariableReturnMappingSolution.h.
◆ _maximum_equivalent_stress
template<ComputeStage compute_stage>
◆ _maximum_gauge_ratio
template<ComputeStage compute_stage>
◆ _minimum_equivalent_stress
template<ComputeStage compute_stage>
◆ _model
template<ComputeStage compute_stage>
◆ _num_resids
template<ComputeStage compute_stage>
◆ _pore_shape
template<ComputeStage compute_stage>
◆ _pore_shape_factor
template<ComputeStage compute_stage>
◆ _porosity_old
template<ComputeStage compute_stage>
◆ _power
template<ComputeStage compute_stage>
◆ _power_factor
template<ComputeStage compute_stage>
◆ _relative_tolerance
template<ComputeStage compute_stage>
◆ _residual
template<ComputeStage compute_stage>
◆ _residual_history
template<ComputeStage compute_stage>
◆ _svrms_name
template<ComputeStage compute_stage>
◆ _total_strain_base_name
template<ComputeStage compute_stage>
◆ _verbose
template<ComputeStage compute_stage>
◆ usingMaterialMembers
template<ComputeStage compute_stage>
◆ usingSingleVariableReturnMappingSolutionMembers
template<ComputeStage compute_stage>
◆ usingStressUpdateBaseMembers
template<ComputeStage compute_stage>
◆ usingViscoplasticityStressUpdateBaseMembers
template<ComputeStage compute_stage>
The documentation for this class was generated from the following files:
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)
virtual ADReal initialGuess(const ADReal &)
Compute an initial guess for the value of the scalar.
const bool _internal_solve_full_iteration_history
Whether to output iteration information all the time (regardless of whether iterations converge)
PoreShapeModel
Enum to choose which pore shape model to use.
enum ADViscoplasticityStressUpdate::ViscoplasticityModel _model
void returnMappingSolve(const ADReal &effective_trial_stress, ADReal &scalar, const ConsoleStream &console)
Perform the return mapping iterations.
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
bool _bracket_solution
Whether to save upper and lower bounds of root for scalar, and set solution to the midpoint between t...
void updateIntermediatePorosity(const ADRankTwoTensor &elastic_strain_increment)
const RankTwoTensor _identity_two
Rank two identity tensor.
virtual void outputIterationStep(std::stringstream *iter_output, const ADReal &effective_trial_stress, const ADReal &scalar, const Real reference_residual)
Output information for a single iteration step to build the convergence history of the model.
void computeInelasticStrainIncrement(ADReal &gauge_stress, ADReal &dpsi_dgauge, ADRankTwoTensor &creep_strain_increment, const ADReal &equiv_stress, const ADRankTwoTensor &dev_stress, const ADRankTwoTensor &stress)
static InputParameters validParams()
unsigned int _iteration
iteration number
const MaterialProperty< Real > & _porosity_old
Material property for the old porosity.
const Real _minimum_equivalent_stress
Minimum value of equivalent stress below which viscoplasticiy is not calculated.
ADReal _initial_residual
Residual values, kept as members to retain solver state for summary outputting.
ADReal _intermediate_porosity
Container for the porosity calculated from all other intelastic models except the current model.
bool convergedAcceptable(const unsigned int it, const Real reference)
Check to see whether the residual is within acceptable convergence limits.
enum ADViscoplasticityStressUpdate::PoreShapeModel _pore_shape
std::vector< Real > _residual_history
History of residuals used to check whether progress is still being made on decreasing the residual.
virtual ADReal maximumPermissibleValue(const ADReal &effective_trial_stress) const
Compute the maximum permissible value of the scalar.
ADReal _hydro_stress
Container for hydrostatic stress.
Real _max_inelastic_increment
Max increment for inelastic strain.
Real _relative_tolerance
Relative convergence tolerance.
virtual void computeStressInitialize(const ADReal &, const ADRankFourTensor &)
Perform any necessary initialization before return mapping iterations.
const MaterialProperty< RankTwoTensor > & _inelastic_strain_old
ViscoplasticityModel
Enum to choose which viscoplastic model to use.
enum ADSingleVariableReturnMappingSolution::InternalSolveOutput _internal_solve_output_on
const RankTwoTensor _dhydro_stress_dsigma
Derivative of hydrostatic stress with respect to the stress tensor.
Real _acceptable_multiplier
Multiplier applied to relative and absolute tolerances for acceptable convergence.
bool converged(const ADReal &residual, const Real reference)
Check to see whether the residual is within the convergence limits.
void updateBounds(const ADReal &scalar, const ADReal &residual, const Real init_resid_sign, ADReal &scalar_upper_bound, ADReal &scalar_lower_bound, std::stringstream *iter_output)
Update the upper and lower bounds of the root for the effective inelastic strain.
const Real _maximum_equivalent_stress
Maximum value of equivalent stress above which an exception is thrown.
const Real _power_factor
Power factor used for LPS model.
const Real _maximum_gauge_ratio
Maximum ratio between the gauge stress and the equilvalent stress.
const Real _power
Exponent on the effective stress.
const unsigned int _max_its
Maximum number of return mapping iterations.
virtual void computeStressFinalize(const ADRankTwoTensor &)
Perform any necessary steps to finalize state after return mapping iterations.
const MaterialProperty< Real > & _effective_inelastic_strain_old
void checkPermissibleRange(ADReal &scalar, ADReal &scalar_increment, const ADReal &scalar_old, const ADReal min_permissible_scalar, const ADReal max_permissible_scalar, std::stringstream *iter_output)
Check to see whether solution is within admissible range, and set it within that range if it is not.
virtual void iterationFinalize(ADReal)
Finalize internal state variables for a model for a given iteration.
virtual ADReal minimumPermissibleValue(const ADReal &effective_trial_stress) const
Compute the minimum permissible value of the scalar.
bool _check_range
Whether to check to see whether iterative solution is within admissible range, and set within that ra...
const std::size_t _num_resids
Number of residuals to be stored in history.
const std::string _base_name
Name used as a prefix for all material properties related to the stress update model.
ADReal _derivative
Container for _derivative.
const bool _verbose
Flag to enable verbose output.
Real _absolute_tolerance
Absolute convergence tolerance.
bool _line_search
Whether to use line searches to improve convergence.
ADRankTwoTensor computeDGaugeDSigma(const ADReal &gauge_stress, const ADReal &equiv_stress, const ADRankTwoTensor &dev_stress, const ADRankTwoTensor &stress)
virtual ADReal computeResidual(const ADReal &effective_trial_stress, const ADReal &scalar)=0
Compute the residual for a predicted value of the scalar.
virtual void outputIterationSummary(std::stringstream *iter_output, const unsigned int total_it)
Output summary information for the convergence history of the model.
SolveState internalSolve(const ADReal effective_trial_stress, ADReal &scalar, std::stringstream *iter_output=nullptr)
Method called from within this class to perform the actual return mappping iterations.
virtual ADReal computeDerivative(const ADReal &effective_trial_stress, const ADReal &scalar)=0
Compute the derivative of the residual as a function of the scalar variable.
virtual Real computeReferenceResidual(const ADReal &effective_trial_stress, const ADReal &scalar)=0
Compute a reference quantity to be used for checking relative convergence.