www.mooseframework.org
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Types | Private Member Functions | Private Attributes | List of all members
RadialReturnStressUpdate Class Referenceabstract

RadialReturnStressUpdate computes the radial return stress increment for an isotropic elastic-viscoplasticity model after interating on the difference between new and old trial stress increments. More...

#include <RadialReturnStressUpdate.h>

Inheritance diagram for RadialReturnStressUpdate:
[legend]

Public Member Functions

 RadialReturnStressUpdate (const InputParameters &parameters)
 
virtual void updateState (RankTwoTensor &strain_increment, RankTwoTensor &inelastic_strain_increment, const RankTwoTensor &rotation_increment, RankTwoTensor &stress_new, const RankTwoTensor &stress_old, const RankFourTensor &elasticity_tensor, const RankTwoTensor &elastic_strain_old, bool compute_full_tangent_operator, RankFourTensor &tangent_operator) override
 A radial return (J2) mapping method is performed with return mapping iterations. More...
 
virtual Real computeReferenceResidual (const Real effective_trial_stress, const Real scalar_effective_inelastic_strain) override
 Compute a reference quantity to be used for checking relative convergence. More...
 
virtual Real minimumPermissibleValue (const Real) const override
 Compute the minimum permissible value of the scalar. More...
 
virtual Real maximumPermissibleValue (const Real effective_trial_stress) const override
 Compute the maximum permissible value of the scalar. 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...
 
bool isIsotropic () override
 Radial return mapped models should be isotropic by default! More...
 
void setQp (unsigned int qp)
 Sets the value of the global variable _qp for inheriting classes. More...
 
virtual void propagateQpStatefulProperties ()
 If updateState is not called during a timestep, this will be. More...
 
virtual TangentCalculationMethod getTangentCalculationMethod ()
 
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
 

Static Public Member Functions

static InputParameters validParams ()
 

Protected Member Functions

virtual void initQpStatefulProperties () override
 
void propagateQpStatefulPropertiesRadialReturn ()
 Propagate the properties pertaining to this intermediate class. More...
 
virtual void computeStressInitialize (const Real, const RankFourTensor &)
 Perform any necessary initialization before return mapping iterations. More...
 
virtual Real computeStressDerivative (const Real, const Real)
 Calculate the derivative of the strain increment with respect to the updated stress. More...
 
virtual void computeStressFinalize (const RankTwoTensor &)
 Perform any necessary steps to finalize state after return mapping iterations. More...
 
void outputIterationSummary (std::stringstream *iter_output, const unsigned int total_it) override
 Output summary information for the convergence history of the model. More...
 
void returnMappingSolve (const Real effective_trial_stress, Real &scalar, const ConsoleStream &console)
 Perform the return mapping iterations. More...
 
virtual Real initialGuess (const Real)
 Compute an initial guess for the value of the scalar. More...
 
virtual Real computeResidual (const Real effective_trial_stress, const Real scalar)=0
 Compute the residual for a predicted value of the scalar. More...
 
virtual Real computeDerivative (const Real effective_trial_stress, const Real scalar)=0
 Compute the derivative of the residual as a function of the scalar variable. More...
 
virtual void iterationFinalize (Real)
 Finalize internal state variables for a model for a given iteration. More...
 
virtual void outputIterationStep (std::stringstream *iter_output, const unsigned int it, const Real effective_trial_stress, const Real scalar, const Real residual, const Real reference_residual)
 Output information for a single iteration step to build the convergence history of the model. More...
 

Protected Attributes

Real _three_shear_modulus
 3 * shear modulus More...
 
MaterialProperty< Real > & _effective_inelastic_strain
 
const MaterialProperty< Real > & _effective_inelastic_strain_old
 
Real _max_inelastic_increment
 
const RankTwoTensor _identity_two
 Rank two identity tensor. More...
 
const RankFourTensor _identity_symmetric_four
 Rank four symmetric identity tensor. More...
 
const RankFourTensor _deviatoric_projection_four
 Rank four deviatoric projection tensor. More...
 
const std::string _base_name
 Name used as a prefix for all material properties related to the stress update model. More...
 
bool _check_range
 Whether to check to see whether iterative solution is within admissible range, and set within that range if outside. More...
 
bool _line_search
 Whether to use line searches to improve convergence. More...
 
bool _bracket_solution
 Whether to save upper and lower bounds of root for scalar, and set solution to the midpoint between those bounds if outside them. More...
 

Private Types

enum  InternalSolveOutput { InternalSolveOutput::NEVER, InternalSolveOutput::ON_ERROR, InternalSolveOutput::ALWAYS }
 
enum  SolveState { SolveState::SUCCESS, SolveState::NAN_INF, SolveState::EXCEEDED_ITERATIONS }
 

Private Member Functions

SolveState internalSolve (const Real effective_trial_stress, Real &scalar, std::stringstream *iter_output=nullptr)
 Method called from within this class to perform the actual return mappping iterations. More...
 
bool converged (const Real residual, const Real reference)
 Check to see whether the residual is within the convergence limits. More...
 
bool convergedAcceptable (const unsigned int it, const Real residual, const Real reference)
 Check to see whether the residual is within acceptable convergence limits. More...
 
void checkPermissibleRange (Real &scalar, Real &scalar_increment, const Real scalar_old, const Real min_permissible_scalar, const Real 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 Real scalar, const Real residual, const Real init_resid_sign, Real &scalar_upper_bound, Real &scalar_lower_bound, std::stringstream *iter_output)
 Update the upper and lower bounds of the root for the effective inelastic strain. More...
 

Private Attributes

enum SingleVariableReturnMappingSolution::InternalSolveOutput _internal_solve_output_on
 
const unsigned int _max_its
 Maximum number of return mapping iterations. More...
 
const bool _internal_solve_full_iteration_history
 Whether to output iteration information all the time (regardless of whether iterations converge) More...
 
Real _relative_tolerance
 Relative convergence tolerance. More...
 
Real _absolute_tolerance
 Absolute convergence tolerance. More...
 
Real _acceptable_multiplier
 Multiplier applied to relative and absolute tolerances for acceptable convergence. More...
 
const std::size_t _num_resids
 Number of residuals to be stored in history. More...
 
std::vector< Real > _residual_history
 History of residuals used to check whether progress is still being made on decreasing the residual. More...
 
unsigned int _iteration
 iteration number More...
 
const std::string _svrms_name
 MOOSE input name of the object performing the solve. More...
 
Real _initial_residual
 Residual values, kept as members to retain solver state for summary outputting. More...
 
Real _residual
 

Detailed Description

RadialReturnStressUpdate computes the radial return stress increment for an isotropic elastic-viscoplasticity model after interating on the difference between new and old trial stress increments.

