Optimist  0.0.0
A C++ library for optimization
Loading...
Searching...
No Matches
Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver > Class Template Reference

Class container for the generic root-finding/optimization problem solver. More...

#include <Solver.hh>

Public Member Functions

 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)
 
template<Integer FunInDim, Integer FunOutDim, typename DerivedFunction>
bool rootfind (Function< Real, FunInDim, FunOutDim, DerivedFunction > const &function, const InputType &x_ini, InputType &x_sol)
 
template<Integer FunInDim, Integer FunOutDim, typename DerivedFunction>
bool optimize (Function< Real, FunInDim, FunOutDim, DerivedFunction > const &function, const InputType &x_ini, InputType &x_sol)
 
std::string name () const
 

Protected Types

using InputType = typename std::conditional_t<SolInDim == 1, Real, Eigen::Vector<Real, SolInDim>>
 
using OutputType = typename std::conditional_t<SolOutDim == 1, Real, Eigen::Vector<Real, SolOutDim>>
 
using TraceType = typename std::vector<InputType>
 
using FirstDerivativeType = std::conditional_t<SolInDim == 1 && SolOutDim == 1, Real, Eigen::Matrix<Real, SolOutDim, SolInDim>>
 
using SecondDerivativeType
 
using FunctionWrapper = typename std::function<void(const InputType &, OutputType &)>
 
using FirstDerivativeWrapper = typename std::function<void(const InputType &, FirstDerivativeType &)>
 
using SecondDerivativeWrapper = typename std::function<void(const InputType &, SecondDerivativeType &)>
 

Protected Member Functions

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)
 
template<Integer FunInDim, Integer FunOutDim, typename DerivedFunction>
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

InputType m_lower_bound
 
InputType m_upper_bound
 
Integer m_function_evaluations {0}
 
Integer m_first_derivative_evaluations {0}
 
Integer m_second_derivative_evaluations {0}
 
Integer m_max_function_evaluations {1000}
 
Integer m_max_first_derivative_evaluations {1000}
 
Integer m_max_second_derivative_evaluations {1000}
 
Integer m_iterations {0}
 
Integer m_max_iterations {100}
 
Real m_alpha {0.8}
 
Integer m_relaxations {0}
 
Integer m_max_relaxations {10}
 
Real m_tolerance {EPSILON_LOW}
 
bool m_verbose {false}
 
bool m_damped {true}
 
std::ostream * m_ostream {&std::cout}
 
std::string m_task {"Undefined"}
 
bool m_converged {false}
 
TraceType m_trace
 

Detailed Description

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
class Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >

Solver

The solvers implemented in Optimist are derived from the Solver class, which provides a common interface for solving nonlinear systems of equations. The Solver class is an abstract base class that defines the methods and attributes required by the derived classes. The derived classes implement the specific algorithms for solving the nonlinear systems of equations.

This class provides a base class for developing solvers that handle a generic function \(\mathbf{f}(\mathbf{x}): \mathbb{R}^{n} \rightarrow \mathbb{R}^{m}\). In general, the derived classes are designed to solve root-finding problems of the form

\[ \mathbf{f}(\mathbf{x}) = \mathbf{0} \text{,} \quad \text{with} \quad \mathbf{f}: \mathbb{R}^{n} \rightarrow \mathbb{R}^{n} \text{,} \]

as well as optimization problems of the form

\[ \min_{\mathbf{x}} f(\mathbf{x}) \text{,} \quad f: \mathbb{R}^{n} \rightarrow \mathbb{R} \text{.} \]

If \(n > 1\) the solver is a multi-dimensional solver, otherwise, it is a scalar solver. The linear algebra operations are performed using the Eigen library, which provides a high-level interface for matrix operations. Otherwise, for \(n = 1\), the base type changes to Real, a scalar type defined in the Optimist library.

Template Parameters
RealScalar number type.
SolInDimThe root-finding/optimization problem input dimension.
SolOutDimThe root-finding/optimization problem output dimension.
DerivedSolverDerived solver class.

Member Typedef Documentation

◆ FirstDerivativeType

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
using Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::FirstDerivativeType = std::conditional_t<SolInDim == 1 && SolOutDim == 1, Real, Eigen::Matrix<Real, SolOutDim, SolInDim>>
protected

First derivative type.

