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

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

#include <Broyden.hh>

Inherits Optimist::RootFinder::QuasiNewton< Real, N, Broyden< Real, N > >.

Public Types

using Method = enum class Method : Integer {GOOD = 0, BAD = 1, COMBINED = 2}
using Vector = typename QuasiNewton<Real, N, Broyden<Real, N>>::Vector
using Matrix = typename QuasiNewton<Real, N, Broyden<Real, N>>::Matrix
using FunctionWrapper = typename QuasiNewton<Real, N, Broyden<Real, N>>::FunctionWrapper
using JacobianWrapper = typename QuasiNewton<Real, N, Broyden<Real, N>>::JacobianWrapper
Public Types inherited from Optimist::RootFinder::QuasiNewton< Real, N, Broyden< Real, N > >
using Method
using Vector
using Matrix
using FunctionWrapper
using JacobianWrapper
Public Types inherited from Optimist::RootFinder::RootFinder< Real, N, Broyden< Real, N >, true >
using Vector
using Matrix
using Tensor
using FunctionWrapper
using JacobianWrapper
using HessianWrapper

Public Member Functions

 Broyden ()
std::string name_impl () const
Method method () const
void method (Method t_method)
void enable_good_method ()
void enable_bad_method ()
void enable_combined_method ()
void set_method (Method t_method)
void update_impl (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 &, Matrix &jacobian_new)
Public Member Functions inherited from Optimist::RootFinder::QuasiNewton< Real, N, Broyden< Real, N > >
 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, Broyden< Real, N >, true >
 RootFinder ()
std::string name () const
Integer jacobian_evaluations () const
Integer max_jacobian_evaluations () const
Integer hessian_evaluations () const
Integer max_hessian_evaluations () const
bool solve (FunctionWrapper function, Vector const &x_ini, Vector &x_sol)
Public Member Functions inherited from Optimist::SolverBase< Real, N, N, Broyden< Real, N >, ForceEigen >
 SolverBase ()
