Optimist  0.0.0
A C++ library for optimization
Loading...
Searching...
No Matches
Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen > Class Template Reference

Class container for the multi-dimensional root finder. More...

#include <RootFinder.hh>

Inherits Optimist::SolverBase< Real, N, N, DerivedSolver, false >.

Public Types

using Vector = typename SolverBase<Real, N, N, DerivedSolver, ForceEigen>::InputType
using Matrix = typename SolverBase<Real, N, N, DerivedSolver, ForceEigen>::FirstDerivativeType
using Tensor = typename SolverBase<Real, N, N, DerivedSolver, ForceEigen>::SecondDerivativeType
Public Types inherited from Optimist::SolverBase< Real, N, N, DerivedSolver, false >
using InputType
using OutputType
using TraceType
using FirstDerivativeType
using SecondDerivativeType

Public Member Functions

 RootFinder ()
std::string name () const
Integer jacobian_evaluations () const
Integer max_jacobian_evaluations () const
void max_jacobian_evaluations (Integer t_jacobian_evaluations)
Integer hessian_evaluations () const
Integer max_hessian_evaluations () const
void max_hessian_evaluations (Integer t_hessian_evaluations)
template<typename FunctionLambda>
bool solve (FunctionLambda &&function, Vector const &x_ini, Vector &x_sol)
template<typename FunctionLambda, typename JacobianLambda>
bool solve (FunctionLambda &&function, JacobianLambda &&jacobian, Vector const &x_ini, Vector &x_sol)
template<typename FunctionLambda, typename JacobianLambda, typename HessianLambda>
bool solve (FunctionLambda &&function, JacobianLambda &&jacobian, HessianLambda &&hessian, Vector const &x_ini, Vector &x_sol)
Public Member Functions inherited from Optimist::SolverBase< Real, N, N, DerivedSolver, false >
 SolverBase ()
InputType const & lower_bound () const
InputType const & upper_bound () const
void bounds (InputType const &t_lower_bound, InputType const &t_upper_bound)
constexpr Integer input_dimension () const
constexpr Integer output_dimension () const
Integer function_evaluations () const
void max_function_evaluations (Integer t_max_function_evaluations)
Integer iterations () const
Integer max_iterations () const
Real alpha () const
Integer relaxations () const
Integer max_relaxations () const
Real tolerance () const
void verbose_mode (bool t_verbose)
void enable_verbose_mode ()
void disable_verbose_mode ()
void damped_mode (bool t_damped)
void enable_damped_mode ()
void disable_damped_mode ()
std::string task () const
bool converged () const
const TraceTypetrace () const
std::ostream & ostream () const
bool solve (FunctionLambda &&function, InputType const &x_ini, InputType &x_sol)
bool rootfind (FunctionBase< Real, FunInDim, FunOutDim, DerivedFunction, ForceEigen &&FunOutDim==1 > const &function, InputType const &x_ini, InputType &x_sol)
bool optimize (FunctionBase< Real, FunInDim, FunOutDim, DerivedFunction, ForceEigen &&FunOutDim==1 > const &function, InputType const &x_ini, InputType &x_sol)
std::string name () const

Static Public Attributes

static constexpr bool is_rootfinder {true}
static constexpr bool is_optimizer {false}

Protected Member Functions

template<typename JacobianLambda>
bool evaluate_jacobian (JacobianLambda &&jacobian, Vector const &x, Matrix &out)
template<typename HessianLambda>
bool evaluate_hessian (HessianLambda &&hessian, Vector const &x, Matrix &out)
Protected Member Functions inherited from Optimist::SolverBase< Real, N, N, DerivedSolver, false >
Integer first_derivative_evaluations () const
Integer max_first_derivative_evaluations () const
Integer second_derivative_evaluations () const
Integer max_second_derivative_evaluations () const
void reset ()
bool evaluate_function (FunctionLambda &&function, InputType const &x, OutputType &out)
bool evaluate_first_derivative (FirstDerivativeLambda &&function, InputType const &x, FirstDerivativeType &out)
bool evaluate_second_derivative (SecondDerivativeLambda &&function, InputType const &x, SecondDerivativeType &out)
void store_trace (InputType const &x)
bool damp (FunctionLambda &&function, InputType const &x_old, InputType const &function_old, InputType const &step_old, InputType &x_new, InputType &function_new, InputType &step_new)
void header ()
void bottom ()
void info (Real residuals, std::string const &notes="-")

