https://mooseframework.inl.gov
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Types | Private Member Functions | Private Attributes | List of all members
SingleVariableReturnMappingSolutionTempl< is_ad > Class Template Referenceabstract

Base class that provides capability for Newton return mapping iterations on a single variable. More...

#include <SingleVariableReturnMappingSolution.h>

Inheritance diagram for SingleVariableReturnMappingSolutionTempl< is_ad >:
[legend]

Public Member Functions

 SingleVariableReturnMappingSolutionTempl (const InputParameters &parameters)
 
virtual ~SingleVariableReturnMappingSolutionTempl ()
 

Static Public Member Functions

static InputParameters validParams ()
 

Protected Member Functions

void returnMappingSolve (const GenericReal< is_ad > &effective_trial_stress, GenericReal< is_ad > &scalar, const ConsoleStream &console)
 Perform the return mapping iterations. More...
 
virtual GenericReal< is_ad > minimumPermissibleValue (const GenericReal< is_ad > &effective_trial_stress) const
 Compute the minimum permissible value of the scalar. More...
 
virtual GenericReal< is_ad > maximumPermissibleValue (const GenericReal< is_ad > &effective_trial_stress) const
 Compute the maximum permissible value of the scalar. More...
 
virtual GenericReal< is_ad > initialGuess (const GenericReal< is_ad > &)
 Compute an initial guess for the value of the scalar. More...
 
virtual GenericReal< is_ad > computeResidual (const GenericReal< is_ad > &, const GenericReal< is_ad > &)=0
 Compute the residual for a predicted value of the scalar. More...
 
virtual GenericReal< is_ad > computeDerivative (const GenericReal< is_ad > &, const GenericReal< is_ad > &)=0
 Compute the derivative of the residual as a function of the scalar variable. More...
 
virtual GenericChainedReal< is_ad > computeResidualAndDerivative (const GenericReal< is_ad > &, const GenericChainedReal< is_ad > &)
 Compute the residual and the derivative for a predicted value of the scalar. More...
 
virtual Real computeReferenceResidual (const GenericReal< is_ad > &effective_trial_stress, const GenericReal< is_ad > &scalar)=0
 Compute a reference quantity to be used for checking relative convergence. More...
 
virtual void preStep (const GenericReal< is_ad > &, const GenericReal< is_ad > &, const GenericReal< is_ad > &)
 This method is called before taking a step in the return mapping algorithm. More...
 
virtual void iterationFinalize (const GenericReal< is_ad > &)
 Finalize internal state variables for a model for a given iteration. More...
 
virtual void outputIterationSummary (std::stringstream *iter_output, const unsigned int total_it)
 Output summary information for the convergence history of the model. More...
 
virtual void outputIterationStep (std::stringstream *iter_output, const GenericReal< is_ad > &effective_trial_stress, const GenericReal< is_ad > &scalar, const Real reference_residual)
 Output information for a single iteration step to build the convergence history of the model. More...
 
bool converged (const GenericReal< is_ad > &residual, const Real reference)
 Check to see whether the residual is within the convergence limits. More...
 

Protected Attributes

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

void computeResidualAndDerivativeHelper (const GenericReal< is_ad > &effective_trial_stress, const GenericReal< is_ad > &scalar)
 Helper function to compute and set the _residual and _derivative. More...
 
