ConicBundle
CMsymdense.hxx
Go to the documentation of this file.
1 
2 
3 #ifndef CONICBUNDLE_CMSYMDENSE_HXX
4 #define CONICBUNDLE_CMSYMDENSE_HXX
5 
15 #include "Coeffmat.hxx"
16 
17 namespace ConicBundle {
18 
22 
27  class CMsymdense: public Coeffmat
28  {
29  private:
31  public:
34  {A=Ain;CM_type=CM_symdense;infop=cip;}
36  virtual ~CMsymdense(){}
37 
39  virtual Coeffmat* clone() const
40  { return new CMsymdense(A,ConicBundle::clone(infop)); }
41 
43  virtual CH_Matrix_Classes::Integer dim() const {return A.rowdim();}
44 
47 
49  virtual void make_symmatrix(CH_Matrix_Classes::Symmatrix& S) const {S=A;}
50 
52  virtual CH_Matrix_Classes::Real norm(void) const {return CH_Matrix_Classes::norm2(A);}
53 
55  virtual Coeffmat* subspace(const CH_Matrix_Classes::Matrix& P) const
57 
60  {A*=d; if (infop) infop->multiply(d);}
61 
64 
67 
70  {
73  const CH_Matrix_Classes::Integer brd=B.rowdim();
74  const CH_Matrix_Classes::Integer pcd=P.coldim();
75  const CH_Matrix_Classes::Integer prd=P.rowdim();
76  const CH_Matrix_Classes::Real *pp=P.get_store()+start_row;
77  const CH_Matrix_Classes::Real *vp=A.get_store();
78  for(CH_Matrix_Classes::Integer i=0;i<brd;i++){
79  CH_Matrix_Classes::mat_xpeya(pcd,bp+i,brd,pp+i,prd,(*vp++));
80  for(CH_Matrix_Classes::Integer j=i+1;j<brd;j++){
81  CH_Matrix_Classes::mat_xpeya(pcd,bp+i,brd,pp+j,prd,(*vp));
82  CH_Matrix_Classes::mat_xpeya(pcd,bp+j,brd,pp+i,prd,(*vp++));
83  }
84  }
85  CH_Matrix_Classes::Real trval=0.;
86  if (Lam==0){
87  for(CH_Matrix_Classes::Integer i=0;i<pcd;i++){
88  trval+=CH_Matrix_Classes::mat_ip(brd,pp,bp);
89  bp+=brd;
90  pp+=prd;
91  }
92  }
93  else {
94  assert(Lam->dim()==P.coldim());
95  for(CH_Matrix_Classes::Integer i=0;i<pcd;i++){
96  trval+=(*Lam)(i)*CH_Matrix_Classes::mat_ip(brd,pp,bp);
97  bp+=brd;
98  pp+=prd;
99  }
100  }
101  return trval;
102  }
103 
106 
109  { B.xpeya(A*C,d); }
110  //B+=d*(*this)*C
111 
114  { CH_Matrix_Classes::genmult(A,C,B,d,1.); }
115  //B+=d*(*this)*C
116 
119  {
120  if (P.coldim()<Q.coldim()){
122  CH_Matrix_Classes::genmult(A,P,tmp1,1.,0.,0);
123  CH_Matrix_Classes::genmult(tmp1,Q,R,1.,0.,1,0);
124  }
125  else {
127  CH_Matrix_Classes::genmult(A,Q,tmp1,1.,0.,0);
128  CH_Matrix_Classes::genmult(P,tmp1,R,1.,0.,1,0);
129  }
130  }
131 
134  { return 2*A.rowdim()*A.rowdim(); }
135 
137  virtual int dense() const
138  {return 1;}
139 
141  virtual int sparse() const
142  { return 0;}
143 
147  CH_Matrix_Classes::Matrix& /* val */,
148  CH_Matrix_Classes::Real /* d=1. */)const
149  { return 0;}
150 
152  virtual int support_in(const CH_Matrix_Classes::Sparsesym& /* S */) const
153  {return 0;}
154 
157  {return CH_Matrix_Classes::ip(A,S);}
158 
161  { S.init(CH_Matrix_Classes::transpose(P)*A*P); }
162 
165  {
168  const CH_Matrix_Classes::Integer brd=B.rowdim();
169  const CH_Matrix_Classes::Integer pcd=P.coldim();
170  const CH_Matrix_Classes::Integer prd=P.rowdim();
171  const CH_Matrix_Classes::Real *pp=P.get_store()+start_row;
172  const CH_Matrix_Classes::Real *vp=A.get_store();
173  for(CH_Matrix_Classes::Integer i=0;i<brd;i++){
174  CH_Matrix_Classes::mat_xpeya(pcd,bp+i,brd,pp+i,prd,(*vp++));
175  for(CH_Matrix_Classes::Integer j=i+1;j<brd;j++){
176  CH_Matrix_Classes::mat_xpeya(pcd,bp+i,brd,pp+j,prd,(*vp));
177  CH_Matrix_Classes::mat_xpeya(pcd,bp+j,brd,pp+i,prd,(*vp++));
178  }
179  }
181  {for(CH_Matrix_Classes::Integer i=0;i<pcd;i++){
182  const CH_Matrix_Classes::Real* bp2=bp;
183  if (alpha==1.){
184  for(CH_Matrix_Classes::Integer j=i;j<pcd;j++){
185  (*sp++)+=CH_Matrix_Classes::mat_ip(brd,pp,bp2);
186  bp2+=brd;
187  }
188  }
189  else {
190  for(CH_Matrix_Classes::Integer j=i;j<pcd;j++){
191  (*sp++)+=alpha*CH_Matrix_Classes::mat_ip(brd,pp,bp2);
192  bp2+=brd;
193  }
194  }
195  bp+=brd;
196  pp+=prd;
197  }}
198  chk_set_init(S,1);
199  }
200 
203  CH_Matrix_Classes::Real alpha=1.,CH_Matrix_Classes::Real beta=0.,int btrans=0) const
204  {
205  return CH_Matrix_Classes::genmult(A,B,C,alpha,beta,btrans);
206  }
207 
210  CH_Matrix_Classes::Real alpha=1.,CH_Matrix_Classes::Real beta=0.,int btrans=0) const
211  {
213  return CH_Matrix_Classes::genmult(B,A,C,alpha,beta,btrans);
214  }
215 
217  virtual int equal(const Coeffmat* p,double tol=1e-6) const
218  {
219  const CMsymdense *pp=dynamic_cast<const CMsymdense *>(p);
220  if (pp==0)
221  return 0;
222  if (A.rowdim()!=(pp->A).rowdim())
223  return 0;
224  return (CH_Matrix_Classes::norm2(A-pp->A)<tol);
225  }
226 
228  virtual std::ostream& display(std::ostream& o) const
229  {o<<"CMsymdense\n"; A.display(o); return o;}
230 
232  virtual std::ostream& out(std::ostream& o) const
233  {return o<<"SYMMETRIC_DENSE\n"<<A;}
234 
236  virtual std::istream& in(std::istream& i)
237  {return i>>A;}
238 
240  CMsymdense(std::istream& is,CoeffmatInfo* cip=0)
241  {CM_type=CM_symdense;infop=cip;in(is);}
242 
243  //--- specific routines
245  const CH_Matrix_Classes::Symmatrix& get_A() const {return A;}
246 
247  };
248 
250 }
251 #endif
252 
#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
virtual CH_Matrix_Classes::Real gramip(const CH_Matrix_Classes::Matrix &P) const
returns ip(*this,PP^T)=trace P^T(*this)P
Definition: CMsymdense.hxx:66
Header declaring the classes ConicBundle::Coeffmat, ConicBundle::CoeffmatPointer, ConicBundle::Coeffm...
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: CMsymdense.hxx:59
Matrix & xpeya(const Matrix &A, Real d=1.)
sets *this+=d*A and returns *this
virtual std::ostream & out(std::ostream &o) const
put entire contents onto outstream with the class type in the beginning so that the derived class can...
Definition: CMsymdense.hxx:232
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: CMsymdense.hxx:164
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: CMsymdense.hxx:152
Integer rowdim() const
returns the row dimension
Definition: matrix.hxx:215
double Real
all real numbers in calculations are of this type
Definition: matop.hxx:50
virtual std::ostream & display(std::ostream &o) const
display constraint information
Definition: CMsymdense.hxx:228
for CMsymdense
Definition: Coeffmat.hxx:31
allows to memorize the scalings applied to a Coeffmat and offers the basis for storing further user d...
Definition: Coeffmat.hxx:52
CMsymdense(const CH_Matrix_Classes::Symmatrix &Ain, CoeffmatInfo *cip=0)
copy Ain and possibly store the user information
Definition: CMsymdense.hxx:33
virtual CH_Matrix_Classes::Real ip(const CH_Matrix_Classes::Sparsesym &S) const
returns the inner product of the constraint matrix with S
Definition: CMsymdense.hxx:156
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: CMsymdense.hxx:46
virtual int dense() const
returns 1 if its structure is as bad as its dense symmetric representation, otherwise 0 ...
Definition: CMsymdense.hxx:137
CH_Matrix_Classes::Symmatrix A
the dense coefficient matrix
Definition: CMsymdense.hxx:30
Real * get_store()
returns the current address of the internal value array; use cautiously, do not use delete! ...
Definition: symmat.hxx:224
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: CMsymdense.hxx:55
defines a base class for coefficient matrices in semidefinite programming, in particular for use with...
Definition: Coeffmat.hxx:125
Symmatrix & init(const Symmatrix &A, double d=1.)
initialize to *this=A*d
Definition: symmat.hxx:753
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
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...
const CH_Matrix_Classes::Symmatrix & get_A() const
returns the const reference to the internal symmetric matrix
Definition: CMsymdense.hxx:245
Symmatrix & xetriu_yza(const Matrix &A, const Matrix &B, Real d=1.)
sets *this(i,j), i<=j to the upper triangle of the matrix product d*transpose(A)*B ...
Real * get_store()
returns the current address of the internal value array; use cautiously, do not use delete! ...
Definition: matrix.hxx:326
virtual void addprodto(CH_Matrix_Classes::Matrix &B, const CH_Matrix_Classes::Sparsemat &C, CH_Matrix_Classes::Real d=1.) const
computes B+=d*(*this)*C
Definition: CMsymdense.hxx:113
Matrix class of symmetric matrices with real values of type Real
Definition: symmat.hxx:43
Integer rowdim() const
returns the row dimension
Definition: symmat.hxx:159
virtual Coeffmat * clone() const
makes an explicit copy of itself and returns a pointer to it
Definition: CMsymdense.hxx:39
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: CMsymdense.hxx:69
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: CMsymdense.hxx:236
virtual void addprodto(CH_Matrix_Classes::Matrix &B, const CH_Matrix_Classes::Matrix &C, CH_Matrix_Classes::Real d=1.) const
computes B+=d*(*this)*C
Definition: CMsymdense.hxx:108
conic bundle method solver for sum of convex functions. See the ConicBundle_Manual for a quick introd...
Definition: CBSolver.hxx:22
Coeffmattype CM_type
in order to enable type identification
Definition: Coeffmat.hxx:134
CoeffmatInfo * infop
allows the user to specify and output additional information
Definition: Coeffmat.hxx:135
virtual const CH_Matrix_Classes::Matrix & postgenmult(const CH_Matrix_Classes::Matrix &B, CH_Matrix_Classes::Matrix &C, CH_Matrix_Classes::Real alpha=1., CH_Matrix_Classes::Real beta=0., int btrans=0) const
computes C= alpha*(*this)*B^(T if btrans) + beta*C, C is also returned
Definition: CMsymdense.hxx:202
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: CMsymdense.hxx:63
implements a general dense symmetric Coeffmat based on CH_Matrix_Classes::Symmatrix (for use with Mat...
Definition: CMsymdense.hxx:27
virtual CH_Matrix_Classes::Integer prodvec_flops() const
returns an estimate of number of flops to compute addprodto for a vector
Definition: CMsymdense.hxx:133
virtual CH_Matrix_Classes::Integer dim() const
returns the order of the represented symmetric matrix
Definition: CMsymdense.hxx:43
Matrix class of symmetric matrices with real values of type Real
Definition: sparssym.hxx:89
void multiply(CH_Matrix_Classes::Real sf)
scales the scale factor
Definition: Coeffmat.hxx:70
virtual void make_symmatrix(CH_Matrix_Classes::Symmatrix &S) const
returns a dense symmetric constraint matrix
Definition: CMsymdense.hxx:49
Integer coldim() const
returns the column dimension
Definition: matrix.hxx:218
Matrix class for real values of type Real
Definition: matrix.hxx:74
Indexmatrix transpose(const Indexmatrix &A)
returns an Indexmatrix that is the transpose of A
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: CMsymdense.hxx:217
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
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
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=...
CMsymdense(std::istream &is, CoeffmatInfo *cip=0)
constructor with istream and possibly additional user information
Definition: CMsymdense.hxx:240
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: CMsymdense.hxx:145
virtual CH_Matrix_Classes::Real norm(void) const
returns the Frobenius norm of the matrix
Definition: CMsymdense.hxx:52
CoeffmatInfo * clone(const CoeffmatInfo *cip)
if cip is not zero, it calls and returns cip->clone() and 0 otherwise
Definition: Coeffmat.hxx:106
virtual int sparse() const
returns 0 if not sparse, otherwise 1
Definition: CMsymdense.hxx:141
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 void project(CH_Matrix_Classes::Symmatrix &S, const CH_Matrix_Classes::Matrix &P) const
computes S=P^T*(*this)*P
Definition: CMsymdense.hxx:160
virtual const CH_Matrix_Classes::Matrix & pregenmult(const CH_Matrix_Classes::Matrix &B, CH_Matrix_Classes::Matrix &C, CH_Matrix_Classes::Real alpha=1., CH_Matrix_Classes::Real beta=0., int btrans=0) const
computes C= alpha*B^(T if btrans)*(*this) + beta*C, C is also returned
Definition: CMsymdense.hxx:209
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: CMsymdense.hxx:118
Symmatrix & xpeya(const Symmatrix &A, Real d=1.)
sets *this+=d*A and returns *this
virtual void addmeto(CH_Matrix_Classes::Symmatrix &S, CH_Matrix_Classes::Real d=1.) const
computes S+=d*(*this);
Definition: CMsymdense.hxx:105