ConicBundle
CMlowranksd.hxx
Go to the documentation of this file.
1 
2 
3 #ifndef CONICBUNDLE_CMLOWRANKSD_HXX
4 #define CONICBUNDLE_CMLOWRANKSD_HXX
5 
14 #include "Coeffmat.hxx"
15 #include "CMlowrankdd.hxx"
16 
17 namespace ConicBundle {
18 
22 
26  class CMlowranksd: public Coeffmat
27  {
28  private:
31  public:
34  {A=Ain;B=Bin;CM_type=CM_lowranksd;infop=cip;}
36  virtual ~CMlowranksd(){}
37 
38 
40  virtual Coeffmat* clone() const
41  {return new CMlowranksd(A,B,ConicBundle::clone(infop));}
42 
44  virtual CH_Matrix_Classes::Integer dim() const { return A.rowdim(); }
45 
48  { return CH_Matrix_Classes::ip(A.row(i),B.row(j))+ CH_Matrix_Classes::ip(B.row(i),A.row(j));}
49 
52  { CH_Matrix_Classes::rank2add(A,B,S,2.);}
53 
55  virtual CH_Matrix_Classes::Real norm(void) const
58  return sqrt(2.*CH_Matrix_Classes::ip(C,D)+d);}
59 
61  virtual Coeffmat* subspace(const CH_Matrix_Classes::Matrix& P) const
63  return new CMlowrankdd(C,D,ConicBundle::clone(infop)); }
64 
67  { A*=d; if (infop) infop->multiply(d);}
68 
72 
76  return 2.*CH_Matrix_Classes::ip(C,D); }
77 
80  {
83  const CH_Matrix_Classes::Integer brd=C.rowdim();
84  const CH_Matrix_Classes::Integer pcd=P.coldim();
85  const CH_Matrix_Classes::Integer prd=P.rowdim();
86  const CH_Matrix_Classes::Real *pp=P.get_store()+start_row;
88  const CH_Matrix_Classes::Real *vp=(A.get_colval()).get_store();
89  const CH_Matrix_Classes::Integer *ip=(A.get_colindex()).get_store();
90  for(CH_Matrix_Classes::Integer j=0;j<colinfo.rowdim();j++){
91  CH_Matrix_Classes::Integer indj=colinfo(j,0);
92  for(CH_Matrix_Classes::Integer i=colinfo(j,1);--i>=0;){
93  CH_Matrix_Classes::mat_xpeya(pcd,bp+indj,brd,pp+(*ip++),prd,(*vp++));
94  }
95  }
98  {for(CH_Matrix_Classes::Integer j=0;j<D.coldim();j++){
99  const CH_Matrix_Classes::Real *pp=P.get_store()+start_row+j*P.rowdim();
100  const CH_Matrix_Classes::Real *ap=B.get_store();
102  for(CH_Matrix_Classes::Integer i=D.rowdim();--i>=0;){
103  (*dp++)=CH_Matrix_Classes::mat_ip(ard,ap,pp);
104  ap+=ard;
105  }
106  }}
107  chk_set_init(D,1);
108  CH_Matrix_Classes::Real trval=0;
109  if (Lam==0) {
110  trval=CH_Matrix_Classes::ip(C,D);
111  }
112  else {
113  assert(Lam->dim()==P.coldim());
114  const CH_Matrix_Classes::Real *cp=C.get_store();
115  const CH_Matrix_Classes::Real *dp=D.get_store();
116  const CH_Matrix_Classes::Real *lp=Lam->get_store();
117  for (CH_Matrix_Classes::Integer i=0;i<C.coldim();i++,cp+=C.rowdim(),dp+=D.rowdim())
118  trval+=(*lp++)*CH_Matrix_Classes::mat_ip(C.rowdim(),cp,dp);
119  }
120  return 2.*trval;
121  }
122 
125  { CH_Matrix_Classes::rank2add(A,B,S,2.*d,1.); }
126 
130  CH_Matrix_Classes::genmult(B,CH_Matrix_Classes::genmult(A,C,E,1.,0.,1),D,d,1.);}
131 
135  CH_Matrix_Classes::genmult(B,CH_Matrix_Classes::genmult(A,C,E,1.,0.,1),D,d,1.);}
136 
139  {
140  CH_Matrix_Classes::Matrix tmp1; CH_Matrix_Classes::genmult(P,A,tmp1,1.,0.,1,0);
141  CH_Matrix_Classes::Matrix tmp2; CH_Matrix_Classes::genmult(B,Q,tmp2,1.,0.,1,0);
142  CH_Matrix_Classes::genmult(tmp1,tmp2,R,1.,0.,0,0);
143  CH_Matrix_Classes::genmult(P,B,tmp1,1.,0.,1,0);
144  CH_Matrix_Classes::genmult(A,Q,tmp2,1.,0.,1,0);
145  CH_Matrix_Classes::genmult(tmp1,tmp2,R,1.,1.,0,0);
146  }
147 
150  { return 4*A.nonzeros()+4*B.rowdim()*B.coldim(); }
151 
153  virtual int dense() const
154  {return 0;}
155 
157  virtual int sparse() const
158  { return 0;}
159 
163  CH_Matrix_Classes::Matrix& /* val */,
164  CH_Matrix_Classes::Real /* d=1. */)const
165  {return 0;}
166 
168  virtual int support_in(const CH_Matrix_Classes::Sparsesym& /* S */) const
169  {return 0;}
170 
173  {return 2.*CH_Matrix_Classes::ip(S*A,B);}
174 
178  CH_Matrix_Classes::rank2add(C,D,S,2.);}
179 
182  {
185  const CH_Matrix_Classes::Integer brd=C.rowdim();
186  const CH_Matrix_Classes::Integer pcd=P.coldim();
187  const CH_Matrix_Classes::Integer prd=P.rowdim();
188  const CH_Matrix_Classes::Real *pp=P.get_store()+start_row;
189  const CH_Matrix_Classes::Indexmatrix &colinfo=A.get_colinfo();
190  const CH_Matrix_Classes::Real *vp=(A.get_colval()).get_store();
191  const CH_Matrix_Classes::Integer *ip=(A.get_colindex()).get_store();
192  for(CH_Matrix_Classes::Integer j=0;j<colinfo.rowdim();j++){
193  CH_Matrix_Classes::Integer indj=colinfo(j,0);
194  for(CH_Matrix_Classes::Integer i=colinfo(j,1);--i>=0;){
195  CH_Matrix_Classes::mat_xpeya(pcd,bp+indj,brd,pp+(*ip++),prd,(*vp++));
196  }
197  }
200  {for(CH_Matrix_Classes::Integer j=0;j<D.coldim();j++){
201  const CH_Matrix_Classes::Real *pp=P.get_store()+start_row+j*P.rowdim();
202  const CH_Matrix_Classes::Real *ap=B.get_store();
204  for(CH_Matrix_Classes::Integer i=D.rowdim();--i>=0;){
205  (*dp++)=CH_Matrix_Classes::mat_ip(ard,ap,pp);
206  ap+=ard;
207  }
208  }}
209  chk_set_init(D,1);
210  CH_Matrix_Classes::rank2add(C,D,S,2.*alpha,1.,1);
211  }
212 
215  CH_Matrix_Classes::Real alpha=1.,CH_Matrix_Classes::Real beta=0.,int dtrans=0) const
216  {
218  CH_Matrix_Classes::genmult(A,CH_Matrix_Classes::genmult(B,D,E,1.,0.,1,dtrans),C,alpha,beta);
219  return CH_Matrix_Classes::genmult(B,CH_Matrix_Classes::genmult(A,D,E,1.,0.,1,dtrans),C,alpha,1.);
220  }
221 
224  CH_Matrix_Classes::Real alpha=1.,CH_Matrix_Classes::Real beta=0.,int dtrans=0) const
225  {
227  CH_Matrix_Classes::genmult(CH_Matrix_Classes::genmult(D,A,E,1.,0.,dtrans),B,C,alpha,beta,0,1);
228  return CH_Matrix_Classes::genmult(CH_Matrix_Classes::genmult(D,B,E,1.,0.,dtrans),A,C,alpha,1.,0,1);
229  }
230 
232  virtual int equal(const Coeffmat* p,double tol=1e-6) const
233  {
234  const CMlowranksd *pp=dynamic_cast<const CMlowranksd *>(p);
235  if (pp==0)
236  return 0;
237  if (CH_Matrix_Classes::equal(A,pp->A,tol) &&
238  (B.rowdim()==(pp->B).rowdim())&&(B.coldim()==(pp->B).coldim()) &&
239  (CH_Matrix_Classes::norm2(B-pp->B)<tol))
240  return 1;
241  return 0;
242  }
243 
245  virtual std::ostream& display(std::ostream& o) const
246  {o<<"CMlowranksd\n";A.display(o);B.display(o);return o;}
247 
249  virtual std::ostream& out(std::ostream& o) const
250  {return o<<"LOWRANK_SPARSE_DENSE\n"<<A<<B;}
251 
253  virtual std::istream& in(std::istream& i)
254  {return i>>A>>B;
255  if((A.rowdim()!=B.rowdim())||(A.coldim()!=B.coldim())){
256  i.clear(std::ios::failbit);
257  if (CH_Matrix_Classes::materrout) (*CH_Matrix_Classes::materrout)<<"*** ERROR: CMlowranksd::in(): dimensions of A and B do not match"<<std::endl;
258  }
259  return i;
260  }
261 
263  CMlowranksd(std::istream& is,CoeffmatInfo* cip=0)
264  {CM_type=CM_lowranksd;infop=cip;in(is);}
265 
266  };
267 
269 
270 }
271 #endif
272 
#define chk_set_init(x, y)
CONICBUNDLE_DEBUG being undefined, the template function is removed. Otherwise it would allow to set ...
Definition: matop.hxx:1774
int Integer
all integer numbers in calculations and indexing are of this type
Definition: matop.hxx:40
Header declaring the classes ConicBundle::Coeffmat, ConicBundle::CoeffmatPointer, ConicBundle::Coeffm...
implements a low rank matrix as Coeffmat with a sparse rectangular CH_Matrix_Classes::Sparsemat and...
Definition: CMlowranksd.hxx:26
Integer rowdim() const
returns the row dimension
Definition: matrix.hxx:215
std::ostream * materrout
if not zero, this is the output stream for runtime error messages, by default it is set to &std::cout...
double Real
all real numbers in calculations are of this type
Definition: matop.hxx:50
allows to memorize the scalings applied to a Coeffmat and offers the basis for storing further user d...
Definition: Coeffmat.hxx:52
const Indexmatrix & get_colinfo() const
returns information on nonzero columns, k by 3, listing: index, %nonzeros, first index in colindex/co...
Definition: sparsmat.hxx:265
defines a base class for coefficient matrices in semidefinite programming, in particular for use with...
Definition: Coeffmat.hxx:125
void display(std::ostream &out, int precision=0, int width=0, int screenwidth=0) const
displays a matrix in a pretty way for bounded screen widths; for variables of value zero default valu...
Matrix class for integral values of type Integer
Definition: indexmat.hxx:195
Real norm2(const Matrix &A)
returns the Frobenius norm of A, i.e., the square root of the sum of A(i,j)*A(i,j) over all i...
Definition: matrix.hxx:1235
virtual Coeffmat * clone() const
makes an explicit copy of itself and returns a pointer to it
Definition: CMlowranksd.hxx:40
virtual void multiply(CH_Matrix_Classes::Real d)
multiply constraint permanentely by d; this is to allow scaling or sign changes in the constraints ...
Definition: CMlowranksd.hxx:66
virtual void addmeto(CH_Matrix_Classes::Symmatrix &S, CH_Matrix_Classes::Real d=1.) const
computes S+=d*(*this);
Definition: CMlowranksd.hxx:124
virtual CH_Matrix_Classes::Real norm(void) const
returns the Frobenius norm of the matrix
Definition: CMlowranksd.hxx:55
virtual int sparse() const
returns 0 if not sparse, otherwise 1
Definition: CMlowranksd.hxx:157
const Matrix & get_colval() const
returns the value vector of the column representation holding the value for each element ...
Definition: sparsmat.hxx:269
virtual void left_right_prod(const CH_Matrix_Classes::Matrix &P, const CH_Matrix_Classes::Matrix &Q, CH_Matrix_Classes::Matrix &R) const
computes R=P^T*(*this)*Q
Definition: CMlowranksd.hxx:138
Real * get_store()
returns the current address of the internal value array; use cautiously, do not use delete! ...
Definition: matrix.hxx:326
Matrix class of symmetric matrices with real values of type Real
Definition: symmat.hxx:43
CH_Matrix_Classes::Matrix B
this is B in A*B^T+B*A^T
Definition: CMlowranksd.hxx:30
virtual CH_Matrix_Classes::Real ip(const CH_Matrix_Classes::Sparsesym &S) const
returns the inner product of the constraint matrix with S
Definition: CMlowranksd.hxx:172
void display(std::ostream &out, int precision=0, int width=0, int screenwidth=0) const
displays a matrix in a pretty way for bounded screen widths; for variables of value zero default valu...
conic bundle method solver for sum of convex functions. See the ConicBundle_Manual for a quick introd...
Definition: CBSolver.hxx:22
virtual std::istream & in(std::istream &i)
counterpart to out(), does not read the class type, though. This is assumed to have been read in orde...
Definition: CMlowranksd.hxx:253
virtual void make_symmatrix(CH_Matrix_Classes::Symmatrix &S) const
returns a dense symmetric constraint matrix
Definition: CMlowranksd.hxx:51
Matrix row(Integer i) const
returns row i copied to a new matrix
Coeffmattype CM_type
in order to enable type identification
Definition: Coeffmat.hxx:134
virtual void project(CH_Matrix_Classes::Symmatrix &S, const CH_Matrix_Classes::Matrix &P) const
computes S=P^T*(*this)*P
Definition: CMlowranksd.hxx:176
CoeffmatInfo * infop
allows the user to specify and output additional information
Definition: Coeffmat.hxx:135
virtual std::ostream & out(std::ostream &o) const
put entire contents onto ostream with the class type in the beginning so that the derived class can b...
Definition: CMlowranksd.hxx:249
virtual void add_projection(CH_Matrix_Classes::Symmatrix &S, const CH_Matrix_Classes::Matrix &P, CH_Matrix_Classes::Real alpha=1., CH_Matrix_Classes::Integer start_row=0) const
computes S+=Q^T(*this)Q for Q=P.rows(start_row,start_row+dim-1)
Definition: CMlowranksd.hxx:181
Integer trace(const Indexmatrix &A)
returns the sum of the diagonal elements A(i,i) over all i
virtual int dense() const
returns 1 if its structure is as bad as its dense symmetric representation, otherwise 0 ...
Definition: CMlowranksd.hxx:153
CMlowranksd(std::istream &is, CoeffmatInfo *cip=0)
constructor with istream and possibly additional user information
Definition: CMlowranksd.hxx:263
virtual CH_Matrix_Classes::Real ip(const CH_Matrix_Classes::Symmatrix &S) const
returns ip(*this,S)=trace(*this*S), the trace inner product
Definition: CMlowranksd.hxx:70
virtual const CH_Matrix_Classes::Matrix & postgenmult(const CH_Matrix_Classes::Matrix &D, CH_Matrix_Classes::Matrix &C, CH_Matrix_Classes::Real alpha=1., CH_Matrix_Classes::Real beta=0., int dtrans=0) const
computes C= alpha*(*this)*D^(T if dtrans) + beta*C, C is also returned
Definition: CMlowranksd.hxx:214
Matrix class of symmetric matrices with real values of type Real
Definition: sparssym.hxx:89
for CMlowranksd
Definition: Coeffmat.hxx:34
void multiply(CH_Matrix_Classes::Real sf)
scales the scale factor
Definition: Coeffmat.hxx:70
implements a low rank matrix as Coeffmat with each a dense rectangular CH_Matrix_Classes::Matrix (f...
Definition: CMlowrankdd.hxx:24
virtual std::ostream & display(std::ostream &o) const
display constraint information
Definition: CMlowranksd.hxx:245
Integer coldim() const
returns the column dimension
Definition: matrix.hxx:218
virtual CH_Matrix_Classes::Real gramip(const CH_Matrix_Classes::Matrix &P, CH_Matrix_Classes::Integer start_row, const CH_Matrix_Classes::Matrix *Lam=0) const
returns ip(*this,QQ^T)=trace Q^T(*this)Q for Q=P.rows(start_row,start_row+dim-1)
Definition: CMlowranksd.hxx:79
virtual int equal(const Coeffmat *p, double tol=1e-6) const
returns 1, if p is the same derived class and entries differ by less than tol, otherwise zero ...
Definition: CMlowranksd.hxx:232
Matrix class for real values of type Real
Definition: matrix.hxx:74
CMlowranksd(const CH_Matrix_Classes::Sparsemat &Ain, const CH_Matrix_Classes::Matrix &Bin, CoeffmatInfo *cip=0)
copy Ain, Bin and store the user information
Definition: CMlowranksd.hxx:33
virtual Coeffmat * subspace(const CH_Matrix_Classes::Matrix &P) const
delivers a new object on the heap corresponding to the matrix P^T(*this)P, the caller is responsible ...
Definition: CMlowranksd.hxx:61
virtual const CH_Matrix_Classes::Matrix & pregenmult(const CH_Matrix_Classes::Matrix &D, CH_Matrix_Classes::Matrix &C, CH_Matrix_Classes::Real alpha=1., CH_Matrix_Classes::Real beta=0., int dtrans=0) const
computes C= alpha*D^(T if dtrans)*(*this) + beta*C, C is also returned
Definition: CMlowranksd.hxx:223
Val mat_ip(Integer len, const Val *x, const Val *y, const Val *d=0)
return sum(x[i]*y[i]) summing over len elements of the arrays x and y.
Definition: matop.hxx:1096
Matrix class of sparse matrices with real values of type Real
Definition: sparsmat.hxx:74
virtual CH_Matrix_Classes::Integer prodvec_flops() const
returns an estimate of number of flops to compute addprodto for a vector
Definition: CMlowranksd.hxx:149
void mat_xpeya(Integer len, Val *x, const Val *y, const Val a)
Set x[i]+=a*y[i] for len elements of the arrays x and y.
Definition: matop.hxx:543
CH_Matrix_Classes::Sparsemat A
this is A in A*B^T+B*A^T
Definition: CMlowranksd.hxx:29
const Indexmatrix & get_colindex() const
returns the index vector of the column representation holding the row index for each element ...
Definition: sparsmat.hxx:267
Integer rowdim() const
returns the row dimension
Definition: indexmat.hxx:321
virtual CH_Matrix_Classes::Real operator()(CH_Matrix_Classes::Integer i, CH_Matrix_Classes::Integer j) const
returns the value of the matrix element (i,j)
Definition: CMlowranksd.hxx:47
virtual CH_Matrix_Classes::Real gramip(const CH_Matrix_Classes::Matrix &P) const
returns ip(*this,PP^T)=trace P^T(*this)P
Definition: CMlowranksd.hxx:74
virtual void addprodto(CH_Matrix_Classes::Matrix &D, const CH_Matrix_Classes::Matrix &C, CH_Matrix_Classes::Real d=1.) const
computes D+=d*(*this)*C
Definition: CMlowranksd.hxx:128
Indexmatrix & genmult(const Indexmatrix &A, const Indexmatrix &B, Indexmatrix &C, Integer alpha=1, Integer beta=0, int atrans=0, int btrans=0)
returns C=beta*C+alpha*A*B, where A and B may be transposed; C must not be equal to A and B; if beta=...
Sparsemat row(Integer i) const
returns row i copied to a new sparse matrix
virtual int support_in(const CH_Matrix_Classes::Sparsesym &) const
returns 0 if the support of the costraint matrix is not contained in the support of the sparse symmet...
Definition: CMlowranksd.hxx:168
Symmatrix & rank2add(const Matrix &A, const Matrix &B, Symmatrix &C, Real alpha=1., Real beta=0., int trans=0)
returns C=beta*C+alpha*(A*B^T+B*A^T)/2 [or for transposed (A^T*B+B^T*A)/2]. If beta==0. then C is initiliazed to the correct size.
virtual int sparse(CH_Matrix_Classes::Indexmatrix &, CH_Matrix_Classes::Indexmatrix &, CH_Matrix_Classes::Matrix &, CH_Matrix_Classes::Real) const
returns 0 if not sparse. If it is sparse it returns 1 and the nonzero structure in I...
Definition: CMlowranksd.hxx:161
bool equal(const Matrix &A, const Matrix &B)
returns true if both matrices have the same size and the same elements
Definition: matrix.hxx:1433
Integer coldim() const
returns the column dimension
Definition: sparsmat.hxx:202
Integer rowdim() const
returns the row dimension
Definition: sparsmat.hxx:199
double sqrt(int a)
return sqrt for int a
Definition: mymath.hxx:121
CoeffmatInfo * clone(const CoeffmatInfo *cip)
if cip is not zero, it calls and returns cip->clone() and 0 otherwise
Definition: Coeffmat.hxx:106
Real ip(const Matrix &A, const Matrix &B)
returns the usual inner product of A and B, i.e., the sum of A(i,j)*B(i,j) over all i...
Definition: matrix.hxx:1165
virtual CH_Matrix_Classes::Integer dim() const
returns the order of the represented symmetric matrix
Definition: CMlowranksd.hxx:44
Header declaring the class ConicBundle::CMlowrankdd (needed for ConicBundle::AffineMatrixFunction) ...
Integer nonzeros() const
returns the number of nonzeros
Definition: sparsmat.hxx:205
virtual void addprodto(CH_Matrix_Classes::Matrix &D, const CH_Matrix_Classes::Sparsemat &C, CH_Matrix_Classes::Real d=1.) const
computes D+=d*(*this)*C
Definition: CMlowranksd.hxx:133