ConicBundle
Public Member Functions | Private Member Functions | Private Attributes | List of all members
ConicBundle::QPDirectKKTSolver Class Reference

implements a direct KKT Solver variant of QPKKTSolverObject More...

#include <QPDirectKKTSolver.hxx>

Inheritance diagram for ConicBundle::QPDirectKKTSolver:
ConicBundle::QPKKTSolverObject ConicBundle::CBout

Public Member Functions

virtual void clear ()
 reset data to empty
 
 QPDirectKKTSolver (bool in_factorize_ABC=false, CBout *cb=0, int cbinc=-1)
 default constructor
 
virtual ~QPDirectKKTSolver ()
 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 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
 
- 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
 
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 compute_DLR_Schur_complement (CH_Matrix_Classes::Symmatrix &sc)
 computes the Schur complement if Hp->is_DLR()==TRUE, i.e., the quadratic term is representable as diagonal plus low rank
 
int compute_dense_Schur_complement (CH_Matrix_Classes::Symmatrix &sc)
 computes the Schur compelement if the quadratic term is only retrievable in dense form
 

Private Attributes

CH_Matrix_Classes::Integer dim
 dimension of the quadratic term
 
CH_Matrix_Classes::Integer Anr
 number of rows in the constrain matrix =(A==0)?0:A->rowdim();
 
CH_Matrix_Classes::Integer bsz
 size of the bundle information in the cutting model
 
CH_Matrix_Classes::Integer csz
 number of constraints within the cutting model
 
const CH_Matrix_Classes::MatrixVp
 used for low rank precondititioning
 
CH_Matrix_Classes::Real Hfactor
 scalar factor for H term
 
CH_Matrix_Classes::Matrix Diag_inv
 inverse of the KKT diagonal in low rank versions
 
CH_Matrix_Classes::Symmatrix Qchol
 Cholesky facotr of KKT Q or of its low rank part.
 
CH_Matrix_Classes::Symmatrix Schur_complement
 Schur complement of ABC with Q if no bounds.
 
CH_Matrix_Classes::Symmatrix AQiAt_inv
 Cholesky or Aasen factor of ABC with Schur complement.
 
CH_Matrix_Classes::Symmatrix CABinvCt_inv
 for forming partial Schur complements
 
CH_Matrix_Classes::Matrix LinvC
 for partial Schur complement w.r.t. C block
 
CH_Matrix_Classes::Indexmatrix piv
 pivot sequence for Aasen
 
bool factorize_ABC
 if true, use Aasen to factor the ABC block
 

Additional Inherited Members

- Protected Attributes inherited from ConicBundle::QPKKTSolverObject
QPSolverProxObjectHp
 points to the quadratic cost representation, may NOT be NULL afer init
 
QPModelBlockObjectmodel
 points to the cutting model information, may be NULL
 
const CH_Matrix_Classes::SparsematA
 points to a possibly present constraint matrix, may be NULL
 
const CH_Matrix_Classes::Indexmatrixeq_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
 

Detailed Description

implements a direct KKT Solver variant of QPKKTSolverObject

See the text to QPKKTSolverObject for the terminology of the primal dual KKT System and the general outline.

For solving the system the current class first forms the Schur complement with respect to $H+D_x$; this is done in dependence of whether the prox term pointed to by Hp is of the form of Diagonal plus low rank (-> low rank inversion) or dense or other (dense Cholesky). The computed Schur complement matrix [A; B; 0]H^{-1}[A^T,B^T,0] is stored in Schur_complement before it is added to ABC. The code then continues in dependence of the value of factorize_ABC

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)

Member Function Documentation

◆ QPinit_KKTdata()

virtual int ConicBundle::QPDirectKKTSolver::QPinit_KKTdata ( QPSolverProxObject Hp,
QPModelBlockObject model,
const CH_Matrix_Classes::Sparsemat A,
const CH_Matrix_Classes::Indexmatrix eq_indices 
)
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

Parameters
Hpmay not be be NULL
modelmay be NULL
Amay be NULL
eq_indicesif not NULL these rows of A correspond to equations

Implements ConicBundle::QPKKTSolverObject.

Referenced by QPDirectKKTSolver().


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