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

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

#include <Greenstadt.hh>

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

Public Types

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

Public Member Functions

 Greenstadt ()
constexpr 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 (const Vector &, const Vector &, 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::QuasiNewton< Vector, Greenstadt< 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, Greenstadt< 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, Greenstadt< 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, Greenstadt< Vector > >
static constexpr bool IsRootFinder
static constexpr bool IsOptimizer

Private Attributes

Method m_method {Method::ONE}

Additional Inherited Members

Protected Member Functions inherited from Optimist::RootFinder::RootFinder< Vector, Greenstadt< 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, Greenstadt< 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, Greenstadt< 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::Greenstadt< Vector >

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

Member Typedef Documentation

◆ Method

template<typename Vector>
using Optimist::RootFinder::Greenstadt< Vector >::Method
Initial value:
enum class Method : Integer {
ONE = 1,
TWO = 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

Greenstadt solver type enumeration.

◆ Scalar

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

◆ VectorTrait

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

Constructor & Destructor Documentation

◆ Greenstadt()

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

Class constructor for the Greenstadt solver.

Member Function Documentation

◆ enable_one_method()

template<typename Vector>
void Optimist::RootFinder::Greenstadt< Vector >::enable_one_method ( )
inline

Enable the Greenstadt1 solver method.

◆ enable_two_method()

template<typename Vector>
void Optimist::RootFinder::Greenstadt< Vector >::enable_two_method ( )
inline

Enable the Greenstadt2 solver method.

◆ method() [1/2]

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

Get the enumeration type of the Greenstadt solver method.

Returns
The Greenstadt solver enumeration type.

◆ method() [2/2]

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

Get the Greenstadt solver name.

Returns
The Greenstadt solver name.

◆ set_method()

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

Set the Greenstadt solver method.

Parameters
[in]t_methodThe Greenstadt solver method enumeration type.

◆ update_impl()

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

Greenstadt solver type.

◆ RequiresFirstDerivative

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

◆ RequiresFunction

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

◆ RequiresSecondDerivative

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

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