gmsh-TingyuanDoc  0.1
An Open-Source Timing-driven Analytical Mixed-size FPGA Placer
fullMatrix.h
Go to the documentation of this file.
1 // Gmsh - Copyright (C) 1997-2022 C. Geuzaine, J.-F. Remacle
2 //
3 // See the LICENSE.txt file in the Gmsh root directory for license information.
4 // Please report all issues on https://gitlab.onelab.info/gmsh/gmsh/issues.
5 
6 #ifndef FULL_MATRIX_H
7 #define FULL_MATRIX_H
8 
9 #include "GmshConfig.h"
10 #include "GmshMessage.h"
11 #include <cmath>
12 #include <cstdio>
13 #include <complex>
14 #include <iostream>
15 
16 #if defined(HAVE_EIGEN)
17 #ifdef Success // in X11 header X.h
18 #undef Success
19 #endif
20 #include <sstream>
21 #include <Eigen/Dense>
22 #endif
23 
24 template <class scalar> class fullMatrix;
25 
26 // An interface for vectors of scalars (real or complex, with simple or double
27 // precision). The first index of a fullVector is 0. fullVectors can own their
28 // scalars, or just be an access point to an other fullVector; such a fullVector
29 // is called a proxy.
30 
31 template <class scalar> class fullVector {
32 private:
33  int _r; // size of the vector
34  scalar *_data; // pointer on the first element
35  bool _ownData;
36  friend class fullMatrix<scalar>;
37 
38 #if defined(HAVE_EIGEN)
39  typedef Eigen::Map<Eigen::Matrix<scalar, Eigen::Dynamic, 1> > EigenVec;
40 #endif
41 
42 public:
43  // Instantiate a zero size fullVector
44  fullVector(void) : _r(0), _data(0), _ownData(1) {}
45  // Instantiate a fullVector of size r filled with zeros
46  fullVector(int r) : _r(r), _ownData(1)
47  {
48  _data = new scalar[_r];
49  setAll(scalar(0.));
50  }
51  // Instantiate a proxy fullVector given by [original[0], original[1], ...,
52  // original[r - 1]]
53  fullVector(scalar *original, int r)
54  {
55  _r = r;
56  _ownData = false;
57  _data = original;
58  }
59  // Instantiate a fullVector, which is a copy (and not a proxy) of other
60  fullVector(const fullVector<scalar> &other) : _r(other._r), _ownData(1)
61  {
62  _data = new scalar[_r];
63  for(int i = 0; i < _r; ++i) _data[i] = other._data[i];
64  }
66  {
67  if(_ownData && _data) delete[] _data;
68  }
69  inline int size() const { return _r; }
70  inline const scalar *getDataPtr() const { return _data; }
71  inline scalar *getDataPtr() { return _data; }
72  inline scalar operator()(int i) const { return _data[i]; }
73  inline scalar &operator()(int i) { return _data[i]; }
75  {
76  if(this != &other) {
77  if((!resize(other._r, false) && _r > 2 * other._r)) {
78  if(_data) delete[] _data;
79  _r = other._r;
80  _data = new scalar[_r];
81  }
82  setAll(other);
83  }
84  return *this;
85  }
86  void copy(const fullVector<scalar> &v, int i0, int ni, int desti0)
87  {
88  for(int i = i0, desti = desti0; i < i0 + ni; i++, desti++)
89  (*this)(desti) = v(i);
90  }
91  inline void set(int r, scalar v)
92  {
93 #ifdef _DEBUG
94  if(r >= _r || r < 0) {
95  Msg::Error("Invalid index in vector: %i (size = %i)", r, _r);
96  return;
97  }
98 #endif
99  (*this)(r) = v;
100  }
101  // Return the L^2 norm
102  scalar norm() const;
103  bool resize(int r, bool resetValue = true)
104  {
105  if(_r < r || !_ownData) {
106  if(_ownData && _data) delete[] _data;
107  _r = r;
108  _data = new scalar[_r];
109  _ownData = true;
110  if(resetValue) setAll(scalar(0.));
111  return true;
112  }
113  _r = r;
114  if(resetValue) setAll(scalar(0.));
115  return false;
116  }
117  // this fullVector becomes a proxy of original, that is: [original(r_start),
118  // ..., original(r_start + r - 1)]. Previous data are lost
119  void setAsProxy(const fullVector<scalar> &original, int r_start, int r)
120  {
121  if(_ownData && _data) delete[] _data;
122  _ownData = false;
123  _r = r;
124  _data = original._data + r_start;
125  }
126  // This fullVector becomes a proxy of original matrix cth row, that is:
127  // [original(0, c), ..., original(original.size1() - 1, c)]. Previous data are
128  // lost
129  void setAsProxy(const fullMatrix<scalar> &original, int c)
130  {
131  if(_ownData && _data) delete[] _data;
132  _ownData = false;
133  _r = original._r;
134  _data = original._data + c * _r;
135  }
136  // This fullVector becomes a proxy of the array. Previous data are lost
137  void setAsProxy(scalar *data, int r)
138  {
139  if(_ownData && _data) delete[] _data;
140  _ownData = false;
141  _r = r;
142  _data = data;
143  }
144  inline void scale(const scalar s)
145  {
146 #if defined(HAVE_EIGEN)
147  EigenVec vv(_data, _r);
148  vv *= s;
149 #else
150  if(s == scalar(0.))
151  for(int i = 0; i < _r; ++i) _data[i] = scalar(0.);
152  else if(s == -1.)
153  for(int i = 0; i < _r; ++i) _data[i] = -_data[i];
154  else
155  for(int i = 0; i < _r; ++i) _data[i] *= s;
156 #endif
157  }
158  inline void setAll(const scalar &m)
159  {
160 #if defined(HAVE_EIGEN)
161  EigenVec vv(_data, _r);
162  vv.setConstant(m);
163 #else
164  for(int i = 0; i < _r; i++) set(i, m);
165 #endif
166  }
167  void setAll(const fullVector<scalar> &m)
168 #if defined(HAVE_EIGEN)
169  {
170  EigenVec vv(_data, _r);
171  EigenVec vm(m._data, m._r);
172  vv = vm;
173  }
174 #elif !defined(HAVE_BLAS)
175  {
176  for(int i = 0; i < _r; i++) _data[i] = m._data[i];
177  }
178 #endif
179  ;
180  // Scalar product
181  inline scalar operator*(const fullVector<scalar> &other)
182  {
183 #if defined(HAVE_EIGEN)
184  EigenVec vv(_data, _r), vother(other._data, other._r);
185  scalar s = vv.dot(vother);
186  return s;
187 #else
188  scalar s = 0.;
189  for(int i = 0; i < _r; ++i) s += _data[i] * other._data[i];
190  return s;
191 #endif
192  }
193  // v(i) = v(i) + alpha * x(i)
194  void axpy(const fullVector<scalar> &x, scalar alpha = 1.)
195 #if defined(HAVE_EIGEN)
196  {
197  EigenVec vv(_data, _r), vx(x._data, x._r);
198  vv.noalias() += alpha * vx;
199  }
200 #elif !defined(HAVE_BLAS)
201  {
202  for(int i = 0; i < _r; i++) _data[i] += alpha * x._data[i];
203  }
204 #endif
205  ;
206  // v(i) = v(i) * x(i)
208  {
209  for(int i = 0; i < _r; i++) _data[i] *= x._data[i];
210  }
211  void print(const std::string name = "", const std::string format = "") const;
212  void binarySave(FILE *f) const { fwrite(_data, sizeof(scalar), _r, f); }
213  void binaryLoad(FILE *f)
214  {
215  if(fread(_data, sizeof(scalar), _r, f) != (size_t)_r) return;
216  }
217  bool getOwnData() const { return _ownData; };
218  void setOwnData(bool ownData) { _ownData = ownData; };
219 };
220 
221 // An interface for dense matrix of scalars (real or complex, with simple or
222 // double precision)
223 
224 template <class scalar> class fullMatrix {
225 private:
226  bool _ownData; // should data be freed on delete?
227  int _r, _c; // size of the matrix
228  scalar *_data; // pointer on the first element
229  friend class fullVector<scalar>;
230 
231 #if defined(HAVE_EIGEN)
232  typedef Eigen::Map<Eigen::Matrix<scalar, Eigen::Dynamic, Eigen::Dynamic> > EigenMat;
233  typedef Eigen::Map<Eigen::Matrix<scalar, Eigen::Dynamic, 1> > EigenVec;
234 #endif
235 
236 public:
237  fullMatrix(scalar *original, int r, int c)
238  {
239  _r = r;
240  _c = c;
241  _ownData = false;
242  _data = original;
243  }
244  fullMatrix(fullMatrix<scalar> &original, int c_start, int c)
245  {
246  _c = c;
247  _r = original._r;
248  _ownData = false;
249  _data = original._data + c_start * _r;
250  }
251  fullMatrix(int r, int c, bool init0 = true) : _r(r), _c(c)
252  {
253  _data = new scalar[_r * _c];
254  _ownData = true;
255  if(init0) setAll(scalar(0.));
256  }
257  fullMatrix(int r, int c, scalar *data)
258  : _r(r), _c(c), _data(data), _ownData(false)
259  {
260  setAll(scalar(0.));
261  }
262  fullMatrix(const fullMatrix<scalar> &other) : _r(other._r), _c(other._c)
263  {
264  _data = new scalar[_r * _c];
265  _ownData = true;
266  for(int i = 0; i < _r * _c; ++i) _data[i] = other._data[i];
267  }
268  fullMatrix() : _ownData(false), _r(0), _c(0), _data(0) {}
270  {
271  if(_data && _ownData) delete[] _data;
272  }
273  // get information (size, value)
274  inline int size1() const { return _r; }
275  inline int size2() const { return _c; }
276  inline scalar get(int r, int c) const
277  {
278 #ifdef _DEBUG
279  if(r >= _r || r < 0 || c >= _c || c < 0) {
280  Msg::Error("Invalid index in dense matrix: %i %i (size = %i %i)", r, c,
281  _r, _c);
282  return 0;
283  }
284 #endif
285  return (*this)(r, c);
286  }
287  inline const scalar *getDataPtr() const { return _data; }
288  inline scalar *getDataPtr() { return _data; }
289  inline void set(int r, int c, scalar v)
290  {
291 #ifdef _DEBUG
292  if(r >= _r || r < 0 || c >= _c || c < 0) {
293  Msg::Error("Invalid index in dense matrix: %i %i (size = %i %i)", r, c,
294  _r, _c);
295  return;
296  }
297 #endif
298  (*this)(r, c) = v;
299  }
300  inline scalar norm() const
301  {
302  scalar n = 0.;
303  for(int i = 0; i < _r; ++i)
304  for(int j = 0; j < _c; ++j) n += (*this)(i, j) * (*this)(i, j);
305  return sqrt(n);
306  }
307  bool resize(int r, int c, bool resetValue = true)
308  {
309  // data will be owned (same as constructor)
310  if((r * c > _r * _c) || !_ownData) {
311  if(_ownData && _data) delete[] _data;
312  _r = r;
313  _c = c;
314  _data = new scalar[_r * _c];
315  _ownData = true;
316  if(resetValue) setAll(scalar(0.));
317  return true;
318  }
319  _r = r;
320  _c = c;
321  if(resetValue) setAll(scalar(0.));
322  return false; // no reallocation
323  }
324  void reshape(int nbRows, int nbColumns)
325  {
326  if(nbRows == -1 && nbColumns != -1) nbRows = _r * _c / nbColumns;
327  if(nbRows != -1 && nbColumns == -1) nbColumns = _r * _c / nbRows;
328  if(nbRows * nbColumns != _r * _c)
329  Msg::Error("Invalid dense matrix reshape: total number of entries must "
330  "be equal (new %i x %i != old %i x %i)",
331  nbRows, nbColumns, _r, _c);
332  _r = nbRows;
333  _c = nbColumns;
334  }
335  void setAsProxy(const fullMatrix<scalar> &original)
336  {
337  if(_data && _ownData) delete[] _data;
338  _c = original._c;
339  _r = original._r;
340  _ownData = false;
341  _data = original._data;
342  }
343  void setAsProxy(const fullMatrix<scalar> &original, int c_start, int c)
344  {
345  if(_data && _ownData) delete[] _data;
346  _c = c;
347  _r = original._r;
348  _ownData = false;
349  _data = original._data + c_start * _r;
350  }
351  void setAsProxy(scalar *data, int r, int c)
352  {
353  if(_data && _ownData) delete[] _data;
354  _c = c;
355  _r = r;
356  _ownData = false;
357  _data = data;
358  }
360  {
361  copy(other);
362  return *this;
363  }
364  void operator+=(const fullMatrix<scalar> &other)
365  {
366  if(_r != other._r || _c != other._c) {
367  Msg::Error("Cannot sum dense matrices of different sizes");
368  return;
369  }
370  for(int i = 0; i < _r * _c; ++i) _data[i] += other._data[i];
371  }
372  inline scalar operator()(int i, int j) const
373  {
374 #ifdef _DEBUG
375  if(i >= _r || i < 0 || j >= _c || j < 0) {
376  Msg::Error("Invalid index to access dense matrix: %i %i (size = %i %i)",
377  i, j, _r, _c);
378  return 0;
379  }
380 #endif
381  return _data[i + _r * j];
382  }
383  inline scalar &operator()(int i, int j)
384  {
385 #ifdef _DEBUG
386  if(i >= _r || i < 0 || j >= _c || j < 0) {
387  Msg::Error("Invalid index to access dense matrix: %i %i (size = %i %i)",
388  i, j, _r, _c);
389  return _data[0];
390  }
391 #endif
392  return _data[i + _r * j];
393  }
394  void copy(const fullMatrix<scalar> &a, int i0, int ni, int j0, int nj,
395  int desti0, int destj0)
396  {
397  for(int i = i0, desti = desti0; i < i0 + ni; i++, desti++)
398  for(int j = j0, destj = destj0; j < j0 + nj; j++, destj++)
399  (*this)(desti, destj) = a(i, j);
400  }
401  void copy(const fullMatrix<scalar> &a)
402  {
403  if(_data && !_ownData) {
404  Msg::Error("Dense matrix copy prohibited for proxies, use setAll "
405  "instead");
406  return;
407  }
408  if(_r != a._r || _c != a._c) {
409  if(_data && _ownData) delete[] _data;
410  _r = a._r;
411  _c = a._c;
412  _data = new scalar[_r * _c];
413  _ownData = true;
414  }
415  setAll(a);
416  }
417  void copyOneColumn(const fullVector<scalar> &x, const int ind) const
418  {
419  int cind = _c * ind;
420  for(int i = 0; i < _r; i++) _data[cind + i] = x(i);
421  }
422  inline void setAll(const scalar &m)
423  {
424 #if defined(HAVE_EIGEN)
425  EigenMat ma(_data, _r, _c);
426  ma.setConstant(m);
427 #else
428  for(int i = 0; i < _r * _c; i++) _data[i] = m;
429 #endif
430  }
431  void setAll(const fullMatrix<scalar> &m)
432 #if defined(HAVE_EIGEN)
433  {
434  EigenMat ma(_data, _r, _c), mm(m._data, m._r, m._c);
435  ma = mm;
436  }
437 #elif !defined(HAVE_BLAS)
438  {
439  if(_r != m._r || _c != m._c) {
440  Msg::Error("Dense matrix sizes do not match in setAll");
441  return;
442  }
443  for(int i = 0; i < _r * _c; i++) _data[i] = m._data[i];
444  }
445 #endif
446  ;
447  void scale(const scalar s)
448 #if defined(HAVE_EIGEN)
449  {
450  EigenMat ma(_data, _r, _c);
451  ma *= s;
452  }
453 #elif !defined(HAVE_BLAS)
454  {
455  if(s == 0.) // this is not really correct nan*0 (or inf*0) is expected to
456  // give nan
457  for(int i = 0; i < _r * _c; ++i) _data[i] = scalar(0.);
458  else
459  for(int i = 0; i < _r * _c; ++i) _data[i] *= s;
460  }
461 #endif
462  ;
463  inline void add(const double &a)
464  {
465  for(int i = 0; i < _r * _c; ++i) _data[i] += a;
466  }
467  inline void add(const fullMatrix<scalar> &m)
468  {
469 #if defined(HAVE_EIGEN)
470  EigenMat ma(_data, _r, _c), mm(m._data, m._r, m._c);
471  ma.noalias() += mm;
472 #else
473  for(int i = 0; i < _r; i++)
474  for(int j = 0; j < _c; j++) (*this)(i, j) += m(i, j);
475 #endif
476  }
477  inline void add(const fullMatrix<scalar> &m, const scalar &a)
478  {
479 #if defined(HAVE_EIGEN)
480  EigenMat ma(_data, _r, _c), mm(m._data, m._r, m._c);
481  ma.noalias() += a * mm;
482 #else
483  for(int i = 0; i < _r; i++)
484  for(int j = 0; j < _c; j++) (*this)(i, j) += a * m(i, j);
485 #endif
486  }
487  void mult(const fullVector<scalar> &x, fullVector<scalar> &y) const
488 #if defined(HAVE_EIGEN)
489  {
490  EigenMat ma(_data, _r, _c);
491  EigenVec vx(x._data, x._r), vy(y._data, y._r);
492  vy = ma * vx;
493  }
494 #elif !defined(HAVE_BLAS)
495  {
496  y.scale(scalar(0.));
497  for(int i = 0; i < _r; i++)
498  for(int j = 0; j < _c; j++) y._data[i] += (*this)(i, j) * x(j);
499  }
500 #endif
501  ;
503 #if defined(HAVE_EIGEN)
504  {
505  EigenMat ma(_data, _r, _c);
506  EigenVec vx(x._data, x._r), vy(y._data, y._r);
507  vy += ma * vx;
508  }
509 #elif !defined(HAVE_BLAS)
510  {
511  for(int i = 0; i < _r; i++)
512  for(int j = 0; j < _c; j++) y._data[i] += (*this)(i, j) * x(j);
513  }
514 #endif
515  ;
517 #if defined(HAVE_EIGEN)
518  {
519  EigenMat ma(_data, _r, _c), mb(b._data, b._r, b._c), mc(c._data, c._r, c._c);
520  mc.noalias() = ma * mb;
521  }
522 #elif !defined(HAVE_BLAS)
523  {
524  c.scale(scalar(0.));
525  for(int i = 0; i < _r; i++)
526  for(int j = 0; j < b._c; j++)
527  for(int k = 0; k < _c; k++)
528  c._data[i + _r * j] += (*this)(i, k) * b(k, j);
529  }
530 #endif
531  ;
533  {
534  for(int i = 0; i < _r * _c; i++) _data[i] *= a._data[i];
535  }
536  void multOnBlock(const fullMatrix<scalar> &b, const int ncol, const int fcol,
537  const int alpha, const int beta, fullVector<scalar> &c) const
538 #if defined(HAVE_EIGEN) || !defined(HAVE_BLAS)
539  {
540  int row = 0;
541  if(beta != 1) c.scale(beta);
542  for(int j = fcol; j < fcol + ncol; j++)
543  for(int k = 0; k < _c; k++)
544  c._data[j] += alpha * (*this)(row, k) * b(k, j);
545  }
546 #endif
547  ;
548  void multWithATranspose(const fullVector<scalar> &x, scalar alpha,
549  scalar beta, fullVector<scalar> &y) const
550 #if defined(HAVE_EIGEN) || !defined(HAVE_BLAS)
551  {
552  y.scale(beta);
553  for(int j = 0; j < _c; j++)
554  for(int i = 0; i < _r; i++) y._data[j] += alpha * (*this)(i, j) * x(i);
555  }
556 #endif
557  ;
559  {
561  for(int i = 0; i < _r; i++)
562  for(int j = 0; j < _c; j++) T(j, i) = (*this)(i, j);
563  return T;
564  }
565  inline void transposeInPlace()
566  {
567  if(_r != _c) {
568  Msg::Error("In-place transposition requires a square matrix "
569  "(size = %d %d)", _r, _c);
570  return;
571  }
572  scalar t;
573  for(int i = 0; i < _r; i++)
574  for(int j = 0; j < i; j++) {
575  t = _data[i + _r * j];
576  _data[i + _r * j] = _data[j + _r * i];
577  _data[j + _r * i] = t;
578  }
579  }
580  void gemm(const fullMatrix<scalar> &a, const fullMatrix<scalar> &b,
581  scalar alpha = 1., scalar beta = 1., bool transposeA = false,
582  bool transposeB = false)
583 #if defined(HAVE_EIGEN) || !defined(HAVE_BLAS)
584  {
585  const fullMatrix<scalar> &A = transposeA ? a.transpose() : a;
586  const fullMatrix<scalar> &B = transposeA ? b.transpose() : b;
587  fullMatrix<scalar> temp(A._r, B._c);
588  A.mult(B, temp);
589  temp.scale(alpha);
590  scale(beta);
591  add(temp);
592  }
593 #endif
594  ;
595  void axpy(const fullMatrix<scalar> &x, scalar alpha = 1.)
596 #if defined(HAVE_EIGEN) || !defined(HAVE_BLAS)
597  {
598  int n = _r * _c;
599  for(int i = 0; i < n; i++) _data[i] += alpha * x._data[i];
600  }
601 #endif
602  ;
604 #if defined(HAVE_EIGEN)
605  {
606  if(_r != _c || _r != rhs._r || _r != result._r) {
607  Msg::Error("Wrong sizes for dense linear system solve (size = %d %d, "
608  "%d, %d)", _r, _c, result._r, rhs._r);
609  return false;
610  }
611  EigenMat ma(_data, _r, _c);
612  EigenVec vb(rhs._data, rhs._r);
613  EigenVec vx(result._data, result._r);
614  vx = ma.colPivHouseholderQr().solve(vb);
615  return true;
616  }
617 #elif !defined(HAVE_LAPACK)
618  {
619  Msg::Error("LU factorization and substitution requires Eigen or LAPACK");
620  return false;
621  }
622 #endif
623  ;
625 #if defined(HAVE_EIGEN) || !defined(HAVE_LAPACK)
626  {
627  Msg::Error("LU factorization requires LAPACK");
628  return false;
629  }
630 #endif
631  ;
633  fullVector<scalar> &result)
634 #if defined(HAVE_EIGEN) || !defined(HAVE_LAPACK)
635  {
636  Msg::Error("LU substitution requires LAPACK");
637  return false;
638  }
639 #endif
640  ;
641  bool invert(fullMatrix<scalar> &result) const
642 #if defined(HAVE_EIGEN)
643  {
644  if(_r != _c) {
645  Msg::Error("Dense matrix inverse requires square matrix (size = %d %d)",
646  _r, _c);
647  return false;
648  }
649  result.resize(_r, _c);
650  EigenMat ma(_data, _r, _c);
651  EigenMat mi(result._data, _r, _c);
652  mi = ma.inverse();
653  return true;
654  }
655 #elif !defined(HAVE_LAPACK)
656  {
657  Msg::Error("Matrix inversion requires Eigen or LAPACK");
658  return false;
659  }
660 #endif
661  ;
663 #if defined(HAVE_EIGEN)
664  {
665  if(_r != _c) {
666  Msg::Error("Dense matrix inversion requires square matrix (size = %d %d)",
667  _r, _c);
668  return false;
669  }
670  EigenMat ma(_data, _r, _c);
671  ma = ma.inverse();
672  return true;
673  }
674 #elif !defined(HAVE_LAPACK)
675  {
676  Msg::Error("Dense matrix inversion requires LAPACK");
677  return false;
678  }
679 #endif
680  ;
681  scalar determinant() const
682 #if defined(HAVE_EIGEN)
683  {
684  EigenMat ma(_data, _r, _c);
685  return ma.determinant();
686  }
687 #elif !defined(HAVE_LAPACK)
688  {
689  Msg::Error("Dense matrix inversion requires Eigen or LAPACK");
690  return false;
691  }
692 #endif
693  ;
694  void swap(scalar *a, int inca, scalar *b, int incb, int n)
695  {
696  scalar tmp;
697  for(int i = 0; i < n; i++, a += inca, b += incb) {
698  tmp = (*a);
699  (*a) = (*b);
700  (*b) = tmp;
701  }
702  }
703  void eigSort(int n, scalar *wr, scalar *wi, scalar *VL, scalar *VR)
704  {
705  // Sort the eigenvalues/vectors in ascending order according to
706  // their real part. Warning: this will screw up the ordering if we
707  // have complex eigenvalues.
708  for(int i = 0; i < n - 1; i++) {
709  int k = i;
710  scalar ek = wr[i];
711  // search for something to swap
712  for(int j = i + 1; j < n; j++) {
713  const scalar ej = wr[j];
714  if(ej < ek) {
715  k = j;
716  ek = ej;
717  }
718  }
719  if(k != i) {
720  swap(&wr[i], 1, &wr[k], 1, 1);
721  swap(&wi[i], 1, &wi[k], 1, 1);
722  swap(&VL[n * i], 1, &VL[n * k], 1, n);
723  swap(&VR[n * i], 1, &VR[n * k], 1, n);
724  }
725  }
726  }
727  bool eig(fullVector<double> &eigenValReal, fullVector<double> &eigenValImag,
728  fullMatrix<scalar> &leftEigenVect, fullMatrix<scalar> &rightEigenVect,
729  bool sortRealPart = false)
730 #if defined(HAVE_EIGEN)
731  {
732  EigenMat ma(_data, _r, _c);
733  Eigen::EigenSolver<Eigen::Matrix<scalar, Eigen::Dynamic, Eigen::Dynamic> > es(ma);
734  if(es.info() != Eigen::Success) {
735  Msg::Warning("Eigen could not compute eigenvalues/eigenvectors");
736  return false;
737  }
738  EigenVec vr(eigenValReal._data, eigenValReal._r);
739  vr = es.eigenvalues().real();
740  EigenVec vi(eigenValImag._data, eigenValImag._r);
741  vi = es.eigenvalues().imag();
742  EigenMat mr(rightEigenVect._data, rightEigenVect._r, rightEigenVect._c);
743  mr = es.eigenvectors().real();
744  EigenMat ml(leftEigenVect._data, leftEigenVect._r, leftEigenVect._c);
745  // FIXME: only true for symmetric matrices - compute the true left eigenvectors!
746  ml = mr;
747 
748  if(sortRealPart)
749  eigSort(_r, eigenValReal._data, eigenValImag._data,
750  leftEigenVect._data, rightEigenVect._data);
751  return true;
752  }
753 #elif !defined(HAVE_LAPACK)
754  {
755  Msg::Error("Eigenvalue computation of dense matrices requires Eigen or "
756  "LAPACK");
757  return false;
758  }
759 #endif
760  ;
762 #if defined(HAVE_EIGEN)
763  {
764  EigenMat ma(_data, _r, _c);
765  EigenMat mv(V._data, V._r, V._c);
766  EigenVec vs(S._data, S._r);
767  Eigen::BDCSVD
768  <Eigen::Matrix<scalar, Eigen::Dynamic, Eigen::Dynamic> >
769  svd(ma, Eigen::ComputeThinU | Eigen::ComputeThinV);
770  ma = svd.matrixU();
771  mv = svd.matrixV();
772  vs = svd.singularValues();
773  return true;
774  }
775 #elif !defined(HAVE_LAPACK)
776  {
777  Msg::Error("Singular value decomposition of dense matrices requires "
778  "Eigen or LAPACK");
779  return false;
780  }
781 #endif
782  ;
783  void print(const std::string &name = "",
784  const std::string &format = "") const;
785  void binarySave(FILE *f) const { fwrite(_data, sizeof(scalar), _r * _c, f); }
786  void binaryLoad(FILE *f)
787  {
788  if(fread(_data, sizeof(scalar), _r * _c, f) != (size_t)_r) return;
789  }
790  bool getOwnData() const { return _ownData; };
791  void setOwnData(bool ownData) { _ownData = ownData; };
792 };
793 
794 #endif
fullMatrix::add
void add(const double &a)
Definition: fullMatrix.h:463
fullVector::_r
int _r
Definition: fullMatrix.h:33
fullVector::setAsProxy
void setAsProxy(const fullVector< scalar > &original, int r_start, int r)
Definition: fullMatrix.h:119
fullVector::set
void set(int r, scalar v)
Definition: fullMatrix.h:91
fullMatrix::set
void set(int r, int c, scalar v)
Definition: fullMatrix.h:289
fullMatrix::copy
void copy(const fullMatrix< scalar > &a)
Definition: fullMatrix.h:401
fullMatrix::add
void add(const fullMatrix< scalar > &m)
Definition: fullMatrix.h:467
fullVector
Definition: MElement.h:26
fullVector::multTByT
void multTByT(const fullVector< scalar > &x)
Definition: fullMatrix.h:207
fullVector::_ownData
bool _ownData
Definition: fullMatrix.h:35
c
static double c(int i, int j, fullMatrix< double > &CA, const std::vector< SPoint3 > &P, const std::vector< SPoint3 > &Q)
Definition: discreteFrechetDistance.cpp:15
fullMatrix::setAsProxy
void setAsProxy(const fullMatrix< scalar > &original, int c_start, int c)
Definition: fullMatrix.h:343
fullVector::fullVector
fullVector(int r)
Definition: fullMatrix.h:46
Msg::Warning
static void Warning(const char *fmt,...)
Definition: GmshMessage.cpp:543
fullMatrix::setAll
void setAll(const scalar &m)
Definition: fullMatrix.h:422
Msg::Error
static void Error(const char *fmt,...)
Definition: GmshMessage.cpp:482
fullVector::norm
scalar norm() const
fullMatrix::_r
int _r
Definition: fullMatrix.h:227
LegendrePolynomials::f
void f(int n, double u, double *val)
Definition: orthogonalBasis.cpp:77
fullMatrix::operator()
scalar & operator()(int i, int j)
Definition: fullMatrix.h:383
fullMatrix::scale
void scale(const scalar s)
Definition: fullMatrix.h:447
fullVector::size
int size() const
Definition: fullMatrix.h:69
fullMatrix::transpose
fullMatrix< scalar > transpose() const
Definition: fullMatrix.h:558
fullVector::scale
void scale(const scalar s)
Definition: fullMatrix.h:144
fullVector::resize
bool resize(int r, bool resetValue=true)
Definition: fullMatrix.h:103
fullVector::binarySave
void binarySave(FILE *f) const
Definition: fullMatrix.h:212
fullMatrix::get
scalar get(int r, int c) const
Definition: fullMatrix.h:276
fullVector::setAll
void setAll(const fullVector< scalar > &m)
Definition: fullMatrix.h:167
fullMatrix::print
void print(const std::string &name="", const std::string &format="") const
fullMatrix::setAsProxy
void setAsProxy(scalar *data, int r, int c)
Definition: fullMatrix.h:351
GmshMessage.h
fullMatrix::luFactor
bool luFactor(fullVector< int > &ipiv)
Definition: fullMatrix.h:624
fullVector::getOwnData
bool getOwnData() const
Definition: fullMatrix.h:217
fullMatrix::norm
scalar norm() const
Definition: fullMatrix.h:300
fullVector::fullVector
fullVector(const fullVector< scalar > &other)
Definition: fullMatrix.h:60
fullMatrix::setOwnData
void setOwnData(bool ownData)
Definition: fullMatrix.h:791
fullVector::_data
scalar * _data
Definition: fullMatrix.h:34
fullVector::fullVector
fullVector(void)
Definition: fullMatrix.h:44
fullMatrix::mult
void mult(const fullMatrix< scalar > &b, fullMatrix< scalar > &c) const
Definition: fullMatrix.h:516
fullMatrix
Definition: MElement.h:27
fullVector::setAsProxy
void setAsProxy(scalar *data, int r)
Definition: fullMatrix.h:137
fullVector::operator=
fullVector< scalar > & operator=(const fullVector< scalar > &other)
Definition: fullMatrix.h:74
fullMatrix::fullMatrix
fullMatrix()
Definition: fullMatrix.h:268
fullMatrix::luSubstitute
bool luSubstitute(const fullVector< scalar > &rhs, fullVector< int > &ipiv, fullVector< scalar > &result)
Definition: fullMatrix.h:632
fullMatrix::multOnBlock
void multOnBlock(const fullMatrix< scalar > &b, const int ncol, const int fcol, const int alpha, const int beta, fullVector< scalar > &c) const
Definition: fullMatrix.h:536
fullMatrix::copyOneColumn
void copyOneColumn(const fullVector< scalar > &x, const int ind) const
Definition: fullMatrix.h:417
fullMatrix::fullMatrix
fullMatrix(fullMatrix< scalar > &original, int c_start, int c)
Definition: fullMatrix.h:244
fullMatrix::reshape
void reshape(int nbRows, int nbColumns)
Definition: fullMatrix.h:324
fullMatrix::multAddy
void multAddy(const fullVector< scalar > &x, fullVector< scalar > &y) const
Definition: fullMatrix.h:502
fullMatrix::getDataPtr
scalar * getDataPtr()
Definition: fullMatrix.h:288
fullVector::fullVector
fullVector(scalar *original, int r)
Definition: fullMatrix.h:53
fullMatrix::binaryLoad
void binaryLoad(FILE *f)
Definition: fullMatrix.h:786
fullMatrix::invertInPlace
bool invertInPlace()
Definition: fullMatrix.h:662
fullMatrix::determinant
scalar determinant() const
Definition: fullMatrix.h:681
fullVector::operator*
scalar operator*(const fullVector< scalar > &other)
Definition: fullMatrix.h:181
fullVector::~fullVector
~fullVector()
Definition: fullMatrix.h:65
fullMatrix::eig
bool eig(fullVector< double > &eigenValReal, fullVector< double > &eigenValImag, fullMatrix< scalar > &leftEigenVect, fullMatrix< scalar > &rightEigenVect, bool sortRealPart=false)
Definition: fullMatrix.h:727
fullMatrix::operator()
scalar operator()(int i, int j) const
Definition: fullMatrix.h:372
fullVector::getDataPtr
scalar * getDataPtr()
Definition: fullMatrix.h:71
fullMatrix::binarySave
void binarySave(FILE *f) const
Definition: fullMatrix.h:785
fullMatrix::gemm
void gemm(const fullMatrix< scalar > &a, const fullMatrix< scalar > &b, scalar alpha=1., scalar beta=1., bool transposeA=false, bool transposeB=false)
Definition: fullMatrix.h:580
S
#define S
Definition: DefaultOptions.h:21
fullMatrix::transposeInPlace
void transposeInPlace()
Definition: fullMatrix.h:565
fullMatrix::size2
int size2() const
Definition: fullMatrix.h:275
fullMatrix::_ownData
bool _ownData
Definition: fullMatrix.h:226
fullMatrix::axpy
void axpy(const fullMatrix< scalar > &x, scalar alpha=1.)
Definition: fullMatrix.h:595
fullMatrix::setAsProxy
void setAsProxy(const fullMatrix< scalar > &original)
Definition: fullMatrix.h:335
fullMatrix::operator=
fullMatrix< scalar > & operator=(const fullMatrix< scalar > &other)
Definition: fullMatrix.h:359
fullMatrix::multWithATranspose
void multWithATranspose(const fullVector< scalar > &x, scalar alpha, scalar beta, fullVector< scalar > &y) const
Definition: fullMatrix.h:548
fullVector::setAsProxy
void setAsProxy(const fullMatrix< scalar > &original, int c)
Definition: fullMatrix.h:129
fullVector::getDataPtr
const scalar * getDataPtr() const
Definition: fullMatrix.h:70
fullMatrix::size1
int size1() const
Definition: fullMatrix.h:274
fullMatrix::_c
int _c
Definition: fullMatrix.h:227
fullMatrix::getOwnData
bool getOwnData() const
Definition: fullMatrix.h:790
fullMatrix::svd
bool svd(fullMatrix< scalar > &V, fullVector< scalar > &S)
Definition: fullMatrix.h:761
fullVector::binaryLoad
void binaryLoad(FILE *f)
Definition: fullMatrix.h:213
fullMatrix::_data
scalar * _data
Definition: fullMatrix.h:228
fullVector::setOwnData
void setOwnData(bool ownData)
Definition: fullMatrix.h:218
fullVector::copy
void copy(const fullVector< scalar > &v, int i0, int ni, int desti0)
Definition: fullMatrix.h:86
fullMatrix::copy
void copy(const fullMatrix< scalar > &a, int i0, int ni, int j0, int nj, int desti0, int destj0)
Definition: fullMatrix.h:394
fullMatrix::getDataPtr
const scalar * getDataPtr() const
Definition: fullMatrix.h:287
fullMatrix::operator+=
void operator+=(const fullMatrix< scalar > &other)
Definition: fullMatrix.h:364
fullMatrix::luSolve
bool luSolve(const fullVector< scalar > &rhs, fullVector< scalar > &result)
Definition: fullMatrix.h:603
fullMatrix::eigSort
void eigSort(int n, scalar *wr, scalar *wi, scalar *VL, scalar *VR)
Definition: fullMatrix.h:703
fullMatrix::mult
void mult(const fullVector< scalar > &x, fullVector< scalar > &y) const
Definition: fullMatrix.h:487
fullMatrix::multTByT
void multTByT(const fullMatrix< scalar > &a)
Definition: fullMatrix.h:532
fullMatrix::add
void add(const fullMatrix< scalar > &m, const scalar &a)
Definition: fullMatrix.h:477
fullMatrix::fullMatrix
fullMatrix(int r, int c, scalar *data)
Definition: fullMatrix.h:257
fullMatrix::resize
bool resize(int r, int c, bool resetValue=true)
Definition: fullMatrix.h:307
fullMatrix::swap
void swap(scalar *a, int inca, scalar *b, int incb, int n)
Definition: fullMatrix.h:694
fullVector::operator()
scalar operator()(int i) const
Definition: fullMatrix.h:72
fullMatrix::~fullMatrix
~fullMatrix()
Definition: fullMatrix.h:269
fullMatrix::fullMatrix
fullMatrix(int r, int c, bool init0=true)
Definition: fullMatrix.h:251
fullVector::setAll
void setAll(const scalar &m)
Definition: fullMatrix.h:158
fullVector::axpy
void axpy(const fullVector< scalar > &x, scalar alpha=1.)
Definition: fullMatrix.h:194
fullMatrix::setAll
void setAll(const fullMatrix< scalar > &m)
Definition: fullMatrix.h:431
fullVector::operator()
scalar & operator()(int i)
Definition: fullMatrix.h:73
fullMatrix::fullMatrix
fullMatrix(const fullMatrix< scalar > &other)
Definition: fullMatrix.h:262
fullMatrix::invert
bool invert(fullMatrix< scalar > &result) const
Definition: fullMatrix.h:641
fullMatrix::fullMatrix
fullMatrix(scalar *original, int r, int c)
Definition: fullMatrix.h:237
fullVector::print
void print(const std::string name="", const std::string format="") const