Friends

class SolverBase< Real, N, N, RootFinder< Real, N, DerivedSolver, ForceEigen > >

Additional Inherited Members

Protected Attributes inherited from Optimist::SolverBase< Real, N, N, DerivedSolver, false >
InputType m_lower_bound
InputType m_upper_bound
Integer m_function_evaluations
Integer m_first_derivative_evaluations
Integer m_second_derivative_evaluations
Integer m_max_function_evaluations
Integer m_max_first_derivative_evaluations
Integer m_max_second_derivative_evaluations
Integer m_iterations
Integer m_max_iterations
Real m_alpha
Integer m_relaxations
Integer m_max_relaxations
Real m_tolerance
bool m_verbose
bool m_damped
std::ostream * m_ostream
std::string m_task
bool m_converged
TraceType m_trace

Detailed Description

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
class Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >

Root finder

This section describes the multi-dimensional root-finders implemented in Optimist. The available root-finding solvers are Newton (derivative) and quasi-Newton (non-derivative) methods. Derivative methods employ the function's derivative to find the root with high accuracy, while non-derivative methods approximate the derivative for improved efficiency in certain scenarios.

Here, the solvers are implemented for solving problems of the form

\[ \mathbf{f}(\mathbf{x}) = \mathbf{0} \text{,} \quad \text{with} \quad \mathbf{f}: \mathbb{R}^{n} \rightarrow \mathbb{R}^{n} \text{,} \]

which consist in finding the root of the function \(\mathbf{f}\) by iteratively updating the current iterate \(\mathbf{x}_k\) until convergence is achieved. The solvers require the function \(\mathbf{f}\) and its Jacobian matrix \(\mathbf{Jf}_{\mathbf{x}}\) to be provided by the user. The Jacobian matrix can be computed analytically or numerically, depending on the problem's complexity and the user's preference. Alternatively, the Jacobian can be approximated numerically using finite differences, depending on the problem's complexity and the user's preference.

Affine-invariant step

By default, the multi-dimensional root-finders in Optimist use optionally an affine-invariant advancing step to ensure convergence. The generic advancing step \(\mathbf{h}_k\) is computed as

\[ \mathbf{x}_{k+1} = \mathbf{x}_k + \alpha_k \mathbf{h}_k \text{,} \]

where \(\alpha_k\) is a damping coefficient that ensures convergence by satisfying

\[ \left\|\mathbf{Jf}_{\mathbf{x}}(\mathbf{x}_k)^{-1} \mathbf{f}(\mathbf{x}_{k+1})\right\| \leq \left(1 - \displaystyle\frac{\alpha_k}{2}\right) \left\|\mathbf{Jf}_{\mathbf{x}}(\mathbf{x}_k)^{-1} \mathbf{f}(\mathbf{x}_k)\right\| = \left(1 - \displaystyle\frac{\alpha_k}{2} \right) \left\|\mathbf{h}_k\right\| \text{.} \]

For detailed information on the affine-invariant Newton's method, refer to this link, or the works by P. Deuflhard.

Template Parameters
RealScalar number type.
NThe dimension of the root-finding problem.
DerivedSolverDerived solver class.
ForceEigenForce the use of Eigen types for input and output.

Member Typedef Documentation

◆ Matrix

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
using Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::Matrix = typename SolverBase<Real, N, N, DerivedSolver, ForceEigen>::FirstDerivativeType

◆ Tensor

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
using Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::Tensor = typename SolverBase<Real, N, N, DerivedSolver, ForceEigen>::SecondDerivativeType

◆ Vector

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
using Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::Vector = typename SolverBase<Real, N, N, DerivedSolver, ForceEigen>::InputType

Constructor & Destructor Documentation

◆ RootFinder()

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::RootFinder ( )
inline

Class constructor for the multi-dimensional root finder.

Member Function Documentation

◆ evaluate_hessian()

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
template<typename HessianLambda>
bool Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::evaluate_hessian ( HessianLambda && hessian,
Vector const & x,
Matrix & out )
inlineprotected

Evaluate the Hessian function.

Template Parameters
HessianLambdaThe Hessian lambda function type.
Parameters
[in]hessianHessian lambda function.
[in]xInput point.
[out]outHessian value.
Returns
The boolean flag for successful evaluation.

