NestedSolveTempl<is_ad> and its instantiations NestedSolve and ADNestedSolve are utility classes that implement a non-linear solver. More...
#include <NestedSolve.h>
Classes | |
struct | CorrespondingJacobianTempl |
Deduce the Jacobian type from the solution type. More... | |
Public Types | |
enum | State { State::NONE, State::CONVERGED_ABS, State::CONVERGED_REL, State::CONVERGED_ACCEPTABLE_ABS, State::CONVERGED_ACCEPTABLE_REL, State::CONVERGED_BOUNDS, State::EXACT_GUESS, State::CONVERGED_XTOL, State::NOT_CONVERGED } |
possible solver states More... | |
template<int N = 0> | |
using | Value = typename std::conditional< N==1, NSReal, typename std::conditional< N==0, NestedSolveTempl< is_ad >::DynamicVector, Eigen::Matrix< NSReal, N, 1 > >::type >::type |
Residual and solution type. More... | |
template<int N = 0> | |
using | Jacobian = typename std::conditional< N==1, NSReal, typename std::conditional< N==0, NestedSolveTempl< is_ad >::DynamicMatrix, Eigen::Matrix< NSReal, N, N > >::type >::type |
Jacobian matrix type. More... | |
using | NSReal = Moose::GenericType< Real, is_ad > |
AD/non-AD switched type shortcuts. More... | |
using | NSRealVectorValue = Moose::GenericType< RealVectorValue, is_ad > |
using | NSRankTwoTensor = Moose::GenericType< RankTwoTensor, is_ad > |
using | DynamicVector = Eigen::Matrix< NSReal, Eigen::Dynamic, 1 > |
Eigen type shortcuts. More... | |
using | DynamicMatrix = Eigen::Matrix< NSReal, Eigen::Dynamic, Eigen::Dynamic > |
template<typename T > | |
using | CorrespondingJacobian = typename NestedSolveInternal::CorrespondingJacobianTempl< T >::type |
Public Member Functions | |
NestedSolveTempl () | |
NestedSolveTempl (const InputParameters ¶ms) | |
template<typename V , typename T > | |
void | nonlinear (V &guess, T &&compute) |
Solve the N*N nonlinear equation system using a built-in Netwon-Raphson loop. More... | |
template<typename V , typename T , typename C > | |
void | nonlinearDamped (V &guess, T &&compute, C &&computeCondition) |
Solve the N*N nonlinear equation system using the damped Netwon-Raphson loop. More... | |
void | setRelativeTolerance (Real rel) |
void | setAbsoluteTolerance (Real abs) |
void | setAcceptableMultiplier (Real am) |
const State & | getState () const |
Get the solver state. More... | |
const std::size_t & | getIterations () |
Get the number of iterations from the last solve. More... | |
template<typename R , typename J > | |
void | nonlinear (NestedSolveTempl< is_ad >::DynamicVector &guess, R &computeResidual, J &computeJacobian) |
template<typename R , typename J > | |
void | nonlinearTight (NestedSolveTempl< is_ad >::DynamicVector &guess, R &computeResidual, J &computeJacobian) |
template<typename R , typename J > | |
void | nonlinear (DynamicVector &guess, R &computeResidual, J &computeJacobian) |
template<typename R , typename J > | |
void | nonlinear (NSReal &guess, R &computeResidual, J &computeJacobian) |
template<typename R , typename J > | |
void | nonlinear (NSRealVectorValue &guess, R &computeResidual, J &computeJacobian) |
template<typename R , typename J > | |
void | nonlinearTight (DynamicVector &guess, R &computeResidual, J &computeJacobian) |
template<typename R , typename J > | |
void | nonlinearTight (NSReal &guess, R &computeResidual, J &computeJacobian) |
template<typename R , typename J > | |
void | nonlinearTight (NSRealVectorValue &guess, R &computeResidual, J &computeJacobian) |
template<typename R , typename J , typename B > | |
void | nonlinearBounded (NSReal &guess, R &computeResidual, J &computeJacobian, B &computeBounds) |
Perform a bounded solve use Eigen::HybridNonLinearSolver. More... | |
Static Public Member Functions | |
static InputParameters | validParams () |
static Real | relativeToleranceDefault () |
default values More... | |
static Real | absoluteToleranceDefault () |
static Real | xToleranceDefault () |
static unsigned int | minIterationsDefault () |
static unsigned int | maxIterationsDefault () |
static Real | acceptableMultiplierDefault () |
static Real | dampingFactorDefault () |
static unsigned int | maxDampingIterationsDefault () |
template<typename V > | |
static Real | normSquare (const V &v) |
Compute squared norm of v (dropping derivatives as this is only for convergence checking) More... | |
static Real | normSquare (const NSReal &v) |
static Real | normSquare (const NSRealVectorValue &v) |
template<typename V > | |
static bool | isRelSmall (const V &a, const V &b, const V &c) |
Check if |a| < |b * c| for all elements in a and b. This checks if 'a' is small relative to. More... | |
static bool | isRelSmall (const NSReal &a, const NSReal &b, const NSReal &c) |
static bool | isRelSmall (const NSRealVectorValue &a, const NSRealVectorValue &b, const NSReal &c) |
static bool | isRelSmall (const DynamicVector &a, const DynamicVector &b, const NSReal &c) |
Public Attributes | |
Real | _relative_tolerance_square |
Real | _absolute_tolerance_square |
Real | _delta_thresh |
Threshold for minimum step size of linear iterations. More... | |
Real | _damping_factor |
Damping factor. More... | |
unsigned int | _max_damping_iterations |
unsigned int | _min_iterations |
unsigned int | _max_iterations |
Real | _acceptable_multiplier |
Protected Member Functions | |
void | sizeItems (const NestedSolveTempl< is_ad >::DynamicVector &guess, NestedSolveTempl< is_ad >::DynamicVector &residual, NestedSolveTempl< is_ad >::DynamicMatrix &jacobian) const |
Size a dynamic Jacobian matrix correctly. More... | |
template<typename V , typename T > | |
void | sizeItems (const V &, V &, T &) const |
Sizing is a no-op for compile time sized types (and scalars) More... | |
bool | isConverged (Real r0_square, Real r_square, bool acceptable) |
Convergence check (updates _status) More... | |
template<typename J , typename V > | |
void | linear (const J &A, V &x, const V &b) const |
Solve A*x=b for x. More... | |
void | linear (const NSReal &A, NSReal &x, const NSReal &b) const |
void | linear (const NSRankTwoTensor &A, NSRealVectorValue &x, const NSRealVectorValue &b) const |
Protected Attributes | |
State | _state |
current solver state More... | |
std::size_t | _n_iterations |
number of nested iterations More... | |
Private Member Functions | |
template<bool store_residual_norm, typename R , typename J > | |
auto | make_adaptor (R &residual, J &jacobian) |
Build a suitable Eigen adaptor functors to interface between moose and Eigen types. More... | |
template<bool store_residual_norm, typename R , typename J > | |
auto | make_real_adaptor (R &residual, J &jacobian) |
template<bool store_residual_norm, typename R , typename J > | |
auto | make_realvector_adaptor (R &residual, J &jacobian) |
NestedSolveTempl<is_ad> and its instantiations NestedSolve and ADNestedSolve are utility classes that implement a non-linear solver.
These classes can be instantiated in any MOOSE object to perform local Newton-Raphson solves.
Definition at line 42 of file NestedSolve.h.
using NestedSolveTempl< is_ad >::CorrespondingJacobian = typename NestedSolveInternal::CorrespondingJacobianTempl<T>::type |
Definition at line 65 of file NestedSolve.h.
using NestedSolveTempl< is_ad >::DynamicMatrix = Eigen::Matrix<NSReal, Eigen::Dynamic, Eigen::Dynamic> |
Definition at line 58 of file NestedSolve.h.
using NestedSolveTempl< is_ad >::DynamicVector = Eigen::Matrix<NSReal, Eigen::Dynamic, 1> |
Eigen type shortcuts.
Definition at line 57 of file NestedSolve.h.
using NestedSolveTempl< is_ad >::Jacobian = typename std::conditional<N == 1, NSReal, typename std::conditional<N == 0, NestedSolveTempl<is_ad>::DynamicMatrix, Eigen::Matrix<NSReal, N, N> >::type>::type |
Jacobian matrix type.
Definition at line 84 of file NestedSolve.h.
using NestedSolveTempl< is_ad >::NSRankTwoTensor = Moose::GenericType<RankTwoTensor, is_ad> |
Definition at line 53 of file NestedSolve.h.
using NestedSolveTempl< is_ad >::NSReal = Moose::GenericType<Real, is_ad> |
AD/non-AD switched type shortcuts.
Definition at line 51 of file NestedSolve.h.
using NestedSolveTempl< is_ad >::NSRealVectorValue = Moose::GenericType<RealVectorValue, is_ad> |
Definition at line 52 of file NestedSolve.h.
using NestedSolveTempl< is_ad >::Value = typename std::conditional<N == 1, NSReal, typename std::conditional<N == 0, NestedSolveTempl<is_ad>::DynamicVector, Eigen::Matrix<NSReal, N, 1> >::type>::type |
Residual and solution type.
Definition at line 75 of file NestedSolve.h.
|
strong |
possible solver states
Enumerator | |
---|---|
NONE | |
CONVERGED_ABS | |
CONVERGED_REL | |
CONVERGED_ACCEPTABLE_ABS | |
CONVERGED_ACCEPTABLE_REL | |
CONVERGED_BOUNDS | |
EXACT_GUESS | |
CONVERGED_XTOL | |
NOT_CONVERGED |
Definition at line 147 of file NestedSolve.h.
NestedSolveTempl< is_ad >::NestedSolveTempl | ( | ) |
Definition at line 52 of file NestedSolve.C.
NestedSolveTempl< is_ad >::NestedSolveTempl | ( | const InputParameters & | params | ) |
Definition at line 67 of file NestedSolve.C.
|
inlinestatic |
Definition at line 121 of file NestedSolve.h.
|
inlinestatic |
Definition at line 125 of file NestedSolve.h.
|
inlinestatic |
Definition at line 126 of file NestedSolve.h.
|
inline |
Get the number of iterations from the last solve.
Definition at line 163 of file NestedSolve.h.
|
inline |
|
protected |
Convergence check (updates _status)
Definition at line 779 of file NestedSolve.h.
|
static |
Check if |a| < |b * c| for all elements in a and b. This checks if 'a' is small relative to.
|
static |
Definition at line 117 of file NestedSolve.C.
|
static |
Definition at line 125 of file NestedSolve.C.
|
static |
Definition at line 140 of file NestedSolve.C.
|
protected |
Solve A*x=b for x.
Definition at line 546 of file NestedSolve.h.
|
inlineprotected |
Definition at line 203 of file NestedSolve.h.
|
protected |
Definition at line 94 of file NestedSolve.C.
|
private |
Build a suitable Eigen adaptor functors to interface between moose and Eigen types.
Definition at line 753 of file NestedSolve.h.
|
private |
Definition at line 761 of file NestedSolve.h.
|
private |
Definition at line 770 of file NestedSolve.h.
|
inlinestatic |
Definition at line 127 of file NestedSolve.h.
|
inlinestatic |
Definition at line 124 of file NestedSolve.h.
|
inlinestatic |
Definition at line 123 of file NestedSolve.h.
void NestedSolveTempl< is_ad >::nonlinear | ( | V & | guess, |
T && | compute | ||
) |
Solve the N*N nonlinear equation system using a built-in Netwon-Raphson loop.
if we exceed the max iterations, we could still be converged (considering the acceptable multiplier)
Definition at line 427 of file NestedSolve.h.
void NestedSolveTempl< is_ad >::nonlinear | ( | DynamicVector & | guess, |
R & | computeResidual, | ||
J & | computeJacobian | ||
) |
The separate residual/Jacobian functor versions use Eigen::HybridNonLinearSolver with a custom backwards compatible convergence check that allows for looser tolerances.
void NestedSolveTempl< is_ad >::nonlinear | ( | NSReal & | guess, |
R & | computeResidual, | ||
J & | computeJacobian | ||
) |
if we exceed the max iterations, we could still be converged (considering the acceptable multiplier)
Definition at line 242 of file NestedSolve.h.
void NestedSolveTempl< is_ad >::nonlinear | ( | NSRealVectorValue & | guess, |
R & | computeResidual, | ||
J & | computeJacobian | ||
) |
Definition at line 288 of file NestedSolve.h.
void NestedSolveTempl< is_ad >::nonlinear | ( | NestedSolveTempl< is_ad >::DynamicVector & | guess, |
R & | computeResidual, | ||
J & | computeJacobian | ||
) |
Definition at line 224 of file NestedSolve.h.
void NestedSolveTempl< is_ad >::nonlinearBounded | ( | NSReal & | guess, |
R & | computeResidual, | ||
J & | computeJacobian, | ||
B & | computeBounds | ||
) |
Perform a bounded solve use Eigen::HybridNonLinearSolver.
Definition at line 377 of file NestedSolve.h.
void NestedSolveTempl< is_ad >::nonlinearDamped | ( | V & | guess, |
T && | compute, | ||
C && | computeCondition | ||
) |
Solve the N*N nonlinear equation system using the damped Netwon-Raphson loop.
Definition at line 482 of file NestedSolve.h.
void NestedSolveTempl< is_ad >::nonlinearTight | ( | DynamicVector & | guess, |
R & | computeResidual, | ||
J & | computeJacobian | ||
) |
The separate residual/Jacobian functor versions use Eigen::HybridNonLinearSolver with the built-in convergence check to tight numerical tolerance.
void NestedSolveTempl< is_ad >::nonlinearTight | ( | NSReal & | guess, |
R & | computeResidual, | ||
J & | computeJacobian | ||
) |
Definition at line 331 of file NestedSolve.h.
void NestedSolveTempl< is_ad >::nonlinearTight | ( | NSRealVectorValue & | guess, |
R & | computeResidual, | ||
J & | computeJacobian | ||
) |
Definition at line 352 of file NestedSolve.h.
void NestedSolveTempl< is_ad >::nonlinearTight | ( | NestedSolveTempl< is_ad >::DynamicVector & | guess, |
R & | computeResidual, | ||
J & | computeJacobian | ||
) |
Definition at line 313 of file NestedSolve.h.
|
static |
Compute squared norm of v (dropping derivatives as this is only for convergence checking)
Definition at line 555 of file NestedSolve.h.
Referenced by NestedSolveInternal::DynamicMatrixEigenAdaptorFunctor< is_ad, R, J, store_residual_norm >::operator()(), NestedSolveInternal::RealEigenAdaptorFunctor< is_ad, R, J, store_residual_norm >::operator()(), and NestedSolveInternal::RealVectorEigenAdaptorFunctor< is_ad, R, J, store_residual_norm >::operator()().
|
static |
Definition at line 103 of file NestedSolve.C.
|
static |
Definition at line 110 of file NestedSolve.C.
|
inlinestatic |
|
inline |
Definition at line 131 of file NestedSolve.h.
|
inline |
Definition at line 132 of file NestedSolve.h.
|
inline |
Definition at line 130 of file NestedSolve.h.
|
protected |
|
inlineprotected |
Sizing is a no-op for compile time sized types (and scalars)
Definition at line 196 of file NestedSolve.h.
|
static |
Definition at line 16 of file NestedSolve.C.
|
inlinestatic |
Definition at line 122 of file NestedSolve.h.
Real NestedSolveTempl< is_ad >::_absolute_tolerance_square |
Definition at line 135 of file NestedSolve.h.
Referenced by NestedSolveTempl< is_ad >::setAbsoluteTolerance().
Real NestedSolveTempl< is_ad >::_acceptable_multiplier |
Definition at line 144 of file NestedSolve.h.
Referenced by NestedSolveTempl< is_ad >::setAcceptableMultiplier().
Real NestedSolveTempl< is_ad >::_damping_factor |
Damping factor.
Definition at line 139 of file NestedSolve.h.
Real NestedSolveTempl< is_ad >::_delta_thresh |
Threshold for minimum step size of linear iterations.
Definition at line 137 of file NestedSolve.h.
unsigned int NestedSolveTempl< is_ad >::_max_damping_iterations |
Definition at line 140 of file NestedSolve.h.
unsigned int NestedSolveTempl< is_ad >::_max_iterations |
Definition at line 143 of file NestedSolve.h.
unsigned int NestedSolveTempl< is_ad >::_min_iterations |
Definition at line 142 of file NestedSolve.h.
|
protected |
number of nested iterations
Definition at line 187 of file NestedSolve.h.
Referenced by NestedSolveTempl< is_ad >::getIterations().
Real NestedSolveTempl< is_ad >::_relative_tolerance_square |
Definition at line 134 of file NestedSolve.h.
Referenced by NestedSolveTempl< is_ad >::setRelativeTolerance().
|
protected |
current solver state
Definition at line 184 of file NestedSolve.h.
Referenced by NestedSolveTempl< is_ad >::getState().