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

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

#include <Broyden.hh>

Inherits Optimist::RootFinder::QuasiNewton< Vector, Broyden< Vector > >.

Public Types

using VectorTrait = TypeTrait<Vector>
using Scalar = typename TypeTrait<Vector>::Scalar
using Method
Public Types inherited from Optimist::RootFinder::QuasiNewton< Vector, Broyden< Vector > >
using VectorTrait
using Scalar
Public Types inherited from Optimist::RootFinder::RootFinder< Vector, Broyden< Vector > >
using Scalar
using Input
using Output
Public Types inherited from Optimist::SolverBase< Vector, Vector, Broyden< Vector > >
using InputTrait
using OutputTrait
using Scalar
using FirstDerivative
using SecondDerivative

Public Member Functions

 Broyden ()
constexpr 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 (const Vector &delta_x_old, const Vector &delta_function_old, const FirstDerivative &jacobian_old, const Vector &delta_x_new, const Vector &delta_function_new, const Vector &, FirstDerivative &jacobian_new)
Public Member Functions inherited from Optimist::RootFinder::QuasiNewton< Vector, Broyden< Vector > >
 QuasiNewton ()
constexpr std::string name_impl () const
bool solve_impl (FunctionLambda &&function, JacobianLambda &&jacobian, const Vector &x_ini, Vector &x_sol)
void update (const Vector &delta_x_old, const Vector &delta_function_old, const FirstDerivative &jacobian_old, const Vector &delta_x_new, const Vector &delta_function_new, const Vector &function_new, FirstDerivative &jacobian_new)
Public Member Functions inherited from Optimist::RootFinder::RootFinder< Vector, Broyden< Vector > >
 RootFinder ()
constexpr std::string name () const
Integer jacobian_evaluations () const
Integer max_jacobian_evaluations () const
Integer hessian_evaluations () const
Integer max_hessian_evaluations () const
bool solve (FunctionLambda &&function, const Input &x_ini, Output &x_sol)
Public Member Functions inherited from Optimist::SolverBase< Vector, Vector, Broyden< Vector > >
 SolverBase ()
void reset_bounds (const Integer n=InputTrait::IsDynamic ? 0 :InputTrait::Dimension)
const Vector & lower_bound () const
const Vector & upper_bound () const
void bounds (const Vector &t_lower_bound, const Vector &t_upper_bound)
constexpr Integer input_dimension () const
constexpr Integer output_dimension () const
Integer function_evaluations () const
void max_function_evaluations (const Integer t_max_function_evaluations)
Integer iterations () const
Integer max_iterations () const
Scalar alpha () const
Integer relaxations () const
Integer max_relaxations () const
Scalar 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
std::ostream & ostream () const
bool solve (FunctionLambda &&function, const Vector &x_ini, Vector &x_sol)
bool rootfind (const FunctionBase< FunctionInput, FunctionOutput, DerivedFunction > &function, const Vector &x_ini, Vector &x_sol)
bool optimize (const FunctionBase< FunctionInput, FunctionOutput, DerivedFunction > &function, const Vector &x_ini, Vector &x_sol)
constexpr std::string name () const

Static Public Attributes

static constexpr bool RequiresFunction {true}
static constexpr bool RequiresFirstDerivative {true}
static constexpr bool RequiresSecondDerivative {false}
Static Public Attributes inherited from Optimist::RootFinder::RootFinder< Vector, Broyden< Vector > >
static constexpr bool IsRootFinder
static constexpr bool IsOptimizer

Private Attributes

Method m_method {Method::COMBINED}

Additional Inherited Members

Protected Member Functions inherited from Optimist::RootFinder::RootFinder< Vector, Broyden< Vector > >
bool evaluate_jacobian (JacobianLambda &&jacobian, const Input &x, FirstDerivative &out)
bool evaluate_hessian (HessianLambda &&hessian, const Input &x, SecondDerivative &out)
Protected Member Functions inherited from Optimist::SolverBase< Vector, Vector, Broyden< Vector > >
Integer first_derivative_evaluations () const
Integer max_first_derivative_evaluations () const
Integer second_derivative_evaluations () const
Integer max_second_derivative_evaluations () const
void reset_counters ()
bool evaluate_function (FunctionLambda &&function, const Vector &x, Vector &out)
bool evaluate_first_derivative (FirstDerivativeLambda &&function, const Vector &x, FirstDerivative &out)
bool evaluate_second_derivative (SecondDerivativeLambda &&function, const Vector &x, SecondDerivative &out)
bool damp (FunctionLambda &&function, const Vector &x_old, const Vector &function_old, const Vector &step_old, Vector &x_new, Vector &function_new, Vector &step_new)
void header ()
void bottom ()
void info (Scalar residuals, const std::string &notes="-")
Protected Attributes inherited from Optimist::SolverBase< Vector, Vector, Broyden< Vector > >
Vector m_lower_bound
Vector 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
Scalar m_alpha
Integer m_relaxations
Integer m_max_relaxations
Scalar m_tolerance
bool m_verbose
bool m_damped
std::ostream * m_ostream
std::string m_task
bool m_converged

