ConicBundle
|
This is a pseudosolver designed for producing comparative statistics on the performance of mainly iterative solvers for the interior point KKT systems of QPSolver. More...
#include <QPKKTSolverComparison.hxx>
Public Member Functions | |
virtual void | clear () |
reset data to empty | |
QPKKTSolverComparison (CBout *cb=0, int cbinc=-1) | |
default constructor | |
virtual | ~QPKKTSolverComparison () |
virtual destructor | |
virtual int | add_solver (QPKKTSolverObject *solver, const char *name) |
the first solver added is the reference solver | |
virtual int | QPinit_KKTdata (QPSolverProxObject *Hp, QPModelBlockObject *model, const CH_Matrix_Classes::Sparsemat *A, const CH_Matrix_Classes::Indexmatrix *eq_indices) |
returns 1 if this class is not applicable in the current data situation, otherwise it stores the data pointers and these need to stay valid throught the use of the other routines but are not deleted here More... | |
virtual int | QPinit_KKTsystem (const CH_Matrix_Classes::Matrix &KKTdiagx, const CH_Matrix_Classes::Matrix &KKTdiagy, CH_Matrix_Classes::Real Hfactor, CH_Matrix_Classes::Real prec, QPSolverParameters *params) |
set up the primal dual KKT system for being solved for predictor and corrector rhs in QPsolve_KKTsystem | |
virtual int | QPsolve_KKTsystem (CH_Matrix_Classes::Matrix &solx, CH_Matrix_Classes::Matrix &soly, const CH_Matrix_Classes::Matrix &primalrhs, const CH_Matrix_Classes::Matrix &dualrhs, CH_Matrix_Classes::Real rhsmu, CH_Matrix_Classes::Real rhscorr, CH_Matrix_Classes::Real prec, QPSolverParameters *params) |
solve the KKTsystem to precision prec for the given right hand sides that have been computed for the value rhsmu of the barrier parameter and in which a rhscorr fraction (out of [0,1] of the corrector term have been included; in iterative solvers solx and soly may be used as starting points | |
virtual CH_Matrix_Classes::Real | QPget_blockH_norm () |
for judging violation this returns (an estimate of) the norm of the H-row in the latest system | |
virtual CH_Matrix_Classes::Real | QPget_blockA_norm () |
for judging violation this returns (an estimate of) the norm of the A-row in the latest system | |
std::vector< std::string > & | get_solvernames () |
returns the names of the solvers | |
std::vector< QPKKT_ProbStats > & | get_probdata () |
returns the stored statistical data | |
int | get_mu_stats (CH_Matrix_Classes::Real lbmu, CH_Matrix_Classes::Real ubmu, CH_Matrix_Classes::Indexmatrix &dims, CH_Matrix_Classes::Matrix &mu, CH_Matrix_Classes::Matrix &prepsecs, CH_Matrix_Classes::Matrix &predsecs, CH_Matrix_Classes::Matrix &corrsecs, CH_Matrix_Classes::Indexmatrix &predcalls, CH_Matrix_Classes::Indexmatrix &corrcalls, CH_Matrix_Classes::Matrix &cond, CH_Matrix_Classes::Indexmatrix &pccols, CH_Matrix_Classes::Matrix &sysviol) |
return those data columns (each a KKT system; columns are more efficient to append than lines) that fall into the given lower and upper bounds on mu More... | |
int | get_prob_stats (CH_Matrix_Classes::Indexmatrix &dims, CH_Matrix_Classes::Indexmatrix &iterations, CH_Matrix_Classes::Matrix &lastmu, CH_Matrix_Classes::Matrix &prepsecs, CH_Matrix_Classes::Matrix &predsecs, CH_Matrix_Classes::Matrix &corrsecs, CH_Matrix_Classes::Indexmatrix &predcalls, CH_Matrix_Classes::Indexmatrix &corrcalls) |
return one data column per subproblem (more efficient to append than lines) with the sum of the time/calls/etc. More... | |
Public Member Functions inherited from ConicBundle::QPKKTSolverObject | |
QPKKTSolverObject (CBout *cb=0, int cbinc=-1) | |
default constructor | |
virtual | ~QPKKTSolverObject () |
virtual destructor | |
virtual CH_Matrix_Classes::Integer | QPget_nmult () const |
for evaluation purposes with iterative solvers, return the number of matrix vector multiplications since the latest call to QPinit_KKTsystem (returns 0 for direct solvers) | |
virtual CH_Matrix_Classes::Real | QPget_condition_number () |
for evaluation purposes with iterative solvers, return an estimate of the condition number (returns -1. for direct solvers or if not supplied) | |
virtual CH_Matrix_Classes::Integer | QPget_precond_rank () |
for evaluation purposes with iterative solvers, return the rank of the precondiontioner used (or the number of n-vector multiplications per call) | |
virtual CH_Matrix_Classes::Integer | QPget_system_size () |
for evaluation purposes with iterative solvers, return the size of the system matrix | |
Public Member Functions inherited from ConicBundle::CBout | |
virtual void | set_out (std::ostream *out=0, int print_level=1) |
Specifies the output level (out==NULL: no output at all, out!=NULL and level=0: errors and warnings, level>0 increasingly detailed information) More... | |
virtual void | set_cbout (const CBout *cb, int incr=-1) |
Specifies the output level relative to the given CBout class. More... | |
void | clear_cbout () |
reset to default settings (out=0,print_level=1) | |
CBout (const CBout *cb=0, int incr=-1) | |
calls set_cbout | |
CBout (std::ostream *outp, int pl=1) | |
initialize correspondingly | |
CBout (const CBout &cb, int incr=0) | |
copy constructor | |
virtual bool | cb_out (int pl=-1) const |
Returns true if out!=0 and (pl<print_level), pl<0 should be used for WARNINGS and ERRORS only, pl==0 for usual output. | |
std::ostream & | get_out () const |
If cb_out() returned true, this returns the output stream, but it will abort if called with out==0. | |
std::ostream * | get_out_ptr () const |
returns the pointer to the output stream | |
int | get_print_level () const |
returns the print_level | |
virtual int | mfile_data (std::ostream &out) const |
writes problem data to the given outstream | |
Private Member Functions | |
int | violation (CH_Matrix_Classes::Matrix &violvec, const CH_Matrix_Classes::Matrix &solx, const CH_Matrix_Classes::Matrix &soly, QPModelBlockObject *inmodel, QPKKT_SolverStats &stats) |
computes and stores the norms of the residualy of the KKT system | |
Private Attributes | |
std::vector< QPKKTSolverObject * > | solver |
solvers add via add_solver() | |
std::vector< std::string > | solvername |
name of the solver for output | |
std::vector< QPKKT_ProbStats > | probdata |
for each bundle subproblem there is a block QPKKT_ProbStats (stored here), which holds for each KKT system a block of QPKKT_KKTStats, which holds for each solver a block QPKKT_SolverStats. | |
std::vector< QPModelBlockObject * > | model |
a seperate clone of the model for each solver in order to avoid side effects | |
QPSolverProxObject * | testHp |
the proximal term/quadratic cost matrix; this points to external object | |
const CH_Matrix_Classes::Sparsemat * | testA |
the constraint matrix; this points to external object | |
CH_Matrix_Classes::Matrix | testKKTdiagx |
stores the input diagonal for the H block | |
CH_Matrix_Classes::Matrix | testKKTdiagy |
stores the input diagonal for the A block | |
CH_Matrix_Classes::Real | testHfactor |
stores a factor for testHp, usually ==1. unless a second order cone variant is used | |
CH_Matrix_Classes::Matrix | testprimalrhs |
stores the input rhs for the A block | |
CH_Matrix_Classes::Matrix | testdualrhs |
stores the input rhs for the H block | |
CH_Tools::Clock | clock |
for taking the time | |
CH_Matrix_Classes::Matrix | insolx |
stores the input solution of the H block (e.g. of predictor) | |
CH_Matrix_Classes::Matrix | insoly |
stores the input solution of the A block (e.g. of predictor) | |
CH_Matrix_Classes::Matrix | outsolx |
stores the output solution of the H block (first solver) | |
CH_Matrix_Classes::Matrix | outsoly |
stores the output solution of the A block (first solver) | |
Friends | |
std::ostream & | operator<< (std::ostream &out, const QPKKTSolverComparison &q) |
(file-) output | |
std::istream & | operator>> (std::istream &in, QPKKTSolverComparison &q) |
(file-) input | |
Additional Inherited Members | |
Protected Attributes inherited from ConicBundle::QPKKTSolverObject | |
QPSolverProxObject * | Hp |
points to the quadratic cost representation, may NOT be NULL afer init | |
QPModelBlockObject * | model |
points to the cutting model information, may be NULL | |
const CH_Matrix_Classes::Sparsemat * | A |
points to a possibly present constraint matrix, may be NULL | |
const CH_Matrix_Classes::Indexmatrix * | eq_indices |
if not NULL, these rows of A correspond to equations; needed for checking applicability of this Object | |
CH_Matrix_Classes::Real | blockH_norm |
for judging violation: (an estimate of) the norm of the H-row in the system | |
CH_Matrix_Classes::Real | blockA_norm |
for judging violation: (an estimate of) the norm of the A-row in the system | |
This is a pseudosolver designed for producing comparative statistics on the performance of mainly iterative solvers for the interior point KKT systems of QPSolver.
Use add_solver() to enter the solvers that are to be compared. The first solver entered is used for actually returning the computations as if this first solver would be the true solver. The others work on exactly the same KKT system and model data via cloning.
The statistics are stored as follows: for each bundle subproblem there is a block QPKKT_ProbStats, which holds for each KKT system a block of QPKKT_KKTStats, which holds for each solver a block QPKKT_SolverStats.
int ConicBundle::QPKKTSolverComparison::get_mu_stats | ( | CH_Matrix_Classes::Real | lbmu, |
CH_Matrix_Classes::Real | ubmu, | ||
CH_Matrix_Classes::Indexmatrix & | dims, | ||
CH_Matrix_Classes::Matrix & | mu, | ||
CH_Matrix_Classes::Matrix & | prepsecs, | ||
CH_Matrix_Classes::Matrix & | predsecs, | ||
CH_Matrix_Classes::Matrix & | corrsecs, | ||
CH_Matrix_Classes::Indexmatrix & | predcalls, | ||
CH_Matrix_Classes::Indexmatrix & | corrcalls, | ||
CH_Matrix_Classes::Matrix & | cond, | ||
CH_Matrix_Classes::Indexmatrix & | pccols, | ||
CH_Matrix_Classes::Matrix & | sysviol | ||
) |
return those data columns (each a KKT system; columns are more efficient to append than lines) that fall into the given lower and upper bounds on mu
lbmu | get KKT systems with last barrier paremeter at least this lower bound |
ubmu | get KKT systems with last barrier paremeter at most this upper bound |
dims | order of H, lowrank of H, # rows of A, # equalities in A, # rows of B, # rows of C |
mu | value of the barrier parameter |
prepsecs | time spent in setting up the subproblem and systems (one row per solver) |
predsecs | time spent in solving the predictors (one row per solver) |
corrsecs | time spent in solving the correctors (one row per solver) |
predcalls | number of matrix vector multiplications in the predictors (one row per solver) |
corrcalls | number of matrix vector multiplications in the correctors (one row per solver) |
cond | condition number of the system (one row per solver) |
pccols | number of columns in the low rank preconditioners |
sysviol | norm of the (corrector) residual of the entire KKT system |
int ConicBundle::QPKKTSolverComparison::get_prob_stats | ( | CH_Matrix_Classes::Indexmatrix & | dims, |
CH_Matrix_Classes::Indexmatrix & | iterations, | ||
CH_Matrix_Classes::Matrix & | lastmu, | ||
CH_Matrix_Classes::Matrix & | prepsecs, | ||
CH_Matrix_Classes::Matrix & | predsecs, | ||
CH_Matrix_Classes::Matrix & | corrsecs, | ||
CH_Matrix_Classes::Indexmatrix & | predcalls, | ||
CH_Matrix_Classes::Indexmatrix & | corrcalls | ||
) |
return one data column per subproblem (more efficient to append than lines) with the sum of the time/calls/etc.
dims | per column/subproblem: order of H, lowrank of H, # rows of A, # equalities in A, # rows of B, # rows of C |
iterations | number of KKTsystem of each subproblem |
lastmu | value of the barrier parameter in the last KKT system of the subproblem |
prepsecs | cummulative time spent in setting up the subproblem and systems (one row per solver) |
predsecs | cummulative time spent in solving the predictors (one row per solver) |
corrsecs | cummulative time spent in solving the correctors (one row per solver) |
predcalls | cummulative number of matrix vector multiplications in the predictors (one row per solver) |
corrcalls | cummulative number of matrix vector multiplications in the correctors (one row per solver) |
|
virtual |
returns 1 if this class is not applicable in the current data situation, otherwise it stores the data pointers and these need to stay valid throught the use of the other routines but are not deleted here
Hp | may not be be NULL |
model | may be NULL |
A | may be NULL |
eq_indices | if not NULL these rows of A correspond to equations |
Implements ConicBundle::QPKKTSolverObject.