fmatvec  0.0.0
symmetric_matrix.h
1/* Copyright (C) 2003-2005 Martin Förg
2
3 * This library is free software; you can redistribute it and/or
4 * modify it under the terms of the GNU Lesser General Public
5 * License as published by the Free Software Foundation; either
6 * version 2.1 of the License, or (at your option) any later version.
7 *
8 * This library is distributed in the hope that it will be useful,
9 * but WITHOUT ANY WARRANTY; without even the implied warranty of
10 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
11 * Lesser General Public License for more details.
12 *
13 * You should have received a copy of the GNU Lesser General Public
14 * License along with this library; if not, write to the Free Software
15 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
16 *
17 * Contact:
18 * martin.o.foerg@googlemail.com
19 *
20 */
21
22#ifndef symmetric_matrix_h
23#define symmetric_matrix_h
24
25#include "types.h"
26#include "indices.h"
27#include "matrix.h"
28#include "_memory.h"
29
30namespace fmatvec {
31
40 template <class AT> class Matrix<Symmetric,Ref,Ref,AT> {
41
42 protected:
43
45
46 Memory<AT> memory;
47 AT *ele;
48 int n{0};
49 int lda{0};
50
51 friend class Matrix<General,Ref,Ref,AT>;
52
53 template <class Type, class Row, class Col> inline Matrix<Symmetric,Ref,Ref,AT>& copy(const Matrix<Type,Row,Col,AT> &A);
56
57 const AT* elePtr(int i, int j) const {
58 return j > i ? ele+i*lda+j : ele+i+j*lda;
59 }
60
61 AT* elePtr(int i, int j) {
62 return j > i ? ele+i*lda+j : ele+i+j*lda;
63 }
64
66
67 public:
68
69 static constexpr bool isVector {false};
70 using value_type = AT;
71 using shape_type = Symmetric;
72
77 explicit Matrix() : memory(), ele(nullptr) {
78 }
79
80 explicit Matrix(int n_, Noinit) : memory(n_*n_), ele((AT*)memory.get()), n(n_), lda(n_) { }
81 explicit Matrix(int n_, Init ini=INIT, const AT &a=AT()) : memory(n_*n_), ele((AT*)memory.get()), n(n_), lda(n_) { init(a); }
82 explicit Matrix(int n_, Eye ini, const AT &a=1) : memory(n_*n_), ele((AT*)memory.get()), n(n_), lda(n_) { init(ini,a); }
83 explicit Matrix(int m_, int n_, Noinit) : memory(n_*n_), ele((AT*)memory.get()), n(n_), lda(n_) { }
84 explicit Matrix(int m_, int n_, Init ini=INIT, const AT &a=AT()) : memory(n_*n_), ele((AT*)memory.get()), n(n_), lda(n_) { init(a); }
85 explicit Matrix(int m_, int n_, Eye ini, const AT &a=1) : memory(n_*n_), ele((AT*)memory.get()), n(n_), lda(n_) { init(ini,a); }
86
92 Matrix(const Matrix<Symmetric,Ref,Ref,AT> &A) : memory(A.n*A.n), ele((AT*)memory.get()) , n(A.n), lda(n) {
93 copy(A);
94 }
95
101 template<class Row>
102 Matrix(const Matrix<Symmetric,Row,Row,AT> &A) : memory(A.size()*A.size()), ele((AT*)memory.get()), n(A.size()), lda(n) {
103 copy(A);
104 }
105
111 template<class Type, class Row, class Col>
112 explicit Matrix(const Matrix<Type,Row,Col,AT> &A) : memory(A.rows()*A.cols()), ele((AT*)memory.get()), n(A.cols()), lda(n) {
113 FMATVEC_ASSERT(A.rows() == A.cols(), AT);
114 copy(A);
115 }
116
124 explicit Matrix(int n_, AT* ele_) : memory(), ele(ele_), n(n_), lda(n_) { }
125
128 ~Matrix() = default;
129
130 Matrix<Symmetric,Ref,Ref,AT>& resize(int n_, Noinit) {
131 n = n_; lda = n;
132 memory.resize(n*n);
133 ele = (AT*)memory.get();
134 return *this;
135 }
136
137 Matrix<Symmetric,Ref,Ref,AT>& resize(int n, Init ini=INIT, const AT &a=AT()) { return resize(n,Noinit()).init(a); }
138
139 Matrix<Symmetric,Ref,Ref,AT>& resize(int n, Eye ini, const AT &a=1) { return resize(n,Noinit()).init(ini,a); }
140
143 void resize(int n, int m) {
144 if(n!=m)
145 throw std::runtime_error("A symmetric matrix cannot have different dimensions for rows and columns.");
146 resize(n);
147 }
148
156 FMATVEC_ASSERT(n == A.rows(), AT);
157 return copy(A);
158 }
159
166 template<class Type, class Row, class Col>
168 FMATVEC_ASSERT(A.rows() == A.cols(), AT);
169 FMATVEC_ASSERT(n == A.rows(), AT);
170 return copy(A);
171 }
172
180 n=A.n;
181 memory = A.memory;
182 ele = A.ele;
183 lda = A.lda;
184 return *this;
185 }
186
193 template<class Type, class Row, class Col>
194 inline Matrix<Symmetric,Ref,Ref,AT>& operator<<=(const Matrix<Type,Row,Col,AT> &A) {
195 FMATVEC_ASSERT(A.rows() == A.cols(), AT);
196 if(n!=A.rows()) resize(A.rows(),NONINIT);
197 return copy(A);
198 }
199
209 AT& operator()(int i, int j) {
210 FMATVEC_ASSERT(i>=0, AT);
211 FMATVEC_ASSERT(j>=0, AT);
212 FMATVEC_ASSERT(i<n, AT);
213 FMATVEC_ASSERT(j<n, AT);
214 return e(i,j);
215 }
216
221 const AT& operator()(int i, int j) const {
222 FMATVEC_ASSERT(i>=0, AT);
223 FMATVEC_ASSERT(j>=0, AT);
224 FMATVEC_ASSERT(i<n, AT);
225 FMATVEC_ASSERT(j<n, AT);
226
227 return e(i,j);
228 }
229
230 AT& ei(int i, int j) {
231 return ele[i+j*lda];
232 }
233
234 const AT& ei(int i, int j) const {
235 return ele[i+j*lda];
236 }
237
238 AT& ej(int i, int j) {
239 return ele[i*lda+j];
240 }
241
242 const AT& ej(int i, int j) const {
243 return ele[i*lda+j];
244 }
245
246 AT& e(int i, int j) {
247 return j > i ? ej(i,j) : ei(i,j);
248 }
249
250 const AT& e(int i, int j) const {
251 return j > i ? ej(i,j) : ei(i,j);
252 }
253
259 AT* operator()() {return ele;}
260
265 const AT* operator()() const {return ele;}
266
271 int size() const {return n;}
272
277 int rows() const {return n;}
278
283 int cols() const {return n;}
284
289 int ldim() const {return lda;}
290
298 CBLAS_ORDER blasOrder() const {
299 return CblasColMajor;
300 }
301
309 CBLAS_UPLO blasUplo() const {
310 return CblasLower;
311 }
312
320 inline Matrix<Symmetric,Ref,Ref,AT>& init(const AT &val);
321 inline Matrix<Symmetric,Ref,Ref,AT>& init(Init, const AT &a=AT()) { return init(a); }
322 inline Matrix<Symmetric,Ref,Ref,AT>& init(Eye, const AT &val=1);
323 inline Matrix<Symmetric,Ref,Ref,AT>& init(Noinit, const AT &a=AT()) { return *this; }
324
327 inline const Matrix<General,Ref,Ref,AT> operator()(const Range<Var,Var> &I, const Range<Var,Var> &J) const;
328
331 inline const Matrix<Symmetric,Ref,Ref,AT> operator()(const Range<Var,Var> &I) const;
332
333 template<class Type, class Row, class Col> inline void set(const Range<Var,Var> &I, const Range<Var,Var> &J, const Matrix<Type,Row,Col,AT> &A);
334 template<class Type, class Row, class Col> inline void add(const Range<Var,Var> &I, const Range<Var,Var> &J, const Matrix<Type,Row,Col,AT> &A);
335
336 template<class Row> inline void set(const Range<Var,Var> &I, const Matrix<Symmetric,Row,Row,AT> &A);
337 template<class Row> inline void add(const Range<Var,Var> &I, const Matrix<Symmetric,Row,Row,AT> &A);
338
339 inline const Matrix<General,Ref,Ref,AT> operator()(const Indices &I, const Indices &J) const;
340
341 inline const Matrix<Symmetric,Ref,Ref,AT> operator()(const Indices &I) const;
342
343 inline void ref(Matrix<Symmetric,Ref,Ref,AT> &A, const Range<Var,Var> &I);
344
349 explicit inline operator std::vector<std::vector<AT>>() const;
350
356 explicit Matrix(const std::vector<std::vector<AT>> &m);
357
358// /*! \brief Cast to AT.
359// *
360// * \return The AT representation of the matrix
361// * */
362// explicit operator AT() const {
363// FMATVEC_ASSERT(n==1, AT);
364// return e(0);
365// }
366//
367// /*! \brief AT Constructor.
368// * Constructs and initializes a matrix with a AT object.
369// * \param x The AT the matrix will be initialized with.
370// * */
371// explicit Matrix(const AT &x) : memory(1), ele((AT*)memory.get()), n(1), lda(1) { ele[0] = x; }
372
373 inline int nonZeroElements() const;
374 };
375
376 template <class AT>
378
379 for(int i=0; i<rows(); i++)
380 for(int j=i; j<cols(); j++)
381 ej(i,j) = val;
382
383 return *this;
384 }
385
386 template <class AT>
388
389 for(int i=0; i<size(); i++) {
390 ej(i,i) = val;
391 for(int j=i+1; j<size(); j++) {
392 ej(i,j) = 0;
393 }
394 }
395 return *this;
396 }
397
398 template <class AT>
400 FMATVEC_ASSERT(I.end()<n, AT);
401 FMATVEC_ASSERT(J.end()<n, AT);
402 Matrix<General,Ref,Ref,AT> A(I.size(),J.size(),NONINIT);
403
404 for(int i=0; i<A.rows(); i++)
405 for(int j=0; j<A.cols(); j++)
406 A.e(i,j) = e(I.start()+i,J.start()+j);
407
408 return A;
409 }
410
411 template <class AT>
413 FMATVEC_ASSERT(I.end()<n, AT);
414
415 Matrix<Symmetric,Ref,Ref,AT> A(I.size(),NONINIT);
416
417 for(int i=0; i<A.rows(); i++)
418 for(int j=i; j<A.cols(); j++)
419 A.e(i,j) = e(I.start()+i,I.start()+j);
420
421 return A;
422 }
423
424 template <class AT> template<class Type, class Row, class Col>
426
427 if(I.start()>=J.start()) FMATVEC_ASSERT(I.start()>=J.end(), AT)
428 else FMATVEC_ASSERT(J.start()>=I.end(), AT);
429 FMATVEC_ASSERT(I.end()<rows(), AT);
430 FMATVEC_ASSERT(J.end()<cols(), AT);
431 FMATVEC_ASSERT(I.size()==A.rows(), AT);
432 FMATVEC_ASSERT(J.size()==A.cols(), AT);
433
434 for(int i=I.start(), ii=0; i<=I.end(); i++, ii++)
435 for(int j=J.start(), jj=0; j<=J.end(); j++, jj++)
436 e(i,j) = A.e(ii,jj);
437 }
438
439 template <class AT> template<class Type, class Row, class Col>
440 inline void Matrix<Symmetric,Ref,Ref,AT>::add(const Range<Var,Var> &I, const Range<Var,Var> &J, const Matrix<Type,Row,Col,AT> &A) {
441
442 if(I.start()>=J.start()) FMATVEC_ASSERT(I.start()>=J.end(), AT)
443 else FMATVEC_ASSERT(J.start()>=I.end(), AT);
444 FMATVEC_ASSERT(I.end()<rows(), AT);
445 FMATVEC_ASSERT(J.end()<cols(), AT);
446 FMATVEC_ASSERT(I.size()==A.rows(), AT);
447 FMATVEC_ASSERT(J.size()==A.cols(), AT);
448
449 for(int i=I.start(), ii=0; i<=I.end(); i++, ii++)
450 for(int j=J.start(), jj=0; j<=J.end(); j++, jj++)
451 e(i,j) += A.e(ii,jj);
452 }
453
454 template <class AT> template<class Row>
455 inline void Matrix<Symmetric,Ref,Ref,AT>::set(const Range<Var,Var> &I, const Matrix<Symmetric,Row,Row,AT> &A) {
456
457 FMATVEC_ASSERT(I.end()<size(), AT);
458 FMATVEC_ASSERT(I.size()==A.size(), AT);
459
460 for(int i=I.start(), ii=0; i<=I.end(); i++, ii++)
461 for(int j=i, jj=ii; j<=I.end(); j++, jj++)
462 e(i,j) = A.e(ii,jj);
463 }
464
465 template <class AT> template<class Row>
466 inline void Matrix<Symmetric,Ref,Ref,AT>::add(const Range<Var,Var> &I, const Matrix<Symmetric,Row,Row,AT> &A) {
467
468 FMATVEC_ASSERT(I.end()<size(), AT);
469 FMATVEC_ASSERT(I.size()==A.size(), AT);
470
471 for(int i=I.start(), ii=0; i<=I.end(); i++, ii++)
472 for(int j=i, jj=ii; j<=I.end(); j++, jj++)
473 e(i,j) += A.e(ii,jj);
474 }
475
476 template <class AT>
477 inline const Matrix<General,Ref,Ref,AT> Matrix<Symmetric,Ref,Ref,AT>::operator()(const Indices &I, const Indices &J) const {
478 FMATVEC_ASSERT(I.max()<size(), AT);
479 FMATVEC_ASSERT(J.max()<size(), AT);
480
481 Matrix<General,Ref,Ref,AT> A(I.size(),J.size(),NONINIT);
482
483 for(int i=0; i<A.rows(); i++)
484 for(int j=0; j<A.cols(); j++)
485 A.e(i,j) = e(I[i],J[j]);
486
487 return A;
488 }
489
490 template <class AT>
491 inline const Matrix<Symmetric,Ref,Ref,AT> Matrix<Symmetric,Ref,Ref,AT>::operator()(const Indices &I) const {
492 FMATVEC_ASSERT(I.max()<size(), AT);
493
494 Matrix<Symmetric,Ref,Ref,AT> A(I.size(),NONINIT);
495
496 for(int i=0; i<A.rows(); i++)
497 for(int j=0; j<A.cols(); j++)
498 A.e(i,j) = e(I[i],I[j]);
499
500 return A;
501 }
502
503 template <class AT>
504 inline void Matrix<Symmetric,Ref,Ref,AT>::ref(Matrix<Symmetric,Ref,Ref,AT> &A, const Range<Var,Var> &I) {
505 FMATVEC_ASSERT(I.end()<A.size(), AT);
506 n=I.size();
507 memory = A.memory;
508 ele = A.elePtr(I.start(),I.start());
509 lda = A.lda;
510 }
511
512 template <class AT>
513 inline Matrix<Symmetric,Ref,Ref,AT>::operator std::vector<std::vector<AT>>() const {
514 std::vector<std::vector<AT>> ret(rows(),std::vector<AT>(cols()));
515 for(int r=0; r<rows(); r++) {
516 for(int c=r; c<cols(); c++) {
517 ret[r][c]=ej(r,c);
518 if(c>r)
519 ret[c][r]=ej(r,c);
520 }
521 }
522 return ret;
523 }
524
525 template <class AT>
526 Matrix<Symmetric,Ref,Ref,AT>::Matrix(const std::vector<std::vector<AT>> &m) : memory(m.size()*m.size()), ele((AT*)memory.get()), n(static_cast<int>(m.size())), lda(static_cast<int>(m.size())) {
527 for(int r=0; r<rows(); r++) {
528 if(static_cast<int>(m[r].size())!=cols())
529 throw std::runtime_error("The rows of the input have different length.");
530 for(int c=r; c<cols(); c++) {
531 ej(r,c)=m[r][c];
532 if(c>r && abs(m[r][c]-m[c][r])>abs(m[r][c]*1e-13+1e-13))
533 throw std::runtime_error("The input is not symmetric.");
534 }
535 }
536 }
537
538 template <class AT>
540 int k = n;
541 for (int j = 0; j < n; j++) {
542 for (int i = j+1; i < n; i++) {
543 if (fabs(e(i, j)) > 1e-16) {
544 k+=2;
545 }
546 }
547 }
548 return k;
549 }
550
552
553 template <class AT>
554 inline Matrix<Symmetric,Ref,Ref,AT>& Matrix<Symmetric,Ref,Ref,AT>::copy(const Matrix<Symmetric,Ref,Ref,AT> &A) {
555 for(int i=0; i<n; i++)
556 for(int j=i; j<n; j++)
557 ej(i,j) = A.ej(i,j);
558 return *this;
559 }
560
561 template <class AT> template <class Type, class Row, class Col>
562 inline Matrix<Symmetric,Ref,Ref,AT>& Matrix<Symmetric,Ref,Ref,AT>::copy(const Matrix<Type,Row,Col,AT> &A) {
563 for(int i=0; i<n; i++)
564 for(int j=i; j<n; j++)
565 ej(i,j) = A.e(i,j);
566 return *this;
567 }
568
569 template <class AT>
570 inline Matrix<Symmetric,Ref,Ref,AT>& Matrix<Symmetric,Ref,Ref,AT>::copy(const Matrix<SymmetricSparse,Ref,Ref,AT> &A) {
571 for(int i=0; i<n; i++) {
572 for(int j=i; j<n; j++)
573 ej(i,j) = 0;
574 for(int j=A.Ip()[i]; j<A.Ip()[i+1]; j++)
575 ej(i,A.Jp()[j]) = A()[j];
576 }
577 return *this;
578 }
579
581
582}
583
584#endif
Definition: types.h:94
Shape class for general matrices.
Definition: types.h:116
Definition: types.h:93
This is a matrix class for general matrices.
Definition: general_matrix.h:42
int cols() const
Number of columns.
Definition: general_matrix.h:313
int rows() const
Number of rows.
Definition: general_matrix.h:307
This is a matrix class for sparse quadratic matrices.
Definition: symmetric_sparse_matrix.h:44
This is a matrix class for symmetric matrices.
Definition: symmetric_matrix.h:40
Matrix(const Matrix< Symmetric, Ref, Ref, AT > &A)
Copy Constructor.
Definition: symmetric_matrix.h:92
Matrix(const Matrix< Symmetric, Row, Row, AT > &A)
Copy Constructor.
Definition: symmetric_matrix.h:102
const AT & operator()(int i, int j) const
Element operator.
Definition: symmetric_matrix.h:221
const AT * operator()() const
Pointer operator.
Definition: symmetric_matrix.h:265
Matrix(const Matrix< Type, Row, Col, AT > &A)
Copy Constructor.
Definition: symmetric_matrix.h:112
int ldim() const
Leading dimension.
Definition: symmetric_matrix.h:289
Matrix< Symmetric, Ref, Ref, AT > & operator=(const Matrix< Symmetric, Ref, Ref, AT > &A)
Assignment operator.
Definition: symmetric_matrix.h:155
Matrix< Symmetric, Ref, Ref, AT > & operator=(const Matrix< Type, Row, Col, AT > &A)
Assignment operator.
Definition: symmetric_matrix.h:167
AT * operator()()
Pointer operator.
Definition: symmetric_matrix.h:259
int size() const
Size.
Definition: symmetric_matrix.h:271
Matrix< Symmetric, Ref, Ref, AT > & operator&=(Matrix< Symmetric, Ref, Ref, AT > &A)
Reference operator.
Definition: symmetric_matrix.h:179
Matrix(int n_, AT *ele_)
Regular Constructor.
Definition: symmetric_matrix.h:124
CBLAS_UPLO blasUplo() const
Symmetry convention.
Definition: symmetric_matrix.h:309
CBLAS_ORDER blasOrder() const
Storage convention.
Definition: symmetric_matrix.h:298
int cols() const
Number of columns.
Definition: symmetric_matrix.h:283
int rows() const
Number of rows.
Definition: symmetric_matrix.h:277
AT & operator()(int i, int j)
Element operator.
Definition: symmetric_matrix.h:209
Matrix< Symmetric, Ref, Ref, AT > & init(const AT &val)
Initialization.
Definition: symmetric_matrix.h:377
Matrix()
Standard constructor.
Definition: symmetric_matrix.h:77
void resize(int n, int m)
Definition: symmetric_matrix.h:143
This is the basic matrix class for arbitrary matrices.
Definition: matrix.h:52
int rows() const
Number of rows.
int cols() const
Number of columns.
AT & operator()(int i, int j)
Standard constructor.
Definition: matrix.h:84
Definition: types.h:92
This is an index class for creating submatrices.
Definition: range.h:44
int end() const
Last element.
Definition: range.h:91
int size() const
Size.
Definition: range.h:97
int start() const
First element.
Definition: range.h:85
Definition: types.h:102
Shape class for symmetric matrices.
Definition: types.h:132
Namespace fmatvec.
Definition: _memory.cc:28