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

Class container for the QuasiNewton's method. More...

#include <QuasiNewton.hh>

Inherits Optimist::RootFinder::RootFinder< Real, N, DerivedSolver >.

Public Types

using Method = enum class Method : Integer {GOOD = 0, BAD = 1, COMBINED = 2}
 
using Vector = typename RootFinder<Real, N, DerivedSolver>::Vector
 
using Matrix = typename RootFinder<Real, N, DerivedSolver>::Matrix
 
using FunctionWrapper = typename RootFinder<Real, N, DerivedSolver>::FunctionWrapper
 
using JacobianWrapper = typename RootFinder<Real, N, DerivedSolver>::JacobianWrapper
 
- Public Types inherited from Optimist::RootFinder::RootFinder< Real, N, DerivedSolver >
using Vector = typename Solver<Real, N, N, DerivedSolver>::InputType
 
using Matrix = typename Solver<Real, N, N, DerivedSolver>::FirstDerivativeType
 
using Tensor = typename Solver<Real, N, N, DerivedSolver>::SecondDerivativeType
 
using FunctionWrapper = typename Solver<Real, N, N, DerivedSolver>::FunctionWrapper
 
using JacobianWrapper = typename Solver<Real, N, N, DerivedSolver>::FirstDerivativeWrapper
 
using HessianWrapper = typename Solver<Real, N, N, DerivedSolver>::SecondDerivativeWrapper
 

Public Member Functions

 QuasiNewton ()
 
std::string name_impl () const
 
bool solve_impl (FunctionWrapper function, JacobianWrapper jacobian, Vector const &x_ini, Vector &x_sol)
 
void update (Vector const &delta_x_old, Vector const &delta_function_old, Matrix const &jacobian_old, Vector const &delta_x_new, Vector const &delta_function_new, Vector const &function_new, Matrix &jacobian_new)
 
- Public Member Functions inherited from Optimist::RootFinder::RootFinder< Real, N, DerivedSolver >
 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)
 
bool solve (FunctionWrapper function, Vector const &x_ini, Vector &x_sol)
 
bool solve (FunctionWrapper function, JacobianWrapper jacobian, Vector const &x_ini, Vector &x_sol)
 
bool solve (FunctionWrapper function, JacobianWrapper jacobian, HessianWrapper hessian, Vector const &x_ini, Vector &x_sol)
 
- Public Member Functions inherited from Optimist::Solver< Real, N, N, DerivedSolver >
 Solver ()
 
 Solver (FunctionWrapper function, const InputType &x_ini, InputType &x_sol)
 
 Solver (FunctionWrapper function, FirstDerivativeWrapper first_derivative, const InputType &x_ini, InputType &x_sol)
 
 Solver (FunctionWrapper function, FirstDerivativeWrapper first_derivative, SecondDerivativeWrapper second_derivative, const InputType &x_ini, InputType &x_sol)
 
const InputTypelower_bound () const
 
void lower_bound (const InputType &t_lower_bound)
 
const InputTypeupper_bound () const
 
void upper_bound (const InputType &t_upper_bound)
 