This radial return mapping class acts as a base class for the radial return creep and plasticity classes / combinations. The stress increment computed by RadialReturnStressUpdate is used by ComputeMultipleInelasticStress which computes the elastic stress for finite strains. This return mapping class is acceptable for finite strains but not total strains. This class is based on the Elasto-viscoplasticity algorithm in F. Dunne and N. Petrinic's Introduction to Computational Plasticity (2004) Oxford University Press.

Definition at line 34 of file RadialReturnStressUpdate.h.

Member Enumeration Documentation

◆ InternalSolveOutput

Enumerator
NEVER 
ON_ERROR 
ALWAYS 

Definition at line 130 of file SingleVariableReturnMappingSolution.h.

131  {
132  NEVER,
133  ON_ERROR,
134  ALWAYS

◆ SolveState

Enumerator
SUCCESS 
NAN_INF 
EXCEEDED_ITERATIONS 

Definition at line 137 of file SingleVariableReturnMappingSolution.h.

138  {
139  SUCCESS,
140  NAN_INF,
141  EXCEEDED_ITERATIONS
142  };

Constructor & Destructor Documentation

◆ RadialReturnStressUpdate()

RadialReturnStressUpdate::RadialReturnStressUpdate ( const InputParameters &  parameters)

Definition at line 37 of file RadialReturnStressUpdate.C.

38  : StressUpdateBase(parameters),
40  _effective_inelastic_strain(declareProperty<Real>(
41  _base_name + getParam<std::string>("effective_inelastic_strain_name"))),
42  _effective_inelastic_strain_old(getMaterialPropertyOld<Real>(
43  _base_name + getParam<std::string>("effective_inelastic_strain_name"))),
44  _max_inelastic_increment(parameters.get<Real>("max_inelastic_increment")),
45  _identity_two(RankTwoTensor::initIdentity),
46  _identity_symmetric_four(RankFourTensor::initIdentitySymmetricFour),
48  _identity_two.outerProduct(_identity_two) / 3.0)
49 {
50 }

Member Function Documentation

◆ checkPermissibleRange()

