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

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

#include <Greenstadt.hh>

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

Public Types

using Method = enum class Method : Integer {ONE = 1, TWO = 2}
using Vector = typename QuasiNewton<Real, N, Greenstadt<Real, N>>::Vector
using Matrix = typename QuasiNewton<Real, N, Greenstadt<Real, N>>::Matrix
using FunctionWrapper = typename QuasiNewton<Real, N, Greenstadt<Real, N>>::FunctionWrapper
using JacobianWrapper = typename QuasiNewton<Real, N, Greenstadt<Real, N>>::JacobianWrapper
Public Types inherited from Optimist::RootFinder::QuasiNewton< Real, N, Greenstadt< Real, N > >
using Method
using Vector
using Matrix
using FunctionWrapper
using JacobianWrapper
Public Types inherited from Optimist::RootFinder::RootFinder< Real, N, Greenstadt< Real, N >, true >
using Vector
using Matrix
using Tensor
using FunctionWrapper
using JacobianWrapper
using HessianWrapper

Public Member Functions

 Greenstadt ()
std::string name_impl () const
Method method () const
void method (Method t_method)
void enable_one_method ()
void enable_two_method ()
void set_method (Method t_method)
void update_impl (Vector const &, Vector const &, 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::QuasiNewton< Real, N, Greenstadt< 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, Greenstadt< 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, Greenstadt< 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, Greenstadt< 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, Greenstadt< 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::ONE}

Additional Inherited Members

Protected Types inherited from Optimist::SolverBase< Real, N, N, Greenstadt< 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, Greenstadt< 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, Greenstadt< 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, Greenstadt< 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::Greenstadt< Real, N >

Greenstadt's method

Similarly to the Broyden's methods, Greenstadt's methods are quasi-Newton methods that approximates the Jacobian matrix \(\mathbf{Jf}_{\mathbf{x}}\) an update rule. The generic Greenstadt'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.

Based on the update rule, Greenstadt's method is classified into two methods: trivially method 1 and method 2. The update of the Jacobian approximation is performed as

\[ \mathbf{H}_{k+1}^{-1} = \mathbf{H}_k^{-1} - \displaystyle\frac{\mathbf{H}_k^{-1} \Delta\mathbf{f}_k - \Delta\mathbf{x}_k}{\mathbf{g} \Delta\mathbf{f}_k} \mathbf{g} \text{,} \]

where \(\Delta\mathbf{x}_k = \mathbf{x}_k - \mathbf{x}_{k-1}\), and \(\Delta\mathbf{f}_k = \mathbf{f}(\mathbf{x}_k) - \mathbf{f}(\mathbf{x}_{k-1})\). The quantity \(\mathbf{g}\) is \(\mathbf{g} = \mathbf{f}_k^\top\) for method 1 and \(\mathbf{g} = \mathbf{H}_{k\top}^{-1} \mathbf{H}_k^{-1} \Delta\mathbf{f}_k^{\top}\) for method 2. The choice of the method is based on the convergence history, switching between the methods to adapt to the problem's behavior.

For more details on the Greenstadt's methods refer to the reference: Spedicato E., Greenstadt J. On some classes of variationally derived Quasi-Newton methods for systems of nonlinear algebraic equations, Numerische Mathematik, 29 (1978), pp. 363-380.

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

Member Typedef Documentation

◆ FunctionWrapper

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

◆ JacobianWrapper

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

◆ Matrix

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

◆ Method

template<typename Real, Integer N>
using Optimist::RootFinder::Greenstadt< Real, N >::Method = enum class Method : Integer {ONE = 1, TWO = 2}

Greenstadt solver type.

◆ Vector

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

Constructor & Destructor Documentation

◆ Greenstadt()

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

Class constructor for the Greenstadt solver.

Member Function Documentation

◆ enable_one_method()

template<typename Real, Integer N>
void Optimist::RootFinder::Greenstadt< Real, N >::enable_one_method ( )
inline

Enable the Greenstadt1 solver method.

◆ enable_two_method()

template<typename Real, Integer N>
void Optimist::RootFinder::Greenstadt< Real, N >::enable_two_method ( )
inline

Enable the Greenstadt2 solver method.

◆ method() [1/2]

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

Get the enumeration type of the Greenstadt solver method.

Returns
The Greenstadt solver enumeration type.

◆ method() [2/2]

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

Set the enumeration type of the Greenstadt solver method.

Parameters
[in]t_methodThe Greenstadt solver method enumeration type.

◆ name_impl()

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

Get the Greenstadt solver name.

Returns
The Greenstadt solver name.

◆ set_method()

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

Set the Greenstadt solver method.

Parameters
[in]t_methodThe Greenstadt solver method enumeration type.

◆ update_impl()

template<typename Real, Integer N>
void Optimist::RootFinder::Greenstadt< Real, N >::update_impl ( Vector const & ,
Vector const & ,
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 Greenstadt'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::Greenstadt< Real, N >::m_method {Method::ONE}
private

Greenstadt solver type.

◆ requires_first_derivative

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

◆ requires_function

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

◆ requires_second_derivative

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

Basic constants.


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