void bounds (const InputType &t_lower_bound, const InputType &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 max_function_evaluations () const
 
Integer iterations () const
 
Integer max_iterations () const
 
void max_iterations (Integer t_max_iterations)
 
Real alpha () const
 
void alpha (Real t_alpha)
 
Integer relaxations () const
 
Integer max_relaxations () const
 
void max_relaxations (Integer t_max_relaxations)
 
Real tolerance () const
 
void tolerance (Real t_tolerance)
 
void verbose_mode (bool t_verbose)
 
bool verbose_mode () const
 
void enable_verbose_mode ()
 
void disable_verbose_mode ()
 
void damped_mode (bool t_damped)
 
bool damped_mode () const
 
void enable_damped_mode ()
 
void disable_damped_mode ()
 
std::string task () const
 
void task (std::string t_task)
 
bool converged () const
 
const TraceTypetrace () const
 
std::ostream & ostream () const
 
void ostream (std::ostream &t_ostream)
 
bool solve (FunctionWrapper function, const InputType &x_ini, InputType &x_sol)
 
bool solve (FunctionWrapper function, FirstDerivativeWrapper first_derivative, const InputType &x_ini, InputType &x_sol)
 
bool solve (FunctionWrapper function, FirstDerivativeWrapper first_derivative, SecondDerivativeWrapper second_derivative, const InputType &x_ini, InputType &x_sol)
 
bool rootfind (Function< Real, FunInDim, FunOutDim, DerivedFunction > const &function, const InputType &x_ini, InputType &x_sol)
 
bool optimize (Function< Real, FunInDim, FunOutDim, DerivedFunction > const &function, const InputType &x_ini, InputType &x_sol)
 
std::string name () const
 

Static Public Attributes

static constexpr bool requires_function {true}
 
static constexpr bool requires_first_derivative {true}
 
static constexpr bool requires_second_derivative {false}
 
- Static Public Attributes inherited from Optimist::RootFinder::RootFinder< Real, N, DerivedSolver >
static constexpr bool is_rootfinder {true}
 
static constexpr bool is_optimizer {false}
 
static constexpr bool requires_function {DerivedSolver::requires_function}
 
static constexpr bool requires_first_derivative {DerivedSolver::requires_first_derivative}
 
static constexpr bool requires_second_derivative {DerivedSolver::requires_second_derivative}
 

Private Attributes

Method m_method {Method::COMBINED}
 

Additional Inherited Members

- Public Attributes inherited from Optimist::RootFinder::RootFinder< Real, N, DerivedSolver >
friend Solver< Real, N, N, RootFinder< Real, N, DerivedSolver > >
 
- Protected Types inherited from Optimist::Solver< Real, N, N, DerivedSolver >
using InputType
 
using OutputType
 
using TraceType
 
using FirstDerivativeType
 
using SecondDerivativeType
 
using FunctionWrapper
 
using FirstDerivativeWrapper
 
using SecondDerivativeWrapper
 
- Protected Member Functions inherited from Optimist::RootFinder::RootFinder< Real, N, DerivedSolver >
void evaluate_jacobian (JacobianWrapper jacobian, const Vector &x, Matrix &out)
 
void evaluate_hessian (HessianWrapper hessian, const Vector &x, Matrix &out)
 
- Protected Member Functions inherited from Optimist::Solver< Real, N, N, DerivedSolver >
Integer first_derivative_evaluations () const
 
Integer max_first_derivative_evaluations () const
 
void max_first_derivative_evaluations (Integer first_derivative_evaluations)
 
Integer second_derivative_evaluations () const
 
Integer max_second_derivative_evaluations () const
 
void max_second_derivative_evaluations (Integer second_derivative_evaluations)
 
bool solve (Function< Real, FunInDim, FunOutDim, DerivedFunction > const &function, const InputType &x_ini, InputType &x_sol, bool is_optimization)
 
void reset ()
 
void evaluate_function (FunctionWrapper function, const InputType &x, OutputType &out)
 
void evaluate_first_derivative (FirstDerivativeWrapper function, const InputType &x, FirstDerivativeType &out)
 
void evaluate_second_derivative (SecondDerivativeWrapper function, const InputType &x, SecondDerivativeType &out)
 
void store_trace (const InputType &x)
 
bool damp (FunctionWrapper 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="-")
 
- Protected Attributes inherited from Optimist::Solver< Real, N, N, DerivedSolver >
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>
class Optimist::RootFinder::QuasiNewton< Real, N, DerivedSolver >

Quasi-Newton method

The quasi-Newton is a family of methods that approximate the Jacobian matrix or its inverse. Such methods are based on the linearization of the function \(\mathbf{f}(\mathbf{x})\) around the current iterate \(\mathbf{x}_k\), which leads to the linear system

\[ \mathbf{Jf}_{\mathbf{x}}(\mathbf{x}_k) \mathbf{h}_k = -\mathbf{f}(\mathbf{x}_k) \text{.} \]

Here, the Jacobian matrix \(\mathbf{Jf}_{\mathbf{x}}\) is supposed not to be explicitly available or computationally expensive to compute. For this reason, quasi-Newton methods approximate it iteratively.

Template Parameters
RealScalar number type.
NDimension of the root-finding problem.

Member Typedef Documentation

◆ FunctionWrapper

template<typename Real, Integer N, typename DerivedSolver>
using Optimist::RootFinder::QuasiNewton< Real, N, DerivedSolver >::FunctionWrapper = typename RootFinder<Real, N, DerivedSolver>::FunctionWrapper

◆ JacobianWrapper

template<typename Real, Integer N, typename DerivedSolver>
using Optimist::RootFinder::QuasiNewton< Real, N, DerivedSolver >::JacobianWrapper = typename RootFinder<Real, N, DerivedSolver>::JacobianWrapper

◆ Matrix

template<typename Real, Integer N, typename DerivedSolver>
using Optimist::RootFinder::QuasiNewton< Real, N, DerivedSolver >::Matrix = typename RootFinder<Real, N, DerivedSolver>::Matrix

◆ Method

template<typename Real, Integer N, typename DerivedSolver>
using Optimist::RootFinder::QuasiNewton< Real, N, DerivedSolver >::Method = enum class Method : Integer {GOOD = 0, BAD = 1, COMBINED = 2}

QuasiNewton solver type.

◆ Vector

template<typename Real, Integer N, typename DerivedSolver>
using Optimist::RootFinder::QuasiNewton< Real, N, DerivedSolver >::Vector = typename RootFinder<Real, N, DerivedSolver>::Vector

Constructor & Destructor Documentation

◆ QuasiNewton()

template<typename Real, Integer N, typename DerivedSolver>
Optimist::RootFinder::QuasiNewton< Real, N, DerivedSolver >::QuasiNewton ( )
inline

Class constructor for the QuasiNewton solver.

Member Function Documentation

◆ name_impl()

template<typename Real, Integer N, typename DerivedSolver>
std::string Optimist::RootFinder::QuasiNewton< Real, N, DerivedSolver >::name_impl ( ) const
inline

Get the QuasiNewton solver name.

Returns
The QuasiNewton solver name.

◆ solve_impl()

template<typename Real, Integer N, typename DerivedSolver>
bool Optimist::RootFinder::QuasiNewton< Real, N, DerivedSolver >::solve_impl ( FunctionWrapper function,
JacobianWrapper jacobian,
Vector const & x_ini,
Vector & x_sol )
inline

Solve the nonlinear system of equations \( \mathbf{f}(\mathbf{x}) = 0 \), with \(\mathbf{f}: \mathbb{R}^n \rightarrow \mathbb{R}^n \).

Parameters
[in]functionFunction wrapper.
[in]jacobianJacobian wrapper.
[in]x_iniInitialization point.
[out]x_solSolution point.
Returns
The convergence boolean flag.

◆ update()

template<typename Real, Integer N, typename DerivedSolver>
void Optimist::RootFinder::QuasiNewton< Real, N, DerivedSolver >::update ( Vector const & delta_x_old,
Vector const & delta_function_old,
Matrix const & jacobian_old,
Vector const & delta_x_new,
Vector const & delta_function_new,
Vector const & function_new,
Matrix & jacobian_new )
inline

Jacobian approximation update rule for the QuasiNewton's method.

Parameters
[in]delta_x_oldOld difference between points.
[in]delta_function_oldOld difference between function values.
[in]jacobian_oldOld jacobian approximation.
[in]delta_x_newNew difference between points.
[in]delta_function_newNew difference between function values.
[in]function_newNew function value.
[out]jacobian_newNew jacobian approximation.

Member Data Documentation

◆ m_method

template<typename Real, Integer N, typename DerivedSolver>
Method Optimist::RootFinder::QuasiNewton< Real, N, DerivedSolver >::m_method {Method::COMBINED}
private

QuasiNewton solver type.

◆ requires_first_derivative

template<typename Real, Integer N, typename DerivedSolver>
bool Optimist::RootFinder::QuasiNewton< Real, N, DerivedSolver >::requires_first_derivative {true}
staticconstexpr

◆ requires_function

template<typename Real, Integer N, typename DerivedSolver>
bool Optimist::RootFinder::QuasiNewton< Real, N, DerivedSolver >::requires_function {true}
staticconstexpr

◆ requires_second_derivative

template<typename Real, Integer N, typename DerivedSolver>
bool Optimist::RootFinder::QuasiNewton< Real, N, DerivedSolver >::requires_second_derivative {false}
staticconstexpr

Basic constants.


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