◆ FirstDerivativeWrapper

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
using Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::FirstDerivativeWrapper = typename std::function<void(const InputType &, FirstDerivativeType &)>
protected

First derivative type.

◆ FunctionWrapper

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
using Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::FunctionWrapper = typename std::function<void(const InputType &, OutputType &)>
protected

Function type.

◆ InputType

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
using Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::InputType = typename std::conditional_t<SolInDim == 1, Real, Eigen::Vector<Real, SolInDim>>
protected

Basic constants. Input type.

◆ OutputType

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
using Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::OutputType = typename std::conditional_t<SolOutDim == 1, Real, Eigen::Vector<Real, SolOutDim>>
protected

Output type.

◆ SecondDerivativeType

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
using Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::SecondDerivativeType
protected
Initial value:
std::conditional_t<SolInDim == 1 && SolOutDim == 1, Real,
std::conditional_t<SolInDim == 1 || SolOutDim == 1, Eigen::Matrix<Real, SolInDim, SolInDim>,
std::vector<Eigen::Matrix<Real, SolInDim, SolInDim>>>>

Second derivative type.

◆ SecondDerivativeWrapper

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
using Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::SecondDerivativeWrapper = typename std::function<void(const InputType &, SecondDerivativeType &)>
protected

Second derivative type.

◆ TraceType

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
using Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::TraceType = typename std::vector<InputType>
protected

Input trace type.

Constructor & Destructor Documentation

◆ Solver() [1/4]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::Solver ( )
inline

Class constructor for the nonlinear solver.

◆ Solver() [2/4]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::Solver ( FunctionWrapper function,
const InputType & x_ini,
InputType & x_sol )
inline

Class constructor for the nonlinear solver.

Parameters
[in]functionFunction wrapper.
[in]x_iniInitialization point.
[out]x_solSolution point.

◆ Solver() [3/4]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::Solver ( FunctionWrapper function,
FirstDerivativeWrapper first_derivative,
const InputType & x_ini,
InputType & x_sol )
inline

Class constructor for the nonlinear solver.

Parameters
[in]functionFunction wrapper.
[in]first_derivativeFirst derivative wrapper.
[in]x_iniInitialization point.
[out]x_solSolution point.

◆ Solver() [4/4]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::Solver ( FunctionWrapper function,
FirstDerivativeWrapper first_derivative,
SecondDerivativeWrapper second_derivative,
const InputType & x_ini,
InputType & x_sol )
inline

Class constructor for the nonlinear solver.

Parameters
[in]functionFunction wrapper.
[in]first_derivativeFirst derivative wrapper.
[in]second_derivativeThe second derivative wrapper.
[in]x_iniInitialization point.
[out]x_solSolution point.

Member Function Documentation