void SingleVariableReturnMappingSolution::checkPermissibleRange ( Real &  scalar,
Real &  scalar_increment,
const Real  scalar_old,
const Real  min_permissible_scalar,
const Real  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
scalarCurrent value of the inelastic strain increment
scalar_incrementIncremental change in scalar from the previous iteration
scalar_oldPrevious value of scalar
min_permissible_scalarMinimum permissible value of scalar
max_permissible_scalarMaximum permissible value of scalar
iter_outputOutput stream

Definition at line 356 of file SingleVariableReturnMappingSolution.C.

362 {
363  if (scalar > max_permissible_scalar)
364  {
365  scalar_increment = (max_permissible_scalar - scalar_old) / 2.0;
366  scalar = scalar_old + scalar_increment;
367  if (iter_output)
368  *iter_output << "Scalar greater than maximum (" << max_permissible_scalar
369  << ") adjusted scalar=" << scalar << " scalar_increment=" << scalar_increment
370  << std::endl;
371  }
372  else if (scalar < min_permissible_scalar)
373  {
374  scalar_increment = (min_permissible_scalar - scalar_old) / 2.0;
375  scalar = scalar_old + scalar_increment;
376  if (iter_output)
377  *iter_output << "Scalar less than minimum (" << min_permissible_scalar
378  << ") adjusted scalar=" << scalar << " scalar_increment=" << scalar_increment
379  << std::endl;
380  }
381 }

Referenced by SingleVariableReturnMappingSolution::internalSolve().

◆ computeDerivative()

virtual Real SingleVariableReturnMappingSolution::computeDerivative ( const Real  effective_trial_stress,
const Real  scalar 
)
protectedpure virtualinherited

Compute the derivative of the residual as a function of the scalar variable.

The residual should be in strain increment units for all models for consistency.

Parameters
effective_trial_stressEffective trial stress
scalarInelastic strain increment magnitude being solved for

Implemented in IsotropicPlasticity, CLSHPlasticModel, PowerLawCreepModel, IsotropicPlasticityStressUpdate, PowerLawCreepStressUpdate, and PowerLawCreepExceptionTest.

Referenced by RadialReturnCreepStressUpdateBase::computeStressDerivative(), and SingleVariableReturnMappingSolution::internalSolve().

◆ computeReferenceResidual()

Real RadialReturnStressUpdate::computeReferenceResidual ( const Real  effective_trial_stress,
const Real  scalar 
)
overridevirtual

Compute a reference quantity to be used for checking relative convergence.

This should be in strain increment units for all models for consistency.

Parameters
effective_trial_stressEffective trial stress
scalarInelastic strain increment magnitude being solved for

Implements SingleVariableReturnMappingSolution.

Definition at line 150 of file RadialReturnStressUpdate.C.

152 {
153  return effective_trial_stress / _three_shear_modulus - scalar_effective_inelastic_strain;
154 }

◆ computeResidual()

virtual Real SingleVariableReturnMappingSolution::computeResidual ( const Real  effective_trial_stress,
const Real  scalar 
)
protectedpure virtualinherited

Compute the residual for a predicted value of the scalar.

This residual should be in strain increment units for all models for consistency.

Parameters
effective_trial_stressEffective trial stress
scalarInelastic strain increment magnitude being solved for

Implemented in IsotropicPlasticity, CLSHPlasticModel, PowerLawCreepModel, IsotropicPlasticityStressUpdate, PowerLawCreepStressUpdate, and PowerLawCreepExceptionTest.

Referenced by SingleVariableReturnMappingSolution::internalSolve().

◆ computeStressDerivative()

virtual Real RadialReturnStressUpdate::computeStressDerivative ( const  Real,
const  Real 
)
inlineprotectedvirtual

Calculate the derivative of the strain increment with respect to the updated stress.

Parameters
effective_trial_stressEffective trial stress
scalarInelastic strain increment magnitude being solved for

Reimplemented in RadialReturnCreepStressUpdateBase.

Definition at line 118 of file RadialReturnStressUpdate.h.

119  {
120  return 0.0;
121  }

Referenced by updateState().

◆ computeStressFinalize()

virtual void RadialReturnStressUpdate::computeStressFinalize ( const RankTwoTensor )
inlineprotectedvirtual

Perform any necessary steps to finalize state after return mapping iterations.

Parameters
inelasticStrainIncrementInelastic strain increment

Reimplemented in IsotropicPlasticityStressUpdate, and RadialReturnCreepStressUpdateBase.

Definition at line 127 of file RadialReturnStressUpdate.h.

127 {}

Referenced by updateState().

◆ computeStressInitialize()

virtual void RadialReturnStressUpdate::computeStressInitialize ( const  Real,
const RankFourTensor  
)
inlineprotectedvirtual

Perform any necessary initialization before return mapping iterations.

Parameters
effective_trial_stressEffective trial stress
elasticityTensorElasticity tensor

Reimplemented in TemperatureDependentHardeningStressUpdate, IsotropicPlasticityStressUpdate, IsotropicPowerLawHardeningStressUpdate, and PowerLawCreepStressUpdate.

Definition at line 108 of file RadialReturnStressUpdate.h.

110  {
111  }

Referenced by updateState().

◆ computeTimeStepLimit()

Real RadialReturnStressUpdate::computeTimeStepLimit ( )
overridevirtual

Compute the limiting value of the time step for this material.

Returns
Limiting time step

Reimplemented from StressUpdateBase.

Definition at line 163 of file RadialReturnStressUpdate.C.

164 {
165  Real scalar_inelastic_strain_incr;
166 
167  scalar_inelastic_strain_incr =
169  if (MooseUtils::absoluteFuzzyEqual(scalar_inelastic_strain_incr, 0.0))
170  return std::numeric_limits<Real>::max();
171 
172  return _dt * _max_inelastic_increment / scalar_inelastic_strain_incr;
173 }

◆ converged()

bool SingleVariableReturnMappingSolution::converged ( const Real  residual,
const Real  reference 
)
privateinherited

Check to see whether the residual is within the convergence limits.

Parameters
residualCurrent value of the residual
referenceCurrent value of the reference quantity
Returns
Whether the model converged

Definition at line 298 of file SingleVariableReturnMappingSolution.C.

299 {
300  return (std::abs(residual) <= _absolute_tolerance ||
301  std::abs(residual / reference) <= _relative_tolerance);
302 }

Referenced by SingleVariableReturnMappingSolution::convergedAcceptable(), and SingleVariableReturnMappingSolution::internalSolve().

◆ convergedAcceptable()

bool SingleVariableReturnMappingSolution::convergedAcceptable ( const unsigned int  it,
const Real  residual,
const Real  reference 
)
privateinherited

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
residualCurrent iteration count
residualCurrent value of the residual
referenceCurrent value of the reference quantity
Returns
Whether the model converged

Definition at line 305 of file SingleVariableReturnMappingSolution.C.

308 {
309  // Require that we have at least done _num_resids evaluations before we allow for
310  // acceptable convergence
311  if (it < _num_resids)
312  return false;
313 
314  // Check to see whether the residual has dropped by convergence_history_factor over
315  // the last _num_resids iterations. If it has (which means it's still making progress),
316  // don't consider it to be converged within the acceptable limits.
317  const Real convergence_history_factor = 10.0;
318  if (std::abs(residual * convergence_history_factor) <
319  std::abs(_residual_history[(it + 1) % _num_resids]))
320  return false;
321 
322  // Now that it's determined that progress is not being made, treat it as converged if
323  // we're within the acceptable convergence limits
324  return converged(residual / _acceptable_multiplier, reference);
325 }

Referenced by SingleVariableReturnMappingSolution::internalSolve().

◆ getTangentCalculationMethod()

virtual TangentCalculationMethod StressUpdateBase::getTangentCalculationMethod ( )
inlinevirtualinherited

Reimplemented in MultiParameterPlasticityStressUpdate, and RadialReturnCreepStressUpdateBase.

Definition at line 116 of file StressUpdateBase.h.

117  {
119  }

Referenced by updateState().

◆ initialGuess()

virtual Real SingleVariableReturnMappingSolution::initialGuess ( const  Real)
inlineprotectedvirtualinherited

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_stressEffective trial stress

Definition at line 64 of file SingleVariableReturnMappingSolution.h.

64 { return 0.0; }

Referenced by SingleVariableReturnMappingSolution::internalSolve().

◆ initQpStatefulProperties()

void RadialReturnStressUpdate::initQpStatefulProperties ( )
overrideprotectedvirtual

◆ internalSolve()

SingleVariableReturnMappingSolution::SolveState SingleVariableReturnMappingSolution::internalSolve ( const Real  effective_trial_stress,
Real &  scalar,
std::stringstream *  iter_output = nullptr 
)
privateinherited

Method called from within this class to perform the actual return mappping iterations.

Parameters
effective_trial_stressEffective trial stress
scalarInelastic strain increment magnitude being solved for
iter_outputOutput stream – if null, no output is produced
Returns
Whether the solution was successful

Definition at line 157 of file SingleVariableReturnMappingSolution.C.

160 {
161  scalar = initialGuess(effective_trial_stress);
162  Real scalar_old = scalar;
163  Real scalar_increment = 0.0;
164  const Real min_permissible_scalar = minimumPermissibleValue(effective_trial_stress);
165  const Real max_permissible_scalar = maximumPermissibleValue(effective_trial_stress);
166  Real scalar_upper_bound = max_permissible_scalar;
167  Real scalar_lower_bound = min_permissible_scalar;
168  _iteration = 0;
169 
170  _initial_residual = _residual = computeResidual(effective_trial_stress, scalar);
171 
172  Real residual_old = _residual;
173  Real init_resid_sign = MathUtils::sign(_residual);
174  Real reference_residual = computeReferenceResidual(effective_trial_stress, scalar);
175 
176  if (converged(_residual, reference_residual))
177  {
178  iterationFinalize(scalar);
180  iter_output, _iteration, effective_trial_stress, scalar, _residual, reference_residual);
181  return SolveState::SUCCESS;
182  }
183 
184  _residual_history.assign(_num_resids, std::numeric_limits<Real>::max());
186 
187  while (_iteration < _max_its && !converged(_residual, reference_residual) &&
188  !convergedAcceptable(_iteration, _residual, reference_residual))
189  {
190  scalar_increment = -_residual / computeDerivative(effective_trial_stress, scalar);
191  scalar = scalar_old + scalar_increment;
192 
193  if (_check_range)
194  checkPermissibleRange(scalar,
195  scalar_increment,
196  scalar_old,
197  min_permissible_scalar,
198  max_permissible_scalar,
199  iter_output);
200 
201  _residual = computeResidual(effective_trial_stress, scalar);
202  reference_residual = computeReferenceResidual(effective_trial_stress, scalar);
203  iterationFinalize(scalar);
204 
205  if (_bracket_solution)
206  updateBounds(
207  scalar, _residual, init_resid_sign, scalar_upper_bound, scalar_lower_bound, iter_output);
208 
209  if (converged(_residual, reference_residual))
210  {
212  iter_output, _iteration, effective_trial_stress, scalar, _residual, reference_residual);
213  break;
214  }
215  else
216  {
217  bool modified_increment = false;
218 
219  // Line Search
220  if (_line_search)
221  {
222  if (residual_old - _residual != 0.0)
223  {
224  Real alpha = residual_old / (residual_old - _residual);
225  alpha = MathUtils::clamp(alpha, 1.0e-2, 1.0);
226 
227  if (alpha != 1.0)
228  {
229  modified_increment = true;
230  scalar_increment *= alpha;
231  if (iter_output)
232  *iter_output << " Line search alpha = " << alpha
233  << " increment = " << scalar_increment << std::endl;
234  }
235  }
236  }
237 
238  if (_bracket_solution)
239  {
240  // Check to see whether trial scalar_increment is outside the bounds, and set it to a point
241  // within the bounds if it is
242  if (scalar_old + scalar_increment >= scalar_upper_bound ||
243  scalar_old + scalar_increment <= scalar_lower_bound)
244  {
245  if (scalar_upper_bound != max_permissible_scalar &&
246  scalar_lower_bound != min_permissible_scalar)
247  {
248  const Real frac = 0.5;
249  scalar_increment =
250  (1.0 - frac) * scalar_lower_bound + frac * scalar_upper_bound - scalar_old;
251  modified_increment = true;
252  if (iter_output)
253  *iter_output << " Trial scalar_increment exceeded bounds. Setting between "
254  "lower/upper bounds. frac: "
255  << frac << std::endl;
256  }
257  }
258  }
259 
260  // Update the trial scalar and recompute residual if the line search or bounds checking
261  // modified the increment
262  if (modified_increment)
263  {
264  scalar = scalar_old + scalar_increment;
265  _residual = computeResidual(effective_trial_stress, scalar);
266  reference_residual = computeReferenceResidual(effective_trial_stress, scalar);
267  iterationFinalize(scalar);
268 
269  if (_bracket_solution)
270  updateBounds(scalar,
271  _residual,
272  init_resid_sign,
273  scalar_upper_bound,
274  scalar_lower_bound,
275  iter_output);
276  }
277  }
278 
280  iter_output, _iteration, effective_trial_stress, scalar, _residual, reference_residual);
281 
282  ++_iteration;
283  residual_old = _residual;
284  scalar_old = scalar;
286  }
287 
288  if (std::isnan(_residual) || std::isinf(_residual))
289  return SolveState::NAN_INF;
290 
291  if (_iteration == _max_its)
293 
294  return SolveState::SUCCESS;
295 }