const InputTypelower_bound () const
const InputTypeupper_bound () const
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 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 (FunctionWrapper function, const InputType &x_ini, InputType &x_sol)
bool rootfind (FunctionBase< Real, FunInDim, FunOutDim, DerivedFunction, ForceEigen &&FunOutDim==1 > const &function, const InputType &x_ini, InputType &x_sol)
bool optimize (FunctionBase< Real, FunInDim, FunOutDim, DerivedFunction, ForceEigen &&FunOutDim==1 > 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::QuasiNewton< Real, N, Broyden< Real, N > >
static constexpr bool requires_function
static constexpr bool requires_first_derivative
static constexpr bool requires_second_derivative
Static Public Attributes inherited from Optimist::RootFinder::RootFinder< Real, N, Broyden< Real, N >, true >
static constexpr bool is_rootfinder
static constexpr bool is_optimizer
static constexpr bool requires_function
static constexpr bool requires_first_derivative
static constexpr bool requires_second_derivative

Private Attributes

Method m_method {Method::COMBINED}

Additional Inherited Members

Protected Types inherited from Optimist::SolverBase< Real, N, N, Broyden< Real, N >, ForceEigen >
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, Broyden< Real, N >, true >
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::SolverBase< Real, N, N, Broyden< Real, N >, ForceEigen >
Integer first_derivative_evaluations () const
Integer max_first_derivative_evaluations () const
Integer second_derivative_evaluations () const
Integer max_second_derivative_evaluations () const
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::SolverBase< Real, N, N, Broyden< Real, N >, ForceEigen >
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>
class Optimist::RootFinder::Broyden< Real, N >

Broyden's method

The Broyden family includes the good, bad, and combined methods. These methods are quasi-Newton approaches that approximate the Jacobian matrix or its inverse. The combined method integrates the updates of the good and bad methods, offering an alternative approach to the classical Newton's method. The Broyden's 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, Broyden's methods approximate the Jacobian matrix iteratively. The Broyden's method is defined as

\[ \mathbf{H}_k(\mathbf{x}_k) \mathbf{h}_k = -\mathbf{f}(\mathbf{x}_k) \text{,} \]

where \(\mathbf{H}_k\) is the (inverse) Jacobian approximation at the \(k\)-th iteration. The generic advancing step is then defined as

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

where \(\alpha_k\) is a damping coefficient that ensures affine-invariant criteria is satisfied.

What distinguishes Broyden's combined method from the generic Broyden's method is the update of the Jacobian approximation. The Broyden's combined method uses the following update rule

\[ \begin{cases} \mathbf{H}_{k+1}^{-1} = \mathbf{H}_k^{-1} - \displaystyle\frac{\mathbf{H}_k^{-1} \Delta\mathbf{f}_k - \Delta\mathbf{x}_k}{\mathbf{C}_g \Delta\mathbf{f}_k} \mathbf{C}_g & \left\| \displaystyle\frac{\Delta\mathbf{x}_k^\top \Delta\mathbf{x}_{k-1}}{\Delta\mathbf{x}_k^\top \mathbf{H}_k^{-1} \Delta\mathbf{f}_k} \right\| < \left\| \displaystyle\frac{\Delta\mathbf{f}_k^\top \Delta\mathbf{f}_{k-1}}{\Delta\mathbf{f}_k^\top\Delta\mathbf{f}_k} \right\| \\ \mathbf{H}_{k+1}^{-1} = \mathbf{H}_k^{-1} - \displaystyle\frac{\mathbf{H}_k^{-1} \Delta\mathbf{f}_k-\Delta\mathbf{x}_k}{\mathbf{C}_b \Delta\mathbf{f}_k} \mathbf{C}_b & \text{otherwise} \end{cases} \text{,} \]

with \(\mathbf{C}_g = \Delta\mathbf{x}_k^\top \mathbf{H}_k^{-1}\), \(\mathbf{C}_b = \Delta\mathbf{f}_k^\top\), \(\Delta\mathbf{x}_k = \mathbf{x}_{k+1} - \mathbf{x}_k\), and \(\Delta\mathbf{f}_k = \mathbf{f}(\mathbf{x}_{k+1}) - \mathbf{f}(\mathbf{x}_k)\). Such a rule allows Broyden's combined method to adapt to the problem's behavior, switching between the good and bad methods based on the convergence history.

For more details on the Broyden's methods refer to the references: Broyden C.G., A class of methods for solving nonlinear simultaneous equations, Mathematics of Computation, 19 (1965), pp. 577-593, 10.1090/s0025-5718-1965-0198670-6, and Martinez J., Ochi L., Sobre dois métodos de broyden, Matemática Aplicada e Computacional, IV Congresso Nacional de Matemática Aplicada e Computacional, Rio de Janeiro, Brasil, setembro de 1981.

NOTE: In Optimist, the implemented Broyden's class for the combined method can be used as a good or bad solver by setting the appropriate parameters.

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

Member Typedef Documentation

◆ FunctionWrapper

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

◆ JacobianWrapper

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

◆ Matrix

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

◆ Method

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

Broyden solver type.

◆ Vector

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

Constructor & Destructor Documentation

◆ Broyden()

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

Class constructor for the Broyden solver.

Member Function Documentation

◆ enable_bad_method()

template<typename Real, Integer N>
void Optimist::RootFinder::Broyden< Real, N >::enable_bad_method ( )
inline

Enable the bad Broyden solver method.

◆ enable_combined_method()

template<typename Real, Integer N>
void Optimist::RootFinder::Broyden< Real, N >::enable_combined_method ( )
inline

Enable the combined Broyden solver method.

◆ enable_good_method()

template<typename Real, Integer N>
void Optimist::RootFinder::Broyden< Real, N >::enable_good_method ( )
inline

Enable the good Broyden solver method.

◆ method() [1/2]

template<typename Real, Integer N>
Method Optimist::RootFinder::Broyden< Real, N >::method ( ) const
inline

Get the enumeration type of the Broyden solver method.

Returns
The Broyden solver enumeration type.

◆ method() [2/2]

template<typename Real, Integer N>
void Optimist::RootFinder::Broyden< Real, N >::method ( Method t_method)
inline

Set the enumeration type of the Broyden solver method.

Parameters
[in]t_methodThe Broyden solver enumeration type.

◆ name_impl()

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

Get the Broyden solver name.

Returns
The Broyden solver name.

◆ set_method()

template<typename Real, Integer N>
void Optimist::RootFinder::Broyden< Real, N >::set_method ( Method t_method)
inline

Set the Broyden solver type.

Parameters
[in]t_methodThe Broyden solver type enumeration.

◆ update_impl()

template<typename Real, Integer N>
void Optimist::RootFinder::Broyden< Real, N >::update_impl ( 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 & ,
Matrix & jacobian_new )
inline

Jacobian approximation update rule for the Broyden'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>
Method Optimist::RootFinder::Broyden< Real, N >::m_method {Method::COMBINED}
private

Broyden solver type.

◆ requires_first_derivative

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

◆ requires_function

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

◆ requires_second_derivative

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

Basic constants.


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