SolveState internalSolve (const GenericReal< is_ad > effective_trial_stress, GenericReal< is_ad > &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 (GenericReal< is_ad > &scalar, GenericReal< is_ad > &scalar_increment, const GenericReal< is_ad > &scalar_old, const GenericReal< is_ad > min_permissible_scalar, const GenericReal< is_ad > 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 GenericReal< is_ad > &scalar, const GenericReal< is_ad > &residual, const Real init_resid_sign, GenericReal< is_ad > &scalar_upper_bound, GenericReal< is_ad > &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 SingleVariableReturnMappingSolutionTempl::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 bool _ad_derivative
 
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...
 
GenericReal< is_ad > _derivative
 Derivative of the residual. More...
 
const std::string _svrms_name
 MOOSE input name of the object performing the solve. More...
 
GenericReal< is_ad > _initial_residual
 Residual values, kept as members to retain solver state for summary outputting. More...
 
GenericReal< is_ad > _residual
 

Detailed Description

template<bool is_ad>
class SingleVariableReturnMappingSolutionTempl< is_ad >

Base class that provides capability for Newton return mapping iterations on a single variable.

Definition at line 24 of file SingleVariableReturnMappingSolution.h.

Member Enumeration Documentation

◆ InternalSolveOutput

template<bool is_ad>
enum SingleVariableReturnMappingSolutionTempl::InternalSolveOutput
strongprivate
Enumerator
NEVER 
ON_ERROR 
ALWAYS 

Definition at line 172 of file SingleVariableReturnMappingSolution.h.

173  {
174  NEVER,
175  ON_ERROR,
176  ALWAYS
enum SingleVariableReturnMappingSolutionTempl::InternalSolveOutput _internal_solve_output_on

◆ SolveState

template<bool is_ad>
enum SingleVariableReturnMappingSolutionTempl::SolveState
strongprivate
Enumerator
SUCCESS 
NAN_INF 
EXCEEDED_ITERATIONS 

Definition at line 179 of file SingleVariableReturnMappingSolution.h.

180  {
181  SUCCESS,
182  NAN_INF,
183  EXCEEDED_ITERATIONS
184  };

Constructor & Destructor Documentation

◆ SingleVariableReturnMappingSolutionTempl()

Definition at line 59 of file SingleVariableReturnMappingSolution.C.

61  : _check_range(false),
62  _line_search(true),
63  _bracket_solution(true),
65  parameters.get<MooseEnum>("internal_solve_output_on").getEnum<InternalSolveOutput>()),
66  _max_its(1000), // Far larger than ever expected to be needed
68  parameters.get<bool>("internal_solve_full_iteration_history")),
69  _relative_tolerance(parameters.get<Real>("relative_tolerance")),
70  _absolute_tolerance(parameters.get<Real>("absolute_tolerance")),
71  _acceptable_multiplier(parameters.get<Real>("acceptable_multiplier")),
72  _ad_derivative(parameters.get<bool>("automatic_differentiation_return_mapping")),
73  _num_resids(30),
74  _residual_history(_num_resids, std::numeric_limits<Real>::max()),
75  _iteration(0),
76  _initial_residual(0.0),
77  _residual(0.0),
78  _svrms_name(parameters.get<std::string>("_object_name"))
79 {
80 }
bool _line_search
Whether to use line searches to improve convergence.
Real _acceptable_multiplier
Multiplier applied to relative and absolute tolerances for acceptable convergence.
std::vector< Real > _residual_history
History of residuals used to check whether progress is still being made on decreasing the residual...
std::vector< std::pair< R1, R2 > > get(const std::string &param1, const std::string &param2) const
Real _relative_tolerance
Relative convergence tolerance.
T getEnum() const
const unsigned int _max_its
Maximum number of return mapping iterations.
const bool _internal_solve_full_iteration_history
Whether to output iteration information all the time (regardless of whether iterations converge) ...
Real _absolute_tolerance
Absolute convergence tolerance.
bool _bracket_solution
Whether to save upper and lower bounds of root for scalar, and set solution to the midpoint between t...
GenericReal< is_ad > _initial_residual
Residual values, kept as members to retain solver state for summary outputting.
enum SingleVariableReturnMappingSolutionTempl::InternalSolveOutput _internal_solve_output_on
const std::string _svrms_name
MOOSE input name of the object performing the solve.
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.

◆ ~SingleVariableReturnMappingSolutionTempl()

template<bool is_ad>
virtual SingleVariableReturnMappingSolutionTempl< is_ad >::~SingleVariableReturnMappingSolutionTempl ( )
inlinevirtual

Definition at line 30 of file SingleVariableReturnMappingSolution.h.

30 {}

Member Function Documentation

◆ checkPermissibleRange()

template<bool is_ad>
void SingleVariableReturnMappingSolutionTempl< is_ad >::checkPermissibleRange ( GenericReal< is_ad > &  scalar,
GenericReal< is_ad > &  scalar_increment,
const GenericReal< is_ad > &  scalar_old,
const GenericReal< is_ad >  min_permissible_scalar,
const GenericReal< is_ad >  max_permissible_scalar,
std::stringstream *  iter_output 
)
private

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 359 of file SingleVariableReturnMappingSolution.C.

366 {
367  if (scalar > max_permissible_scalar)
368  {
369  scalar_increment = (max_permissible_scalar - scalar_old) / 2.0;
370  scalar = scalar_old + scalar_increment;
371  if (iter_output)
372  *iter_output << "Scalar greater than maximum ("
373  << MetaPhysicL::raw_value(max_permissible_scalar)
374  << ") adjusted scalar=" << MetaPhysicL::raw_value(scalar)
375  << " scalar_increment=" << MetaPhysicL::raw_value(scalar_increment) << std::endl;
376  }
377  else if (scalar < min_permissible_scalar)
378  {
379  scalar_increment = (min_permissible_scalar - scalar_old) / 2.0;
380  scalar = scalar_old + scalar_increment;
381  if (iter_output)
382  *iter_output << "Scalar less than minimum (" << MetaPhysicL::raw_value(min_permissible_scalar)
383  << ") adjusted scalar=" << MetaPhysicL::raw_value(scalar)
384  << " scalar_increment=" << MetaPhysicL::raw_value(scalar_increment) << std::endl;
385  }
386 }
auto raw_value(const Eigen::Map< T > &in)

◆ computeDerivative()

template<bool is_ad>
virtual GenericReal<is_ad> SingleVariableReturnMappingSolutionTempl< is_ad >::computeDerivative ( const GenericReal< is_ad > &  ,
const GenericReal< is_ad > &   
)
protectedpure virtual

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 LAROMANCEStressUpdateBaseTempl< is_ad >, CompositePowerLawCreepStressUpdateTempl< is_ad >, IsotropicPlasticityStressUpdateTempl< is_ad >, PowerLawCreepStressUpdateTempl< is_ad >, PowerLawCreepTestTempl< is_ad >, and CombinedNonlinearHardeningPlasticityTempl< is_ad >.

◆ computeReferenceResidual()

template<bool is_ad>
virtual Real SingleVariableReturnMappingSolutionTempl< is_ad >::computeReferenceResidual ( const GenericReal< is_ad > &  effective_trial_stress,
const GenericReal< is_ad > &  scalar 
)
protectedpure virtual

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

Implemented in RadialReturnStressUpdateTempl< is_ad >.

◆ computeResidual()

template<bool is_ad>
virtual GenericReal<is_ad> SingleVariableReturnMappingSolutionTempl< is_ad >::computeResidual ( const GenericReal< is_ad > &  ,
const GenericReal< is_ad > &   
)
protectedpure virtual

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 LAROMANCEStressUpdateBaseTempl< is_ad >, IsotropicPlasticityStressUpdateTempl< is_ad >, CompositePowerLawCreepStressUpdateTempl< is_ad >, PowerLawCreepStressUpdateTempl< is_ad >, PowerLawCreepTestTempl< is_ad >, and CombinedNonlinearHardeningPlasticityTempl< is_ad >.

◆ computeResidualAndDerivative()

template<bool is_ad>
virtual GenericChainedReal<is_ad> SingleVariableReturnMappingSolutionTempl< is_ad >::computeResidualAndDerivative ( const GenericReal< is_ad > &  ,
const GenericChainedReal< is_ad > &   
)
inlineprotectedvirtual

Compute the residual and the derivative 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

Reimplemented in CompositePowerLawCreepStressUpdateTempl< is_ad >, and PowerLawCreepStressUpdateTempl< is_ad >.

Definition at line 97 of file SingleVariableReturnMappingSolution.h.

99  {
100  mooseError("computeResidualAndDerivative has to be implemented if "
101  "automatic_differentiation_return_mapping = true.");
102  return 0;
103  };
void mooseError(Args &&... args)

◆ computeResidualAndDerivativeHelper()

template<bool is_ad>
void SingleVariableReturnMappingSolutionTempl< is_ad >::computeResidualAndDerivativeHelper ( const GenericReal< is_ad > &  effective_trial_stress,
const GenericReal< is_ad > &  scalar 
)
private

Helper function to compute and set the _residual and _derivative.

Definition at line 307 of file SingleVariableReturnMappingSolution.C.

309 {
310  if (_ad_derivative)
311  {
312  GenericChainedReal<is_ad> residual_and_derivative =
313  computeResidualAndDerivative(effective_trial_stress, GenericChainedReal<is_ad>(scalar, 1));
314  _residual = residual_and_derivative.value();
315  _derivative = residual_and_derivative.derivatives();
316  }
317  else
318  {
319  _residual = computeResidual(effective_trial_stress, scalar);
320  _derivative = computeDerivative(effective_trial_stress, scalar);
321  }
322 }
virtual GenericReal< is_ad > computeDerivative(const GenericReal< is_ad > &, const GenericReal< is_ad > &)=0
Compute the derivative of the residual as a function of the scalar variable.
Moose::GenericType< ChainedReal, is_ad > GenericChainedReal
virtual GenericChainedReal< is_ad > computeResidualAndDerivative(const GenericReal< is_ad > &, const GenericChainedReal< is_ad > &)
Compute the residual and the derivative for a predicted value of the scalar.
virtual GenericReal< is_ad > computeResidual(const GenericReal< is_ad > &, const GenericReal< is_ad > &)=0
Compute the residual for a predicted value of the scalar.
GenericReal< is_ad > _derivative
Derivative of the residual.

◆ converged()

template<bool is_ad>
bool SingleVariableReturnMappingSolutionTempl< is_ad >::converged ( const GenericReal< is_ad > &  residual,
const Real  reference 
)
protected

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 326 of file SingleVariableReturnMappingSolution.C.

328 {
329  const Real residual = MetaPhysicL::raw_value(ad_residual);
330  return (std::abs(residual) <= _absolute_tolerance ||
331  std::abs(residual / reference) <= _relative_tolerance);
332 }
auto raw_value(const Eigen::Map< T > &in)
Real _relative_tolerance
Relative convergence tolerance.
Real _absolute_tolerance
Absolute convergence tolerance.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ convergedAcceptable()

template<bool is_ad>
bool SingleVariableReturnMappingSolutionTempl< is_ad >::convergedAcceptable ( const unsigned int  it,
const Real  reference 
)
private

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 336 of file SingleVariableReturnMappingSolution.C.

338 {
339  // Require that we have at least done _num_resids evaluations before we allow for
340  // acceptable convergence
341  if (it < _num_resids)
342  return false;
343 
344  // Check to see whether the residual has dropped by convergence_history_factor over
345  // the last _num_resids iterations. If it has (which means it's still making progress),
346  // don't consider it to be converged within the acceptable limits.
347  const Real convergence_history_factor = 10.0;
348  if (std::abs(_residual * convergence_history_factor) <
349  std::abs(_residual_history[(it + 1) % _num_resids]))
350  return false;
351 
352  // Now that it's determined that progress is not being made, treat it as converged if
353  // we're within the acceptable convergence limits
354  return converged(_residual / _acceptable_multiplier, reference);
355 }
Real _acceptable_multiplier
Multiplier applied to relative and absolute tolerances for acceptable convergence.
std::vector< Real > _residual_history
History of residuals used to check whether progress is still being made on decreasing the residual...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool converged(const GenericReal< is_ad > &residual, const Real reference)
Check to see whether the residual is within the convergence limits.
const std::size_t _num_resids
Number of residuals to be stored in history.

◆ initialGuess()

template<bool is_ad>
virtual GenericReal<is_ad> SingleVariableReturnMappingSolutionTempl< is_ad >::initialGuess ( const GenericReal< is_ad > &  )
inlineprotectedvirtual

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

Reimplemented in PowerLawCreepTestTempl< is_ad >.

Definition at line 66 of file SingleVariableReturnMappingSolution.h.

67  {
68  return 0.0;
69  }

◆ internalSolve()

template<bool is_ad>
SingleVariableReturnMappingSolutionTempl< is_ad >::SolveState SingleVariableReturnMappingSolutionTempl< is_ad >::internalSolve ( const GenericReal< is_ad >  effective_trial_stress,
GenericReal< is_ad > &  scalar,
std::stringstream *  iter_output = nullptr 
)
private

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 163 of file SingleVariableReturnMappingSolution.C.

167 {
168  scalar = initialGuess(effective_trial_stress);
169  GenericReal<is_ad> scalar_old = scalar;
170  GenericReal<is_ad> scalar_increment = 0.0;
171  const GenericReal<is_ad> min_permissible_scalar = minimumPermissibleValue(effective_trial_stress);
172  const GenericReal<is_ad> max_permissible_scalar = maximumPermissibleValue(effective_trial_stress);
173  GenericReal<is_ad> scalar_upper_bound = max_permissible_scalar;
174  GenericReal<is_ad> scalar_lower_bound = min_permissible_scalar;
175  _iteration = 0;
176 
177  computeResidualAndDerivativeHelper(effective_trial_stress, scalar);
179 
180  GenericReal<is_ad> residual_old = _residual;
182  Real reference_residual = computeReferenceResidual(effective_trial_stress, scalar);
183 
184  if (converged(_residual, reference_residual))
185  {
186  iterationFinalize(scalar);
187  outputIterationStep(iter_output, effective_trial_stress, scalar, reference_residual);
188  return SolveState::SUCCESS;
189  }
190 
191  _residual_history.assign(_num_resids, std::numeric_limits<Real>::max());
193 
194  while (_iteration < _max_its && !converged(_residual, reference_residual) &&
195  !convergedAcceptable(_iteration, reference_residual))
196  {
197  preStep(scalar_old, _residual, _derivative);
198 
199  scalar_increment = -_residual / _derivative;
200  scalar = scalar_old + scalar_increment;
201 
202  if (_check_range)
203  checkPermissibleRange(scalar,
204  scalar_increment,
205  scalar_old,
206  min_permissible_scalar,
207  max_permissible_scalar,
208  iter_output);
209 
210  computeResidualAndDerivativeHelper(effective_trial_stress, scalar);
211  reference_residual = computeReferenceResidual(effective_trial_stress, scalar);
212  iterationFinalize(scalar);
213 
214  if (_bracket_solution)
215  updateBounds(
216  scalar, _residual, init_resid_sign, scalar_upper_bound, scalar_lower_bound, iter_output);
217 
218  if (converged(_residual, reference_residual))
219  {
220  outputIterationStep(iter_output, effective_trial_stress, scalar, reference_residual);
221  break;
222  }
223  else
224  {
225  bool modified_increment = false;
226 
227  // Line Search
228  if (_line_search)
229  {
230  if (residual_old - _residual != 0.0)
231  {
232  GenericReal<is_ad> alpha = residual_old / (residual_old - _residual);
233  alpha = MathUtils::clamp(alpha, 1.0e-2, 1.0);
234 
235  if (alpha != 1.0)
236  {
237  modified_increment = true;
238  scalar_increment *= alpha;
239  if (iter_output)
240  *iter_output << " Line search alpha = " << MetaPhysicL::raw_value(alpha)
241  << " increment = " << MetaPhysicL::raw_value(scalar_increment)
242  << std::endl;
243  }
244  }
245  }
246 
247  if (_bracket_solution)
248  {
249  // Check to see whether trial scalar_increment is outside the bounds, and set it to a point
250  // within the bounds if it is
251  if (scalar_old + scalar_increment >= scalar_upper_bound ||
252  scalar_old + scalar_increment <= scalar_lower_bound)
253  {
254  if (scalar_upper_bound != max_permissible_scalar &&
255  scalar_lower_bound != min_permissible_scalar)
256  {
257  const Real frac = 0.5;
258  scalar_increment =
259  (1.0 - frac) * scalar_lower_bound + frac * scalar_upper_bound - scalar_old;
260  modified_increment = true;
261  if (iter_output)
262  *iter_output << " Trial scalar_increment exceeded bounds. Setting between "
263  "lower/upper bounds. frac: "
264  << frac << std::endl;
265  }
266  }
267  }
268 
269  // Update the trial scalar and recompute residual if the line search or bounds checking
270  // modified the increment
271  if (modified_increment)
272  {
273  scalar = scalar_old + scalar_increment;
274  computeResidualAndDerivativeHelper(effective_trial_stress, scalar);
275  reference_residual = computeReferenceResidual(effective_trial_stress, scalar);
276  iterationFinalize(scalar);
277 
278  if (_bracket_solution)
279  updateBounds(scalar,
280  _residual,
281  init_resid_sign,
282  scalar_upper_bound,
283  scalar_lower_bound,
284  iter_output);
285  }
286  }
287 
288  outputIterationStep(iter_output, effective_trial_stress, scalar, reference_residual);
289 
290  ++_iteration;
291  residual_old = _residual;
292  scalar_old = scalar;
294  }
295 
296  if (std::isnan(_residual) || std::isinf(MetaPhysicL::raw_value(_residual)))
297  return SolveState::NAN_INF;
298 
299  if (_iteration == _max_its)
301 
302  return SolveState::SUCCESS;
303 }
Moose::GenericType< Real, is_ad > GenericReal
bool _line_search
Whether to use line searches to improve convergence.
virtual void outputIterationStep(std::stringstream *iter_output, const GenericReal< is_ad > &effective_trial_stress, const GenericReal< is_ad > &scalar, const Real reference_residual)
Output information for a single iteration step to build the convergence history of the model...
std::vector< Real > _residual_history
History of residuals used to check whether progress is still being made on decreasing the residual...
virtual GenericReal< is_ad > maximumPermissibleValue(const GenericReal< is_ad > &effective_trial_stress) const
Compute the maximum permissible value of the scalar.
bool convergedAcceptable(const unsigned int it, const Real reference)
Check to see whether the residual is within acceptable convergence limits.
auto raw_value(const Eigen::Map< T > &in)
virtual void preStep(const GenericReal< is_ad > &, const GenericReal< is_ad > &, const GenericReal< is_ad > &)
This method is called before taking a step in the return mapping algorithm.
const unsigned int _max_its
Maximum number of return mapping iterations.
void checkPermissibleRange(GenericReal< is_ad > &scalar, GenericReal< is_ad > &scalar_increment, const GenericReal< is_ad > &scalar_old, const GenericReal< is_ad > min_permissible_scalar, const GenericReal< is_ad > 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...
void updateBounds(const GenericReal< is_ad > &scalar, const GenericReal< is_ad > &residual, const Real init_resid_sign, GenericReal< is_ad > &scalar_upper_bound, GenericReal< is_ad > &scalar_lower_bound, std::stringstream *iter_output)
Update the upper and lower bounds of the root for the effective inelastic strain. ...
virtual GenericReal< is_ad > minimumPermissibleValue(const GenericReal< is_ad > &effective_trial_stress) const
Compute the minimum permissible value of the scalar.
virtual Real computeReferenceResidual(const GenericReal< is_ad > &effective_trial_stress, const GenericReal< is_ad > &scalar)=0
Compute a reference quantity to be used for checking relative convergence.
T sign(T x)
void computeResidualAndDerivativeHelper(const GenericReal< is_ad > &effective_trial_stress, const GenericReal< is_ad > &scalar)
Helper function to compute and set the _residual and _derivative.
bool _bracket_solution
Whether to save upper and lower bounds of root for scalar, and set solution to the midpoint between t...
GenericReal< is_ad > _initial_residual
Residual values, kept as members to retain solver state for summary outputting.
GenericReal< is_ad > _derivative
Derivative of the residual.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static const std::string alpha
Definition: NS.h:134
bool converged(const GenericReal< is_ad > &residual, const Real reference)
Check to see whether the residual is within the convergence limits.
virtual void iterationFinalize(const GenericReal< is_ad > &)
Finalize internal state variables for a model for a given iteration.
T clamp(const T &x, T2 lowerlimit, T2 upperlimit)
bool _check_range
Whether to check to see whether iterative solution is within admissible range, and set within that ra...
virtual GenericReal< is_ad > initialGuess(const GenericReal< is_ad > &)
Compute an initial guess for the value of the scalar.
const std::size_t _num_resids
Number of residuals to be stored in history.

◆ iterationFinalize()

template<bool is_ad>
virtual void SingleVariableReturnMappingSolutionTempl< is_ad >::iterationFinalize ( const GenericReal< is_ad > &  )
inlineprotectedvirtual

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

Parameters
scalarInelastic strain increment magnitude being solved for

Reimplemented in IsotropicPlasticityStressUpdateTempl< is_ad >, and CombinedNonlinearHardeningPlasticityTempl< is_ad >.

Definition at line 128 of file SingleVariableReturnMappingSolution.h.

128 {}

◆ maximumPermissibleValue()

template<bool is_ad>
GenericReal< is_ad > SingleVariableReturnMappingSolutionTempl< is_ad >::maximumPermissibleValue ( const GenericReal< is_ad > &  effective_trial_stress) const
protectedvirtual

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 in RadialReturnStressUpdateTempl< is_ad >, LAROMANCEStressUpdateBaseTempl< is_ad >, and PowerLawCreepTestTempl< is_ad >.

Definition at line 92 of file SingleVariableReturnMappingSolution.C.

94 {
95  return std::numeric_limits<Real>::max();
96 }

◆ minimumPermissibleValue()

template<bool is_ad>
GenericReal< is_ad > SingleVariableReturnMappingSolutionTempl< is_ad >::minimumPermissibleValue ( const GenericReal< is_ad > &  effective_trial_stress) const
protectedvirtual

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 in RadialReturnStressUpdateTempl< is_ad >, and PowerLawCreepTestTempl< is_ad >.

Definition at line 84 of file SingleVariableReturnMappingSolution.C.

86 {
87  return std::numeric_limits<Real>::lowest();
88 }

◆ outputIterationStep()

template<bool is_ad>
void SingleVariableReturnMappingSolutionTempl< is_ad >::outputIterationStep ( std::stringstream *  iter_output,
const GenericReal< is_ad > &  effective_trial_stress,
const GenericReal< is_ad > &  scalar,
const Real  reference_residual 
)
protectedvirtual

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

Parameters
iter_outputOutput stream
effective_trial_stressEffective trial stress
residualCurrent value of the residual
referenceCurrent value of the reference quantity

Reimplemented in LAROMANCEStressUpdateBaseTempl< is_ad >.

Definition at line 419 of file SingleVariableReturnMappingSolution.C.

Referenced by LAROMANCEStressUpdateBaseTempl< is_ad >::outputIterationStep().

424 {
425  if (iter_output)
426  {
427  const unsigned int it = _iteration;
428  const Real residual = MetaPhysicL::raw_value(_residual);
429 
430  *iter_output << " iteration=" << it
431  << " trial_stress=" << MetaPhysicL::raw_value(effective_trial_stress)
432  << " scalar=" << MetaPhysicL::raw_value(scalar) << " residual=" << residual
433  << " ref_res=" << reference_residual
434  << " rel_res=" << std::abs(residual) / reference_residual
435  << " rel_tol=" << _relative_tolerance << " abs_res=" << std::abs(residual)
436  << " abs_tol=" << _absolute_tolerance << '\n';
437  }
438 }
auto raw_value(const Eigen::Map< T > &in)
Real _relative_tolerance
Relative convergence tolerance.
Real _absolute_tolerance
Absolute convergence tolerance.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ outputIterationSummary()

template<bool is_ad>
void SingleVariableReturnMappingSolutionTempl< is_ad >::outputIterationSummary ( std::stringstream *  iter_output,
const unsigned int  total_it 
)
protectedvirtual

Output summary information for the convergence history of the model.

Parameters
iter_outputOutput stream
total_itTotal iteration count

Reimplemented in RadialReturnStressUpdateTempl< is_ad >, LAROMANCEStressUpdateBaseTempl< is_ad >, and ADViscoplasticityStressUpdate.

Definition at line 442 of file SingleVariableReturnMappingSolution.C.

Referenced by ADViscoplasticityStressUpdate::outputIterationSummary(), LAROMANCEStressUpdateBaseTempl< is_ad >::outputIterationSummary(), and RadialReturnStressUpdateTempl< is_ad >::outputIterationSummary().

444 {
445  if (iter_output)
446  *iter_output << "In " << total_it << " iterations the residual went from "
448  << MetaPhysicL::raw_value(_residual) << " in '" << _svrms_name << "'."
449  << std::endl;
450 }
auto raw_value(const Eigen::Map< T > &in)
GenericReal< is_ad > _initial_residual
Residual values, kept as members to retain solver state for summary outputting.
const std::string _svrms_name
MOOSE input name of the object performing the solve.

◆ preStep()

template<bool is_ad>
virtual void SingleVariableReturnMappingSolutionTempl< is_ad >::preStep ( const GenericReal< is_ad > &  ,
const GenericReal< is_ad > &  ,
const GenericReal< is_ad > &   
)
inlineprotectedvirtual

This method is called before taking a step in the return mapping algorithm.

A typical use case is to accumulate the exact algorithmic tangent during return mapping.

Definition at line 118 of file SingleVariableReturnMappingSolution.h.

121  {
122  }

◆ returnMappingSolve()

template<bool is_ad>
void SingleVariableReturnMappingSolutionTempl< is_ad >::returnMappingSolve ( const GenericReal< is_ad > &  effective_trial_stress,
GenericReal< is_ad > &  scalar,
const ConsoleStream console 
)
protected

Perform the return mapping iterations.

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

Definition at line 100 of file SingleVariableReturnMappingSolution.C.

Referenced by ADViscoplasticityStressUpdate::computeInelasticStrainIncrement(), and ComputeSimoHughesJ2PlasticityStress::computeQpPK1Stress().

104 {
105  // construct the stringstream here only if the debug level is set to ALL
106  std::unique_ptr<std::stringstream> iter_output =
108  ? std::make_unique<std::stringstream>()
109  : nullptr;
110 
111  // do the internal solve and capture iteration info during the first round
112  // iff full history output is requested regardless of whether the solve failed or succeeded
113  auto solve_state =
114  internalSolve(effective_trial_stress,
115  scalar,
116  _internal_solve_full_iteration_history ? iter_output.get() : nullptr);
117  if (solve_state != SolveState::SUCCESS &&
119  {
120  // output suppressed by user, throw immediately
122  mooseException("");
123 
124  // user expects some kind of output, if necessary setup output stream now
125  if (!iter_output)
126  iter_output = std::make_unique<std::stringstream>();
127 
128  // add the appropriate error message to the output
129  switch (solve_state)
130  {
131  case SolveState::NAN_INF:
132  *iter_output << "Encountered inf or nan in material return mapping iterations.\n";
133  break;
134 
136  *iter_output << "Exceeded maximum iterations in material return mapping iterations.\n";
137  break;
138 
139  default:
140  mooseError("Unhandled solver state");
141  }
142 
143  // if full history output is only requested for failed solves we have to repeat
144  // the solve a second time
146  internalSolve(effective_trial_stress, scalar, iter_output.get());
147 
148  // Append summary and throw exception
149  outputIterationSummary(iter_output.get(), _iteration);
150  mooseException(iter_output->str());
151  }
152 
154  {
155  // the solve did not fail but the user requested debug output anyways
156  outputIterationSummary(iter_output.get(), _iteration);
157  console << iter_output->str() << std::flush;
158  }
159 }
void mooseError(Args &&... args)
SolveState internalSolve(const GenericReal< is_ad > effective_trial_stress, GenericReal< is_ad > &scalar, std::stringstream *iter_output=nullptr)
Method called from within this class to perform the actual return mappping iterations.
const bool _internal_solve_full_iteration_history
Whether to output iteration information all the time (regardless of whether iterations converge) ...
virtual void outputIterationSummary(std::stringstream *iter_output, const unsigned int total_it)
Output summary information for the convergence history of the model.
enum SingleVariableReturnMappingSolutionTempl::InternalSolveOutput _internal_solve_output_on

◆ updateBounds()

template<bool is_ad>
void SingleVariableReturnMappingSolutionTempl< is_ad >::updateBounds ( const GenericReal< is_ad > &  scalar,
const GenericReal< is_ad > &  residual,
const Real  init_resid_sign,
GenericReal< is_ad > &  scalar_upper_bound,
GenericReal< is_ad > &  scalar_lower_bound,
std::stringstream *  iter_output 
)
private

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 390 of file SingleVariableReturnMappingSolution.C.

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

◆ validParams()

template<bool is_ad>
InputParameters SingleVariableReturnMappingSolutionTempl< is_ad >::validParams ( )
static

Definition at line 26 of file SingleVariableReturnMappingSolution.C.

Referenced by ADViscoplasticityStressUpdate::validParams(), ComputeSimoHughesJ2PlasticityStress::validParams(), and RadialReturnStressUpdateTempl< is_ad >::validParams().

27 {
29  params.addParam<Real>(
30  "relative_tolerance", 1e-8, "Relative convergence tolerance for Newton iteration");
31  params.addParam<Real>(
32  "absolute_tolerance", 1e-11, "Absolute convergence tolerance for Newton iteration");
33  params.addParam<Real>("acceptable_multiplier",
34  10,
35  "Factor applied to relative and absolute "
36  "tolerance for acceptable convergence if "
37  "iterations are no longer making progress");
38 
39  // diagnostic output parameters
40  MooseEnum internal_solve_output_on_enum("never on_error always", "on_error");
41  params.addParam<MooseEnum>("internal_solve_output_on",
42  internal_solve_output_on_enum,
43  "When to output internal Newton solve information");
44  params.addParam<bool>("internal_solve_full_iteration_history",
45  false,
46  "Set true to output full internal Newton iteration history at times "
47  "determined by `internal_solve_output_on`. If false, only a summary is "
48  "output.");
49  params.addParamNamesToGroup("internal_solve_output_on internal_solve_full_iteration_history",
50  "Debug");
51 
52  params.addParam<bool>("automatic_differentiation_return_mapping",
53  false,
54  "Whether to use automatic differentiation to compute the derivative.");
55  return params;
56 }
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
InputParameters emptyInputParameters()
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)

Member Data Documentation

◆ _absolute_tolerance

template<bool is_ad>
Real SingleVariableReturnMappingSolutionTempl< is_ad >::_absolute_tolerance
private

Absolute convergence tolerance.

Definition at line 197 of file SingleVariableReturnMappingSolution.h.

◆ _acceptable_multiplier

template<bool is_ad>
Real SingleVariableReturnMappingSolutionTempl< is_ad >::_acceptable_multiplier
private

Multiplier applied to relative and absolute tolerances for acceptable convergence.

Definition at line 200 of file SingleVariableReturnMappingSolution.h.

◆ _ad_derivative

template<bool is_ad>
const bool SingleVariableReturnMappingSolutionTempl< is_ad >::_ad_derivative
private

Definition at line 205 of file SingleVariableReturnMappingSolution.h.

◆ _bracket_solution

template<bool is_ad>
bool SingleVariableReturnMappingSolutionTempl< is_ad >::_bracket_solution
protected

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

◆ _check_range

template<bool is_ad>
bool SingleVariableReturnMappingSolutionTempl< is_ad >::_check_range
protected

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

Definition at line 138 of file SingleVariableReturnMappingSolution.h.

Referenced by ADViscoplasticityStressUpdate::ADViscoplasticityStressUpdate(), LAROMANCEStressUpdateBaseTempl< is_ad >::LAROMANCEStressUpdateBaseTempl(), and PowerLawCreepTestTempl< is_ad >::PowerLawCreepTestTempl().

◆ _derivative

template<bool is_ad>
GenericReal<is_ad> SingleVariableReturnMappingSolutionTempl< is_ad >::_derivative
private

Derivative of the residual.

Definition at line 222 of file SingleVariableReturnMappingSolution.h.

◆ _initial_residual

template<bool is_ad>
GenericReal<is_ad> SingleVariableReturnMappingSolutionTempl< is_ad >::_initial_residual
private

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

Definition at line 217 of file SingleVariableReturnMappingSolution.h.

◆ _internal_solve_full_iteration_history

template<bool is_ad>
const bool SingleVariableReturnMappingSolutionTempl< is_ad >::_internal_solve_full_iteration_history
private

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

Definition at line 191 of file SingleVariableReturnMappingSolution.h.

◆ _internal_solve_output_on

template<bool is_ad>
enum SingleVariableReturnMappingSolutionTempl::InternalSolveOutput SingleVariableReturnMappingSolutionTempl< is_ad >::_internal_solve_output_on
private

◆ _iteration

template<bool is_ad>
unsigned int SingleVariableReturnMappingSolutionTempl< is_ad >::_iteration
private

iteration number

Definition at line 214 of file SingleVariableReturnMappingSolution.h.

◆ _line_search

template<bool is_ad>
bool SingleVariableReturnMappingSolutionTempl< is_ad >::_line_search
protected

Whether to use line searches to improve convergence.

Definition at line 141 of file SingleVariableReturnMappingSolution.h.

◆ _max_its

template<bool is_ad>
const unsigned int SingleVariableReturnMappingSolutionTempl< is_ad >::_max_its
private

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

◆ _num_resids

template<bool is_ad>
const std::size_t SingleVariableReturnMappingSolutionTempl< is_ad >::_num_resids
private

Number of residuals to be stored in history.

Definition at line 208 of file SingleVariableReturnMappingSolution.h.

◆ _relative_tolerance

template<bool is_ad>
Real SingleVariableReturnMappingSolutionTempl< is_ad >::_relative_tolerance
private

Relative convergence tolerance.

Definition at line 194 of file SingleVariableReturnMappingSolution.h.

◆ _residual

template<bool is_ad>
GenericReal<is_ad> SingleVariableReturnMappingSolutionTempl< is_ad >::_residual
private

Definition at line 218 of file SingleVariableReturnMappingSolution.h.

◆ _residual_history

template<bool is_ad>
std::vector<Real> SingleVariableReturnMappingSolutionTempl< is_ad >::_residual_history
private

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

Definition at line 211 of file SingleVariableReturnMappingSolution.h.

◆ _svrms_name

template<bool is_ad>
const std::string SingleVariableReturnMappingSolutionTempl< is_ad >::_svrms_name
private

MOOSE input name of the object performing the solve.

Definition at line 225 of file SingleVariableReturnMappingSolution.h.


The documentation for this class was generated from the following files: