Select Git revision
optimization.h
optimization.h 135.48 KiB
/*************************************************************************
Copyright (c) Sergey Bochkanov (ALGLIB project).
>>> SOURCE LICENSE >>>
This program is free software; you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation (www.fsf.org); either version 2 of the
License, or (at your option) any later version.
This program is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
A copy of the GNU General Public License is available at
http://www.fsf.org/licensing/licenses
>>> END OF LICENSE >>>
*************************************************************************/
#ifndef _optimization_pkg_h
#define _optimization_pkg_h
#include "ap.h"
#include "alglibinternal.h"
#include "alglibmisc.h"
#include "linalg.h"
#include "solvers.h"
/////////////////////////////////////////////////////////////////////////
//
// THIS SECTION CONTAINS COMPUTATIONAL CORE DECLARATIONS (DATATYPES)
//
/////////////////////////////////////////////////////////////////////////
namespace alglib_impl
{
typedef struct
{
ae_int_t n;
double epsg;
double epsf;
double epsx;
ae_int_t maxits;
double stpmax;
double suggestedstep;
ae_bool xrep;
ae_bool drep;
ae_int_t cgtype;
ae_int_t prectype;
ae_vector diagh;
ae_vector diaghl2;
ae_matrix vcorr;
ae_int_t vcnt;
ae_vector s;
double diffstep;
ae_int_t nfev;
ae_int_t mcstage;
ae_int_t k;
ae_vector xk;
ae_vector dk;
ae_vector xn;
ae_vector dn;
ae_vector d;
double fold;
double stp;
double curstpmax;
ae_vector yk;
double laststep;
double lastscaledstep;
ae_int_t mcinfo;
ae_bool innerresetneeded;
ae_bool terminationneeded;
double trimthreshold;
ae_int_t rstimer;
ae_vector x;
double f;
ae_vector g;
ae_bool needf;
ae_bool needfg;
ae_bool xupdated;
ae_bool algpowerup;
ae_bool lsstart;
ae_bool lsend;
rcommstate rstate;
ae_int_t repiterationscount;
ae_int_t repnfev;
ae_int_t repterminationtype;
ae_int_t debugrestartscount;
linminstate lstate;
double fbase;
double fm2;
double fm1;
double fp1;
double fp2;
double betahs;
double betady;
ae_vector work0;
ae_vector work1;
} mincgstate;
typedef struct
{
ae_int_t iterationscount;
ae_int_t nfev;
ae_int_t terminationtype;
} mincgreport;
typedef struct
{
ae_int_t nmain;
ae_int_t nslack;
double innerepsg;
double innerepsf;
double innerepsx;
double outerepsx;
double outerepsi;
ae_int_t maxits;
ae_bool xrep;
double stpmax;
double diffstep;
ae_int_t prectype;
ae_vector diaghoriginal;
ae_vector diagh;
ae_vector x;
double f;
ae_vector g;
ae_bool needf;
ae_bool needfg;
ae_bool xupdated;
rcommstate rstate;
ae_int_t repinneriterationscount;
ae_int_t repouteriterationscount;
ae_int_t repnfev;
ae_int_t repterminationtype;
double repdebugeqerr;
double repdebugfs;
double repdebugff;
double repdebugdx;
ae_vector xcur;
ae_vector xprev;
ae_vector xstart;
ae_int_t itsleft;
ae_vector xend;
ae_vector lastg;
double trimthreshold;
ae_matrix ceoriginal;
ae_matrix ceeffective;
ae_matrix cecurrent;
ae_vector ct;
ae_int_t cecnt;
ae_int_t cedim;
ae_vector xe;
ae_vector hasbndl;
ae_vector hasbndu;
ae_vector bndloriginal;
ae_vector bnduoriginal;
ae_vector bndleffective;
ae_vector bndueffective;
ae_vector activeconstraints;
ae_vector constrainedvalues;
ae_vector transforms;
ae_vector seffective;
ae_vector soriginal;
ae_vector w;
ae_vector tmp0;
ae_vector tmp1;
ae_vector tmp2;
ae_vector r;
ae_matrix lmmatrix;
double v0;
double v1;
double v2;
double t;
double errfeas;
double gnorm;
double mpgnorm;
double mba;
ae_int_t variabletofreeze;
double valuetofreeze;
double fbase;
double fm2;
double fm1;
double fp1;
double fp2;
double xm1;
double xp1;
mincgstate cgstate;
mincgreport cgrep;
ae_int_t optdim;
} minbleicstate;
typedef struct
{
ae_int_t inneriterationscount;
ae_int_t outeriterationscount;
ae_int_t nfev;
ae_int_t terminationtype;
double debugeqerr;
double debugfs;
double debugff;
double debugdx;
} minbleicreport;
typedef struct
{
ae_int_t n;
ae_int_t m;
double epsg;
double epsf;
double epsx;
ae_int_t maxits;
ae_bool xrep;
double stpmax;
ae_vector s;
double diffstep;
ae_int_t nfev;
ae_int_t mcstage;
ae_int_t k;
ae_int_t q;
ae_int_t p;
ae_vector rho;
ae_matrix yk;
ae_matrix sk;
ae_vector theta;
ae_vector d;
double stp;
ae_vector work;
double fold;
double trimthreshold;
ae_int_t prectype;
double gammak;
ae_matrix denseh;
ae_vector diagh;
double fbase;
double fm2;
double fm1;
double fp1;
double fp2;
ae_vector autobuf;
ae_vector x;
double f;
ae_vector g;
ae_bool needf;
ae_bool needfg;
ae_bool xupdated;
rcommstate rstate;
ae_int_t repiterationscount;
ae_int_t repnfev;
ae_int_t repterminationtype;
linminstate lstate;
} minlbfgsstate;
typedef struct
{
ae_int_t iterationscount;
ae_int_t nfev;
ae_int_t terminationtype;
} minlbfgsreport;
typedef struct
{
ae_int_t n;
ae_int_t algokind;
ae_int_t akind;
ae_matrix densea;
ae_vector diaga;
ae_vector b;
ae_vector bndl;
ae_vector bndu;
ae_vector havebndl;
ae_vector havebndu;
ae_vector xorigin;
ae_vector startx;
ae_bool havex;
ae_vector xc;
ae_vector gc;
ae_vector activeconstraints;
ae_vector prevactiveconstraints;
double constterm;
ae_vector workbndl;
ae_vector workbndu;
ae_int_t repinneriterationscount;
ae_int_t repouteriterationscount;
ae_int_t repncholesky;
ae_int_t repnmv;
ae_int_t repterminationtype;
ae_vector tmp0;
ae_vector tmp1;
ae_vector itmp0;
ae_vector p2;
ae_matrix bufa;
ae_vector bufb;
ae_vector bufx;
apbuffers buf;
} minqpstate;
typedef struct
{
ae_int_t inneriterationscount;
ae_int_t outeriterationscount;
ae_int_t nmv;
ae_int_t ncholesky;
ae_int_t terminationtype;
} minqpreport;
typedef struct
{
ae_int_t n;
ae_int_t m;
double diffstep;
double epsg;
double epsf;
double epsx;
ae_int_t maxits;
ae_bool xrep;
double stpmax;
ae_int_t maxmodelage;
ae_bool makeadditers;
ae_vector x;
double f;
ae_vector fi;
ae_matrix j;
ae_matrix h;
ae_vector g;
ae_bool needf;
ae_bool needfg;
ae_bool needfgh;
ae_bool needfij;
ae_bool needfi;
ae_bool xupdated;
ae_int_t algomode;
ae_bool hasf;
ae_bool hasfi;
ae_bool hasg;
ae_vector xbase;
double fbase;
ae_vector fibase;
ae_vector gbase;
ae_matrix quadraticmodel;
ae_vector bndl;
ae_vector bndu;
ae_vector havebndl;
ae_vector havebndu;
ae_vector s;
double lambdav;
double nu;
ae_int_t modelage;
ae_vector xdir;
ae_vector deltax;
ae_vector deltaf;
ae_bool deltaxready;
ae_bool deltafready;
ae_int_t repiterationscount;
ae_int_t repterminationtype;
ae_int_t repnfunc;
ae_int_t repnjac;
ae_int_t repngrad;
ae_int_t repnhess;
ae_int_t repncholesky;
rcommstate rstate;
ae_vector choleskybuf;
ae_vector tmp0;
double actualdecrease;
double predicteddecrease;
double xm1;
double xp1;
ae_vector fm1;
ae_vector fp1;
minlbfgsstate internalstate;
minlbfgsreport internalrep;
minqpstate qpstate;
minqpreport qprep;
} minlmstate;
typedef struct
{
ae_int_t iterationscount;
ae_int_t terminationtype;
ae_int_t nfunc;
ae_int_t njac;
ae_int_t ngrad;
ae_int_t nhess;
ae_int_t ncholesky;
} minlmreport;
typedef struct
{
ae_int_t n;
double epsg;
double epsf;
double epsx;
ae_int_t maxits;
ae_bool xrep;
double stpmax;
ae_int_t cgtype;
ae_int_t k;
ae_int_t nfev;
ae_int_t mcstage;
ae_vector bndl;
ae_vector bndu;
ae_int_t curalgo;
ae_int_t acount;
double mu;
double finit;
double dginit;
ae_vector ak;
ae_vector xk;
ae_vector dk;
ae_vector an;
ae_vector xn;
ae_vector dn;
ae_vector d;
double fold;
double stp;
ae_vector work;
ae_vector yk;
ae_vector gc;
double laststep;
ae_vector x;
double f;
ae_vector g;
ae_bool needfg;
ae_bool xupdated;
rcommstate rstate;
ae_int_t repiterationscount;
ae_int_t repnfev;
ae_int_t repterminationtype;
ae_int_t debugrestartscount;
linminstate lstate;
double betahs;
double betady;
} minasastate;
typedef struct
{
ae_int_t iterationscount;
ae_int_t nfev;
ae_int_t terminationtype;
ae_int_t activeconstraints;
} minasareport;
}
/////////////////////////////////////////////////////////////////////////
//
// THIS SECTION CONTAINS C++ INTERFACE
//
/////////////////////////////////////////////////////////////////////////
namespace alglib
{
/*************************************************************************
This object stores state of the nonlinear CG optimizer.
You should use ALGLIB functions to work with this object.
*************************************************************************/
class _mincgstate_owner
{
public:
_mincgstate_owner();
_mincgstate_owner(const _mincgstate_owner &rhs);
_mincgstate_owner& operator=(const _mincgstate_owner &rhs);
virtual ~_mincgstate_owner();
alglib_impl::mincgstate* c_ptr();
alglib_impl::mincgstate* c_ptr() const;
protected:
alglib_impl::mincgstate *p_struct;
};
class mincgstate : public _mincgstate_owner
{
public:
mincgstate();
mincgstate(const mincgstate &rhs);
mincgstate& operator=(const mincgstate &rhs);
virtual ~mincgstate();
ae_bool &needf;
ae_bool &needfg;
ae_bool &xupdated;
double &f;
real_1d_array g;
real_1d_array x;
};
/*************************************************************************
*************************************************************************/
class _mincgreport_owner
{
public:
_mincgreport_owner();
_mincgreport_owner(const _mincgreport_owner &rhs);
_mincgreport_owner& operator=(const _mincgreport_owner &rhs);
virtual ~_mincgreport_owner();
alglib_impl::mincgreport* c_ptr();
alglib_impl::mincgreport* c_ptr() const;
protected:
alglib_impl::mincgreport *p_struct;
};
class mincgreport : public _mincgreport_owner
{
public:
mincgreport();
mincgreport(const mincgreport &rhs);
mincgreport& operator=(const mincgreport &rhs);
virtual ~mincgreport();
ae_int_t &iterationscount;
ae_int_t &nfev;
ae_int_t &terminationtype;
};
/*************************************************************************
This object stores nonlinear optimizer state.
You should use functions provided by MinBLEIC subpackage to work with this
object
*************************************************************************/
class _minbleicstate_owner
{
public:
_minbleicstate_owner();
_minbleicstate_owner(const _minbleicstate_owner &rhs);
_minbleicstate_owner& operator=(const _minbleicstate_owner &rhs);
virtual ~_minbleicstate_owner();
alglib_impl::minbleicstate* c_ptr();
alglib_impl::minbleicstate* c_ptr() const;
protected:
alglib_impl::minbleicstate *p_struct;
};
class minbleicstate : public _minbleicstate_owner
{
public:
minbleicstate();
minbleicstate(const minbleicstate &rhs);
minbleicstate& operator=(const minbleicstate &rhs);
virtual ~minbleicstate();
ae_bool &needf;
ae_bool &needfg;
ae_bool &xupdated;
double &f;
real_1d_array g;
real_1d_array x;
};
/*************************************************************************
This structure stores optimization report:
* InnerIterationsCount number of inner iterations
* OuterIterationsCount number of outer iterations
* NFEV number of gradient evaluations
* TerminationType termination type (see below)
TERMINATION CODES
TerminationType field contains completion code, which can be:
-10 unsupported combination of algorithm settings:
1) StpMax is set to non-zero value,
AND 2) non-default preconditioner is used.
You can't use both features at the same moment,
so you have to choose one of them (and to turn
off another one).
-3 inconsistent constraints. Feasible point is
either nonexistent or too hard to find. Try to
restart optimizer with better initial
approximation
4 conditions on constraints are fulfilled
with error less than or equal to EpsC
5 MaxIts steps was taken
7 stopping conditions are too stringent,
further improvement is impossible,
X contains best point found so far.
ADDITIONAL FIELDS
There are additional fields which can be used for debugging:
* DebugEqErr error in the equality constraints (2-norm)
* DebugFS f, calculated at projection of initial point
to the feasible set
* DebugFF f, calculated at the final point
* DebugDX |X_start-X_final|
*************************************************************************/
class _minbleicreport_owner
{
public:
_minbleicreport_owner();
_minbleicreport_owner(const _minbleicreport_owner &rhs);
_minbleicreport_owner& operator=(const _minbleicreport_owner &rhs);
virtual ~_minbleicreport_owner();
alglib_impl::minbleicreport* c_ptr();
alglib_impl::minbleicreport* c_ptr() const;
protected:
alglib_impl::minbleicreport *p_struct;
};
class minbleicreport : public _minbleicreport_owner
{
public:
minbleicreport();
minbleicreport(const minbleicreport &rhs);
minbleicreport& operator=(const minbleicreport &rhs);
virtual ~minbleicreport();
ae_int_t &inneriterationscount;
ae_int_t &outeriterationscount;
ae_int_t &nfev;
ae_int_t &terminationtype;
double &debugeqerr;
double &debugfs;
double &debugff;
double &debugdx;
};
/*************************************************************************
*************************************************************************/
class _minlbfgsstate_owner
{
public:
_minlbfgsstate_owner();
_minlbfgsstate_owner(const _minlbfgsstate_owner &rhs);
_minlbfgsstate_owner& operator=(const _minlbfgsstate_owner &rhs);
virtual ~_minlbfgsstate_owner();
alglib_impl::minlbfgsstate* c_ptr();
alglib_impl::minlbfgsstate* c_ptr() const;
protected:
alglib_impl::minlbfgsstate *p_struct;
};
class minlbfgsstate : public _minlbfgsstate_owner
{
public:
minlbfgsstate();
minlbfgsstate(const minlbfgsstate &rhs);
minlbfgsstate& operator=(const minlbfgsstate &rhs);
virtual ~minlbfgsstate();
ae_bool &needf;
ae_bool &needfg;
ae_bool &xupdated;
double &f;
real_1d_array g;
real_1d_array x;
};
/*************************************************************************
*************************************************************************/
class _minlbfgsreport_owner
{
public:
_minlbfgsreport_owner();
_minlbfgsreport_owner(const _minlbfgsreport_owner &rhs);
_minlbfgsreport_owner& operator=(const _minlbfgsreport_owner &rhs);
virtual ~_minlbfgsreport_owner();
alglib_impl::minlbfgsreport* c_ptr();
alglib_impl::minlbfgsreport* c_ptr() const;
protected:
alglib_impl::minlbfgsreport *p_struct;
};
class minlbfgsreport : public _minlbfgsreport_owner
{
public:
minlbfgsreport();
minlbfgsreport(const minlbfgsreport &rhs);
minlbfgsreport& operator=(const minlbfgsreport &rhs);
virtual ~minlbfgsreport();
ae_int_t &iterationscount;
ae_int_t &nfev;
ae_int_t &terminationtype;
};
/*************************************************************************
This object stores nonlinear optimizer state.
You should use functions provided by MinQP subpackage to work with this
object
*************************************************************************/
class _minqpstate_owner
{
public:
_minqpstate_owner();
_minqpstate_owner(const _minqpstate_owner &rhs);
_minqpstate_owner& operator=(const _minqpstate_owner &rhs);
virtual ~_minqpstate_owner();
alglib_impl::minqpstate* c_ptr();
alglib_impl::minqpstate* c_ptr() const;
protected:
alglib_impl::minqpstate *p_struct;
};
class minqpstate : public _minqpstate_owner
{
public:
minqpstate();
minqpstate(const minqpstate &rhs);
minqpstate& operator=(const minqpstate &rhs);
virtual ~minqpstate();
};
/*************************************************************************
This structure stores optimization report:
* InnerIterationsCount number of inner iterations
* OuterIterationsCount number of outer iterations
* NCholesky number of Cholesky decomposition
* NMV number of matrix-vector products
(only products calculated as part of iterative
process are counted)
* TerminationType completion code (see below)
Completion codes:
* -5 inappropriate solver was used:
* Cholesky solver for semidefinite or indefinite problems
* Cholesky solver for problems with non-boundary constraints
* -3 inconsistent constraints (or, maybe, feasible point is
too hard to find). If you are sure that constraints are feasible,
try to restart optimizer with better initial approximation.
* 4 successful completion
* 5 MaxIts steps was taken
* 7 stopping conditions are too stringent,
further improvement is impossible,
X contains best point found so far.
*************************************************************************/
class _minqpreport_owner
{
public:
_minqpreport_owner();
_minqpreport_owner(const _minqpreport_owner &rhs);
_minqpreport_owner& operator=(const _minqpreport_owner &rhs);
virtual ~_minqpreport_owner();
alglib_impl::minqpreport* c_ptr();
alglib_impl::minqpreport* c_ptr() const;
protected:
alglib_impl::minqpreport *p_struct;
};
class minqpreport : public _minqpreport_owner
{
public:
minqpreport();
minqpreport(const minqpreport &rhs);
minqpreport& operator=(const minqpreport &rhs);
virtual ~minqpreport();
ae_int_t &inneriterationscount;
ae_int_t &outeriterationscount;
ae_int_t &nmv;
ae_int_t &ncholesky;
ae_int_t &terminationtype;
};
/*************************************************************************
Levenberg-Marquardt optimizer.
This structure should be created using one of the MinLMCreate???()
functions. You should not access its fields directly; use ALGLIB functions
to work with it.
*************************************************************************/
class _minlmstate_owner
{
public:
_minlmstate_owner();
_minlmstate_owner(const _minlmstate_owner &rhs);
_minlmstate_owner& operator=(const _minlmstate_owner &rhs);
virtual ~_minlmstate_owner();
alglib_impl::minlmstate* c_ptr();
alglib_impl::minlmstate* c_ptr() const;
protected:
alglib_impl::minlmstate *p_struct;
};
class minlmstate : public _minlmstate_owner
{
public:
minlmstate();
minlmstate(const minlmstate &rhs);
minlmstate& operator=(const minlmstate &rhs);
virtual ~minlmstate();
ae_bool &needf;
ae_bool &needfg;
ae_bool &needfgh;
ae_bool &needfi;
ae_bool &needfij;
ae_bool &xupdated;
double &f;
real_1d_array fi;
real_1d_array g;
real_2d_array h;
real_2d_array j;
real_1d_array x;
};
/*************************************************************************
Optimization report, filled by MinLMResults() function
FIELDS:
* TerminationType, completetion code:
* -9 derivative correctness check failed;
see Rep.WrongNum, Rep.WrongI, Rep.WrongJ for
more information.
* 1 relative function improvement is no more than
EpsF.
* 2 relative step is no more than EpsX.
* 4 gradient is no more than EpsG.
* 5 MaxIts steps was taken
* 7 stopping conditions are too stringent,
further improvement is impossible
* IterationsCount, contains iterations count
* NFunc, number of function calculations
* NJac, number of Jacobi matrix calculations
* NGrad, number of gradient calculations
* NHess, number of Hessian calculations
* NCholesky, number of Cholesky decomposition calculations
*************************************************************************/
class _minlmreport_owner
{
public:
_minlmreport_owner();
_minlmreport_owner(const _minlmreport_owner &rhs);
_minlmreport_owner& operator=(const _minlmreport_owner &rhs);
virtual ~_minlmreport_owner();
alglib_impl::minlmreport* c_ptr();
alglib_impl::minlmreport* c_ptr() const;
protected:
alglib_impl::minlmreport *p_struct;
};
class minlmreport : public _minlmreport_owner
{
public:
minlmreport();
minlmreport(const minlmreport &rhs);
minlmreport& operator=(const minlmreport &rhs);
virtual ~minlmreport();
ae_int_t &iterationscount;
ae_int_t &terminationtype;
ae_int_t &nfunc;
ae_int_t &njac;
ae_int_t &ngrad;
ae_int_t &nhess;
ae_int_t &ncholesky;
};
/*************************************************************************
*************************************************************************/
class _minasastate_owner
{
public:
_minasastate_owner();
_minasastate_owner(const _minasastate_owner &rhs);
_minasastate_owner& operator=(const _minasastate_owner &rhs);
virtual ~_minasastate_owner();
alglib_impl::minasastate* c_ptr();
alglib_impl::minasastate* c_ptr() const;
protected:
alglib_impl::minasastate *p_struct;
};
class minasastate : public _minasastate_owner
{
public:
minasastate();
minasastate(const minasastate &rhs);
minasastate& operator=(const minasastate &rhs);
virtual ~minasastate();
ae_bool &needfg;
ae_bool &xupdated;
double &f;
real_1d_array g;
real_1d_array x;
};
/*************************************************************************
*************************************************************************/
class _minasareport_owner
{
public:
_minasareport_owner();
_minasareport_owner(const _minasareport_owner &rhs);
_minasareport_owner& operator=(const _minasareport_owner &rhs);
virtual ~_minasareport_owner();
alglib_impl::minasareport* c_ptr();
alglib_impl::minasareport* c_ptr() const;
protected:
alglib_impl::minasareport *p_struct;
};
class minasareport : public _minasareport_owner
{
public:
minasareport();
minasareport(const minasareport &rhs);
minasareport& operator=(const minasareport &rhs);
virtual ~minasareport();
ae_int_t &iterationscount;
ae_int_t &nfev;
ae_int_t &terminationtype;
ae_int_t &activeconstraints;
};
/*************************************************************************
NONLINEAR CONJUGATE GRADIENT METHOD
DESCRIPTION:
The subroutine minimizes function F(x) of N arguments by using one of the
nonlinear conjugate gradient methods.
These CG methods are globally convergent (even on non-convex functions) as
long as grad(f) is Lipschitz continuous in a some neighborhood of the
L = { x : f(x)<=f(x0) }.
REQUIREMENTS:
Algorithm will request following information during its operation:
* function value F and its gradient G (simultaneously) at given point X
USAGE:
1. User initializes algorithm state with MinCGCreate() call
2. User tunes solver parameters with MinCGSetCond(), MinCGSetStpMax() and
other functions
3. User calls MinCGOptimize() function which takes algorithm state and
pointer (delegate, etc.) to callback function which calculates F/G.
4. User calls MinCGResults() to get solution
5. Optionally, user may call MinCGRestartFrom() to solve another problem
with same N but another starting point and/or another function.
MinCGRestartFrom() allows to reuse already initialized structure.
INPUT PARAMETERS:
N - problem dimension, N>0:
* if given, only leading N elements of X are used
* if not given, automatically determined from size of X
X - starting point, array[0..N-1].
OUTPUT PARAMETERS:
State - structure which stores algorithm state
-- ALGLIB --
Copyright 25.03.2010 by Bochkanov Sergey
*************************************************************************/
void mincgcreate(const ae_int_t n, const real_1d_array &x, mincgstate &state);
void mincgcreate(const real_1d_array &x, mincgstate &state);
/*************************************************************************
The subroutine is finite difference variant of MinCGCreate(). It uses
finite differences in order to differentiate target function.
Description below contains information which is specific to this function
only. We recommend to read comments on MinCGCreate() in order to get more
information about creation of CG optimizer.
INPUT PARAMETERS:
N - problem dimension, N>0:
* if given, only leading N elements of X are used
* if not given, automatically determined from size of X
X - starting point, array[0..N-1].
DiffStep- differentiation step, >0
OUTPUT PARAMETERS:
State - structure which stores algorithm state
NOTES:
1. algorithm uses 4-point central formula for differentiation.
2. differentiation step along I-th axis is equal to DiffStep*S[I] where
S[] is scaling vector which can be set by MinCGSetScale() call.
3. we recommend you to use moderate values of differentiation step. Too
large step will result in too large truncation errors, while too small
step will result in too large numerical errors. 1.0E-6 can be good
value to start with.
4. Numerical differentiation is very inefficient - one gradient
calculation needs 4*N function evaluations. This function will work for
any N - either small (1...10), moderate (10...100) or large (100...).
However, performance penalty will be too severe for any N's except for
small ones.
We should also say that code which relies on numerical differentiation
is less robust and precise. L-BFGS needs exact gradient values.
Imprecise gradient may slow down convergence, especially on highly
nonlinear problems.
Thus we recommend to use this function for fast prototyping on small-
dimensional problems only, and to implement analytical gradient as soon
as possible.
-- ALGLIB --
Copyright 16.05.2011 by Bochkanov Sergey
*************************************************************************/
void mincgcreatef(const ae_int_t n, const real_1d_array &x, const double diffstep, mincgstate &state);
void mincgcreatef(const real_1d_array &x, const double diffstep, mincgstate &state);
/*************************************************************************
This function sets stopping conditions for CG optimization algorithm.
INPUT PARAMETERS:
State - structure which stores algorithm state
EpsG - >=0
The subroutine finishes its work if the condition
|v|<EpsG is satisfied, where:
* |.| means Euclidian norm
* v - scaled gradient vector, v[i]=g[i]*s[i]
* g - gradient
* s - scaling coefficients set by MinCGSetScale()
EpsF - >=0
The subroutine finishes its work if on k+1-th iteration
the condition |F(k+1)-F(k)|<=EpsF*max{|F(k)|,|F(k+1)|,1}
is satisfied.
EpsX - >=0
The subroutine finishes its work if on k+1-th iteration
the condition |v|<=EpsX is fulfilled, where:
* |.| means Euclidian norm
* v - scaled step vector, v[i]=dx[i]/s[i]
* dx - ste pvector, dx=X(k+1)-X(k)
* s - scaling coefficients set by MinCGSetScale()
MaxIts - maximum number of iterations. If MaxIts=0, the number of
iterations is unlimited.
Passing EpsG=0, EpsF=0, EpsX=0 and MaxIts=0 (simultaneously) will lead to
automatic stopping criterion selection (small EpsX).
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void mincgsetcond(const mincgstate &state, const double epsg, const double epsf, const double epsx, const ae_int_t maxits);
/*************************************************************************
This function sets scaling coefficients for CG optimizer.
ALGLIB optimizers use scaling matrices to test stopping conditions (step
size and gradient are scaled before comparison with tolerances). Scale of
the I-th variable is a translation invariant measure of:
a) "how large" the variable is
b) how large the step should be to make significant changes in the function
Scaling is also used by finite difference variant of CG optimizer - step
along I-th axis is equal to DiffStep*S[I].
In most optimizers (and in the CG too) scaling is NOT a form of
preconditioning. It just affects stopping conditions. You should set
preconditioner by separate call to one of the MinCGSetPrec...() functions.
There is special preconditioning mode, however, which uses scaling
coefficients to form diagonal preconditioning matrix. You can turn this
mode on, if you want. But you should understand that scaling is not the
same thing as preconditioning - these are two different, although related
forms of tuning solver.
INPUT PARAMETERS:
State - structure stores algorithm state
S - array[N], non-zero scaling coefficients
S[i] may be negative, sign doesn't matter.
-- ALGLIB --
Copyright 14.01.2011 by Bochkanov Sergey
*************************************************************************/
void mincgsetscale(const mincgstate &state, const real_1d_array &s);
/*************************************************************************
This function turns on/off reporting.
INPUT PARAMETERS:
State - structure which stores algorithm state
NeedXRep- whether iteration reports are needed or not
If NeedXRep is True, algorithm will call rep() callback function if it is
provided to MinCGOptimize().
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void mincgsetxrep(const mincgstate &state, const bool needxrep);
/*************************************************************************
This function sets CG algorithm.
INPUT PARAMETERS:
State - structure which stores algorithm state
CGType - algorithm type:
* -1 automatic selection of the best algorithm
* 0 DY (Dai and Yuan) algorithm
* 1 Hybrid DY-HS algorithm
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void mincgsetcgtype(const mincgstate &state, const ae_int_t cgtype);
/*************************************************************************
This function sets maximum step length
INPUT PARAMETERS:
State - structure which stores algorithm state
StpMax - maximum step length, >=0. Set StpMax to 0.0, if you don't
want to limit step length.
Use this subroutine when you optimize target function which contains exp()
or other fast growing functions, and optimization algorithm makes too
large steps which leads to overflow. This function allows us to reject
steps that are too large (and therefore expose us to the possible
overflow) without actually calculating function value at the x+stp*d.
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void mincgsetstpmax(const mincgstate &state, const double stpmax);
/*************************************************************************
This function allows to suggest initial step length to the CG algorithm.
Suggested step length is used as starting point for the line search. It
can be useful when you have badly scaled problem, i.e. when ||grad||
(which is used as initial estimate for the first step) is many orders of
magnitude different from the desired step.
Line search may fail on such problems without good estimate of initial
step length. Imagine, for example, problem with ||grad||=10^50 and desired
step equal to 0.1 Line search function will use 10^50 as initial step,
then it will decrease step length by 2 (up to 20 attempts) and will get
10^44, which is still too large.
This function allows us to tell than line search should be started from
some moderate step length, like 1.0, so algorithm will be able to detect
desired step length in a several searches.
Default behavior (when no step is suggested) is to use preconditioner, if
it is available, to generate initial estimate of step length.
This function influences only first iteration of algorithm. It should be
called between MinCGCreate/MinCGRestartFrom() call and MinCGOptimize call.
Suggested step is ignored if you have preconditioner.
INPUT PARAMETERS:
State - structure used to store algorithm state.
Stp - initial estimate of the step length.
Can be zero (no estimate).
-- ALGLIB --
Copyright 30.07.2010 by Bochkanov Sergey
*************************************************************************/
void mincgsuggeststep(const mincgstate &state, const double stp);
/*************************************************************************
Modification of the preconditioner: preconditioning is turned off.
INPUT PARAMETERS:
State - structure which stores algorithm state
NOTE: you can change preconditioner "on the fly", during algorithm
iterations.
-- ALGLIB --
Copyright 13.10.2010 by Bochkanov Sergey
*************************************************************************/
void mincgsetprecdefault(const mincgstate &state);
/*************************************************************************
Modification of the preconditioner: diagonal of approximate Hessian is
used.
INPUT PARAMETERS:
State - structure which stores algorithm state
D - diagonal of the approximate Hessian, array[0..N-1],
(if larger, only leading N elements are used).
NOTE: you can change preconditioner "on the fly", during algorithm
iterations.
NOTE 2: D[i] should be positive. Exception will be thrown otherwise.
NOTE 3: you should pass diagonal of approximate Hessian - NOT ITS INVERSE.
-- ALGLIB --
Copyright 13.10.2010 by Bochkanov Sergey
*************************************************************************/
void mincgsetprecdiag(const mincgstate &state, const real_1d_array &d);
/*************************************************************************
Modification of the preconditioner: scale-based diagonal preconditioning.
This preconditioning mode can be useful when you don't have approximate
diagonal of Hessian, but you know that your variables are badly scaled
(for example, one variable is in [1,10], and another in [1000,100000]),
and most part of the ill-conditioning comes from different scales of vars.
In this case simple scale-based preconditioner, with H[i] = 1/(s[i]^2),
can greatly improve convergence.
IMPRTANT: you should set scale of your variables with MinCGSetScale() call
(before or after MinCGSetPrecScale() call). Without knowledge of the scale
of your variables scale-based preconditioner will be just unit matrix.
INPUT PARAMETERS:
State - structure which stores algorithm state
NOTE: you can change preconditioner "on the fly", during algorithm
iterations.
-- ALGLIB --
Copyright 13.10.2010 by Bochkanov Sergey
*************************************************************************/
void mincgsetprecscale(const mincgstate &state);
/*************************************************************************
This function provides reverse communication interface
Reverse communication interface is not documented or recommended to use.
See below for functions which provide better documented API
*************************************************************************/
bool mincgiteration(const mincgstate &state);
/*************************************************************************
This family of functions is used to launcn iterations of nonlinear optimizer
These functions accept following parameters:
state - algorithm state
func - callback which calculates function (or merit function)
value func at given point x
grad - callback which calculates function (or merit function)
value func and gradient grad at given point x
rep - optional callback which is called after each iteration
can be NULL
ptr - optional pointer which is passed to func/grad/hess/jac/rep
can be NULL
NOTES:
1. This function has two different implementations: one which uses exact
(analytical) user-supplied gradient, and one which uses function value
only and numerically differentiates function in order to obtain
gradient.
Depending on the specific function used to create optimizer object
(either MinCGCreate() for analytical gradient or MinCGCreateF() for
numerical differentiation) you should choose appropriate variant of
MinCGOptimize() - one which accepts function AND gradient or one which
accepts function ONLY.
Be careful to choose variant of MinCGOptimize() which corresponds to
your optimization scheme! Table below lists different combinations of
callback (function/gradient) passed to MinCGOptimize() and specific
function used to create optimizer.
| USER PASSED TO MinCGOptimize()
CREATED WITH | function only | function and gradient
------------------------------------------------------------
MinCGCreateF() | work FAIL
MinCGCreate() | FAIL work
Here "FAIL" denotes inappropriate combinations of optimizer creation
function and MinCGOptimize() version. Attemps to use such combination
(for example, to create optimizer with MinCGCreateF() and to pass
gradient information to MinCGOptimize()) will lead to exception being
thrown. Either you did not pass gradient when it WAS needed or you
passed gradient when it was NOT needed.
-- ALGLIB --
Copyright 20.04.2009 by Bochkanov Sergey
*************************************************************************/
void mincgoptimize(mincgstate &state,
void (*func)(const real_1d_array &x, double &func, void *ptr),
void (*rep)(const real_1d_array &x, double func, void *ptr) = NULL,
void *ptr = NULL);
void mincgoptimize(mincgstate &state,
void (*grad)(const real_1d_array &x, double &func, real_1d_array &grad, void *ptr),
void (*rep)(const real_1d_array &x, double func, void *ptr) = NULL,
void *ptr = NULL);
/*************************************************************************
Conjugate gradient results
INPUT PARAMETERS:
State - algorithm state
OUTPUT PARAMETERS:
X - array[0..N-1], solution
Rep - optimization report:
* Rep.TerminationType completetion code:
* 1 relative function improvement is no more than
EpsF.
* 2 relative step is no more than EpsX.
* 4 gradient norm is no more than EpsG
* 5 MaxIts steps was taken
* 7 stopping conditions are too stringent,
further improvement is impossible,
we return best X found so far
* 8 terminated by user
* Rep.IterationsCount contains iterations count
* NFEV countains number of function calculations
-- ALGLIB --
Copyright 20.04.2009 by Bochkanov Sergey
*************************************************************************/
void mincgresults(const mincgstate &state, real_1d_array &x, mincgreport &rep);
/*************************************************************************
Conjugate gradient results
Buffered implementation of MinCGResults(), which uses pre-allocated buffer
to store X[]. If buffer size is too small, it resizes buffer. It is
intended to be used in the inner cycles of performance critical algorithms
where array reallocation penalty is too large to be ignored.
-- ALGLIB --
Copyright 20.04.2009 by Bochkanov Sergey
*************************************************************************/
void mincgresultsbuf(const mincgstate &state, real_1d_array &x, mincgreport &rep);
/*************************************************************************
This subroutine restarts CG algorithm from new point. All optimization
parameters are left unchanged.
This function allows to solve multiple optimization problems (which
must have same number of dimensions) without object reallocation penalty.
INPUT PARAMETERS:
State - structure used to store algorithm state.
X - new starting point.
-- ALGLIB --
Copyright 30.07.2010 by Bochkanov Sergey
*************************************************************************/
void mincgrestartfrom(const mincgstate &state, const real_1d_array &x);
/*************************************************************************
BOUND CONSTRAINED OPTIMIZATION
WITH ADDITIONAL LINEAR EQUALITY AND INEQUALITY CONSTRAINTS
DESCRIPTION:
The subroutine minimizes function F(x) of N arguments subject to any
combination of:
* bound constraints
* linear inequality constraints
* linear equality constraints
REQUIREMENTS:
* user must provide function value and gradient
* starting point X0 must be feasible or
not too far away from the feasible set
* grad(f) must be Lipschitz continuous on a level set:
L = { x : f(x)<=f(x0) }
* function must be defined everywhere on the feasible set F
USAGE:
Constrained optimization if far more complex than the unconstrained one.
Here we give very brief outline of the BLEIC optimizer. We strongly recommend
you to read examples in the ALGLIB Reference Manual and to read ALGLIB User Guide
on optimization, which is available at http://www.alglib.net/optimization/
1. User initializes algorithm state with MinBLEICCreate() call
2. USer adds boundary and/or linear constraints by calling
MinBLEICSetBC() and MinBLEICSetLC() functions.
3. User sets stopping conditions for underlying unconstrained solver
with MinBLEICSetInnerCond() call.
This function controls accuracy of underlying optimization algorithm.
4. User sets stopping conditions for outer iteration by calling
MinBLEICSetOuterCond() function.
This function controls handling of boundary and inequality constraints.
5. Additionally, user may set limit on number of internal iterations
by MinBLEICSetMaxIts() call.
This function allows to prevent algorithm from looping forever.
6. User calls MinBLEICOptimize() function which takes algorithm state and
pointer (delegate, etc.) to callback function which calculates F/G.
7. User calls MinBLEICResults() to get solution
8. Optionally user may call MinBLEICRestartFrom() to solve another problem
with same N but another starting point.
MinBLEICRestartFrom() allows to reuse already initialized structure.
INPUT PARAMETERS:
N - problem dimension, N>0:
* if given, only leading N elements of X are used
* if not given, automatically determined from size ofX
X - starting point, array[N]:
* it is better to set X to a feasible point
* but X can be infeasible, in which case algorithm will try
to find feasible point first, using X as initial
approximation.
OUTPUT PARAMETERS:
State - structure stores algorithm state
-- ALGLIB --
Copyright 28.11.2010 by Bochkanov Sergey
*************************************************************************/
void minbleiccreate(const ae_int_t n, const real_1d_array &x, minbleicstate &state);
void minbleiccreate(const real_1d_array &x, minbleicstate &state);
/*************************************************************************
The subroutine is finite difference variant of MinBLEICCreate(). It uses
finite differences in order to differentiate target function.
Description below contains information which is specific to this function
only. We recommend to read comments on MinBLEICCreate() in order to get
more information about creation of BLEIC optimizer.
INPUT PARAMETERS:
N - problem dimension, N>0:
* if given, only leading N elements of X are used
* if not given, automatically determined from size of X
X - starting point, array[0..N-1].
DiffStep- differentiation step, >0
OUTPUT PARAMETERS:
State - structure which stores algorithm state
NOTES:
1. algorithm uses 4-point central formula for differentiation.
2. differentiation step along I-th axis is equal to DiffStep*S[I] where
S[] is scaling vector which can be set by MinBLEICSetScale() call.
3. we recommend you to use moderate values of differentiation step. Too
large step will result in too large truncation errors, while too small
step will result in too large numerical errors. 1.0E-6 can be good
value to start with.
4. Numerical differentiation is very inefficient - one gradient
calculation needs 4*N function evaluations. This function will work for
any N - either small (1...10), moderate (10...100) or large (100...).
However, performance penalty will be too severe for any N's except for
small ones.
We should also say that code which relies on numerical differentiation
is less robust and precise. CG needs exact gradient values. Imprecise
gradient may slow down convergence, especially on highly nonlinear
problems.
Thus we recommend to use this function for fast prototyping on small-
dimensional problems only, and to implement analytical gradient as soon
as possible.
-- ALGLIB --
Copyright 16.05.2011 by Bochkanov Sergey
*************************************************************************/
void minbleiccreatef(const ae_int_t n, const real_1d_array &x, const double diffstep, minbleicstate &state);
void minbleiccreatef(const real_1d_array &x, const double diffstep, minbleicstate &state);
/*************************************************************************
This function sets boundary constraints for BLEIC optimizer.
Boundary constraints are inactive by default (after initial creation).
They are preserved after algorithm restart with MinBLEICRestartFrom().
INPUT PARAMETERS:
State - structure stores algorithm state
BndL - lower bounds, array[N].
If some (all) variables are unbounded, you may specify
very small number or -INF.
BndU - upper bounds, array[N].
If some (all) variables are unbounded, you may specify
very large number or +INF.
NOTE 1: it is possible to specify BndL[i]=BndU[i]. In this case I-th
variable will be "frozen" at X[i]=BndL[i]=BndU[i].
NOTE 2: this solver has following useful properties:
* bound constraints are always satisfied exactly
* function is evaluated only INSIDE area specified by bound constraints,
even when numerical differentiation is used (algorithm adjusts nodes
according to boundary constraints)
-- ALGLIB --
Copyright 28.11.2010 by Bochkanov Sergey
*************************************************************************/
void minbleicsetbc(const minbleicstate &state, const real_1d_array &bndl, const real_1d_array &bndu);
/*************************************************************************
This function sets linear constraints for BLEIC optimizer.
Linear constraints are inactive by default (after initial creation).
They are preserved after algorithm restart with MinBLEICRestartFrom().
INPUT PARAMETERS:
State - structure previously allocated with MinBLEICCreate call.
C - linear constraints, array[K,N+1].
Each row of C represents one constraint, either equality
or inequality (see below):
* first N elements correspond to coefficients,
* last element corresponds to the right part.
All elements of C (including right part) must be finite.
CT - type of constraints, array[K]:
* if CT[i]>0, then I-th constraint is C[i,*]*x >= C[i,n+1]
* if CT[i]=0, then I-th constraint is C[i,*]*x = C[i,n+1]
* if CT[i]<0, then I-th constraint is C[i,*]*x <= C[i,n+1]
K - number of equality/inequality constraints, K>=0:
* if given, only leading K elements of C/CT are used
* if not given, automatically determined from sizes of C/CT
NOTE 1: linear (non-bound) constraints are satisfied only approximately:
* there always exists some minor violation (about Epsilon in magnitude)
due to rounding errors
* numerical differentiation, if used, may lead to function evaluations
outside of the feasible area, because algorithm does NOT change
numerical differentiation formula according to linear constraints.
If you want constraints to be satisfied exactly, try to reformulate your
problem in such manner that all constraints will become boundary ones
(this kind of constraints is always satisfied exactly, both in the final
solution and in all intermediate points).
-- ALGLIB --
Copyright 28.11.2010 by Bochkanov Sergey
*************************************************************************/
void minbleicsetlc(const minbleicstate &state, const real_2d_array &c, const integer_1d_array &ct, const ae_int_t k);
void minbleicsetlc(const minbleicstate &state, const real_2d_array &c, const integer_1d_array &ct);
/*************************************************************************
This function sets stopping conditions for the underlying nonlinear CG
optimizer. It controls overall accuracy of solution. These conditions
should be strict enough in order for algorithm to converge.
INPUT PARAMETERS:
State - structure which stores algorithm state
EpsG - >=0
The subroutine finishes its work if the condition
|v|<EpsG is satisfied, where:
* |.| means Euclidian norm
* v - scaled gradient vector, v[i]=g[i]*s[i]
* g - gradient
* s - scaling coefficients set by MinBLEICSetScale()
EpsF - >=0
The subroutine finishes its work if on k+1-th iteration
the condition |F(k+1)-F(k)|<=EpsF*max{|F(k)|,|F(k+1)|,1}
is satisfied.
EpsX - >=0
The subroutine finishes its work if on k+1-th iteration
the condition |v|<=EpsX is fulfilled, where:
* |.| means Euclidian norm
* v - scaled step vector, v[i]=dx[i]/s[i]
* dx - ste pvector, dx=X(k+1)-X(k)
* s - scaling coefficients set by MinBLEICSetScale()
Passing EpsG=0, EpsF=0 and EpsX=0 (simultaneously) will lead to
automatic stopping criterion selection.
These conditions are used to terminate inner iterations. However, you
need to tune termination conditions for outer iterations too.
-- ALGLIB --
Copyright 28.11.2010 by Bochkanov Sergey
*************************************************************************/
void minbleicsetinnercond(const minbleicstate &state, const double epsg, const double epsf, const double epsx);
/*************************************************************************
This function sets stopping conditions for outer iteration of BLEIC algo.
These conditions control accuracy of constraint handling and amount of
infeasibility allowed in the solution.
INPUT PARAMETERS:
State - structure which stores algorithm state
EpsX - >0, stopping condition on outer iteration step length
EpsI - >0, stopping condition on infeasibility
Both EpsX and EpsI must be non-zero.
MEANING OF EpsX
EpsX is a stopping condition for outer iterations. Algorithm will stop
when solution of the current modified subproblem will be within EpsX
(using 2-norm) of the previous solution.
MEANING OF EpsI
EpsI controls feasibility properties - algorithm won't stop until all
inequality constraints will be satisfied with error (distance from current
point to the feasible area) at most EpsI.
-- ALGLIB --
Copyright 28.11.2010 by Bochkanov Sergey
*************************************************************************/
void minbleicsetoutercond(const minbleicstate &state, const double epsx, const double epsi);
/*************************************************************************
This function sets scaling coefficients for BLEIC optimizer.
ALGLIB optimizers use scaling matrices to test stopping conditions (step
size and gradient are scaled before comparison with tolerances). Scale of
the I-th variable is a translation invariant measure of:
a) "how large" the variable is
b) how large the step should be to make significant changes in the function
Scaling is also used by finite difference variant of the optimizer - step
along I-th axis is equal to DiffStep*S[I].
In most optimizers (and in the BLEIC too) scaling is NOT a form of
preconditioning. It just affects stopping conditions. You should set
preconditioner by separate call to one of the MinBLEICSetPrec...()
functions.
There is a special preconditioning mode, however, which uses scaling
coefficients to form diagonal preconditioning matrix. You can turn this
mode on, if you want. But you should understand that scaling is not the
same thing as preconditioning - these are two different, although related
forms of tuning solver.
INPUT PARAMETERS:
State - structure stores algorithm state
S - array[N], non-zero scaling coefficients
S[i] may be negative, sign doesn't matter.
-- ALGLIB --
Copyright 14.01.2011 by Bochkanov Sergey
*************************************************************************/
void minbleicsetscale(const minbleicstate &state, const real_1d_array &s);
/*************************************************************************
Modification of the preconditioner: preconditioning is turned off.
INPUT PARAMETERS:
State - structure which stores algorithm state
-- ALGLIB --
Copyright 13.10.2010 by Bochkanov Sergey
*************************************************************************/
void minbleicsetprecdefault(const minbleicstate &state);
/*************************************************************************
Modification of the preconditioner: diagonal of approximate Hessian is
used.
INPUT PARAMETERS:
State - structure which stores algorithm state
D - diagonal of the approximate Hessian, array[0..N-1],
(if larger, only leading N elements are used).
NOTE 1: D[i] should be positive. Exception will be thrown otherwise.
NOTE 2: you should pass diagonal of approximate Hessian - NOT ITS INVERSE.
-- ALGLIB --
Copyright 13.10.2010 by Bochkanov Sergey
*************************************************************************/
void minbleicsetprecdiag(const minbleicstate &state, const real_1d_array &d);
/*************************************************************************
Modification of the preconditioner: scale-based diagonal preconditioning.
This preconditioning mode can be useful when you don't have approximate
diagonal of Hessian, but you know that your variables are badly scaled
(for example, one variable is in [1,10], and another in [1000,100000]),
and most part of the ill-conditioning comes from different scales of vars.
In this case simple scale-based preconditioner, with H[i] = 1/(s[i]^2),
can greatly improve convergence.
IMPRTANT: you should set scale of your variables with MinBLEICSetScale()
call (before or after MinBLEICSetPrecScale() call). Without knowledge of
the scale of your variables scale-based preconditioner will be just unit
matrix.
INPUT PARAMETERS:
State - structure which stores algorithm state
-- ALGLIB --
Copyright 13.10.2010 by Bochkanov Sergey
*************************************************************************/
void minbleicsetprecscale(const minbleicstate &state);
/*************************************************************************
This function allows to stop algorithm after specified number of inner
iterations.
INPUT PARAMETERS:
State - structure which stores algorithm state
MaxIts - maximum number of inner iterations.
If MaxIts=0, the number of iterations is unlimited.
-- ALGLIB --
Copyright 28.11.2010 by Bochkanov Sergey
*************************************************************************/
void minbleicsetmaxits(const minbleicstate &state, const ae_int_t maxits);
/*************************************************************************
This function turns on/off reporting.
INPUT PARAMETERS:
State - structure which stores algorithm state
NeedXRep- whether iteration reports are needed or not
If NeedXRep is True, algorithm will call rep() callback function if it is
provided to MinBLEICOptimize().
-- ALGLIB --
Copyright 28.11.2010 by Bochkanov Sergey
*************************************************************************/
void minbleicsetxrep(const minbleicstate &state, const bool needxrep);
/*************************************************************************
This function sets maximum step length
IMPORTANT: this feature is hard to combine with preconditioning. You can't
set upper limit on step length, when you solve optimization problem with
linear (non-boundary) constraints AND preconditioner turned on.
When non-boundary constraints are present, you have to either a) use
preconditioner, or b) use upper limit on step length. YOU CAN'T USE BOTH!
In this case algorithm will terminate with appropriate error code.
INPUT PARAMETERS:
State - structure which stores algorithm state
StpMax - maximum step length, >=0. Set StpMax to 0.0, if you don't
want to limit step length.
Use this subroutine when you optimize target function which contains exp()
or other fast growing functions, and optimization algorithm makes too
large steps which lead to overflow. This function allows us to reject
steps that are too large (and therefore expose us to the possible
overflow) without actually calculating function value at the x+stp*d.
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void minbleicsetstpmax(const minbleicstate &state, const double stpmax);
/*************************************************************************
This function provides reverse communication interface
Reverse communication interface is not documented or recommended to use.
See below for functions which provide better documented API
*************************************************************************/
bool minbleiciteration(const minbleicstate &state);
/*************************************************************************
This family of functions is used to launcn iterations of nonlinear optimizer
These functions accept following parameters:
state - algorithm state
func - callback which calculates function (or merit function)
value func at given point x
grad - callback which calculates function (or merit function)
value func and gradient grad at given point x
rep - optional callback which is called after each iteration
can be NULL
ptr - optional pointer which is passed to func/grad/hess/jac/rep
can be NULL
NOTES:
1. This function has two different implementations: one which uses exact
(analytical) user-supplied gradient, and one which uses function value
only and numerically differentiates function in order to obtain
gradient.
Depending on the specific function used to create optimizer object
(either MinBLEICCreate() for analytical gradient or MinBLEICCreateF()
for numerical differentiation) you should choose appropriate variant of
MinBLEICOptimize() - one which accepts function AND gradient or one
which accepts function ONLY.
Be careful to choose variant of MinBLEICOptimize() which corresponds to
your optimization scheme! Table below lists different combinations of
callback (function/gradient) passed to MinBLEICOptimize() and specific
function used to create optimizer.
| USER PASSED TO MinBLEICOptimize()
CREATED WITH | function only | function and gradient
------------------------------------------------------------
MinBLEICCreateF() | work FAIL
MinBLEICCreate() | FAIL work
Here "FAIL" denotes inappropriate combinations of optimizer creation
function and MinBLEICOptimize() version. Attemps to use such
combination (for example, to create optimizer with MinBLEICCreateF()
and to pass gradient information to MinCGOptimize()) will lead to
exception being thrown. Either you did not pass gradient when it WAS
needed or you passed gradient when it was NOT needed.
-- ALGLIB --
Copyright 28.11.2010 by Bochkanov Sergey
*************************************************************************/
void minbleicoptimize(minbleicstate &state,
void (*func)(const real_1d_array &x, double &func, void *ptr),
void (*rep)(const real_1d_array &x, double func, void *ptr) = NULL,
void *ptr = NULL);
void minbleicoptimize(minbleicstate &state,
void (*grad)(const real_1d_array &x, double &func, real_1d_array &grad, void *ptr),
void (*rep)(const real_1d_array &x, double func, void *ptr) = NULL,
void *ptr = NULL);
/*************************************************************************
BLEIC results
INPUT PARAMETERS:
State - algorithm state
OUTPUT PARAMETERS:
X - array[0..N-1], solution
Rep - optimization report. You should check Rep.TerminationType
in order to distinguish successful termination from
unsuccessful one.
More information about fields of this structure can be
found in the comments on MinBLEICReport datatype.
-- ALGLIB --
Copyright 28.11.2010 by Bochkanov Sergey
*************************************************************************/
void minbleicresults(const minbleicstate &state, real_1d_array &x, minbleicreport &rep);
/*************************************************************************
BLEIC results
Buffered implementation of MinBLEICResults() which uses pre-allocated buffer
to store X[]. If buffer size is too small, it resizes buffer. It is
intended to be used in the inner cycles of performance critical algorithms
where array reallocation penalty is too large to be ignored.
-- ALGLIB --
Copyright 28.11.2010 by Bochkanov Sergey
*************************************************************************/
void minbleicresultsbuf(const minbleicstate &state, real_1d_array &x, minbleicreport &rep);
/*************************************************************************
This subroutine restarts algorithm from new point.
All optimization parameters (including constraints) are left unchanged.
This function allows to solve multiple optimization problems (which
must have same number of dimensions) without object reallocation penalty.
INPUT PARAMETERS:
State - structure previously allocated with MinBLEICCreate call.
X - new starting point.
-- ALGLIB --
Copyright 28.11.2010 by Bochkanov Sergey
*************************************************************************/
void minbleicrestartfrom(const minbleicstate &state, const real_1d_array &x);
/*************************************************************************
LIMITED MEMORY BFGS METHOD FOR LARGE SCALE OPTIMIZATION
DESCRIPTION:
The subroutine minimizes function F(x) of N arguments by using a quasi-
Newton method (LBFGS scheme) which is optimized to use a minimum amount
of memory.
The subroutine generates the approximation of an inverse Hessian matrix by
using information about the last M steps of the algorithm (instead of N).
It lessens a required amount of memory from a value of order N^2 to a
value of order 2*N*M.
REQUIREMENTS:
Algorithm will request following information during its operation:
* function value F and its gradient G (simultaneously) at given point X
USAGE:
1. User initializes algorithm state with MinLBFGSCreate() call
2. User tunes solver parameters with MinLBFGSSetCond() MinLBFGSSetStpMax()
and other functions
3. User calls MinLBFGSOptimize() function which takes algorithm state and
pointer (delegate, etc.) to callback function which calculates F/G.
4. User calls MinLBFGSResults() to get solution
5. Optionally user may call MinLBFGSRestartFrom() to solve another problem
with same N/M but another starting point and/or another function.
MinLBFGSRestartFrom() allows to reuse already initialized structure.
INPUT PARAMETERS:
N - problem dimension. N>0
M - number of corrections in the BFGS scheme of Hessian
approximation update. Recommended value: 3<=M<=7. The smaller
value causes worse convergence, the bigger will not cause a
considerably better convergence, but will cause a fall in the
performance. M<=N.
X - initial solution approximation, array[0..N-1].
OUTPUT PARAMETERS:
State - structure which stores algorithm state
NOTES:
1. you may tune stopping conditions with MinLBFGSSetCond() function
2. if target function contains exp() or other fast growing functions, and
optimization algorithm makes too large steps which leads to overflow,
use MinLBFGSSetStpMax() function to bound algorithm's steps. However,
L-BFGS rarely needs such a tuning.
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void minlbfgscreate(const ae_int_t n, const ae_int_t m, const real_1d_array &x, minlbfgsstate &state);
void minlbfgscreate(const ae_int_t m, const real_1d_array &x, minlbfgsstate &state);
/*************************************************************************
The subroutine is finite difference variant of MinLBFGSCreate(). It uses
finite differences in order to differentiate target function.
Description below contains information which is specific to this function
only. We recommend to read comments on MinLBFGSCreate() in order to get
more information about creation of LBFGS optimizer.
INPUT PARAMETERS:
N - problem dimension, N>0:
* if given, only leading N elements of X are used
* if not given, automatically determined from size of X
M - number of corrections in the BFGS scheme of Hessian
approximation update. Recommended value: 3<=M<=7. The smaller
value causes worse convergence, the bigger will not cause a
considerably better convergence, but will cause a fall in the
performance. M<=N.
X - starting point, array[0..N-1].
DiffStep- differentiation step, >0
OUTPUT PARAMETERS:
State - structure which stores algorithm state
NOTES:
1. algorithm uses 4-point central formula for differentiation.
2. differentiation step along I-th axis is equal to DiffStep*S[I] where
S[] is scaling vector which can be set by MinLBFGSSetScale() call.
3. we recommend you to use moderate values of differentiation step. Too
large step will result in too large truncation errors, while too small
step will result in too large numerical errors. 1.0E-6 can be good
value to start with.
4. Numerical differentiation is very inefficient - one gradient
calculation needs 4*N function evaluations. This function will work for
any N - either small (1...10), moderate (10...100) or large (100...).
However, performance penalty will be too severe for any N's except for
small ones.
We should also say that code which relies on numerical differentiation
is less robust and precise. LBFGS needs exact gradient values.
Imprecise gradient may slow down convergence, especially on highly
nonlinear problems.
Thus we recommend to use this function for fast prototyping on small-
dimensional problems only, and to implement analytical gradient as soon
as possible.
-- ALGLIB --
Copyright 16.05.2011 by Bochkanov Sergey
*************************************************************************/
void minlbfgscreatef(const ae_int_t n, const ae_int_t m, const real_1d_array &x, const double diffstep, minlbfgsstate &state);
void minlbfgscreatef(const ae_int_t m, const real_1d_array &x, const double diffstep, minlbfgsstate &state);
/*************************************************************************
This function sets stopping conditions for L-BFGS optimization algorithm.
INPUT PARAMETERS:
State - structure which stores algorithm state
EpsG - >=0
The subroutine finishes its work if the condition
|v|<EpsG is satisfied, where:
* |.| means Euclidian norm
* v - scaled gradient vector, v[i]=g[i]*s[i]
* g - gradient
* s - scaling coefficients set by MinLBFGSSetScale()
EpsF - >=0
The subroutine finishes its work if on k+1-th iteration
the condition |F(k+1)-F(k)|<=EpsF*max{|F(k)|,|F(k+1)|,1}
is satisfied.
EpsX - >=0
The subroutine finishes its work if on k+1-th iteration
the condition |v|<=EpsX is fulfilled, where:
* |.| means Euclidian norm
* v - scaled step vector, v[i]=dx[i]/s[i]
* dx - ste pvector, dx=X(k+1)-X(k)
* s - scaling coefficients set by MinLBFGSSetScale()
MaxIts - maximum number of iterations. If MaxIts=0, the number of
iterations is unlimited.
Passing EpsG=0, EpsF=0, EpsX=0 and MaxIts=0 (simultaneously) will lead to
automatic stopping criterion selection (small EpsX).
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void minlbfgssetcond(const minlbfgsstate &state, const double epsg, const double epsf, const double epsx, const ae_int_t maxits);
/*************************************************************************
This function turns on/off reporting.
INPUT PARAMETERS:
State - structure which stores algorithm state
NeedXRep- whether iteration reports are needed or not
If NeedXRep is True, algorithm will call rep() callback function if it is
provided to MinLBFGSOptimize().
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void minlbfgssetxrep(const minlbfgsstate &state, const bool needxrep);
/*************************************************************************
This function sets maximum step length
INPUT PARAMETERS:
State - structure which stores algorithm state
StpMax - maximum step length, >=0. Set StpMax to 0.0 (default), if
you don't want to limit step length.
Use this subroutine when you optimize target function which contains exp()
or other fast growing functions, and optimization algorithm makes too
large steps which leads to overflow. This function allows us to reject
steps that are too large (and therefore expose us to the possible
overflow) without actually calculating function value at the x+stp*d.
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void minlbfgssetstpmax(const minlbfgsstate &state, const double stpmax);
/*************************************************************************
This function sets scaling coefficients for LBFGS optimizer.
ALGLIB optimizers use scaling matrices to test stopping conditions (step
size and gradient are scaled before comparison with tolerances). Scale of
the I-th variable is a translation invariant measure of:
a) "how large" the variable is
b) how large the step should be to make significant changes in the function
Scaling is also used by finite difference variant of the optimizer - step
along I-th axis is equal to DiffStep*S[I].
In most optimizers (and in the LBFGS too) scaling is NOT a form of
preconditioning. It just affects stopping conditions. You should set
preconditioner by separate call to one of the MinLBFGSSetPrec...()
functions.
There is special preconditioning mode, however, which uses scaling
coefficients to form diagonal preconditioning matrix. You can turn this
mode on, if you want. But you should understand that scaling is not the
same thing as preconditioning - these are two different, although related
forms of tuning solver.
INPUT PARAMETERS:
State - structure stores algorithm state
S - array[N], non-zero scaling coefficients
S[i] may be negative, sign doesn't matter.
-- ALGLIB --
Copyright 14.01.2011 by Bochkanov Sergey
*************************************************************************/
void minlbfgssetscale(const minlbfgsstate &state, const real_1d_array &s);
/*************************************************************************
Modification of the preconditioner: default preconditioner (simple
scaling, same for all elements of X) is used.
INPUT PARAMETERS:
State - structure which stores algorithm state
NOTE: you can change preconditioner "on the fly", during algorithm
iterations.
-- ALGLIB --
Copyright 13.10.2010 by Bochkanov Sergey
*************************************************************************/
void minlbfgssetprecdefault(const minlbfgsstate &state);
/*************************************************************************
Modification of the preconditioner: Cholesky factorization of approximate
Hessian is used.
INPUT PARAMETERS:
State - structure which stores algorithm state
P - triangular preconditioner, Cholesky factorization of
the approximate Hessian. array[0..N-1,0..N-1],
(if larger, only leading N elements are used).
IsUpper - whether upper or lower triangle of P is given
(other triangle is not referenced)
After call to this function preconditioner is changed to P (P is copied
into the internal buffer).
NOTE: you can change preconditioner "on the fly", during algorithm
iterations.
NOTE 2: P should be nonsingular. Exception will be thrown otherwise.
-- ALGLIB --
Copyright 13.10.2010 by Bochkanov Sergey
*************************************************************************/
void minlbfgssetpreccholesky(const minlbfgsstate &state, const real_2d_array &p, const bool isupper);
/*************************************************************************
Modification of the preconditioner: diagonal of approximate Hessian is
used.
INPUT PARAMETERS:
State - structure which stores algorithm state
D - diagonal of the approximate Hessian, array[0..N-1],
(if larger, only leading N elements are used).
NOTE: you can change preconditioner "on the fly", during algorithm
iterations.
NOTE 2: D[i] should be positive. Exception will be thrown otherwise.
NOTE 3: you should pass diagonal of approximate Hessian - NOT ITS INVERSE.
-- ALGLIB --
Copyright 13.10.2010 by Bochkanov Sergey
*************************************************************************/
void minlbfgssetprecdiag(const minlbfgsstate &state, const real_1d_array &d);
/*************************************************************************
Modification of the preconditioner: scale-based diagonal preconditioning.
This preconditioning mode can be useful when you don't have approximate
diagonal of Hessian, but you know that your variables are badly scaled
(for example, one variable is in [1,10], and another in [1000,100000]),
and most part of the ill-conditioning comes from different scales of vars.
In this case simple scale-based preconditioner, with H[i] = 1/(s[i]^2),
can greatly improve convergence.
IMPRTANT: you should set scale of your variables with MinLBFGSSetScale()
call (before or after MinLBFGSSetPrecScale() call). Without knowledge of
the scale of your variables scale-based preconditioner will be just unit
matrix.
INPUT PARAMETERS:
State - structure which stores algorithm state
-- ALGLIB --
Copyright 13.10.2010 by Bochkanov Sergey
*************************************************************************/
void minlbfgssetprecscale(const minlbfgsstate &state);
/*************************************************************************
This function provides reverse communication interface
Reverse communication interface is not documented or recommended to use.
See below for functions which provide better documented API
*************************************************************************/
bool minlbfgsiteration(const minlbfgsstate &state);
/*************************************************************************
This family of functions is used to launcn iterations of nonlinear optimizer
These functions accept following parameters:
state - algorithm state
func - callback which calculates function (or merit function)
value func at given point x
grad - callback which calculates function (or merit function)
value func and gradient grad at given point x
rep - optional callback which is called after each iteration
can be NULL
ptr - optional pointer which is passed to func/grad/hess/jac/rep
can be NULL
NOTES:
1. This function has two different implementations: one which uses exact
(analytical) user-supplied gradient, and one which uses function value
only and numerically differentiates function in order to obtain
gradient.
Depending on the specific function used to create optimizer object
(either MinLBFGSCreate() for analytical gradient or MinLBFGSCreateF()
for numerical differentiation) you should choose appropriate variant of
MinLBFGSOptimize() - one which accepts function AND gradient or one
which accepts function ONLY.
Be careful to choose variant of MinLBFGSOptimize() which corresponds to
your optimization scheme! Table below lists different combinations of
callback (function/gradient) passed to MinLBFGSOptimize() and specific
function used to create optimizer.
| USER PASSED TO MinLBFGSOptimize()
CREATED WITH | function only | function and gradient
------------------------------------------------------------
MinLBFGSCreateF() | work FAIL
MinLBFGSCreate() | FAIL work
Here "FAIL" denotes inappropriate combinations of optimizer creation
function and MinLBFGSOptimize() version. Attemps to use such
combination (for example, to create optimizer with MinLBFGSCreateF() and
to pass gradient information to MinCGOptimize()) will lead to exception
being thrown. Either you did not pass gradient when it WAS needed or
you passed gradient when it was NOT needed.
-- ALGLIB --
Copyright 20.03.2009 by Bochkanov Sergey
*************************************************************************/
void minlbfgsoptimize(minlbfgsstate &state,
void (*func)(const real_1d_array &x, double &func, void *ptr),
void (*rep)(const real_1d_array &x, double func, void *ptr) = NULL,
void *ptr = NULL);
void minlbfgsoptimize(minlbfgsstate &state,
void (*grad)(const real_1d_array &x, double &func, real_1d_array &grad, void *ptr),
void (*rep)(const real_1d_array &x, double func, void *ptr) = NULL,
void *ptr = NULL);
/*************************************************************************
L-BFGS algorithm results
INPUT PARAMETERS:
State - algorithm state
OUTPUT PARAMETERS:
X - array[0..N-1], solution
Rep - optimization report:
* Rep.TerminationType completetion code:
* -2 rounding errors prevent further improvement.
X contains best point found.
* -1 incorrect parameters were specified
* 1 relative function improvement is no more than
EpsF.
* 2 relative step is no more than EpsX.
* 4 gradient norm is no more than EpsG
* 5 MaxIts steps was taken
* 7 stopping conditions are too stringent,
further improvement is impossible
* Rep.IterationsCount contains iterations count
* NFEV countains number of function calculations
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void minlbfgsresults(const minlbfgsstate &state, real_1d_array &x, minlbfgsreport &rep);
/*************************************************************************
L-BFGS algorithm results
Buffered implementation of MinLBFGSResults which uses pre-allocated buffer
to store X[]. If buffer size is too small, it resizes buffer. It is
intended to be used in the inner cycles of performance critical algorithms
where array reallocation penalty is too large to be ignored.
-- ALGLIB --
Copyright 20.08.2010 by Bochkanov Sergey
*************************************************************************/
void minlbfgsresultsbuf(const minlbfgsstate &state, real_1d_array &x, minlbfgsreport &rep);
/*************************************************************************
This subroutine restarts LBFGS algorithm from new point. All optimization
parameters are left unchanged.
This function allows to solve multiple optimization problems (which
must have same number of dimensions) without object reallocation penalty.
INPUT PARAMETERS:
State - structure used to store algorithm state
X - new starting point.
-- ALGLIB --
Copyright 30.07.2010 by Bochkanov Sergey
*************************************************************************/
void minlbfgsrestartfrom(const minlbfgsstate &state, const real_1d_array &x);
/*************************************************************************
CONSTRAINED QUADRATIC PROGRAMMING
The subroutine creates QP optimizer. After initial creation, it contains
default optimization problem with zero quadratic and linear terms and no
constraints. You should set quadratic/linear terms with calls to functions
provided by MinQP subpackage.
INPUT PARAMETERS:
N - problem size
OUTPUT PARAMETERS:
State - optimizer with zero quadratic/linear terms
and no constraints
-- ALGLIB --
Copyright 11.01.2011 by Bochkanov Sergey
*************************************************************************/
void minqpcreate(const ae_int_t n, minqpstate &state);
/*************************************************************************
This function sets linear term for QP solver.
By default, linear term is zero.
INPUT PARAMETERS:
State - structure which stores algorithm state
B - linear term, array[N].
-- ALGLIB --
Copyright 11.01.2011 by Bochkanov Sergey
*************************************************************************/
void minqpsetlinearterm(const minqpstate &state, const real_1d_array &b);
/*************************************************************************
This function sets quadratic term for QP solver.
By default quadratic term is zero.
IMPORTANT: this solver minimizes following function:
f(x) = 0.5*x'*A*x + b'*x.
Note that quadratic term has 0.5 before it. So if you want to minimize
f(x) = x^2 + x
you should rewrite your problem as follows:
f(x) = 0.5*(2*x^2) + x
and your matrix A will be equal to [[2.0]], not to [[1.0]]
INPUT PARAMETERS:
State - structure which stores algorithm state
A - matrix, array[N,N]
IsUpper - (optional) storage type:
* if True, symmetric matrix A is given by its upper
triangle, and the lower triangle isnt used
* if False, symmetric matrix A is given by its lower
triangle, and the upper triangle isnt used
* if not given, both lower and upper triangles must be
filled.
-- ALGLIB --
Copyright 11.01.2011 by Bochkanov Sergey
*************************************************************************/
void minqpsetquadraticterm(const minqpstate &state, const real_2d_array &a, const bool isupper);
void minqpsetquadraticterm(const minqpstate &state, const real_2d_array &a);
/*************************************************************************
This function sets starting point for QP solver. It is useful to have
good initial approximation to the solution, because it will increase
speed of convergence and identification of active constraints.
INPUT PARAMETERS:
State - structure which stores algorithm state
X - starting point, array[N].
-- ALGLIB --
Copyright 11.01.2011 by Bochkanov Sergey
*************************************************************************/
void minqpsetstartingpoint(const minqpstate &state, const real_1d_array &x);
/*************************************************************************
This function sets origin for QP solver. By default, following QP program
is solved:
min(0.5*x'*A*x+b'*x)
This function allows to solve different problem:
min(0.5*(x-x_origin)'*A*(x-x_origin)+b'*(x-x_origin))
INPUT PARAMETERS:
State - structure which stores algorithm state
XOrigin - origin, array[N].
-- ALGLIB --
Copyright 11.01.2011 by Bochkanov Sergey
*************************************************************************/
void minqpsetorigin(const minqpstate &state, const real_1d_array &xorigin);
/*************************************************************************
This function tells solver to use Cholesky-based algorithm.
Cholesky-based algorithm can be used when:
* problem is convex
* there is no constraints or only boundary constraints are present
This algorithm has O(N^3) complexity for unconstrained problem and is up
to several times slower on bound constrained problems (these additional
iterations are needed to identify active constraints).
INPUT PARAMETERS:
State - structure which stores algorithm state
-- ALGLIB --
Copyright 11.01.2011 by Bochkanov Sergey
*************************************************************************/
void minqpsetalgocholesky(const minqpstate &state);
/*************************************************************************
This function sets boundary constraints for QP solver
Boundary constraints are inactive by default (after initial creation).
After being set, they are preserved until explicitly turned off with
another SetBC() call.
INPUT PARAMETERS:
State - structure stores algorithm state
BndL - lower bounds, array[N].
If some (all) variables are unbounded, you may specify
very small number or -INF (latter is recommended because
it will allow solver to use better algorithm).
BndU - upper bounds, array[N].
If some (all) variables are unbounded, you may specify
very large number or +INF (latter is recommended because
it will allow solver to use better algorithm).
NOTE: it is possible to specify BndL[i]=BndU[i]. In this case I-th
variable will be "frozen" at X[i]=BndL[i]=BndU[i].
-- ALGLIB --
Copyright 11.01.2011 by Bochkanov Sergey
*************************************************************************/
void minqpsetbc(const minqpstate &state, const real_1d_array &bndl, const real_1d_array &bndu);
/*************************************************************************
This function solves quadratic programming problem.
You should call it after setting solver options with MinQPSet...() calls.
INPUT PARAMETERS:
State - algorithm state
You should use MinQPResults() function to access results after calls
to this function.
-- ALGLIB --
Copyright 11.01.2011 by Bochkanov Sergey
*************************************************************************/
void minqpoptimize(const minqpstate &state);
/*************************************************************************
QP solver results
INPUT PARAMETERS:
State - algorithm state
OUTPUT PARAMETERS:
X - array[0..N-1], solution
Rep - optimization report. You should check Rep.TerminationType,
which contains completion code, and you may check another
fields which contain another information about algorithm
functioning.
-- ALGLIB --
Copyright 11.01.2011 by Bochkanov Sergey
*************************************************************************/
void minqpresults(const minqpstate &state, real_1d_array &x, minqpreport &rep);
/*************************************************************************
QP results
Buffered implementation of MinQPResults() which uses pre-allocated buffer
to store X[]. If buffer size is too small, it resizes buffer. It is
intended to be used in the inner cycles of performance critical algorithms
where array reallocation penalty is too large to be ignored.
-- ALGLIB --
Copyright 11.01.2011 by Bochkanov Sergey
*************************************************************************/
void minqpresultsbuf(const minqpstate &state, real_1d_array &x, minqpreport &rep);
/*************************************************************************
IMPROVED LEVENBERG-MARQUARDT METHOD FOR
NON-LINEAR LEAST SQUARES OPTIMIZATION
DESCRIPTION:
This function is used to find minimum of function which is represented as
sum of squares:
F(x) = f[0]^2(x[0],...,x[n-1]) + ... + f[m-1]^2(x[0],...,x[n-1])
using value of function vector f[] and Jacobian of f[].
REQUIREMENTS:
This algorithm will request following information during its operation:
* function vector f[] at given point X
* function vector f[] and Jacobian of f[] (simultaneously) at given point
There are several overloaded versions of MinLMOptimize() function which
correspond to different LM-like optimization algorithms provided by this
unit. You should choose version which accepts fvec() and jac() callbacks.
First one is used to calculate f[] at given point, second one calculates
f[] and Jacobian df[i]/dx[j].
You can try to initialize MinLMState structure with VJ function and then
use incorrect version of MinLMOptimize() (for example, version which
works with general form function and does not provide Jacobian), but it
will lead to exception being thrown after first attempt to calculate
Jacobian.
USAGE:
1. User initializes algorithm state with MinLMCreateVJ() call
2. User tunes solver parameters with MinLMSetCond(), MinLMSetStpMax() and
other functions
3. User calls MinLMOptimize() function which takes algorithm state and
callback functions.
4. User calls MinLMResults() to get solution
5. Optionally, user may call MinLMRestartFrom() to solve another problem
with same N/M but another starting point and/or another function.
MinLMRestartFrom() allows to reuse already initialized structure.
INPUT PARAMETERS:
N - dimension, N>1
* if given, only leading N elements of X are used
* if not given, automatically determined from size of X
M - number of functions f[i]
X - initial solution, array[0..N-1]
OUTPUT PARAMETERS:
State - structure which stores algorithm state
NOTES:
1. you may tune stopping conditions with MinLMSetCond() function
2. if target function contains exp() or other fast growing functions, and
optimization algorithm makes too large steps which leads to overflow,
use MinLMSetStpMax() function to bound algorithm's steps.
-- ALGLIB --
Copyright 30.03.2009 by Bochkanov Sergey
*************************************************************************/
void minlmcreatevj(const ae_int_t n, const ae_int_t m, const real_1d_array &x, minlmstate &state);
void minlmcreatevj(const ae_int_t m, const real_1d_array &x, minlmstate &state);
/*************************************************************************
IMPROVED LEVENBERG-MARQUARDT METHOD FOR
NON-LINEAR LEAST SQUARES OPTIMIZATION
DESCRIPTION:
This function is used to find minimum of function which is represented as
sum of squares:
F(x) = f[0]^2(x[0],...,x[n-1]) + ... + f[m-1]^2(x[0],...,x[n-1])
using value of function vector f[] only. Finite differences are used to
calculate Jacobian.
REQUIREMENTS:
This algorithm will request following information during its operation:
* function vector f[] at given point X
There are several overloaded versions of MinLMOptimize() function which
correspond to different LM-like optimization algorithms provided by this
unit. You should choose version which accepts fvec() callback.
You can try to initialize MinLMState structure with VJ function and then
use incorrect version of MinLMOptimize() (for example, version which
works with general form function and does not accept function vector), but
it will lead to exception being thrown after first attempt to calculate
Jacobian.
USAGE:
1. User initializes algorithm state with MinLMCreateV() call
2. User tunes solver parameters with MinLMSetCond(), MinLMSetStpMax() and
other functions
3. User calls MinLMOptimize() function which takes algorithm state and
callback functions.
4. User calls MinLMResults() to get solution
5. Optionally, user may call MinLMRestartFrom() to solve another problem
with same N/M but another starting point and/or another function.
MinLMRestartFrom() allows to reuse already initialized structure.
INPUT PARAMETERS:
N - dimension, N>1
* if given, only leading N elements of X are used
* if not given, automatically determined from size of X
M - number of functions f[i]
X - initial solution, array[0..N-1]
DiffStep- differentiation step, >0
OUTPUT PARAMETERS:
State - structure which stores algorithm state
See also MinLMIteration, MinLMResults.
NOTES:
1. you may tune stopping conditions with MinLMSetCond() function
2. if target function contains exp() or other fast growing functions, and
optimization algorithm makes too large steps which leads to overflow,
use MinLMSetStpMax() function to bound algorithm's steps.
-- ALGLIB --
Copyright 30.03.2009 by Bochkanov Sergey
*************************************************************************/
void minlmcreatev(const ae_int_t n, const ae_int_t m, const real_1d_array &x, const double diffstep, minlmstate &state);
void minlmcreatev(const ae_int_t m, const real_1d_array &x, const double diffstep, minlmstate &state);
/*************************************************************************
LEVENBERG-MARQUARDT-LIKE METHOD FOR NON-LINEAR OPTIMIZATION
DESCRIPTION:
This function is used to find minimum of general form (not "sum-of-
-squares") function
F = F(x[0], ..., x[n-1])
using its gradient and Hessian. Levenberg-Marquardt modification with
L-BFGS pre-optimization and internal pre-conditioned L-BFGS optimization
after each Levenberg-Marquardt step is used.
REQUIREMENTS:
This algorithm will request following information during its operation:
* function value F at given point X
* F and gradient G (simultaneously) at given point X
* F, G and Hessian H (simultaneously) at given point X
There are several overloaded versions of MinLMOptimize() function which
correspond to different LM-like optimization algorithms provided by this
unit. You should choose version which accepts func(), grad() and hess()
function pointers. First pointer is used to calculate F at given point,
second one calculates F(x) and grad F(x), third one calculates F(x),
grad F(x), hess F(x).
You can try to initialize MinLMState structure with FGH-function and then
use incorrect version of MinLMOptimize() (for example, version which does
not provide Hessian matrix), but it will lead to exception being thrown
after first attempt to calculate Hessian.
USAGE:
1. User initializes algorithm state with MinLMCreateFGH() call
2. User tunes solver parameters with MinLMSetCond(), MinLMSetStpMax() and
other functions
3. User calls MinLMOptimize() function which takes algorithm state and
pointers (delegates, etc.) to callback functions.
4. User calls MinLMResults() to get solution
5. Optionally, user may call MinLMRestartFrom() to solve another problem
with same N but another starting point and/or another function.
MinLMRestartFrom() allows to reuse already initialized structure.
INPUT PARAMETERS:
N - dimension, N>1
* if given, only leading N elements of X are used
* if not given, automatically determined from size of X
X - initial solution, array[0..N-1]
OUTPUT PARAMETERS:
State - structure which stores algorithm state
NOTES:
1. you may tune stopping conditions with MinLMSetCond() function
2. if target function contains exp() or other fast growing functions, and
optimization algorithm makes too large steps which leads to overflow,
use MinLMSetStpMax() function to bound algorithm's steps.
-- ALGLIB --
Copyright 30.03.2009 by Bochkanov Sergey
*************************************************************************/
void minlmcreatefgh(const ae_int_t n, const real_1d_array &x, minlmstate &state);
void minlmcreatefgh(const real_1d_array &x, minlmstate &state);
/*************************************************************************
This function sets stopping conditions for Levenberg-Marquardt optimization
algorithm.
INPUT PARAMETERS:
State - structure which stores algorithm state
EpsG - >=0
The subroutine finishes its work if the condition
|v|<EpsG is satisfied, where:
* |.| means Euclidian norm
* v - scaled gradient vector, v[i]=g[i]*s[i]
* g - gradient
* s - scaling coefficients set by MinLMSetScale()
EpsF - >=0
The subroutine finishes its work if on k+1-th iteration
the condition |F(k+1)-F(k)|<=EpsF*max{|F(k)|,|F(k+1)|,1}
is satisfied.
EpsX - >=0
The subroutine finishes its work if on k+1-th iteration
the condition |v|<=EpsX is fulfilled, where:
* |.| means Euclidian norm
* v - scaled step vector, v[i]=dx[i]/s[i]
* dx - ste pvector, dx=X(k+1)-X(k)
* s - scaling coefficients set by MinLMSetScale()
MaxIts - maximum number of iterations. If MaxIts=0, the number of
iterations is unlimited. Only Levenberg-Marquardt
iterations are counted (L-BFGS/CG iterations are NOT
counted because their cost is very low compared to that of
LM).
Passing EpsG=0, EpsF=0, EpsX=0 and MaxIts=0 (simultaneously) will lead to
automatic stopping criterion selection (small EpsX).
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void minlmsetcond(const minlmstate &state, const double epsg, const double epsf, const double epsx, const ae_int_t maxits);
/*************************************************************************
This function turns on/off reporting.
INPUT PARAMETERS:
State - structure which stores algorithm state
NeedXRep- whether iteration reports are needed or not
If NeedXRep is True, algorithm will call rep() callback function if it is
provided to MinLMOptimize(). Both Levenberg-Marquardt and internal L-BFGS
iterations are reported.
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void minlmsetxrep(const minlmstate &state, const bool needxrep);
/*************************************************************************
This function sets maximum step length
INPUT PARAMETERS:
State - structure which stores algorithm state
StpMax - maximum step length, >=0. Set StpMax to 0.0, if you don't
want to limit step length.
Use this subroutine when you optimize target function which contains exp()
or other fast growing functions, and optimization algorithm makes too
large steps which leads to overflow. This function allows us to reject
steps that are too large (and therefore expose us to the possible
overflow) without actually calculating function value at the x+stp*d.
NOTE: non-zero StpMax leads to moderate performance degradation because
intermediate step of preconditioned L-BFGS optimization is incompatible
with limits on step size.
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void minlmsetstpmax(const minlmstate &state, const double stpmax);
/*************************************************************************
This function sets scaling coefficients for LM optimizer.
ALGLIB optimizers use scaling matrices to test stopping conditions (step
size and gradient are scaled before comparison with tolerances). Scale of
the I-th variable is a translation invariant measure of:
a) "how large" the variable is
b) how large the step should be to make significant changes in the function
Generally, scale is NOT considered to be a form of preconditioner. But LM
optimizer is unique in that it uses scaling matrix both in the stopping
condition tests and as Marquardt damping factor.
Proper scaling is very important for the algorithm performance. It is less
important for the quality of results, but still has some influence (it is
easier to converge when variables are properly scaled, so premature
stopping is possible when very badly scalled variables are combined with
relaxed stopping conditions).
INPUT PARAMETERS:
State - structure stores algorithm state
S - array[N], non-zero scaling coefficients
S[i] may be negative, sign doesn't matter.
-- ALGLIB --
Copyright 14.01.2011 by Bochkanov Sergey
*************************************************************************/
void minlmsetscale(const minlmstate &state, const real_1d_array &s);
/*************************************************************************
This function sets boundary constraints for LM optimizer
Boundary constraints are inactive by default (after initial creation).
They are preserved until explicitly turned off with another SetBC() call.
INPUT PARAMETERS:
State - structure stores algorithm state
BndL - lower bounds, array[N].
If some (all) variables are unbounded, you may specify
very small number or -INF (latter is recommended because
it will allow solver to use better algorithm).
BndU - upper bounds, array[N].
If some (all) variables are unbounded, you may specify
very large number or +INF (latter is recommended because
it will allow solver to use better algorithm).
NOTE 1: it is possible to specify BndL[i]=BndU[i]. In this case I-th
variable will be "frozen" at X[i]=BndL[i]=BndU[i].
NOTE 2: this solver has following useful properties:
* bound constraints are always satisfied exactly
* function is evaluated only INSIDE area specified by bound constraints
or at its boundary
-- ALGLIB --
Copyright 14.01.2011 by Bochkanov Sergey
*************************************************************************/
void minlmsetbc(const minlmstate &state, const real_1d_array &bndl, const real_1d_array &bndu);
/*************************************************************************
This function is used to change acceleration settings
You can choose between three acceleration strategies:
* AccType=0, no acceleration.
* AccType=1, secant updates are used to update quadratic model after each
iteration. After fixed number of iterations (or after model breakdown)
we recalculate quadratic model using analytic Jacobian or finite
differences. Number of secant-based iterations depends on optimization
settings: about 3 iterations - when we have analytic Jacobian, up to 2*N
iterations - when we use finite differences to calculate Jacobian.
AccType=1 is recommended when Jacobian calculation cost is prohibitive
high (several Mx1 function vector calculations followed by several NxN
Cholesky factorizations are faster than calculation of one M*N Jacobian).
It should also be used when we have no Jacobian, because finite difference
approximation takes too much time to compute.
Table below list optimization protocols (XYZ protocol corresponds to
MinLMCreateXYZ) and acceleration types they support (and use by default).
ACCELERATION TYPES SUPPORTED BY OPTIMIZATION PROTOCOLS:
protocol 0 1 comment
V + +
VJ + +
FGH +
DAFAULT VALUES:
protocol 0 1 comment
V x without acceleration it is so slooooooooow
VJ x
FGH x
NOTE: this function should be called before optimization. Attempt to call
it during algorithm iterations may result in unexpected behavior.
NOTE: attempt to call this function with unsupported protocol/acceleration
combination will result in exception being thrown.
-- ALGLIB --
Copyright 14.10.2010 by Bochkanov Sergey
*************************************************************************/
void minlmsetacctype(const minlmstate &state, const ae_int_t acctype);
/*************************************************************************
This function provides reverse communication interface
Reverse communication interface is not documented or recommended to use.
See below for functions which provide better documented API
*************************************************************************/
bool minlmiteration(const minlmstate &state);
/*************************************************************************
This family of functions is used to launcn iterations of nonlinear optimizer
These functions accept following parameters:
state - algorithm state
func - callback which calculates function (or merit function)
value func at given point x
grad - callback which calculates function (or merit function)
value func and gradient grad at given point x
hess - callback which calculates function (or merit function)
value func, gradient grad and Hessian hess at given point x
fvec - callback which calculates function vector fi[]
at given point x
jac - callback which calculates function vector fi[]
and Jacobian jac at given point x
rep - optional callback which is called after each iteration
can be NULL
ptr - optional pointer which is passed to func/grad/hess/jac/rep
can be NULL
NOTES:
1. Depending on function used to create state structure, this algorithm
may accept Jacobian and/or Hessian and/or gradient. According to the
said above, there ase several versions of this function, which accept
different sets of callbacks.
This flexibility opens way to subtle errors - you may create state with
MinLMCreateFGH() (optimization using Hessian), but call function which
does not accept Hessian. So when algorithm will request Hessian, there
will be no callback to call. In this case exception will be thrown.
Be careful to avoid such errors because there is no way to find them at
compile time - you can see them at runtime only.
-- ALGLIB --
Copyright 10.03.2009 by Bochkanov Sergey
*************************************************************************/
void minlmoptimize(minlmstate &state,
void (*fvec)(const real_1d_array &x, real_1d_array &fi, void *ptr),
void (*rep)(const real_1d_array &x, double func, void *ptr) = NULL,
void *ptr = NULL);
void minlmoptimize(minlmstate &state,
void (*fvec)(const real_1d_array &x, real_1d_array &fi, void *ptr),
void (*jac)(const real_1d_array &x, real_1d_array &fi, real_2d_array &jac, void *ptr),
void (*rep)(const real_1d_array &x, double func, void *ptr) = NULL,
void *ptr = NULL);
void minlmoptimize(minlmstate &state,
void (*func)(const real_1d_array &x, double &func, void *ptr),
void (*grad)(const real_1d_array &x, double &func, real_1d_array &grad, void *ptr),
void (*hess)(const real_1d_array &x, double &func, real_1d_array &grad, real_2d_array &hess, void *ptr),
void (*rep)(const real_1d_array &x, double func, void *ptr) = NULL,
void *ptr = NULL);
void minlmoptimize(minlmstate &state,
void (*func)(const real_1d_array &x, double &func, void *ptr),
void (*jac)(const real_1d_array &x, real_1d_array &fi, real_2d_array &jac, void *ptr),
void (*rep)(const real_1d_array &x, double func, void *ptr) = NULL,
void *ptr = NULL);
void minlmoptimize(minlmstate &state,
void (*func)(const real_1d_array &x, double &func, void *ptr),
void (*grad)(const real_1d_array &x, double &func, real_1d_array &grad, void *ptr),
void (*jac)(const real_1d_array &x, real_1d_array &fi, real_2d_array &jac, void *ptr),
void (*rep)(const real_1d_array &x, double func, void *ptr) = NULL,
void *ptr = NULL);
/*************************************************************************
Levenberg-Marquardt algorithm results
INPUT PARAMETERS:
State - algorithm state
OUTPUT PARAMETERS:
X - array[0..N-1], solution
Rep - optimization report;
see comments for this structure for more info.
-- ALGLIB --
Copyright 10.03.2009 by Bochkanov Sergey
*************************************************************************/
void minlmresults(const minlmstate &state, real_1d_array &x, minlmreport &rep);
/*************************************************************************
Levenberg-Marquardt algorithm results
Buffered implementation of MinLMResults(), which uses pre-allocated buffer
to store X[]. If buffer size is too small, it resizes buffer. It is
intended to be used in the inner cycles of performance critical algorithms
where array reallocation penalty is too large to be ignored.
-- ALGLIB --
Copyright 10.03.2009 by Bochkanov Sergey
*************************************************************************/
void minlmresultsbuf(const minlmstate &state, real_1d_array &x, minlmreport &rep);
/*************************************************************************
This subroutine restarts LM algorithm from new point. All optimization
parameters are left unchanged.
This function allows to solve multiple optimization problems (which
must have same number of dimensions) without object reallocation penalty.
INPUT PARAMETERS:
State - structure used for reverse communication previously
allocated with MinLMCreateXXX call.
X - new starting point.
-- ALGLIB --
Copyright 30.07.2010 by Bochkanov Sergey
*************************************************************************/
void minlmrestartfrom(const minlmstate &state, const real_1d_array &x);
/*************************************************************************
This is obsolete function.
Since ALGLIB 3.3 it is equivalent to MinLMCreateVJ().
-- ALGLIB --
Copyright 30.03.2009 by Bochkanov Sergey
*************************************************************************/
void minlmcreatevgj(const ae_int_t n, const ae_int_t m, const real_1d_array &x, minlmstate &state);
void minlmcreatevgj(const ae_int_t m, const real_1d_array &x, minlmstate &state);
/*************************************************************************
This is obsolete function.
Since ALGLIB 3.3 it is equivalent to MinLMCreateFJ().
-- ALGLIB --
Copyright 30.03.2009 by Bochkanov Sergey
*************************************************************************/
void minlmcreatefgj(const ae_int_t n, const ae_int_t m, const real_1d_array &x, minlmstate &state);
void minlmcreatefgj(const ae_int_t m, const real_1d_array &x, minlmstate &state);
/*************************************************************************
This function is considered obsolete since ALGLIB 3.1.0 and is present for
backward compatibility only. We recommend to use MinLMCreateVJ, which
provides similar, but more consistent and feature-rich interface.
-- ALGLIB --
Copyright 30.03.2009 by Bochkanov Sergey
*************************************************************************/
void minlmcreatefj(const ae_int_t n, const ae_int_t m, const real_1d_array &x, minlmstate &state);
void minlmcreatefj(const ae_int_t m, const real_1d_array &x, minlmstate &state);
/*************************************************************************
Obsolete function, use MinLBFGSSetPrecDefault() instead.
-- ALGLIB --
Copyright 13.10.2010 by Bochkanov Sergey
*************************************************************************/
void minlbfgssetdefaultpreconditioner(const minlbfgsstate &state);
/*************************************************************************
Obsolete function, use MinLBFGSSetCholeskyPreconditioner() instead.
-- ALGLIB --
Copyright 13.10.2010 by Bochkanov Sergey
*************************************************************************/
void minlbfgssetcholeskypreconditioner(const minlbfgsstate &state, const real_2d_array &p, const bool isupper);
/*************************************************************************
This is obsolete function which was used by previous version of the BLEIC
optimizer. It does nothing in the current version of BLEIC.
-- ALGLIB --
Copyright 28.11.2010 by Bochkanov Sergey
*************************************************************************/
void minbleicsetbarrierwidth(const minbleicstate &state, const double mu);
/*************************************************************************
This is obsolete function which was used by previous version of the BLEIC
optimizer. It does nothing in the current version of BLEIC.
-- ALGLIB --
Copyright 28.11.2010 by Bochkanov Sergey
*************************************************************************/
void minbleicsetbarrierdecay(const minbleicstate &state, const double mudecay);
/*************************************************************************
Obsolete optimization algorithm.
Was replaced by MinBLEIC subpackage.
-- ALGLIB --
Copyright 25.03.2010 by Bochkanov Sergey
*************************************************************************/
void minasacreate(const ae_int_t n, const real_1d_array &x, const real_1d_array &bndl, const real_1d_array &bndu, minasastate &state);
void minasacreate(const real_1d_array &x, const real_1d_array &bndl, const real_1d_array &bndu, minasastate &state);
/*************************************************************************
Obsolete optimization algorithm.
Was replaced by MinBLEIC subpackage.
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void minasasetcond(const minasastate &state, const double epsg, const double epsf, const double epsx, const ae_int_t maxits);
/*************************************************************************
Obsolete optimization algorithm.
Was replaced by MinBLEIC subpackage.
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void minasasetxrep(const minasastate &state, const bool needxrep);
/*************************************************************************
Obsolete optimization algorithm.
Was replaced by MinBLEIC subpackage.
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void minasasetalgorithm(const minasastate &state, const ae_int_t algotype);
/*************************************************************************
Obsolete optimization algorithm.
Was replaced by MinBLEIC subpackage.
-- ALGLIB --
Copyright 02.04.2010 by Bochkanov Sergey
*************************************************************************/
void minasasetstpmax(const minasastate &state, const double stpmax);
/*************************************************************************
This function provides reverse communication interface
Reverse communication interface is not documented or recommended to use.
See below for functions which provide better documented API
*************************************************************************/
bool minasaiteration(const minasastate &state);
/*************************************************************************
This family of functions is used to launcn iterations of nonlinear optimizer
These functions accept following parameters:
state - algorithm state
grad - callback which calculates function (or merit function)
value func and gradient grad at given point x
rep - optional callback which is called after each iteration
can be NULL
ptr - optional pointer which is passed to func/grad/hess/jac/rep
can be NULL
-- ALGLIB --
Copyright 20.03.2009 by Bochkanov Sergey
*************************************************************************/
void minasaoptimize(minasastate &state,
void (*grad)(const real_1d_array &x, double &func, real_1d_array &grad, void *ptr),
void (*rep)(const real_1d_array &x, double func, void *ptr) = NULL,
void *ptr = NULL);
/*************************************************************************
Obsolete optimization algorithm.
Was replaced by MinBLEIC subpackage.
-- ALGLIB --
Copyright 20.03.2009 by Bochkanov Sergey
*************************************************************************/
void minasaresults(const minasastate &state, real_1d_array &x, minasareport &rep);
/*************************************************************************
Obsolete optimization algorithm.
Was replaced by MinBLEIC subpackage.
-- ALGLIB --
Copyright 20.03.2009 by Bochkanov Sergey
*************************************************************************/
void minasaresultsbuf(const minasastate &state, real_1d_array &x, minasareport &rep);
/*************************************************************************
Obsolete optimization algorithm.
Was replaced by MinBLEIC subpackage.
-- ALGLIB --
Copyright 30.07.2010 by Bochkanov Sergey
*************************************************************************/
void minasarestartfrom(const minasastate &state, const real_1d_array &x, const real_1d_array &bndl, const real_1d_array &bndu);
}
/////////////////////////////////////////////////////////////////////////
//
// THIS SECTION CONTAINS COMPUTATIONAL CORE DECLARATIONS (FUNCTIONS)
//
/////////////////////////////////////////////////////////////////////////
namespace alglib_impl
{
void mincgcreate(ae_int_t n,
/* Real */ ae_vector* x,
mincgstate* state,
ae_state *_state);
void mincgcreatef(ae_int_t n,
/* Real */ ae_vector* x,
double diffstep,
mincgstate* state,
ae_state *_state);
void mincgsetcond(mincgstate* state,
double epsg,
double epsf,
double epsx,
ae_int_t maxits,
ae_state *_state);
void mincgsetscale(mincgstate* state,
/* Real */ ae_vector* s,
ae_state *_state);
void mincgsetxrep(mincgstate* state, ae_bool needxrep, ae_state *_state);
void mincgsetdrep(mincgstate* state, ae_bool needdrep, ae_state *_state);
void mincgsetcgtype(mincgstate* state, ae_int_t cgtype, ae_state *_state);
void mincgsetstpmax(mincgstate* state, double stpmax, ae_state *_state);
void mincgsuggeststep(mincgstate* state, double stp, ae_state *_state);
void mincgsetprecdefault(mincgstate* state, ae_state *_state);
void mincgsetprecdiag(mincgstate* state,
/* Real */ ae_vector* d,
ae_state *_state);
void mincgsetprecscale(mincgstate* state, ae_state *_state);
ae_bool mincgiteration(mincgstate* state, ae_state *_state);
void mincgresults(mincgstate* state,
/* Real */ ae_vector* x,
mincgreport* rep,
ae_state *_state);
void mincgresultsbuf(mincgstate* state,
/* Real */ ae_vector* x,
mincgreport* rep,
ae_state *_state);
void mincgrestartfrom(mincgstate* state,
/* Real */ ae_vector* x,
ae_state *_state);
void mincgsetprecdiagfast(mincgstate* state,
/* Real */ ae_vector* d,
ae_state *_state);
void mincgsetpreclowrankfast(mincgstate* state,
/* Real */ ae_vector* d1,
/* Real */ ae_vector* c,
/* Real */ ae_matrix* v,
ae_int_t vcnt,
ae_state *_state);
void mincgsetprecvarpart(mincgstate* state,
/* Real */ ae_vector* d2,
ae_state *_state);
ae_bool _mincgstate_init(mincgstate* p, ae_state *_state, ae_bool make_automatic);
ae_bool _mincgstate_init_copy(mincgstate* dst, mincgstate* src, ae_state *_state, ae_bool make_automatic);
void _mincgstate_clear(mincgstate* p);
ae_bool _mincgreport_init(mincgreport* p, ae_state *_state, ae_bool make_automatic);
ae_bool _mincgreport_init_copy(mincgreport* dst, mincgreport* src, ae_state *_state, ae_bool make_automatic);
void _mincgreport_clear(mincgreport* p);
void minbleiccreate(ae_int_t n,
/* Real */ ae_vector* x,
minbleicstate* state,
ae_state *_state);
void minbleiccreatef(ae_int_t n,
/* Real */ ae_vector* x,
double diffstep,
minbleicstate* state,
ae_state *_state);
void minbleicsetbc(minbleicstate* state,
/* Real */ ae_vector* bndl,
/* Real */ ae_vector* bndu,
ae_state *_state);
void minbleicsetlc(minbleicstate* state,
/* Real */ ae_matrix* c,
/* Integer */ ae_vector* ct,
ae_int_t k,
ae_state *_state);
void minbleicsetinnercond(minbleicstate* state,
double epsg,
double epsf,
double epsx,
ae_state *_state);
void minbleicsetoutercond(minbleicstate* state,
double epsx,
double epsi,
ae_state *_state);
void minbleicsetscale(minbleicstate* state,
/* Real */ ae_vector* s,
ae_state *_state);
void minbleicsetprecdefault(minbleicstate* state, ae_state *_state);
void minbleicsetprecdiag(minbleicstate* state,
/* Real */ ae_vector* d,
ae_state *_state);
void minbleicsetprecscale(minbleicstate* state, ae_state *_state);
void minbleicsetmaxits(minbleicstate* state,
ae_int_t maxits,
ae_state *_state);
void minbleicsetxrep(minbleicstate* state,
ae_bool needxrep,
ae_state *_state);
void minbleicsetstpmax(minbleicstate* state,
double stpmax,
ae_state *_state);
ae_bool minbleiciteration(minbleicstate* state, ae_state *_state);
void minbleicresults(minbleicstate* state,
/* Real */ ae_vector* x,
minbleicreport* rep,
ae_state *_state);
void minbleicresultsbuf(minbleicstate* state,
/* Real */ ae_vector* x,
minbleicreport* rep,
ae_state *_state);
void minbleicrestartfrom(minbleicstate* state,
/* Real */ ae_vector* x,
ae_state *_state);
ae_bool _minbleicstate_init(minbleicstate* p, ae_state *_state, ae_bool make_automatic);
ae_bool _minbleicstate_init_copy(minbleicstate* dst, minbleicstate* src, ae_state *_state, ae_bool make_automatic);
void _minbleicstate_clear(minbleicstate* p);
ae_bool _minbleicreport_init(minbleicreport* p, ae_state *_state, ae_bool make_automatic);
ae_bool _minbleicreport_init_copy(minbleicreport* dst, minbleicreport* src, ae_state *_state, ae_bool make_automatic);
void _minbleicreport_clear(minbleicreport* p);
void minlbfgscreate(ae_int_t n,
ae_int_t m,
/* Real */ ae_vector* x,
minlbfgsstate* state,
ae_state *_state);
void minlbfgscreatef(ae_int_t n,
ae_int_t m,
/* Real */ ae_vector* x,
double diffstep,
minlbfgsstate* state,
ae_state *_state);
void minlbfgssetcond(minlbfgsstate* state,
double epsg,
double epsf,
double epsx,
ae_int_t maxits,
ae_state *_state);
void minlbfgssetxrep(minlbfgsstate* state,
ae_bool needxrep,
ae_state *_state);
void minlbfgssetstpmax(minlbfgsstate* state,
double stpmax,
ae_state *_state);
void minlbfgssetscale(minlbfgsstate* state,
/* Real */ ae_vector* s,
ae_state *_state);
void minlbfgscreatex(ae_int_t n,
ae_int_t m,
/* Real */ ae_vector* x,
ae_int_t flags,
double diffstep,
minlbfgsstate* state,
ae_state *_state);
void minlbfgssetprecdefault(minlbfgsstate* state, ae_state *_state);
void minlbfgssetpreccholesky(minlbfgsstate* state,
/* Real */ ae_matrix* p,
ae_bool isupper,
ae_state *_state);
void minlbfgssetprecdiag(minlbfgsstate* state,
/* Real */ ae_vector* d,
ae_state *_state);
void minlbfgssetprecscale(minlbfgsstate* state, ae_state *_state);
ae_bool minlbfgsiteration(minlbfgsstate* state, ae_state *_state);
void minlbfgsresults(minlbfgsstate* state,
/* Real */ ae_vector* x,
minlbfgsreport* rep,
ae_state *_state);
void minlbfgsresultsbuf(minlbfgsstate* state,
/* Real */ ae_vector* x,
minlbfgsreport* rep,
ae_state *_state);
void minlbfgsrestartfrom(minlbfgsstate* state,
/* Real */ ae_vector* x,
ae_state *_state);
ae_bool _minlbfgsstate_init(minlbfgsstate* p, ae_state *_state, ae_bool make_automatic);
ae_bool _minlbfgsstate_init_copy(minlbfgsstate* dst, minlbfgsstate* src, ae_state *_state, ae_bool make_automatic);
void _minlbfgsstate_clear(minlbfgsstate* p);
ae_bool _minlbfgsreport_init(minlbfgsreport* p, ae_state *_state, ae_bool make_automatic);
ae_bool _minlbfgsreport_init_copy(minlbfgsreport* dst, minlbfgsreport* src, ae_state *_state, ae_bool make_automatic);
void _minlbfgsreport_clear(minlbfgsreport* p);
void minqpcreate(ae_int_t n, minqpstate* state, ae_state *_state);
void minqpsetlinearterm(minqpstate* state,
/* Real */ ae_vector* b,
ae_state *_state);
void minqpsetquadraticterm(minqpstate* state,
/* Real */ ae_matrix* a,
ae_bool isupper,
ae_state *_state);
void minqpsetstartingpoint(minqpstate* state,
/* Real */ ae_vector* x,
ae_state *_state);
void minqpsetorigin(minqpstate* state,
/* Real */ ae_vector* xorigin,
ae_state *_state);
void minqpsetalgocholesky(minqpstate* state, ae_state *_state);
void minqpsetbc(minqpstate* state,
/* Real */ ae_vector* bndl,
/* Real */ ae_vector* bndu,
ae_state *_state);
void minqpoptimize(minqpstate* state, ae_state *_state);
void minqpresults(minqpstate* state,
/* Real */ ae_vector* x,
minqpreport* rep,
ae_state *_state);
void minqpresultsbuf(minqpstate* state,
/* Real */ ae_vector* x,
minqpreport* rep,
ae_state *_state);
void minqpsetlineartermfast(minqpstate* state,
/* Real */ ae_vector* b,
ae_state *_state);
void minqpsetquadratictermfast(minqpstate* state,
/* Real */ ae_matrix* a,
ae_bool isupper,
double s,
ae_state *_state);
void minqprewritediagonal(minqpstate* state,
/* Real */ ae_vector* s,
ae_state *_state);
void minqpsetstartingpointfast(minqpstate* state,
/* Real */ ae_vector* x,
ae_state *_state);
void minqpsetoriginfast(minqpstate* state,
/* Real */ ae_vector* xorigin,
ae_state *_state);
ae_bool _minqpstate_init(minqpstate* p, ae_state *_state, ae_bool make_automatic);
ae_bool _minqpstate_init_copy(minqpstate* dst, minqpstate* src, ae_state *_state, ae_bool make_automatic);
void _minqpstate_clear(minqpstate* p);
ae_bool _minqpreport_init(minqpreport* p, ae_state *_state, ae_bool make_automatic);
ae_bool _minqpreport_init_copy(minqpreport* dst, minqpreport* src, ae_state *_state, ae_bool make_automatic);
void _minqpreport_clear(minqpreport* p);
void minlmcreatevj(ae_int_t n,
ae_int_t m,
/* Real */ ae_vector* x,
minlmstate* state,
ae_state *_state);
void minlmcreatev(ae_int_t n,
ae_int_t m,
/* Real */ ae_vector* x,
double diffstep,
minlmstate* state,
ae_state *_state);
void minlmcreatefgh(ae_int_t n,
/* Real */ ae_vector* x,
minlmstate* state,
ae_state *_state);
void minlmsetcond(minlmstate* state,
double epsg,
double epsf,
double epsx,
ae_int_t maxits,
ae_state *_state);
void minlmsetxrep(minlmstate* state, ae_bool needxrep, ae_state *_state);
void minlmsetstpmax(minlmstate* state, double stpmax, ae_state *_state);
void minlmsetscale(minlmstate* state,
/* Real */ ae_vector* s,
ae_state *_state);
void minlmsetbc(minlmstate* state,
/* Real */ ae_vector* bndl,
/* Real */ ae_vector* bndu,
ae_state *_state);
void minlmsetacctype(minlmstate* state,
ae_int_t acctype,
ae_state *_state);
ae_bool minlmiteration(minlmstate* state, ae_state *_state);
void minlmresults(minlmstate* state,
/* Real */ ae_vector* x,
minlmreport* rep,
ae_state *_state);
void minlmresultsbuf(minlmstate* state,
/* Real */ ae_vector* x,
minlmreport* rep,
ae_state *_state);
void minlmrestartfrom(minlmstate* state,
/* Real */ ae_vector* x,
ae_state *_state);
void minlmcreatevgj(ae_int_t n,
ae_int_t m,
/* Real */ ae_vector* x,
minlmstate* state,
ae_state *_state);
void minlmcreatefgj(ae_int_t n,
ae_int_t m,
/* Real */ ae_vector* x,
minlmstate* state,
ae_state *_state);
void minlmcreatefj(ae_int_t n,
ae_int_t m,
/* Real */ ae_vector* x,
minlmstate* state,
ae_state *_state);
ae_bool _minlmstate_init(minlmstate* p, ae_state *_state, ae_bool make_automatic);
ae_bool _minlmstate_init_copy(minlmstate* dst, minlmstate* src, ae_state *_state, ae_bool make_automatic);
void _minlmstate_clear(minlmstate* p);
ae_bool _minlmreport_init(minlmreport* p, ae_state *_state, ae_bool make_automatic);
ae_bool _minlmreport_init_copy(minlmreport* dst, minlmreport* src, ae_state *_state, ae_bool make_automatic);
void _minlmreport_clear(minlmreport* p);
void minlbfgssetdefaultpreconditioner(minlbfgsstate* state,
ae_state *_state);
void minlbfgssetcholeskypreconditioner(minlbfgsstate* state,
/* Real */ ae_matrix* p,
ae_bool isupper,
ae_state *_state);
void minbleicsetbarrierwidth(minbleicstate* state,
double mu,
ae_state *_state);
void minbleicsetbarrierdecay(minbleicstate* state,
double mudecay,
ae_state *_state);
void minasacreate(ae_int_t n,
/* Real */ ae_vector* x,
/* Real */ ae_vector* bndl,
/* Real */ ae_vector* bndu,
minasastate* state,
ae_state *_state);
void minasasetcond(minasastate* state,
double epsg,
double epsf,
double epsx,
ae_int_t maxits,
ae_state *_state);
void minasasetxrep(minasastate* state, ae_bool needxrep, ae_state *_state);
void minasasetalgorithm(minasastate* state,
ae_int_t algotype,
ae_state *_state);
void minasasetstpmax(minasastate* state, double stpmax, ae_state *_state);
ae_bool minasaiteration(minasastate* state, ae_state *_state);
void minasaresults(minasastate* state,
/* Real */ ae_vector* x,
minasareport* rep,
ae_state *_state);
void minasaresultsbuf(minasastate* state,
/* Real */ ae_vector* x,
minasareport* rep,
ae_state *_state);
void minasarestartfrom(minasastate* state,
/* Real */ ae_vector* x,
/* Real */ ae_vector* bndl,
/* Real */ ae_vector* bndu,
ae_state *_state);
ae_bool _minasastate_init(minasastate* p, ae_state *_state, ae_bool make_automatic);
ae_bool _minasastate_init_copy(minasastate* dst, minasastate* src, ae_state *_state, ae_bool make_automatic);
void _minasastate_clear(minasastate* p);
ae_bool _minasareport_init(minasareport* p, ae_state *_state, ae_bool make_automatic);
ae_bool _minasareport_init_copy(minasareport* dst, minasareport* src, ae_state *_state, ae_bool make_automatic);
void _minasareport_clear(minasareport* p);
}
#endif