◆ alpha() [1/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Real Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::alpha ( ) const
inline

Get relaxation factor \( \alpha \).

Returns
The relaxation factor \( \alpha \).

◆ alpha() [2/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::alpha ( Real t_alpha)
inline

Set relaxation factor \( \alpha \).

Parameters
[in]t_alphaThe relaxation factor \( \alpha \).

◆ bottom()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::bottom ( )
inlineprotected

Print the table bottom solver information.

Note
This has to be properly placed in the derived classes.

◆ bounds()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::bounds ( const InputType & t_lower_bound,
const InputType & t_upper_bound )
inline

Set the bounds.

Parameters
[in]t_lower_boundThe lower bound.
[in]t_upper_boundThe upper bound.

◆ converged()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
bool Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::converged ( ) const
inline

Get the convergence boolean flag.

Returns
The convergence boolean flag.

◆ damp()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
bool Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::damp ( FunctionWrapper function,
InputType const & x_old,
InputType const & function_old,
InputType const & step_old,
InputType & x_new,
InputType & function_new,
InputType & step_new )
inlineprotected

Damp the step using the affine invariant criterion.

Parameters
[in]functionFunction wrapper.
[in]x_oldOld point.
[in]function_oldOld function value.
[in]step_oldOld step.
[in]x_newNew point.
[in]function_newNew function value.
[out]step_newNew step.
Returns
The damping boolean flag, true if the damping is successful, false otherwise.

◆ damped_mode() [1/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
bool Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::damped_mode ( ) const
inline

Get the damped mode boolean flag.

Returns
The damped mode boolean flag.

◆ damped_mode() [2/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::damped_mode ( bool t_damped)
inline

Set the damped mode boolean flag.

Parameters
[in]t_dampedThe damped mode boolean flag.

◆ disable_damped_mode()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::disable_damped_mode ( )
inline

Disable solver's damped mode.

◆ disable_verbose_mode()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::disable_verbose_mode ( )
inline

Disable solver's verbose mode.

◆ enable_damped_mode()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::enable_damped_mode ( )
inline

Enable solver's damped mode.

◆ enable_verbose_mode()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::enable_verbose_mode ( )
inline

Enable solver's verbose mode.

◆ evaluate_first_derivative()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::evaluate_first_derivative ( FirstDerivativeWrapper function,
const InputType & x,
FirstDerivativeType & out )
inlineprotected

Evaluate the first derivative.

Parameters
[in]functionFirst derivative wrapper.
[in]xInput point.
[out]outFirst derivative value.

◆ evaluate_function()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::evaluate_function ( FunctionWrapper function,
const InputType & x,
OutputType & out )
inlineprotected

Evaluate the function.

Parameters
[in]functionFunction wrapper.
[in]xInput point.
[out]outFunction value.

◆ evaluate_second_derivative()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::evaluate_second_derivative ( SecondDerivativeWrapper function,
const InputType & x,
SecondDerivativeType & out )
inlineprotected

Evaluate the second derivative.

Parameters
[in]functionSecond derivative wrapper.
[in]xInput point.
[out]outSecond derivative value.

◆ first_derivative_evaluations()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::first_derivative_evaluations ( ) const
inlineprotected

Get the number of function derivative evaluations.

Returns
The number of function derivative evaluations.

◆ function_evaluations()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::function_evaluations ( ) const
inline

Get the number of function evaluations.

Returns
The number of function evaluations.

◆ header()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::header ( )
inlineprotected

Print the table header solver information.

Note
This has to be properly placed in the derived classes.

◆ info()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::info ( Real residuals,
std::string const & notes = "-" )
inlineprotected

Print the solver information during the algorithm iterations.

Note
This has to be properly placed in the derived classes.

◆ input_dimension()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::input_dimension ( ) const
inlineconstexpr

Get the input dimension of the function.

Returns
The input dimension of the function.

◆ iterations()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::iterations ( ) const
inline

Get the number of algorithm iterations.

◆ lower_bound() [1/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
const InputType & Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::lower_bound ( ) const
inline

Get the lower bound.

Returns
The lower bound.

◆ lower_bound() [2/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::lower_bound ( const InputType & t_lower_bound)
inline

Set the lower bound.

Parameters
[in]t_lower_boundThe lower bound.

◆ max_first_derivative_evaluations() [1/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::max_first_derivative_evaluations ( ) const
inlineprotected

Get the number of maximum allowed first derivative evaluations.

Returns
The number of maximum allowed first derivative evaluations.

◆ max_first_derivative_evaluations() [2/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::max_first_derivative_evaluations ( Integer first_derivative_evaluations)
inlineprotected

Set the number of maximum allowed first derivative evaluations.

Parameters
[in]first_derivative_evaluationsThe number of maximum allowed first derivative evaluations.

◆ max_function_evaluations() [1/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::max_function_evaluations ( ) const
inline

Get the number of maximum allowed function evaluations.

Returns
The number of maximum allowed function evaluations.

◆ max_function_evaluations() [2/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::max_function_evaluations ( Integer t_max_function_evaluations)
inline

Set the number of maximum allowed function evaluations.

Parameters
[in]t_max_function_evaluationsThe number of maximum allowed function evaluations.

◆ max_iterations() [1/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::max_iterations ( ) const
inline

Get the number of maximum allowed iterations.

Returns
The number of maximum allowed iterations.

◆ max_iterations() [2/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::max_iterations ( Integer t_max_iterations)
inline

Set the number of maximum allowed iterations.

Parameters
[in]t_max_iterationsThe number of maximum allowed iterations.

◆ max_relaxations() [1/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::max_relaxations ( ) const
inline

Get the number of maximum allowed relaxations.

Returns
The number of maximum allowed relaxations.

◆ max_relaxations() [2/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::max_relaxations ( Integer t_max_relaxations)
inline

Set the number of maximum allowed relaxations.

Parameters
[in]t_max_relaxationsThe number of maximum allowed relaxations.

◆ max_second_derivative_evaluations() [1/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::max_second_derivative_evaluations ( ) const
inlineprotected

Get the number of maximum allowed second derivative evaluations.

Returns
The number of maximum allowed second derivative evaluations.

◆ max_second_derivative_evaluations() [2/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::max_second_derivative_evaluations ( Integer second_derivative_evaluations)
inlineprotected

Set the number of maximum allowed second derivative evaluations.

Parameters
[in]second_derivative_evaluationsThe number of maximum allowed second derivative evaluations.

◆ name()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
std::string Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::name ( ) const
inline

Get the solver name.

Returns
The solver name.

◆ optimize()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
template<Integer FunInDim, Integer FunOutDim, typename DerivedFunction>
bool Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::optimize ( Function< Real, FunInDim, FunOutDim, DerivedFunction > const & function,
const InputType & x_ini,
InputType & x_sol )
inline

Solve the optimization problem given a function class.

Parameters
[in]functionFunction class.
[in]x_iniInitialization point.
[out]x_solSolution point.
Template Parameters
FunInDimThe function output dimension.
FunOutDimThe function output dimension.
DerivedFunctionDerived function class.

◆ ostream() [1/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
std::ostream & Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::ostream ( ) const
inline

Get the output stream for verbose mode.

Returns
The output stream for verbose mode.

◆ ostream() [2/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::ostream ( std::ostream & t_ostream)
inline

Set the output stream for verbose mode.

Parameters
[in]t_ostreamThe output stream for verbose mode.

◆ output_dimension()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::output_dimension ( ) const
inlineconstexpr

Get the output dimension of the function.

Returns
The output dimension of the function.

◆ relaxations()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::relaxations ( ) const
inline

Get the number of algorithm relaxations.

Returns
The number of algorithm relaxations.

◆ reset()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::reset ( )
inlineprotected

Reset solver internal counters and variables.

◆ rootfind()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
template<Integer FunInDim, Integer FunOutDim, typename DerivedFunction>
bool Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::rootfind ( Function< Real, FunInDim, FunOutDim, DerivedFunction > const & function,
const InputType & x_ini,
InputType & x_sol )
inline

Solve the root-finding problem given a function class.

Parameters
[in]functionFunction class.
[in]x_iniInitialization point.
[out]x_solSolution point.
Template Parameters
FunInDimThe function output dimension.
FunOutDimThe function output dimension.
DerivedFunctionDerived function class.

◆ second_derivative_evaluations()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::second_derivative_evaluations ( ) const
inlineprotected

Get the number of second derivative evaluations.

Returns
The number of second derivative evaluations.

◆ solve() [1/4]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
template<Integer FunInDim, Integer FunOutDim, typename DerivedFunction>
bool Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::solve ( Function< Real, FunInDim, FunOutDim, DerivedFunction > const & function,
const InputType & x_ini,
InputType & x_sol,
bool is_optimization )
inlineprotected

Solve the root-finding/optimization problem given a function class.

Parameters
[in]functionFunction class.
[in]x_iniInitialization point.
[out]x_solSolution point.
[in]is_optimizationBoolean flag for optimization.
Template Parameters
FunInDimThe function output dimension.
FunOutDimThe function output dimension.
DerivedFunctionDerived function class.

◆ solve() [2/4]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
bool Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::solve ( FunctionWrapper function,
const InputType & x_ini,
InputType & x_sol )
inline

Solve the root-finding/optimization problem given the function, and without derivatives.

Parameters
[in]functionFunction wrapper.
[in]x_iniInitialization point.
[out]x_solSolution point.

◆ solve() [3/4]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
bool Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::solve ( FunctionWrapper function,
FirstDerivativeWrapper first_derivative,
const InputType & x_ini,
InputType & x_sol )
inline

Solve the root-finding/optimization problem given the function, and its first derivative.

Parameters
[in]functionFunction wrapper.
[in]first_derivativeFirst derivative wrapper
[in]x_iniInitialization point.
[out]x_solSolution point.

◆ solve() [4/4]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
bool Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::solve ( FunctionWrapper function,
FirstDerivativeWrapper first_derivative,
SecondDerivativeWrapper second_derivative,
const InputType & x_ini,
InputType & x_sol )
inline

Solve the root-finding/optimization problem given the function, and its first and second derivatives.

Parameters
[in]functionFunction wrapper.
[in]first_derivativeFirst derivative wrapper.
[in]second_derivativeThe second derivative wrapper.
[in]x_iniInitialization point.
[out]x_solSolution point.

◆ store_trace()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::store_trace ( const InputType & x)
inlineprotected

Update the history of the solver with the current point and function value.

Parameters
[in]xThe point \( \mathbf{x} \).

◆ task() [1/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
std::string Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::task ( ) const
inline

Get the task name.

Returns
The task name.

◆ task() [2/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::task ( std::string t_task)
inline

Set the task name.

Parameters
[in]t_taskThe task name.

◆ tolerance() [1/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Real Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::tolerance ( ) const
inline

Get the tolerance \( \epsilon \).

Returns
The tolerance \( \epsilon \).

◆ tolerance() [2/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::tolerance ( Real t_tolerance)
inline

Set the tolerance \( \epsilon \) for which the nonlinear solver stops, i.e., \( \left\| \mathbf{F}(\mathbf{x}) \right\|_{2} < \epsilon \).

Parameters
[in]t_toleranceThe tolerance \( \epsilon \).

◆ trace()

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
const TraceType & Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::trace ( ) const
inline

Get the trace of input values during the algorithm iterations.

Returns
The trace of input values.

◆ upper_bound() [1/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
const InputType & Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::upper_bound ( ) const
inline

Get the upper bound.

Returns
The upper bound.

◆ upper_bound() [2/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::upper_bound ( const InputType & t_upper_bound)
inline

Set the upper bound.

Parameters
[in]t_upper_boundThe upper bound.

◆ verbose_mode() [1/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
bool Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::verbose_mode ( ) const
inline

Get the verbose mode boolean flag.

Returns
The verbose mode boolean flag.

◆ verbose_mode() [2/2]

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
void Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::verbose_mode ( bool t_verbose)
inline

Set the verbose mode boolean flag.

Parameters
[in]t_verboseThe verbose mode boolean flag.

Member Data Documentation

◆ m_alpha

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Real Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_alpha {0.8}
protected

Relaxation factor \( \alpha \).

◆ m_converged

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
bool Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_converged {false}
protected

Convergence boolean flag.

◆ m_damped

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
bool Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_damped {true}
protected

Damped mode boolean flag.

◆ m_first_derivative_evaluations

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_first_derivative_evaluations {0}
protected

First derivative evaluations.

◆ m_function_evaluations

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_function_evaluations {0}
protected

Function evaluations.

◆ m_iterations

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_iterations {0}
protected

Algorithm iterations.

◆ m_lower_bound

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
InputType Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_lower_bound
protected

Lower bound.

◆ m_max_first_derivative_evaluations

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_max_first_derivative_evaluations {1000}
protected

Maximum allowed first derivative evaluations.

◆ m_max_function_evaluations

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_max_function_evaluations {1000}
protected

Maximum allowed function evaluations.

◆ m_max_iterations

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_max_iterations {100}
protected

Maximum allowed algorithm iterations.

◆ m_max_relaxations

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_max_relaxations {10}
protected

Maximum allowed algorithm relaxations.

◆ m_max_second_derivative_evaluations

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_max_second_derivative_evaluations {1000}
protected

Maximum allowed second derivative evaluations.

◆ m_ostream

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
std::ostream* Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_ostream {&std::cout}
protected

Output stream for verbose mode.

◆ m_relaxations

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_relaxations {0}
protected

Algorithm relaxations.

◆ m_second_derivative_evaluations

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Integer Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_second_derivative_evaluations {0}
protected

Second derivative evaluations.

◆ m_task

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
std::string Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_task {"Undefined"}
protected

Task name.

◆ m_tolerance

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
Real Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_tolerance {EPSILON_LOW}
protected

Solver tolerance \( \epsilon \) for convergence.

◆ m_trace

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
TraceType Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_trace
protected

Trace for points \( \mathbf{x} \) values.

◆ m_upper_bound

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
InputType Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_upper_bound
protected

Upper bound.

◆ m_verbose

template<typename Real, Integer SolInDim, Integer SolOutDim, typename DerivedSolver>
bool Optimist::Solver< Real, SolInDim, SolOutDim, DerivedSolver >::m_verbose {false}
protected

Verbose mode boolean flag.


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