ConicBundle
|
Iterative solver for the reduced symmetric, in general indefinite primal dual KKT System within QPSolverBasicStructures, where only H and A blocks remain, B and C are removed by Schur complements. If there is no A, PCG may used. More...
#include <QPIterativeKKTHASolver.hxx>
Public Member Functions | |
QPIterativeKKTHASolver (CH_Matrix_Classes::IterativeSolverObject *insolver, QPKKTPrecondObject *inprecond=0, CBout *cb=0, int cbinc=-1) | |
default constructor | |
virtual | ~QPIterativeKKTHASolver () |
virtual destructor | |
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 int | ItSys_mult (const CH_Matrix_Classes::Matrix &in_vec, CH_Matrix_Classes::Matrix &out_vec) |
returns out_vec=(system matrix)*in_vec | |
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::QPIterativeKKTSolver | |
virtual void | clear () |
reset data to empty but keep solver and preconditioner | |
QPIterativeKKTSolver (CH_Matrix_Classes::IterativeSolverObject *insolver, QPKKTPrecondObject *inprecond=0, CBout *cb=0, int cbinc=-1) | |
default constructor | |
virtual | ~QPIterativeKKTSolver () |
virtual destructor | |
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 const CH_Matrix_Classes::Matrix & | ItSys_rhs () |
returns the right hand side vector (dense) | |
virtual int | ItSys_precondM1 (CH_Matrix_Classes::Matrix &vec) |
returns M1^{-1}*vec; default: M1=I | |
virtual int | ItSys_precondM2 (CH_Matrix_Classes::Matrix &vec) |
returns M2^{-1}vec; default: M2=I | |
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) | |
Public Member Functions inherited from ConicBundle::QPKKTSolverObject | |
QPKKTSolverObject (CBout *cb=0, int cbinc=-1) | |
default constructor | |
virtual | ~QPKKTSolverObject () |
virtual destructor | |
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 | |
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 | |
Additional Inherited Members | |
Protected Attributes inherited from ConicBundle::QPIterativeKKTSolver | |
CH_Matrix_Classes::IterativeSolverObject * | solver |
the iterative solution method to be used, may not be NULL, deleted on destruction or replacemnt | |
QPKKTPrecondObject * | precond |
a preconditioning routine compatible with this system AND the solver, may be NULL (no preconditioning); consistency cannot be checked here and must be guaranteed externally; deleted on destruction or replacement | |
CH_Matrix_Classes::Real | Hfactor |
factor for H, usually ==1. unless a SOC model is used for the quadratic term | |
CH_Matrix_Classes::Matrix | KKTdiagx |
diagonal term to be added to the H-block due to bounds on the qudratic variables | |
CH_Matrix_Classes::Matrix | KKTdiagy |
diagonal term to be added to the A-block due to bounds on the constraints | |
CH_Matrix_Classes::Matrix | sysrhs |
the right hand side of the system | |
CH_Matrix_Classes::Matrix | sol |
solution vector of the entire system | |
CH_Matrix_Classes::Matrix | solmod |
temporary solution part of the model variables | |
CH_Matrix_Classes::Matrix | solcstr |
temporary solution part of the model constraint variables | |
CH_Matrix_Classes::Matrix | in_vecx |
temporary storage for x part of in_vec in ItSys_mult | |
CH_Matrix_Classes::Matrix | in_vecy |
temporary storage for y part of in_vec in ItSys_mult | |
CH_Matrix_Classes::Matrix | in_model |
temporary storage for model part of in_vec in ItSys_mult | |
CH_Matrix_Classes::Integer | maxit_bnd |
dynamic bound for maximum number of iterations (bnd relative to maximum encountered so far), reset in QPinitKKTdata | |
CH_Matrix_Classes::Integer | nmult |
holds number of calls to ItSys_mult, reset to 0 in QPinitKKTdata | |
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 | |
Iterative solver for the reduced symmetric, in general indefinite primal dual KKT System within QPSolverBasicStructures, where only H and A blocks remain, B and C are removed by Schur complements. If there is no A, PCG may used.
See the text to QPKKTSolverObject for the terminology of the primal dual KKT System and the general outline.
For this class the preconditioner offered by QPKKTSubspaceHPrecond may be worth to try.
The most important routines of the model described in the QPModelBlockObject that are required here are (besides sizes and multiplications with the bundle matrix B)