Referenced by SingleVariableReturnMappingSolution::returnMappingSolve().

◆ isIsotropic()

bool RadialReturnStressUpdate::isIsotropic ( )
inlineoverridevirtual

Radial return mapped models should be isotropic by default!

Reimplemented from StressUpdateBase.

Definition at line 88 of file RadialReturnStressUpdate.h.

88 { return true; };

◆ iterationFinalize()

virtual void SingleVariableReturnMappingSolution::iterationFinalize ( Real  )
inlineprotectedvirtualinherited

Finalize internal state variables for a model for a given iteration.

Parameters
scalarInelastic strain increment magnitude being solved for

Reimplemented in IsotropicPlasticityStressUpdate, IsotropicPlasticity, and CLSHPlasticModel.

Definition at line 94 of file SingleVariableReturnMappingSolution.h.

94 {}

Referenced by SingleVariableReturnMappingSolution::internalSolve().

◆ maximumPermissibleValue()

Real RadialReturnStressUpdate::maximumPermissibleValue ( const Real  effective_trial_stress) const
overridevirtual

Compute the maximum permissible value of the scalar.

For some models, the magnitude of this may be known.

Parameters
effective_trial_stressEffective trial stress

Reimplemented from SingleVariableReturnMappingSolution.

Definition at line 157 of file RadialReturnStressUpdate.C.

158 {
159  return effective_trial_stress / _three_shear_modulus;
160 }

◆ minimumPermissibleValue()

virtual Real RadialReturnStressUpdate::minimumPermissibleValue ( const  effective_trial_stress) const
inlineoverridevirtual

Compute the minimum permissible value of the scalar.

For some models, the magnitude of this may be known.

Parameters
effective_trial_stressEffective trial stress

Reimplemented from SingleVariableReturnMappingSolution.

Definition at line 67 of file RadialReturnStressUpdate.h.

68  {
69  return 0.0;
70  }

◆ outputIterationStep()

void SingleVariableReturnMappingSolution::outputIterationStep ( std::stringstream *  iter_output,
const unsigned int  it,
const Real  effective_trial_stress,
const Real  scalar,
const Real  residual,
const Real  reference_residual 
)
protectedvirtualinherited

Output information for a single iteration step to build the convergence history of the model.

Parameters
iter_outputOutput stream
itCurrent iteration count
effective_trial_stressEffective trial stress
scalarInelastic strain increment magnitude being solved for
residualCurrent value of the residual
referenceCurrent value of the reference quantity

Definition at line 328 of file SingleVariableReturnMappingSolution.C.

334 {
335  if (iter_output)
336  {
337  *iter_output << " iteration=" << it << " trial_stress=" << effective_trial_stress
338  << " scalar=" << scalar << " residual=" << residual
339  << " ref_res=" << reference_residual
340  << " rel_res=" << std::abs(residual) / reference_residual
341  << " rel_tol=" << _relative_tolerance << " abs_res=" << std::abs(residual)
342  << " abs_tol=" << _absolute_tolerance << '\n';
343  }
344 }

Referenced by SingleVariableReturnMappingSolution::internalSolve().

◆ outputIterationSummary()

void RadialReturnStressUpdate::outputIterationSummary ( std::stringstream *  iter_output,
const unsigned int  total_it 
)
overrideprotectedvirtual

Output summary information for the convergence history of the model.

Parameters
iter_outputOutput stream
total_itTotal iteration count

Reimplemented from SingleVariableReturnMappingSolution.

Definition at line 176 of file RadialReturnStressUpdate.C.

178 {
179  if (iter_output)
180  {
181  *iter_output << "At element " << _current_elem->id() << " _qp=" << _qp << " Coordinates "
182  << _q_point[_qp] << " block=" << _current_elem->subdomain_id() << '\n';
183  }
185 }

◆ propagateQpStatefulProperties()

void StressUpdateBase::propagateQpStatefulProperties ( )
virtualinherited

If updateState is not called during a timestep, this will be.

This method allows derived classes to set internal parameters from their Old values, for instance

Reimplemented in MultiParameterPlasticityStressUpdate, LinearViscoelasticStressUpdate, IsotropicPlasticityStressUpdate, and RadialReturnCreepStressUpdateBase.

