![]() |
SG++-Doxygen-Documentation
|
Levenberg-Marquardt algorithm for least squares optimization. More...
#include <LevenbergMarquardt.hpp>
Public Member Functions | |
void | clone (std::unique_ptr< LeastSquaresOptimizer > &clone) const override |
double | getAcceptanceThreshold () const |
double | getEffectivenessThreshold () const |
double | getInitialDamping () const |
base::VectorFunctionGradient & | getPhiGradient () const |
double | getTolerance () const |
LevenbergMarquardt (const base::VectorFunction &phi, const base::VectorFunctionGradient &phiGradient, size_t maxItCount, double tolerance, double initialDamping, double acceptanceThreshold, double effectivenessThreshold, const base::sle_solver::SLESolver &sleSolver) | |
Constructor. | |
LevenbergMarquardt (const base::VectorFunction &phi, const base::VectorFunctionGradient &phiGradient, size_t maxItCount=DEFAULT_N, double tolerance=DEFAULT_TOLERANCE, double initialDamping=DEFAULT_INITIAL_DAMPING, double acceptanceThreshold=DEFAULT_ACCEPTANCE_THRESHOLD, double effectivenessThreshold=DEFAULT_EFFECTIVENESS_THRESHOLD) | |
Constructor. | |
LevenbergMarquardt (const LevenbergMarquardt &other) | |
Copy constructor. | |
void | optimize () override |
Pure virtual method for optimization of the objective function. | |
void | setAcceptanceThreshold (double acceptanceThreshold) |
void | setEffectivenessThreshold (double effectivenessThreshold) |
void | setInitialDamping (double initialDamping) |
void | setTolerance (double tolerance) |
~LevenbergMarquardt () override | |
Destructor. | |
![]() | |
const base::DataMatrix & | getHistoryOfOptimalPoints () const |
const base::DataVector & | getHistoryOfOptimalValues () const |
size_t | getN () const |
const base::DataVector & | getOptimalPoint () const |
double | getOptimalValue () const |
base::VectorFunction & | getPhiFunction () const |
const base::DataVector & | getStartingPoint () const |
LeastSquaresOptimizer (const base::VectorFunction &phi, size_t N=DEFAULT_N) | |
Constructor. | |
LeastSquaresOptimizer (const LeastSquaresOptimizer &other) | |
Copy constructor. | |
void | setN (size_t N) |
void | setPhiFunction (const base::VectorFunction &phi) |
void | setStartingPoint (const base::DataVector &startingPoint) |
virtual | ~LeastSquaresOptimizer () |
Destructor. | |
Static Public Attributes | |
static constexpr double | DEFAULT_ACCEPTANCE_THRESHOLD = 0.3 |
default acceptance threshold | |
static constexpr double | DEFAULT_EFFECTIVENESS_THRESHOLD = 0.9 |
default effectiveness threshold | |
static constexpr double | DEFAULT_INITIAL_DAMPING = 1.0 |
default initial damping | |
static constexpr double | DEFAULT_TOLERANCE = 1e-6 |
default tolerance | |
![]() | |
static const size_t | DEFAULT_N = 1000 |
default maximal number of iterations or function evaluations | |
Protected Attributes | |
double | beta0 |
acceptance threshold | |
double | beta1 |
effectiveness threshold | |
const base::sle_solver::GaussianElimination | defaultSleSolver |
default linear solver | |
double | mu0 |
initial damping | |
std::unique_ptr< base::VectorFunctionGradient > | phiGradient |
phi gradient | |
const base::sle_solver::SLESolver & | sleSolver |
linear solver | |
double | tol |
tolerance | |
![]() | |
base::DataVector | fHist |
search history vector (optimal values) | |
double | fOpt |
result of optimization (optimal function value) | |
size_t | N |
maximal number of iterations or function evaluations | |
std::unique_ptr< base::VectorFunction > | phi |
phi function | |
base::DataVector | x0 |
starting point | |
base::DataMatrix | xHist |
search history matrix (optimal points) | |
base::DataVector | xOpt |
result of optimization (location of optimum) | |
Levenberg-Marquardt algorithm for least squares optimization.
sgpp::optimization::optimizer::LevenbergMarquardt::LevenbergMarquardt | ( | const base::VectorFunction & | phi, |
const base::VectorFunctionGradient & | phiGradient, | ||
size_t | maxItCount = DEFAULT_N , |
||
double | tolerance = DEFAULT_TOLERANCE , |
||
double | initialDamping = DEFAULT_INITIAL_DAMPING , |
||
double | acceptanceThreshold = DEFAULT_ACCEPTANCE_THRESHOLD , |
||
double | effectivenessThreshold = DEFAULT_EFFECTIVENESS_THRESHOLD |
||
) |
Constructor.
By default, GaussianElimination is used to solve the linear systems.
phi | base function |
phiGradient | Jacobian of phi |
maxItCount | maximal number of function evaluations |
tolerance | tolerance |
initialDamping | initial damping |
acceptanceThreshold | acceptance threshold |
effectivenessThreshold | effectiveness threshold |
References sgpp::optimization::optimizer::LeastSquaresOptimizer::clone(), and phiGradient.
sgpp::optimization::optimizer::LevenbergMarquardt::LevenbergMarquardt | ( | const base::VectorFunction & | phi, |
const base::VectorFunctionGradient & | phiGradient, | ||
size_t | maxItCount, | ||
double | tolerance, | ||
double | initialDamping, | ||
double | acceptanceThreshold, | ||
double | effectivenessThreshold, | ||
const base::sle_solver::SLESolver & | sleSolver | ||
) |
Constructor.
Do not destruct the solver before this object!
phi | phi function |
phiGradient | Jacobian of phi |
maxItCount | maximal number of function evaluations |
tolerance | tolerance |
initialDamping | initial damping |
acceptanceThreshold | acceptance threshold |
effectivenessThreshold | effectiveness threshold |
sleSolver | reference to linear solver for solving the linear systems |
References phiGradient.
sgpp::optimization::optimizer::LevenbergMarquardt::LevenbergMarquardt | ( | const LevenbergMarquardt & | other | ) |
|
override |
Destructor.
|
overridevirtual |
[out] | clone | pointer to cloned object |
Implements sgpp::optimization::optimizer::LeastSquaresOptimizer.
References clone().
Referenced by clone().
double sgpp::optimization::optimizer::LevenbergMarquardt::getAcceptanceThreshold | ( | ) | const |
References beta0.
double sgpp::optimization::optimizer::LevenbergMarquardt::getEffectivenessThreshold | ( | ) | const |
References beta1.
double sgpp::optimization::optimizer::LevenbergMarquardt::getInitialDamping | ( | ) | const |
References mu0.
base::VectorFunctionGradient & sgpp::optimization::optimizer::LevenbergMarquardt::getPhiGradient | ( | ) | const |
References phiGradient.
double sgpp::optimization::optimizer::LevenbergMarquardt::getTolerance | ( | ) | const |
References tol.
|
overridevirtual |
Pure virtual method for optimization of the objective function.
The result of the optimization process can be obtained by member functions, e.g., getOptimalPoint() and getOptimalValue().
Implements sgpp::optimization::optimizer::LeastSquaresOptimizer.
References sgpp::base::DataVector::append(), sgpp::base::DataMatrix::appendRow(), beta0, beta1, sgpp::base::Printer::disableStatusPrinting(), sgpp::base::Printer::enableStatusPrinting(), sgpp::optimization::optimizer::LeastSquaresOptimizer::fHist, sgpp::optimization::optimizer::LeastSquaresOptimizer::fOpt, sgpp::base::Printer::getInstance(), sgpp::base::Printer::isStatusPrintingEnabled(), m, mu, mu0, sgpp::base::DataMatrix::mult(), sgpp::optimization::optimizer::LeastSquaresOptimizer::N, sgpp::optimization::optimizer::LeastSquaresOptimizer::phi, phiGradient, sgpp::base::Printer::printStatusBegin(), sgpp::base::Printer::printStatusEnd(), sgpp::base::Printer::printStatusUpdate(), sgpp::base::DataMatrix::resize(), sleSolver, sgpp::base::sle_solver::SLESolver::solve(), tol, sgpp::base::DataVector::toString(), sgpp::optimization::optimizer::LeastSquaresOptimizer::x0, sgpp::optimization::optimizer::LeastSquaresOptimizer::xHist, and sgpp::optimization::optimizer::LeastSquaresOptimizer::xOpt.
void sgpp::optimization::optimizer::LevenbergMarquardt::setAcceptanceThreshold | ( | double | acceptanceThreshold | ) |
acceptanceThreshold | acceptance threshold |
References beta0.
void sgpp::optimization::optimizer::LevenbergMarquardt::setEffectivenessThreshold | ( | double | effectivenessThreshold | ) |
effectivenessThreshold | effectiveness threshold |
References beta1.
void sgpp::optimization::optimizer::LevenbergMarquardt::setInitialDamping | ( | double | initialDamping | ) |
initialDamping | initial damping |
References mu0.
void sgpp::optimization::optimizer::LevenbergMarquardt::setTolerance | ( | double | tolerance | ) |
tolerance | tolerance |
References tol.
|
protected |
acceptance threshold
Referenced by getAcceptanceThreshold(), optimize(), and setAcceptanceThreshold().
|
protected |
effectiveness threshold
Referenced by getEffectivenessThreshold(), optimize(), and setEffectivenessThreshold().
|
staticconstexpr |
default acceptance threshold
|
staticconstexpr |
default effectiveness threshold
|
staticconstexpr |
default initial damping
|
staticconstexpr |
default tolerance
|
protected |
default linear solver
|
protected |
initial damping
Referenced by getInitialDamping(), optimize(), and setInitialDamping().
|
protected |
phi gradient
Referenced by getPhiGradient(), LevenbergMarquardt(), LevenbergMarquardt(), LevenbergMarquardt(), and optimize().
|
protected |
linear solver
Referenced by optimize().
|
protected |
tolerance
Referenced by getTolerance(), optimize(), and setTolerance().