◆ evaluate_jacobian()

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
template<typename JacobianLambda>
bool Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::evaluate_jacobian ( JacobianLambda && jacobian,
Vector const & x,
Matrix & out )
inlineprotected

Evaluate the Jacobian function.

Template Parameters
JacobianLambdaThe Jacobian lambda function type.
Parameters
[in]jacobianJacobian lambda function.
[in]xInput point.
[out]outJacobian value.
Returns
The boolean flag for successful evaluation.

◆ hessian_evaluations()

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
Integer Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::hessian_evaluations ( ) const
inline

Get the number of Hessian evaluations.

Returns
The number of Hessian evaluations.

◆ jacobian_evaluations()

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
Integer Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::jacobian_evaluations ( ) const
inline

Get the number of Jacobian evaluations.

Returns
The number of Jacobian evaluations.

◆ max_hessian_evaluations() [1/2]

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
Integer Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::max_hessian_evaluations ( ) const
inline

Get the number of maximum allowed Hessian evaluations.

Returns
The number of maximum allowed Hessian evaluations.

◆ max_hessian_evaluations() [2/2]

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
void Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::max_hessian_evaluations ( Integer t_hessian_evaluations)
inline

Set the number of maximum allowed Hessian evaluations.

Parameters
[in]t_hessian_evaluationsThe number of maximum allowed Hessian evaluations.

◆ max_jacobian_evaluations() [1/2]

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
Integer Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::max_jacobian_evaluations ( ) const
inline

Get the number of maximum allowed Jacobian evaluations.

Returns
The number of maximum allowed Jacobian evaluations.

◆ max_jacobian_evaluations() [2/2]

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
void Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::max_jacobian_evaluations ( Integer t_jacobian_evaluations)
inline

Set the number of maximum allowed Jacobian evaluations.

Parameters
[in]t_jacobian_evaluationsThe number of maximum allowed Jacobian evaluations.

◆ name()

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
std::string Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::name ( ) const
inline

Get the solver name.

Returns
The solver name.

◆ solve() [1/3]

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
template<typename FunctionLambda, typename JacobianLambda, typename HessianLambda>
bool Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::solve ( FunctionLambda && function,
JacobianLambda && jacobian,
HessianLambda && hessian,
Vector const & x_ini,
Vector & x_sol )
inline

Solve the root-finding problem given the function, and its Jacobian and Hessian.

Template Parameters
FunctionLambdaThe lambda function type.
JacobianLambdaThe Jacobian lambda function type.
HessianLambdaThe Hessian lambda function type.
Parameters
[in]functionFunction lambda.
[in]jacobianThe Jacobian lambda function.
[in]hessianThe Hessian lambda function.
[in]x_iniInitialization point.
[out]x_solSolution point.
Returns
The convergence boolean flag.

◆ solve() [2/3]

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
template<typename FunctionLambda, typename JacobianLambda>
bool Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::solve ( FunctionLambda && function,
JacobianLambda && jacobian,
Vector const & x_ini,
Vector & x_sol )
inline

Solve the root-finding problem given the function, and its Jacobian.

Template Parameters
FunctionLambdaThe lambda function type.
JacobianLambdaThe Jacobian lambda function type.
Parameters
[in]functionFunction lambda.
[in]jacobianThe Jacobian lambda function.
[in]x_iniInitialization point.
[out]x_solSolution point.
Returns
The convergence boolean flag.

◆ solve() [3/3]

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
template<typename FunctionLambda>
bool Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::solve ( FunctionLambda && function,
Vector const & x_ini,
Vector & x_sol )
inline

Solve the root-finding problem given the function, and without derivatives.

Template Parameters
FunctionLambdaThe lambda function type.
Parameters
[in]functionFunction lambda.
[in]x_iniInitialization point.
[out]x_solSolution point.
Returns
The convergence boolean flag.

◆ SolverBase< Real, N, N, RootFinder< Real, N, DerivedSolver, ForceEigen > >

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
friend class SolverBase< Real, N, N, RootFinder< Real, N, DerivedSolver, ForceEigen > >
friend

Member Data Documentation

◆ is_optimizer

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
bool Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::is_optimizer {false}
staticconstexpr

◆ is_rootfinder

template<typename Real, Integer N, typename DerivedSolver, bool ForceEigen = false>
bool Optimist::RootFinder::RootFinder< Real, N, DerivedSolver, ForceEigen >::is_rootfinder {true}
staticconstexpr

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