Definition at line 49 of file StressUpdateBase.C.

50 {
51  mooseError(
52  "propagateQpStatefulProperties called: it needs to be implemented by your inelastic model");
53 }

◆ propagateQpStatefulPropertiesRadialReturn()

void RadialReturnStressUpdate::propagateQpStatefulPropertiesRadialReturn ( )
protected

Propagate the properties pertaining to this intermediate class.

This is intended to be called from propagateQpStatefulProperties() in classes that inherit from this one. This is intentionally named uniquely because almost all models that derive from this class have their own stateful properties, and this forces them to define their own implementations of propagateQpStatefulProperties().

Definition at line 59 of file RadialReturnStressUpdate.C.

Referenced by RadialReturnCreepStressUpdateBase::propagateQpStatefulProperties(), and IsotropicPlasticityStressUpdate::propagateQpStatefulProperties().

◆ requiresIsotropicTensor()

bool RadialReturnStressUpdate::requiresIsotropicTensor ( )
inlineoverridevirtual

Does the model require the elasticity tensor to be isotropic?

Implements StressUpdateBase.

Definition at line 83 of file RadialReturnStressUpdate.h.

83 { return true; }

◆ resetProperties()

void StressUpdateBase::resetProperties ( )
inlinefinalinherited

Definition at line 123 of file StressUpdateBase.h.

123 {}

◆ resetQpProperties()

void StressUpdateBase::resetQpProperties ( )
inlinefinalinherited

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 122 of file StressUpdateBase.h.

122 {}

◆ returnMappingSolve()

void SingleVariableReturnMappingSolution::returnMappingSolve ( const Real  effective_trial_stress,
Real &  scalar,
const ConsoleStream &  console 
)
protectedinherited

Perform the return mapping iterations.

Parameters
effective_trial_stressEffective trial stress
scalarInelastic strain increment magnitude being solved for
consoleConsole output

Definition at line 96 of file SingleVariableReturnMappingSolution.C.

