ConicBundle
QPKKTSubspaceHPrecond.hxx
Go to the documentation of this file.
1 
2 
3 #ifndef CONICBUNDLE_QPKKTSUBSPACEHPRECOND_HXX
4 #define CONICBUNDLE_QPKKTSUBSPACEHPRECOND_HXX
5 
14 #include "QPKKTPrecondObject.hxx"
15 
16 namespace ConicBundle {
17 
18 class QPSolverParameters;
19 
20 
24 
72 {
73 protected:
74 
75  //--- data describing the KKT system
76  // QPSolverProxObject* Hp; ///< points to the quadratic cost representation, may NOT be NULL afer init
77  // QPModelBlockObject* model; ///< points to the cutting model information, may be NULL
78  // const CH_Matrix_Classes::Sparsemat* A; ///< points to a possibly present constraint matrix, may be NULL
79  // const CH_Matrix_Classes::Indexmatrix* eq_indices; ///< if not NULL, these rows of A correspond to equations; needed for checking applicability of this Object
80 
82 
85 
88 
91 
97 
104 
106 
112 
113 
114 public:
115  // reset data to empty
116  //virtual void clear()
117  //{Hp=0;model=0;A=0;eq_indices=0;}
118 
120  QPKKTSubspaceHPrecond(CH_Matrix_Classes::Integer inmethod=0,CBout* cb=0,int cbinc=-1);
121 
123  virtual ~QPKKTSubspaceHPrecond();
124 
126  virtual int init_data(QPSolverProxObject* Hp,
130  bool SchurComplAineq
131  );
132 
134  virtual int init_system(const CH_Matrix_Classes::Matrix& KKTdiagx,
135  const CH_Matrix_Classes::Matrix& KKTdiagy,
138  QPSolverParameters* params);
139 
141  virtual int precondM1(CH_Matrix_Classes::Matrix& vec);
142 
143  //returns M2^{-1}vec; default: M2=I
144  //virtual int precondM2(CH_Matrix_Classes::Matrix& /* vec */) {return 0;}
145 
147  virtual int set_subspace(const CH_Matrix_Classes::Matrix& insubspace)
148  { if (method/10==3) subspace=insubspace; return 0;}
149 
152 
154  virtual int precond_invG1(CH_Matrix_Classes::Matrix& vec);
155 
157  virtual int precond_invG1tran(CH_Matrix_Classes::Matrix& vec);
158 
161  {return diagH.rowdim();}
162 
165  const CH_Matrix_Classes::Matrix& KKTdiagx,
166  const CH_Matrix_Classes::Matrix& KKTdiagy);
169  {return eigvals.rowdim();}
170 
173  {return t_precond_mult;}
174 
176  virtual void reset_t_precond_mult()
177  { t_precond_mult=0;}
178 };
179 
180 
181 
182 
184 
185 }
186 
187 #endif
188 
int Integer
all integer numbers in calculations and indexing are of this type
Definition: matop.hxx:40
CH_Matrix_Classes::Indexmatrix pivlowrank
the pivot sequence used
Definition: QPKKTSubspaceHPrecond.hxx:94
parameters for steering the termination criteria and solution method of the solver ...
Definition: QPSolverParameters.hxx:26
virtual int init_data(QPSolverProxObject *Hp, QPModelBlockObject *model, const CH_Matrix_Classes::Sparsemat *A, const CH_Matrix_Classes::Indexmatrix *eq_indices, bool SchurComplAineq)
returns 1 if this class is not applicable in the current data situation, otherwise it stores the data...
CH_Matrix_Classes::Symmatrix tmpsym
temporary symmetric matrix
Definition: QPKKTSubspaceHPrecond.hxx:100
virtual int precond_invG1tran(CH_Matrix_Classes::Matrix &vec)
for estimating the condition number with M1=G*G^T this returns G^{-T}*vec; default: G=I ...
Integer rowdim() const
returns the row dimension
Definition: matrix.hxx:215
virtual CH_Matrix_Classes::Integer get_precond_rank()
for evaluation purposes with iterative solvers, return the rank of the precondiontioner used (or the ...
Definition: QPKKTSubspaceHPrecond.hxx:168
double Real
all real numbers in calculations are of this type
Definition: matop.hxx:50
CH_Matrix_Classes::Matrix eigvals
the eigenvalues of the low rank approximation on the subspace
Definition: QPKKTSubspaceHPrecond.hxx:95
Matrix class for integral values of type Integer
Definition: indexmat.hxx:195
virtual CH_Matrix_Classes::Real get_lmin_invM1()
return (an estimate of) the minimum eigenvalue of the preconditioner M1^{-1}; this is used...
CH_Matrix_Classes::Matrix Q
[only used in testing] QR factorization of D^(-.5)*lowrank*eigvecs
Definition: QPKKTSubspaceHPrecond.hxx:105
CH_Tools::Microseconds t_comp_lowrank
time spent in generating the lowrank matrix by multiplying with subspace
Definition: QPKKTSubspaceHPrecond.hxx:109
CH_Matrix_Classes::Matrix subspace
the subspace used for projection
Definition: QPKKTSubspaceHPrecond.hxx:92
CH_Tools::Clock clock
for taking the time spent in various parts
Definition: QPKKTSubspaceHPrecond.hxx:107
QPSolverProxObject * Hp
points to the quadratic cost representation, may NOT be NULL afer init
Definition: QPKKTPrecondObject.hxx:57
CH_Tools::Microseconds t_comp_svd
time spent in computing the singular value decomposition
Definition: QPKKTSubspaceHPrecond.hxx:110
CH_Matrix_Classes::Matrix diagH
diagonal of the prox term
Definition: QPKKTSubspaceHPrecond.hxx:83
Subspace projection preconditioner for the H-block of the KKT-System assuming that B and C have been ...
Definition: QPKKTSubspaceHPrecond.hxx:71
CH_Tools::Microseconds t_gen_subspace
time spent in generating the subspace
Definition: QPKKTSubspaceHPrecond.hxx:108
CH_Matrix_Classes::Real Hfactor
the prox term of Hp is multiplied by this
Definition: QPKKTPrecondObject.hxx:63
Matrix class of symmetric matrices with real values of type Real
Definition: symmat.hxx:43
CH_Matrix_Classes::Real max_sigma
if >0 it gives the last maximum singular value found
Definition: QPKKTSubspaceHPrecond.hxx:87
allows measuring time difference to its initialization time in Microseconds
Definition: clock.hxx:282
CH_Matrix_Classes::Matrix rotmat
rotation to update the subspace
Definition: QPKKTSubspaceHPrecond.hxx:103
CH_Tools::Microseconds t_precond_mult
time spent in multiplying with the preconditioner
Definition: QPKKTSubspaceHPrecond.hxx:111
conic bundle method solver for sum of convex functions. See the ConicBundle_Manual for a quick introd...
Definition: CBSolver.hxx:22
QPModelBlockObject * model
points to the cutting model information, may be NULL
Definition: QPKKTPrecondObject.hxx:58
base class for uniform use of WARNINGS and ERRORS (at some point in time)
Definition: CBout.hxx:30
extra long integer number for expressing and computing time measurements in microseconds.
Definition: clock.hxx:46
const CH_Matrix_Classes::Sparsemat * A
points to a possibly present constraint matrix, may be NULL
Definition: QPKKTPrecondObject.hxx:59
virtual void reset_t_precond_mult()
for evaluation purposes with iterative solvers, reset the time spent in the multiplication with the p...
Definition: QPKKTSubspaceHPrecond.hxx:176
CH_Matrix_Classes::Matrix keepvecs
eignevectors used to update the subspace
Definition: QPKKTSubspaceHPrecond.hxx:102
CH_Matrix_Classes::Matrix keepeigs
eigenvalues used to update the subspace
Definition: QPKKTSubspaceHPrecond.hxx:101
virtual int precondM1(CH_Matrix_Classes::Matrix &vec)
returns M1^{-1}*vec; default: M1=I
const CH_Matrix_Classes::Matrix * Vp
lowrank part of the prox term
Definition: QPKKTSubspaceHPrecond.hxx:84
CH_Matrix_Classes::Matrix lowrank
the selected lowrank representation
Definition: QPKKTSubspaceHPrecond.hxx:93
CH_Matrix_Classes::Integer last_nmult
last number of multplications in iterative solver
Definition: QPKKTSubspaceHPrecond.hxx:86
Header declaring the abstract class ConicBundle::QPKKTPrecondObject.
virtual int set_subspace(const CH_Matrix_Classes::Matrix &insubspace)
if the method admits this, let the subspace be chosen externally
Definition: QPKKTSubspaceHPrecond.hxx:147
Matrix class for real values of type Real
Definition: matrix.hxx:74
virtual ~QPKKTSubspaceHPrecond()
virtual destructor
virtual CH_Tools::Microseconds get_t_precond_mult()
for evaluation purposes with iterative solvers, return the time spent in the multiplication with ...
Definition: QPKKTSubspaceHPrecond.hxx:172
Matrix class of sparse matrices with real values of type Real
Definition: sparsmat.hxx:74
QPKKTSubspaceHPrecond(CH_Matrix_Classes::Integer inmethod=0, CBout *cb=0, int cbinc=-1)
default constructor
abstract interface for model blocks in the constrained QPSolver
Definition: QPModelBlockObject.hxx:89
Abstract Interface for preconditioners to be used with a QPIterativeKKTSolver and a CH_Matrix_Classes...
Definition: QPKKTPrecondObject.hxx:52
CH_Matrix_Classes::Matrix tmpvec
temporary matrix
Definition: QPKKTSubspaceHPrecond.hxx:99
virtual int cond_number_mult(CH_Matrix_Classes::Matrix &vec, const CH_Matrix_Classes::Matrix &KKTdiagx, const CH_Matrix_Classes::Matrix &KKTdiagy)
for estimating the condition number directly for the preconditioned part only
bool SchurComplAineq
if true, the inequalities of A are Schur complemented into the H block
Definition: QPKKTPrecondObject.hxx:61
CH_Matrix_Classes::Matrix tmpmat
temporary matrix
Definition: QPKKTSubspaceHPrecond.hxx:98
virtual CH_Matrix_Classes::Integer precond_size()
for estimating the condition number directly for the preconditioned part only; negative numbers indic...
Definition: QPKKTSubspaceHPrecond.hxx:160
CH_Matrix_Classes::Matrix Diag_inv
diagonal with KKT diagonal part, inverted, possibly the sqrt of it (if there is a low rank part) ...
Definition: QPKKTSubspaceHPrecond.hxx:89
CH_Matrix_Classes::Real diaginvval
if the value is > 0 then Diag_inv must be this values times the all ones vector
Definition: QPKKTSubspaceHPrecond.hxx:90
CH_Matrix_Classes::Integer method
selects generation method for preconditioner
Definition: QPKKTSubspaceHPrecond.hxx:81
const CH_Matrix_Classes::Indexmatrix * eq_indices
if not NULL, these rows of A correspond to equations; needed for checking applicability of this Objec...
Definition: QPKKTPrecondObject.hxx:60
virtual int init_system(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; the input objects...
virtual int precond_invG1(CH_Matrix_Classes::Matrix &vec)
for estimating the condition number with M1=G*G^T this returns G^{-1}*vec; default: G=I ...
CH_Matrix_Classes::Matrix eigvecs
the eigenvectors within the subspace
Definition: QPKKTSubspaceHPrecond.hxx:96
in order to pass a ConicBundle::BundleProxObject, see Quadratic Proximal Terms, to a custzomized QPSo...
Definition: QPSolverObject.hxx:55