Detailed Description

template<typename Vector>
requires TypeTrait<Vector>::IsEigen && (!TypeTrait<Vector>::IsFixed || TypeTrait<Vector>::Dimension > 0)
class Optimist::RootFinder::Broyden< Vector >

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
VectorEigen vector type.

Member Typedef Documentation

◆ Method

template<typename Vector>
using Optimist::RootFinder::Broyden< Vector >::Method
Initial value:
enum class Method : Integer {
GOOD = 0,
BAD = 1,
COMBINED = 2
}
enum class Method :Integer { STANDARD=0, SPENDLEY=1 } Method
Definition NelderMead.hh:101
OPTIMIST_DEFAULT_INTEGER_TYPE Integer
The Integer type as used for the API.
Definition Optimist.hh:97

Broyden solver type enumeration.

◆ Scalar

template<typename Vector>
using Optimist::RootFinder::Broyden< Vector >::Scalar = typename TypeTrait<Vector>::Scalar

◆ VectorTrait

template<typename Vector>
using Optimist::RootFinder::Broyden< Vector >::VectorTrait = TypeTrait<Vector>

Constructor & Destructor Documentation

◆ Broyden()

template<typename Vector>
Optimist::RootFinder::Broyden< Vector >::Broyden ( )
inline

Class constructor for the Broyden solver.

Member Function Documentation

◆ enable_bad_method()

template<typename Vector>
void Optimist::RootFinder::Broyden< Vector >::enable_bad_method ( )
inline

Enable the bad Broyden solver method.

◆ enable_combined_method()

template<typename Vector>
void Optimist::RootFinder::Broyden< Vector >::enable_combined_method ( )
inline

Enable the combined Broyden solver method.

◆ enable_good_method()

template<typename Vector>
void Optimist::RootFinder::Broyden< Vector >::enable_good_method ( )
inline

Enable the good Broyden solver method.

◆ method() [1/2]

template<typename Vector>
Method Optimist::RootFinder::Broyden< Vector >::method ( ) const
inline

Get the enumeration type of the Broyden solver method.

Returns
The Broyden solver enumeration type.

◆ method() [2/2]

template<typename Vector>
void Optimist::RootFinder::Broyden< Vector >::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 Vector>
std::string Optimist::RootFinder::Broyden< Vector >::name_impl ( ) const
inlineconstexpr

Get the Broyden solver name.

Returns
The Broyden solver name.

◆ set_method()

template<typename Vector>
void Optimist::RootFinder::Broyden< Vector >::set_method ( Method t_method)
inline

Set the Broyden solver type.

Parameters
[in]t_methodThe Broyden solver type enumeration.

◆ update_impl()

template<typename Vector>
void Optimist::RootFinder::Broyden< Vector >::update_impl ( const Vector & delta_x_old,
const Vector & delta_function_old,
const FirstDerivative & jacobian_old,
const Vector & delta_x_new,
const Vector & delta_function_new,
const Vector & ,
FirstDerivative & 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 Vector>
Method Optimist::RootFinder::Broyden< Vector >::m_method {Method::COMBINED}
private

Broyden solver method.

◆ RequiresFirstDerivative

template<typename Vector>
bool Optimist::RootFinder::Broyden< Vector >::RequiresFirstDerivative {true}
staticconstexpr

◆ RequiresFunction

template<typename Vector>
bool Optimist::RootFinder::Broyden< Vector >::RequiresFunction {true}
staticconstexpr

◆ RequiresSecondDerivative

template<typename Vector>
bool Optimist::RootFinder::Broyden< Vector >::RequiresSecondDerivative {false}
staticconstexpr

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