99 {
100  // construct the stringstream here only if the debug level is set to ALL
101  std::unique_ptr<std::stringstream> iter_output =
103  ? libmesh_make_unique<std::stringstream>()
104  : nullptr;
105 
106  // do the internal solve and capture iteration info during the first round
107  // iff full history output is requested regardless of whether the solve failed or succeeded
108  auto solve_state =
109  internalSolve(effective_trial_stress,
110  scalar,
111  _internal_solve_full_iteration_history ? iter_output.get() : nullptr);
112  if (solve_state != SolveState::SUCCESS &&
114  {
115  // output suppressed by user, throw immediately
117  throw MooseException("");
118 
119  // user expects some kind of output, if necessary setup output stream now
120  if (!iter_output)
121  iter_output = libmesh_make_unique<std::stringstream>();
122 
123  // add the appropriate error message to the output
124  switch (solve_state)
125  {
126  case SolveState::NAN_INF:
127  *iter_output << "Encountered inf or nan in material return mapping iterations.\n";
128  break;
129 
131  *iter_output << "Exceeded maximum iterations in material return mapping iterations.\n";
132  break;
133 
134  default:
135  mooseError("Unhandled solver state");
136  }
137 
138  // if full history output is only requested for failed solves we have to repeat
139  // the solve a second time
141  internalSolve(effective_trial_stress, scalar, iter_output.get());
142 
143  // Append summary and throw exception
144  outputIterationSummary(iter_output.get(), _iteration);
145  throw MooseException(iter_output->str());
146  }
147 
149  {
150  // the solve did not fail but the user requested debug output anyways
151  outputIterationSummary(iter_output.get(), _iteration);
152  console << iter_output->str();
153  }
154 }

Referenced by ReturnMappingModel::computeStress(), and updateState().

◆ setQp()

void StressUpdateBase::setQp ( unsigned int  qp)
inherited

Sets the value of the global variable _qp for inheriting classes.

Definition at line 43 of file StressUpdateBase.C.

44 {
45  _qp = qp;
46 }

◆ updateBounds()

void SingleVariableReturnMappingSolution::updateBounds ( const Real  scalar,
const Real  residual,
const Real  init_resid_sign,
Real &  scalar_upper_bound,
Real &  scalar_lower_bound,
std::stringstream *  iter_output 
)
privateinherited

Update the upper and lower bounds of the root for the effective inelastic strain.

Parameters
scalarCurrent value of the inelastic strain increment
residualCurrent value of the residual
init_resid_signSign of the initial value of the residual
scalar_upper_boundUpper bound value of scalar
scalar_lower_boundLower bound value of scalar
iter_outputOutput stream

Definition at line 384 of file SingleVariableReturnMappingSolution.C.

390 {
391  // Update upper/lower bounds as applicable
392  if (residual * init_resid_sign < 0.0 && scalar < scalar_upper_bound)
393  {
394  scalar_upper_bound = scalar;
395  if (scalar_upper_bound < scalar_lower_bound)
396  {
397  scalar_upper_bound = scalar_lower_bound;
398  scalar_lower_bound = 0.0;
399  if (iter_output)
400  *iter_output << " Corrected for scalar_upper_bound < scalar_lower_bound" << std::endl;
401  }
402  }
403  // Don't permit setting scalar_lower_bound > scalar_upper_bound (but do permit the reverse).
404  // This ensures that if we encounter multiple roots, we pick the lowest one.
405  else if (residual * init_resid_sign > 0.0 && scalar > scalar_lower_bound &&
406  scalar < scalar_upper_bound)
407  scalar_lower_bound = scalar;
408 }

Referenced by SingleVariableReturnMappingSolution::internalSolve().

◆ updateState()

void RadialReturnStressUpdate::updateState ( RankTwoTensor strain_increment,
RankTwoTensor inelastic_strain_increment,
const RankTwoTensor rotation_increment,
RankTwoTensor stress_new,
const RankTwoTensor stress_old,
const RankFourTensor elasticity_tensor,
const RankTwoTensor elastic_strain_old,
bool  compute_full_tangent_operator,
RankFourTensor tangent_operator 
)
overridevirtual

A radial return (J2) mapping method is performed with return mapping iterations.

Parameters
strain_incrementSum of elastic and inelastic strain increments
inelastic_strain_incrementInelastic strain increment calculated by this class
rotationincrement Not used by this class
stress_newNew trial stress from pure elastic calculation
stress_oldOld state of stress
elasticity_tensorRank 4 C_{ijkl}, must be isotropic
elastic_strain_oldOld state of total elastic strain
compute_full_tangent_operatorFlag currently unused by this class
tangent_operatorCurrently a copy of the elasticity tensor in this class

Implements StressUpdateBase.

Definition at line 65 of file RadialReturnStressUpdate.C.

74 {
75  // compute the deviatoric trial stress and trial strain from the current intermediate
76  // configuration
77  RankTwoTensor deviatoric_trial_stress = stress_new.deviatoric();
78 
79  // compute the effective trial stress
80  Real dev_trial_stress_squared =
81  deviatoric_trial_stress.doubleContraction(deviatoric_trial_stress);
82  Real effective_trial_stress = std::sqrt(3.0 / 2.0 * dev_trial_stress_squared);
83 
84  // Set the value of 3 * shear modulus for use as a reference residual value
86 
87  computeStressInitialize(effective_trial_stress, elasticity_tensor);
88 
89  // Use Newton iteration to determine the scalar effective inelastic strain increment
90  Real scalar_effective_inelastic_strain = 0.0;
91  if (!MooseUtils::absoluteFuzzyEqual(effective_trial_stress, 0.0))
92  {
93  returnMappingSolve(effective_trial_stress, scalar_effective_inelastic_strain, _console);
94  if (scalar_effective_inelastic_strain != 0.0)
95  inelastic_strain_increment =
96  deviatoric_trial_stress *
97  (1.5 * scalar_effective_inelastic_strain / effective_trial_stress);
98  else
99  inelastic_strain_increment.zero();
100  }
101  else
102  inelastic_strain_increment.zero();
103 
104  strain_increment -= inelastic_strain_increment;
106  _effective_inelastic_strain_old[_qp] + scalar_effective_inelastic_strain;
107 
108  // Use the old elastic strain here because we require tensors used by this class
109  // to be isotropic and this method natively allows for changing in time
110  // elasticity tensors
111  stress_new = elasticity_tensor * (strain_increment + elastic_strain_old);
112 
113  computeStressFinalize(inelastic_strain_increment);
114 
115  if (compute_full_tangent_operator &&
117  {
118  if (MooseUtils::absoluteFuzzyEqual(scalar_effective_inelastic_strain, 0.0))
119  tangent_operator.zero();
120  else
121  {
122  // mu = _three_shear_modulus / 3.0;
123  // norm_dev_stress = ||s_n+1||
124  // effective_trial_stress = von mises trial stress = std::sqrt(3.0 / 2.0) * ||s_n+1^trial||
125  // scalar_effective_inelastic_strain = Delta epsilon^cr_n+1
126  // deriv = derivative of scalar_effective_inelastic_strain w.r.t. von mises stress
127  // deriv = std::sqrt(3.0 / 2.0) partial Delta epsilon^cr_n+1n over partial ||s_n+1^trial||
128 
129  mooseAssert(_three_shear_modulus != 0.0, "Shear modulus is zero");
130 
131  const RankTwoTensor deviatoric_stress = stress_new.deviatoric();
132  const Real norm_dev_stress =
133  std::sqrt(deviatoric_stress.doubleContraction(deviatoric_stress));
134  mooseAssert(norm_dev_stress != 0.0, "Norm of the deviatoric is zero");
135 
136  const RankTwoTensor flow_direction = deviatoric_stress / norm_dev_stress;
137  const RankFourTensor flow_direction_dyad = flow_direction.outerProduct(flow_direction);
138  const Real deriv =
139  computeStressDerivative(effective_trial_stress, scalar_effective_inelastic_strain);
140  const Real scalar_one = _three_shear_modulus * scalar_effective_inelastic_strain /
141  std::sqrt(1.5) / norm_dev_stress;
142 
143  tangent_operator = scalar_one * _deviatoric_projection_four +
144  (_three_shear_modulus * deriv - scalar_one) * flow_direction_dyad;
145  }
146  }
147 }

◆ validParams()

InputParameters RadialReturnStressUpdate::validParams ( )
static

Definition at line 18 of file RadialReturnStressUpdate.C.

19 {
20  InputParameters params = StressUpdateBase::validParams();
21  params.addClassDescription("Calculates the effective inelastic strain increment required to "
22  "return the isotropic stress state to a J2 yield surface. This class "
23  "is intended to be a parent class for classes with specific "
24  "constitutive models.");
26  params.addParam<Real>("max_inelastic_increment",
27  1e-4,
28  "The maximum inelastic strain increment allowed in a time step");
29  params.addRequiredParam<std::string>(
30  "effective_inelastic_strain_name",
31  "Name of the material property that stores the effective inelastic strain");
32  params.addParamNamesToGroup("effective_inelastic_strain_name", "Advanced");
33 
34  return params;
35 }

Referenced by RadialReturnCreepStressUpdateBase::validParams(), and IsotropicPlasticityStressUpdate::validParams().

Member Data Documentation

◆ _absolute_tolerance

Real SingleVariableReturnMappingSolution::_absolute_tolerance
privateinherited

◆ _acceptable_multiplier

Real SingleVariableReturnMappingSolution::_acceptable_multiplier
privateinherited

Multiplier applied to relative and absolute tolerances for acceptable convergence.

Definition at line 158 of file SingleVariableReturnMappingSolution.h.

Referenced by SingleVariableReturnMappingSolution::convergedAcceptable().

◆ _base_name

const std::string StressUpdateBase::_base_name
protectedinherited

Name used as a prefix for all material properties related to the stress update model.

Definition at line 128 of file StressUpdateBase.h.

◆ _bracket_solution

bool SingleVariableReturnMappingSolution::_bracket_solution
protectedinherited

Whether to save upper and lower bounds of root for scalar, and set solution to the midpoint between those bounds if outside them.

Definition at line 127 of file SingleVariableReturnMappingSolution.h.

Referenced by SingleVariableReturnMappingSolution::internalSolve().

◆ _check_range

bool SingleVariableReturnMappingSolution::_check_range
protectedinherited

Whether to check to see whether iterative solution is within admissible range, and set within that range if outside.

Definition at line 120 of file SingleVariableReturnMappingSolution.h.

Referenced by SingleVariableReturnMappingSolution::internalSolve().

◆ _deviatoric_projection_four

const RankFourTensor RadialReturnStressUpdate::_deviatoric_projection_four
protected

Rank four deviatoric projection tensor.

Definition at line 152 of file RadialReturnStressUpdate.h.

Referenced by updateState().

◆ _effective_inelastic_strain

MaterialProperty<Real>& RadialReturnStressUpdate::_effective_inelastic_strain
protected

◆ _effective_inelastic_strain_old

const MaterialProperty<Real>& RadialReturnStressUpdate::_effective_inelastic_strain_old
protected

◆ _identity_symmetric_four

const RankFourTensor RadialReturnStressUpdate::_identity_symmetric_four
protected

Rank four symmetric identity tensor.

Definition at line 147 of file RadialReturnStressUpdate.h.

◆ _identity_two

const RankTwoTensor RadialReturnStressUpdate::_identity_two
protected

Rank two identity tensor.

Definition at line 142 of file RadialReturnStressUpdate.h.

◆ _initial_residual

Real SingleVariableReturnMappingSolution::_initial_residual
privateinherited

Residual values, kept as members to retain solver state for summary outputting.

Definition at line 170 of file SingleVariableReturnMappingSolution.h.

Referenced by SingleVariableReturnMappingSolution::internalSolve(), and SingleVariableReturnMappingSolution::outputIterationSummary().

◆ _internal_solve_full_iteration_history

const bool SingleVariableReturnMappingSolution::_internal_solve_full_iteration_history
privateinherited

Whether to output iteration information all the time (regardless of whether iterations converge)

Definition at line 149 of file SingleVariableReturnMappingSolution.h.

Referenced by SingleVariableReturnMappingSolution::returnMappingSolve().

◆ _internal_solve_output_on

enum SingleVariableReturnMappingSolution::InternalSolveOutput SingleVariableReturnMappingSolution::_internal_solve_output_on
privateinherited

◆ _iteration

unsigned int SingleVariableReturnMappingSolution::_iteration
privateinherited

◆ _line_search

bool SingleVariableReturnMappingSolution::_line_search
protectedinherited

Whether to use line searches to improve convergence.

Definition at line 123 of file SingleVariableReturnMappingSolution.h.

Referenced by SingleVariableReturnMappingSolution::internalSolve().

◆ _max_inelastic_increment

Real RadialReturnStressUpdate::_max_inelastic_increment
protected

Definition at line 137 of file RadialReturnStressUpdate.h.

Referenced by computeTimeStepLimit().

◆ _max_its

const unsigned int SingleVariableReturnMappingSolution::_max_its
privateinherited

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 146 of file SingleVariableReturnMappingSolution.h.

Referenced by SingleVariableReturnMappingSolution::internalSolve().

◆ _num_resids

const std::size_t SingleVariableReturnMappingSolution::_num_resids
privateinherited

◆ _relative_tolerance

Real SingleVariableReturnMappingSolution::_relative_tolerance
privateinherited

◆ _residual

Real SingleVariableReturnMappingSolution::_residual
privateinherited

◆ _residual_history

std::vector<Real> SingleVariableReturnMappingSolution::_residual_history
privateinherited

History of residuals used to check whether progress is still being made on decreasing the residual.

Definition at line 164 of file SingleVariableReturnMappingSolution.h.

Referenced by SingleVariableReturnMappingSolution::convergedAcceptable(), and SingleVariableReturnMappingSolution::internalSolve().

◆ _svrms_name

const std::string SingleVariableReturnMappingSolution::_svrms_name
privateinherited

MOOSE input name of the object performing the solve.

Definition at line 175 of file SingleVariableReturnMappingSolution.h.

Referenced by SingleVariableReturnMappingSolution::outputIterationSummary().

◆ _three_shear_modulus

Real RadialReturnStressUpdate::_three_shear_modulus
protected

The documentation for this class was generated from the following files:
SingleVariableReturnMappingSolution::maximumPermissibleValue
virtual Real maximumPermissibleValue(const Real effective_trial_stress) const
Compute the maximum permissible value of the scalar.
Definition: SingleVariableReturnMappingSolution.C:89
SingleVariableReturnMappingSolution::_check_range
bool _check_range
Whether to check to see whether iterative solution is within admissible range, and set within that ra...
Definition: SingleVariableReturnMappingSolution.h:120
SingleVariableReturnMappingSolution::_initial_residual
Real _initial_residual
Residual values, kept as members to retain solver state for summary outputting.
Definition: SingleVariableReturnMappingSolution.h:170
SingleVariableReturnMappingSolution::_bracket_solution
bool _bracket_solution
Whether to save upper and lower bounds of root for scalar, and set solution to the midpoint between t...
Definition: SingleVariableReturnMappingSolution.h:127
SingleVariableReturnMappingSolution::validParams
static InputParameters validParams()
Definition: SingleVariableReturnMappingSolution.C:28
SingleVariableReturnMappingSolution::SolveState::NAN_INF
SingleVariableReturnMappingSolution::_acceptable_multiplier
Real _acceptable_multiplier
Multiplier applied to relative and absolute tolerances for acceptable convergence.
Definition: SingleVariableReturnMappingSolution.h:158
SingleVariableReturnMappingSolution::minimumPermissibleValue
virtual Real minimumPermissibleValue(const Real effective_trial_stress) const
Compute the minimum permissible value of the scalar.
Definition: SingleVariableReturnMappingSolution.C:82
SingleVariableReturnMappingSolution::outputIterationSummary
virtual void outputIterationSummary(std::stringstream *iter_output, const unsigned int total_it)
Output summary information for the convergence history of the model.
Definition: SingleVariableReturnMappingSolution.C:347
RadialReturnStressUpdate::_identity_two
const RankTwoTensor _identity_two
Rank two identity tensor.
Definition: RadialReturnStressUpdate.h:142
RadialReturnStressUpdate::_deviatoric_projection_four
const RankFourTensor _deviatoric_projection_four
Rank four deviatoric projection tensor.
Definition: RadialReturnStressUpdate.h:152
SingleVariableReturnMappingSolution::SolveState::SUCCESS
RadialReturnStressUpdate::_max_inelastic_increment
Real _max_inelastic_increment
Definition: RadialReturnStressUpdate.h:137
TangentCalculationMethod::ELASTIC
StressUpdateBase::_base_name
const std::string _base_name
Name used as a prefix for all material properties related to the stress update model.
Definition: StressUpdateBase.h:128
SingleVariableReturnMappingSolution::_iteration
unsigned int _iteration
iteration number
Definition: SingleVariableReturnMappingSolution.h:167
RadialReturnStressUpdate::computeStressDerivative
virtual Real computeStressDerivative(const Real, const Real)
Calculate the derivative of the strain increment with respect to the updated stress.
Definition: RadialReturnStressUpdate.h:118
SingleVariableReturnMappingSolution::_max_its
const unsigned int _max_its
Maximum number of return mapping iterations.
Definition: SingleVariableReturnMappingSolution.h:146
SingleVariableReturnMappingSolution::updateBounds
void updateBounds(const Real scalar, const Real residual, const Real init_resid_sign, Real &scalar_upper_bound, Real &scalar_lower_bound, std::stringstream *iter_output)
Update the upper and lower bounds of the root for the effective inelastic strain.
Definition: SingleVariableReturnMappingSolution.C:384
SingleVariableReturnMappingSolution::converged
bool converged(const Real residual, const Real reference)
Check to see whether the residual is within the convergence limits.
Definition: SingleVariableReturnMappingSolution.C:298
SingleVariableReturnMappingSolution::initialGuess
virtual Real initialGuess(const Real)
Compute an initial guess for the value of the scalar.
Definition: SingleVariableReturnMappingSolution.h:64
StressUpdateBase::StressUpdateBase
StressUpdateBase(const InputParameters &parameters)
Definition: StressUpdateBase.C:36
SingleVariableReturnMappingSolution::InternalSolveOutput::ALWAYS
SingleVariableReturnMappingSolution::internalSolve
SolveState internalSolve(const Real effective_trial_stress, Real &scalar, std::stringstream *iter_output=nullptr)
Method called from within this class to perform the actual return mappping iterations.
Definition: SingleVariableReturnMappingSolution.C:157
SingleVariableReturnMappingSolution::InternalSolveOutput::NEVER
RadialReturnStressUpdate::_effective_inelastic_strain
MaterialProperty< Real > & _effective_inelastic_strain
Definition: RadialReturnStressUpdate.h:135
SingleVariableReturnMappingSolution::iterationFinalize
virtual void iterationFinalize(Real)
Finalize internal state variables for a model for a given iteration.
Definition: SingleVariableReturnMappingSolution.h:94
SingleVariableReturnMappingSolution::_internal_solve_output_on
enum SingleVariableReturnMappingSolution::InternalSolveOutput _internal_solve_output_on
RadialReturnStressUpdate::computeStressInitialize
virtual void computeStressInitialize(const Real, const RankFourTensor &)
Perform any necessary initialization before return mapping iterations.
Definition: RadialReturnStressUpdate.h:108
RadialReturnStressUpdate::_effective_inelastic_strain_old
const MaterialProperty< Real > & _effective_inelastic_strain_old
Definition: RadialReturnStressUpdate.h:136
SingleVariableReturnMappingSolution::_residual
Real _residual
Definition: SingleVariableReturnMappingSolution.h:171
StressUpdateBase::getTangentCalculationMethod
virtual TangentCalculationMethod getTangentCalculationMethod()
Definition: StressUpdateBase.h:116
SingleVariableReturnMappingSolution::outputIterationStep
virtual void outputIterationStep(std::stringstream *iter_output, const unsigned int it, const Real effective_trial_stress, const Real scalar, const Real residual, const Real reference_residual)
Output information for a single iteration step to build the convergence history of the model.
Definition: SingleVariableReturnMappingSolution.C:328
SingleVariableReturnMappingSolution::computeDerivative
virtual Real computeDerivative(const Real effective_trial_stress, const Real scalar)=0
Compute the derivative of the residual as a function of the scalar variable.
SingleVariableReturnMappingSolution::_line_search
bool _line_search
Whether to use line searches to improve convergence.
Definition: SingleVariableReturnMappingSolution.h:123
RadialReturnStressUpdate::_three_shear_modulus
Real _three_shear_modulus
3 * shear modulus
Definition: RadialReturnStressUpdate.h:133
SingleVariableReturnMappingSolution::checkPermissibleRange
void checkPermissibleRange(Real &scalar, Real &scalar_increment, const Real scalar_old, const Real min_permissible_scalar, const Real 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.
Definition: SingleVariableReturnMappingSolution.C:356
SingleVariableReturnMappingSolution::returnMappingSolve
void returnMappingSolve(const Real effective_trial_stress, Real &scalar, const ConsoleStream &console)
Perform the return mapping iterations.
Definition: SingleVariableReturnMappingSolution.C:96
SingleVariableReturnMappingSolution::computeResidual
virtual Real computeResidual(const Real effective_trial_stress, const Real scalar)=0
Compute the residual for a predicted value of the scalar.
SingleVariableReturnMappingSolution::_internal_solve_full_iteration_history
const bool _internal_solve_full_iteration_history
Whether to output iteration information all the time (regardless of whether iterations converge)
Definition: SingleVariableReturnMappingSolution.h:149
RadialReturnStressUpdate::computeStressFinalize
virtual void computeStressFinalize(const RankTwoTensor &)
Perform any necessary steps to finalize state after return mapping iterations.
Definition: RadialReturnStressUpdate.h:127
SingleVariableReturnMappingSolution::SolveState::EXCEEDED_ITERATIONS
ElasticityTensorTools::getIsotropicShearModulus
T getIsotropicShearModulus(const RankFourTensorTempl< T > &elasticity_tensor)
Get the shear modulus for an isotropic elasticity tensor param elasticity_tensor the tensor (must be ...
Definition: ElasticityTensorTools.h:71
SingleVariableReturnMappingSolution::_residual_history
std::vector< Real > _residual_history
History of residuals used to check whether progress is still being made on decreasing the residual.
Definition: SingleVariableReturnMappingSolution.h:164
RankFourTensorTempl
Definition: ACGrGrElasticDrivingForce.h:20
SingleVariableReturnMappingSolution::convergedAcceptable
bool convergedAcceptable(const unsigned int it, const Real residual, const Real reference)
Check to see whether the residual is within acceptable convergence limits.
Definition: SingleVariableReturnMappingSolution.C:305
SingleVariableReturnMappingSolution::_num_resids
const std::size_t _num_resids
Number of residuals to be stored in history.
Definition: SingleVariableReturnMappingSolution.h:161
TangentCalculationMethod::PARTIAL
SingleVariableReturnMappingSolution::computeReferenceResidual
virtual Real computeReferenceResidual(const Real effective_trial_stress, const Real scalar)=0
Compute a reference quantity to be used for checking relative convergence.
RankTwoTensorTempl< Real >
RadialReturnStressUpdate::_identity_symmetric_four
const RankFourTensor _identity_symmetric_four
Rank four symmetric identity tensor.
Definition: RadialReturnStressUpdate.h:147
SingleVariableReturnMappingSolution::_relative_tolerance
Real _relative_tolerance
Relative convergence tolerance.
Definition: SingleVariableReturnMappingSolution.h:152
SingleVariableReturnMappingSolution::_absolute_tolerance
Real _absolute_tolerance
Absolute convergence tolerance.
Definition: SingleVariableReturnMappingSolution.h:155
SingleVariableReturnMappingSolution::SingleVariableReturnMappingSolution
SingleVariableReturnMappingSolution(const InputParameters &parameters)
Definition: SingleVariableReturnMappingSolution.C:59
StressUpdateBase::validParams
static InputParameters validParams()
Definition: StressUpdateBase.C:17