IFEM  90A354
matrix.h
Go to the documentation of this file.
1 // $Id$
2 //==============================================================================
17 //==============================================================================
18 
19 #ifndef UTL_MATRIX_H
20 #define UTL_MATRIX_H
21 
22 #include <vector>
23 #include <iostream>
24 #include <algorithm>
25 #include <cstring>
26 #include <cctype>
27 #include <cmath>
28 #include "BLAS.h"
29 #include "print_tol.h"
30 
31 #ifdef INDEX_CHECK
32 #if INDEX_CHECK > 1
33 #define ABORT_ON_INDEX_CHECK abort()
34 #else
35 #define ABORT_ON_INDEX_CHECK
36 #endif
37 #define CHECK_INDEX(label,i,n) if (i < 1 || i > n) { \
38  std::cerr << label << i <<" is out of range [1,"<< n <<"]"<< std::endl; \
39  ABORT_ON_INDEX_CHECK; }
40 #else
41 #define CHECK_INDEX(label,i,n)
42 #define ABORT_ON_INDEX_CHECK
43 #endif
44 
45 #ifdef SING_CHECK
46 #define ABORT_ON_SINGULARITY abort()
47 #else
48 #define ABORT_ON_SINGULARITY
49 #endif
50 
51 
52 namespace utl
53 {
55  const char RETAIN = 2;
56 
63  template<class T> class vector
64  {
65  public:
67  vector() {}
69  explicit vector(size_t n) { this->resize(n); }
71  vector(const T* values, size_t n) { this->fill(values,n); }
73  vector(const std::vector<T>& X) : myVec(X) {}
74 
76  vector<T>& operator=(const std::vector<T>& X)
77  {
78  myVec = X;
79  return *this;
80  }
81 
83  T* ptr() { return myVec.empty() ? nullptr : myVec.data(); }
85  const T* ptr() const { return myVec.empty() ? nullptr : myVec.data(); }
86 
88  size_t size() const { return myVec.size(); }
90  bool empty() const { return myVec.empty(); }
92  bool zero(T tol = T(0)) const
93  {
94  return std::all_of(myVec.begin(), myVec.end(),
95  [tol](T v) { return std::fabs(v) <= tol; });
96  }
97 
99  using ConstVecIter = typename std::vector<T>::const_iterator;
101  using VecIter = typename std::vector<T>::iterator;
102 
104  ConstVecIter begin() const { return myVec.begin(); }
106  ConstVecIter end() const { return myVec.end(); }
108  VecIter begin(){ return myVec.begin(); }
110  VecIter end() { return myVec.end(); }
111 
113  operator const std::vector<T>&() const { return myVec; }
115  operator std::vector<T>&() { return myVec; }
116 
118  T& operator[](size_t i) { return myVec[i]; }
120  const T& operator[](size_t i) const { return myVec[i]; }
121 
123  T& operator()(size_t i)
124  {
125  CHECK_INDEX("vector::operator(): Index ",i,myVec.size());
126  return myVec[i-1];
127  }
128 
130  const T& operator()(size_t i) const
131  {
132  CHECK_INDEX("vector::operator(): Index ",i,myVec.size());
133  return myVec[i-1];
134  }
135 
137  void fill(T s) { std::fill(myVec.begin(),myVec.end(),s); }
139  void fill(const T* values, size_t n = 0)
140  {
141  if (n > myVec.size())
142  myVec.resize(n);
143  memcpy(myVec.data(),values,myVec.size()*sizeof(T));
144  }
145 
147  void push_back(T c) { myVec.push_back(c); }
148 
151  {
152  myVec.insert(myVec.end(),i1,i2);
153  }
154 
156  void push_back(const T* p, const T* q) { myVec.insert(myVec.end(),p,q); }
157 
159  void swap(vector<T>& vec) { myVec.swap(vec.myVec); }
160 
162  vector<T>& operator*=(T c);
164  vector<T>& operator/=(T d) { return this->operator*=(T(1)/d); }
165 
167  vector<T>& operator*=(const std::vector<T>& X)
168  {
169  for (size_t i = 0; i < myVec.size() && i < X.size(); i++)
170  myVec[i] *= X[i];
171  return *this;
172  }
174  vector<T>& operator/=(const std::vector<T>& X)
175  {
176  for (size_t i = 0; i < myVec.size() && i < X.size(); i++)
177  myVec[i] *= (X[i] == T(0) ? T(0) : T(1)/X[i]);
178  return *this;
179  }
180 
182  vector<T>& operator+=(const vector<T>& X) { return this->add(X); }
184  vector<T>& operator-=(const vector<T>& X) { return this->add(X,T(-1)); }
186  vector<T>& add(const std::vector<T>& X, const T& alfa = T(1),
187  unsigned int ofsx = 0, int stridex = 1,
188  unsigned int ofsy = 0, int stridey = 1);
189 
192  vector<T>& relax(T alfa, const std::vector<T>& X)
193  {
194  if (alfa != T(1))
195  {
196  this->operator*=(alfa);
197  this->add(X,T(1)-alfa);
198  }
199  return *this;
200  }
203  vector<T>& relax(T alfa, const std::vector<T>& X, const std::vector<T>& Y)
204  {
205  return this->operator=(Y).relax(alfa,X);
206  }
207 
215  T dot(const T* v, size_t nv,
216  size_t off1 = 0, int inc1 = 1,
217  size_t off2 = 0, int inc2 = 1) const;
218 
225  T dot(const std::vector<T>& v,
226  size_t off1 = 0, int inc1 = 1,
227  size_t off2 = 0, int inc2 = 1) const
228  {
229  return this->dot(v.data(),v.size(),off1,inc1,off2,inc2);
230  }
231 
235  T norm2(size_t off = 0, int inc = 1) const;
241  T normInf(size_t& off, int inc = 1, bool sign = false) const;
244  T normInf(int inc = 1) const { size_t o = 0; return this->normInf(o,inc); }
245 
247  T max() const { return *std::max_element(myVec.begin(),myVec.end()); }
249  T min() const { return *std::min_element(myVec.begin(),myVec.end()); }
250 
254  T asum(size_t off = 0, int inc = 1) const;
255 
260  T sum(size_t off = 0, int inc = 1, size_t max = 0) const
261  {
262  T xsum = T(0);
263  if (inc < 1 || myVec.empty())
264  return xsum;
265 
266  if (max == 0 || max > myVec.size())
267  max = myVec.size();
268  for (size_t i = off; i < max; i += inc)
269  xsum += myVec[i];
270  return xsum;
271  }
272 
277  bool resize(size_t n, char forceClear = 0)
278  {
279  if (n == myVec.size())
280  {
281  if (forceClear == 1)
282  this->fill(T(0)); // Erase previous content
283  return false; // Size is not changed
284  }
285 
286  if (forceClear < RETAIN)
287  myVec.clear();
288  myVec.resize(n,T(0));
289  return true;
290  }
291 
293  void reserve(size_t n) { myVec.reserve(n); }
295  void clear() { myVec.clear(); }
296 
297  private:
298  std::vector<T> myVec;
299  };
300 
301 
307  template<class T> class matrixBase
308  {
309  protected:
311  matrixBase() : n{0,0,0,0}, elem(myElem) {}
313  explicit matrixBase(vector<T>& vec) : n{0,0,0,0}, elem(vec) {}
316  matrixBase(size_t n_1, size_t n_2, size_t n_3 = 1, size_t n_4 = 1)
317  : n{n_1,n_2,n_3,n_4}, elem(myElem), myElem(n_1*n_2*n_3*n_4) {}
318 
322  matrixBase(const matrixBase<T>& mat, bool copyContent = true)
323  : elem(myElem), myElem(mat.size())
324  {
325  memcpy(n,mat.n,sizeof(n));
326  if (copyContent)
327  elem = mat.elem;
328  }
329 
338  void redim(size_t n_1, size_t n_2, size_t n_3, size_t n_4, bool forceClear)
339  {
340  if (forceClear)
341  {
342  // Erase previous content
343  if (this->size() == n_1*n_2*n_3*n_4)
344  this->fill(T(0));
345  else
346  this->clear();
347  }
348 
349  if (n[0] == n_1 && n[1] == n_2 && n[2] == n_3 && n[3] == n_4)
350  return; // nothing to do
351 
352  size_t oldn1 = n[0];
353  size_t oldn2 = n[1];
354  size_t oldn3 = n[2];
355  size_t oldSize = this->size();
356  n[0] = n_1;
357  n[1] = n_2;
358  n[2] = n_3;
359  n[3] = n_4;
360  if (this->size() == oldSize)
361  return; // no more to do, size is unchanged
362 
363  // If the size in any of the matrix dimensions, except for the last one,
364  // are changed the previous matrix content must be cleared
365  if (!forceClear)
366  this->clearIfNrowChanged(oldn1,oldn2,oldn3);
367 
368  elem.resize(n[0]*n[1]*n[2]*n[3],RETAIN);
369  }
370 
372  virtual void clearIfNrowChanged(size_t n1, size_t n2, size_t n3) = 0;
373 
374  public:
376  size_t dim(short int d = 1) const { return d > 0 && d <= 4 ? n[d-1] : 0; }
378  size_t size() const { return n[0]*n[1]*n[2]*n[3]; }
380  bool empty() const { return elem.empty(); }
382  bool zero(T tol = T(0)) const { return elem.zero(tol); }
383 
385  const vector<T>& toVec() const { return elem; }
387  operator const std::vector<T>&() const { return elem; }
389  operator std::vector<T>&() { return elem; }
390 
392  T* ptr(size_t c = 0)
393  {
394  return n[0]*c < elem.size() ? elem.ptr() + n[0]*c : nullptr;
395  }
397  const T* ptr(size_t c = 0) const
398  {
399  return n[0]*c < elem.size() ? elem.ptr() + n[0]*c : nullptr;
400  }
401 
403  typename std::vector<T>::iterator begin() { return elem.begin(); }
405  typename std::vector<T>::iterator end() { return elem.end(); }
406 
408  void clear() { n[0] = n[1] = n[2] = n[3] = 0; elem.clear(); }
409 
411  void fill(T s) { std::fill(elem.begin(),elem.end(),s); }
413  void fill(const T* values, size_t n = 0) { elem.fill(values,n); }
414 
416  matrixBase<T>& add(const matrixBase<T>& A, const T& alfa);
418  matrixBase<T>& multiply(const T& c);
419 
422  T norm2(int inc = 1) const { return elem.norm2(0,inc); }
425  T asum(int inc = 1) const { return elem.asum(0,inc); }
429  T sum(int inc = 1) const
430  {
431  if (inc > 0)
432  return elem.sum(0,inc);
433  else if (inc == 0 || (inc *= -1) > static_cast<int>(n[1]))
434  return T(0);
435  else
436  return elem.sum((inc-1)*n[0],1,inc*n[0]);
437  }
438 
439  protected:
440  size_t n[4];
442 
443  private:
445  };
446 
447 
455  template<class T> class matrix : public matrixBase<T>
456  {
457  public:
459  matrix() : nrow(this->n[0]), ncol(this->n[1]) {}
461  explicit matrix(vector<T>& vec)
462  : matrixBase<T>(vec), nrow(this->n[0]), ncol(this->n[1]) {}
464  matrix(size_t r, size_t c)
465  : matrixBase<T>(r,c), nrow(this->n[0]), ncol(this->n[1]) {}
467  matrix(const matrix<T>& mat, bool transposed = false)
468  : matrixBase<T>(mat,false), nrow(this->n[0]), ncol(this->n[1])
469  {
470  nrow = transposed ? mat.ncol : mat.nrow;
471  ncol = transposed ? mat.nrow : mat.ncol;
472  if (transposed)
473  for (size_t r = 0; r < ncol; r++)
474  for (size_t c = 0; c < nrow; c++)
475  this->elem[c+nrow*r] = mat.elem[r+ncol*c];
476  else if (!mat.elem.empty())
477  this->elem.fill(mat.elem.ptr());
478  }
480  virtual ~matrix() {}
481 
488  void resize(size_t r, size_t c, bool forceClear = false)
489  {
490  this->redim(r,c,1,1,forceClear);
491  }
492 
494  matrix<T>& expandRows(int incRows)
495  {
496  int newRows = nrow + incRows;
497  if (newRows < 1 || ncol < 1)
498  // The matrix is empty
499  this->clear();
500  else if (incRows < 0)
501  {
502  // The matrix size is reduced
503  T* newMat = this->ptr() + newRows;
504  for (size_t c = 1; c < ncol; c++, newMat += newRows)
505  memmove(newMat,this->ptr(c),newRows*sizeof(T));
506  nrow = newRows;
507  this->elem.resize(nrow*ncol,RETAIN);
508  }
509  else if (incRows > 0)
510  {
511  // The matrix size is increased
512  size_t oldRows = nrow;
513  nrow = newRows;
514  this->elem.resize(nrow*ncol,RETAIN);
515  T* oldMat = this->ptr() + oldRows*(ncol-1);
516  for (size_t c = ncol-1; c > 0; c--, oldMat -= oldRows)
517  {
518  memmove(this->ptr(c),oldMat,oldRows*sizeof(T));
519  for (size_t r = nrow-1; r >= oldRows; r--)
520  this->elem[r+nrow*(c-1)] = T(0);
521  }
522  }
523  return *this;
524  }
525 
527  bool augmentRows(const matrix<T>& B)
528  {
529  if (B.ncol != ncol)
530  return false;
531 
532  size_t oldRows = nrow;
533  nrow += B.nrow;
534  this->elem.resize(nrow*ncol,RETAIN);
535  T* oldMat = this->ptr() + oldRows*(ncol-1);
536  for (size_t c = ncol; c > 0; c--, oldMat -= oldRows)
537  {
538  if (c > 1)
539  memmove(this->ptr(c-1),oldMat,oldRows*sizeof(T));
540  for (size_t r = nrow; r > oldRows; r--)
541  this->elem[r-1+nrow*(c-1)] = B(r-oldRows,c);
542  }
543  return true;
544  }
545 
547  bool augmentCols(const matrix<T>& B)
548  {
549  if (B.nrow != nrow)
550  return false;
551 
552  this->elem.push_back(B.elem.begin(),B.elem.end());
553  ncol += B.ncol;
554  return true;
555  }
556 
558  size_t rows() const { return nrow; }
560  size_t cols() const { return ncol; }
561 
564  {
565  if (&A == this)
566  return *this;
567 
568  memcpy(this->n,A.n,sizeof(A.n));
569  this->elem = A.elem;
570  return *this;
571  }
572 
574  matrix<T>& operator=(const std::vector<T>& X)
575  {
576  // Do not use vector<T>::operator= because we don't want to alter size
577  size_t nval = X.size() < this->elem.size() ? X.size() : this->elem.size();
578  std::copy(X.begin(),X.begin()+nval,this->elem.begin());
579  std::fill(this->elem.begin()+nval,this->elem.end(),T(0));
580  return *this;
581  }
582 
585  T& operator()(size_t r, size_t c)
586  {
587  CHECK_INDEX("matrix::operator(): Row-index ",r,nrow);
588  CHECK_INDEX("matrix::operator(): Column-index ",c,ncol);
589  return this->elem[r-1+nrow*(c-1)];
590  }
591 
594  const T& operator()(size_t r, size_t c) const
595  {
596  CHECK_INDEX("matrix::operator(): Row-index ",r,nrow);
597  CHECK_INDEX("matrix::operator(): Column-index ",c,ncol);
598  return this->elem[r-1+nrow*(c-1)];
599  }
600 
602  vector<T> getRow(size_t r) const
603  {
604  CHECK_INDEX("matrix::getRow: Row-index ",r,nrow);
605  if (nrow < 2)
606  return this->elem;
607 
608  vector<T> row(ncol);
609  for (size_t i = 0; i < ncol; i++)
610  row[i] = this->elem[r-1+nrow*i];
611  return row;
612  }
613 
615  std::vector<T> getColumn(size_t c) const
616  {
617  CHECK_INDEX("matrix::getColumn: Column-index ",c,ncol);
618  if (ncol < 2)
619  return this->elem;
620 
621  std::vector<T> col(nrow);
622  memcpy(col.data(),this->ptr(c-1),nrow*sizeof(T));
623  return col;
624  }
625 
626  using matrixBase<T>::fill;
628  void fill(const std::vector<T>& v, size_t n, size_t m = 0)
629  {
630  if (n == 0 || v.size() < n)
631  return;
632  if (m == 0) m = v.size()/n;
633  this->resize(n,m,true);
634  if (n*m == v.size())
635  this->elem.fill(v.data());
636  else if ((n = v.size()/m) > nrow)
637  for (size_t c = 0; c < ncol; c++)
638  this->fillColumn(c+1,v.data()+c*n);
639  else // n < nrow
640  for (size_t c = 0; c < ncol; c++)
641  for (size_t r = 0; r < n; r++)
642  this->elem[r+c*nrow] = v[r+c*n];
643  }
644 
646  void fillColumn(size_t c, const std::vector<T>& data)
647  {
648  CHECK_INDEX("matrix::fillColumn: Column-index ",c,ncol);
649  size_t ndata = nrow > data.size() ? data.size() : nrow;
650  memcpy(this->ptr(c-1),data.data(),ndata*sizeof(T));
651  }
652 
654  void fillColumn(size_t c, const T* data)
655  {
656  CHECK_INDEX("matrix::fillColumn: Column-index ",c,ncol);
657  memcpy(this->ptr(c-1),data,nrow*sizeof(T));
658  }
659 
661  void fillRow(size_t r, const T* data)
662  {
663  CHECK_INDEX("matrix::fillRow: Row-index ",r,nrow);
664  if (nrow < 2)
665  this->elem.fill(data);
666  else for (size_t i = 0; i < ncol; i++)
667  this->elem[r-1+nrow*i] = data[i];
668  }
669 
671  void fillBlock(const matrix<T>& block, size_t r, size_t c,
672  bool transposed = false)
673  {
674  size_t nr = transposed ? block.cols() : block.rows();
675  size_t nc = transposed ? block.rows() : block.cols();
676  for (size_t i = 1; i <= nr && i+r-1 <= nrow; i++)
677  {
678  size_t ip = i+r-2 + nrow*(c-1);
679  for (size_t j = 1; j <= nc && j+c-1 <= ncol; j++, ip += nrow)
680  this->elem[ip] = transposed ? block(j,i) : block(i,j);
681  }
682  }
683 
685  void addBlock(const matrix<T>& block, T s, size_t r, size_t c,
686  bool transposed = false)
687  {
688  size_t nr = transposed ? block.cols() : block.rows();
689  size_t nc = transposed ? block.rows() : block.cols();
690  for (size_t i = 1; i <= nr && i+r-1 <= nrow; i++)
691  {
692  size_t ip = i+r-2 + nrow*(c-1);
693  for (size_t j = 1; j <= nc && j+c-1 <= ncol; j++, ip += nrow)
694  this->elem[ip] += s*(transposed ? block(j,i) : block(i,j));
695  }
696  }
697 
699  void extractBlock(matrix<T>& block, size_t r, size_t c,
700  bool addTo = false, bool transposed = false) const
701  {
702  size_t nr = transposed ? block.cols() : block.rows();
703  size_t nc = transposed ? block.rows() : block.cols();
704  for (size_t i = 1; i <= nr && i+r-1 <= nrow; i++)
705  {
706  size_t ip = i+r-2 + nrow*(c-1);
707  for (size_t j = 1; j <= nc && j+c-1 <= ncol; j++, ip += nrow)
708  if (addTo)
709  (transposed ? block(j,i) : block(i,j)) += this->elem[ip];
710  else
711  (transposed ? block(j,i) : block(i,j)) = this->elem[ip];
712  }
713  }
714 
716  matrix<T>& diag(T d, size_t dim = 0)
717  {
718  if (dim > 0)
719  this->resize(dim,dim,true);
720  else
721  this->resize(nrow,ncol,true);
722  for (size_t r = 0; r < nrow && r < ncol; r++)
723  this->elem[r+nrow*r] = d;
724  return *this;
725  }
726 
729  {
730  matrix<T> tmp(*this);
731  for (size_t r = 0; r < nrow; r++)
732  for (size_t c = 0; c < ncol; c++)
733  this->elem[c+ncol*r] = tmp.elem[r+nrow*c];
734 
735  nrow = tmp.ncol;
736  ncol = tmp.nrow;
737  return *this;
738  }
739 
741  T trace() const { return this->elem.sum(0,nrow+1); }
743  T rowsum(size_t r) const { return this->elem.sum(r-1,nrow); }
745  T colsum(size_t c) const { return this->elem.sum(nrow*(c-1),1,nrow*c); }
746 
747 #define THIS(i,j) this->operator()(i,j)
748 
750  T det() const
751  {
752  if (ncol == 1 && nrow >= 1)
753  return THIS(1,1);
754  else if (ncol == 2 && nrow >= 2)
755  return THIS(1,1)*THIS(2,2) - THIS(2,1)*THIS(1,2);
756  else if (ncol == 3 && nrow >= 3)
757  return THIS(1,1)*(THIS(2,2)*THIS(3,3) - THIS(3,2)*THIS(2,3))
758  - THIS(1,2)*(THIS(2,1)*THIS(3,3) - THIS(3,1)*THIS(2,3))
759  + THIS(1,3)*(THIS(2,1)*THIS(3,2) - THIS(3,1)*THIS(2,2));
760  else if (ncol > 0 && nrow > 0) {
761  std::cerr <<"matrix::det: Not available for "
762  << nrow <<"x"<< ncol <<" matrices"<< std::endl;
763  ABORT_ON_SINGULARITY;
764  return T(-999);
765  }
766  else
767  return T(0);
768  }
769 
773  T inverse(T tol = T(0))
774  {
775  T Det = this->det();
776  if (Det == T(-999))
777  return Det;
778  else if (Det <= tol && Det >= -tol) {
779  std::cerr <<"matrix::inverse: Singular matrix |A|="<< Det << std::endl;
780  ABORT_ON_SINGULARITY;
781  return T(0);
782  }
783 
784  if (ncol == 1)
785  THIS(1,1) = T(1) / Det;
786  else if (ncol == 2) {
787  matrix<T> B(2,2);
788  B(1,1) = THIS(2,2) / Det;
789  B(2,1) = -THIS(2,1) / Det;
790  B(1,2) = -THIS(1,2) / Det;
791  B(2,2) = THIS(1,1) / Det;
792  *this = B;
793  }
794  else if (ncol == 3) {
795  matrix<T> B(3,3);
796  B(1,1) = (THIS(2,2)*THIS(3,3) - THIS(3,2)*THIS(2,3)) / Det;
797  B(2,1) = -(THIS(2,1)*THIS(3,3) - THIS(3,1)*THIS(2,3)) / Det;
798  B(3,1) = (THIS(2,1)*THIS(3,2) - THIS(3,1)*THIS(2,2)) / Det;
799  B(1,2) = -(THIS(1,2)*THIS(3,3) - THIS(3,2)*THIS(1,3)) / Det;
800  B(2,2) = (THIS(1,1)*THIS(3,3) - THIS(3,1)*THIS(1,3)) / Det;
801  B(3,2) = -(THIS(1,1)*THIS(3,2) - THIS(3,1)*THIS(1,2)) / Det;
802  B(1,3) = (THIS(1,2)*THIS(2,3) - THIS(2,2)*THIS(1,3)) / Det;
803  B(2,3) = -(THIS(1,1)*THIS(2,3) - THIS(2,1)*THIS(1,3)) / Det;
804  B(3,3) = (THIS(1,1)*THIS(2,2) - THIS(2,1)*THIS(1,2)) / Det;
805  *this = B;
806  }
807 
808  return Det;
809  }
810 
813  bool isSymmetric(T tol = T(0)) const
814  {
815  if (nrow != ncol)
816  return false;
817 
818  for (size_t r = 0; r < nrow; r++)
819  for (size_t c = 0; c < r; c++)
820  {
821  T diff = this->elem[r+nrow*c] - this->elem[c+nrow*r];
822  if (diff < -tol || diff > tol)
823  return false;
824  }
825 
826  return true;
827  }
828 
830  matrix<T>& operator+=(const matrix<T>& A) { return this->add(A); }
832  matrix<T>& operator-=(const matrix<T>& A) { return this->add(A,T(-1)); }
834  matrix<T>& add(const matrix<T>& A, T alfa = T(1))
835  {
836  return static_cast<matrix<T>&>(this->matrixBase<T>::add(A,alfa));
837  }
838 
840  matrix<T>& operator*=(T c) { return this->multiply(c); }
842  matrix<T>& operator/=(T d) { return this->multiply(T(1)/d); }
845  {
846  return static_cast<matrix<T>&>(this->matrixBase<T>::multiply(c));
847  }
848 
860  matrix<T>& multiply(const matrix<T>& A, const matrix<T>& B,
861  bool transA = false, bool transB = false,
862  bool addTo = false, const T& alpha = T(1));
863 
876  bool multiplyMat(const matrix<T>& A, const std::vector<T>& B,
877  bool transA = false, bool addTo = false);
878 
891  bool multiplyMat(const std::vector<T>& A, const matrix<T>& B,
892  bool transB = false, bool addTo = false);
893 
903  bool multiply(const std::vector<T>& X, std::vector<T>& Y,
904  bool transA = false, char addTo = 0) const;
905 
911  bool multiply(const std::vector<T>& X, std::vector<T>& Y,
912  const T& alpha, const T& beta = T(0),
913  bool transA = false, int stridex = 1, int stridey = 1,
914  unsigned int ofsx = 0, unsigned int ofsy = 0) const;
915 
917  bool outer_product(const std::vector<T>& X, const std::vector<T>& Y,
918  bool addTo = false, T alpha = T(1));
919 
921  T normInf() const
922  {
923  if (nrow == 0)
924  return T(0);
925 
926  // Compute row sums
927  vector<T> sums(nrow);
928  for (size_t i = 0; i < nrow; i++)
929  sums[i] = this->elem.asum(i,nrow);
930  return *std::max_element(sums.begin(),sums.end());
931  }
932 
933  private:
935  bool compatible(const std::vector<T>& X, bool transA) const
936  {
937  if (nrow > 0 && ncol > 0)
938  if ((transA ? nrow : ncol) == X.size())
939  return true;
940 
941  std::cerr <<"matrix::multiply: Incompatible matrices: A("
942  << nrow <<','<< ncol <<"), X("<< X.size() <<")\n"
943  <<" when computing Y = "
944  << (transA ? "A^t":"A") <<" * X"<< std::endl;
945  ABORT_ON_INDEX_CHECK;
946  return false;
947  }
948 
950  bool compatible(const matrix<T>& A, const matrix<T>& B,
951  bool transA, bool transB, size_t& M, size_t& N, size_t& K)
952  {
953  M = transA ? A.ncol : A.nrow;
954  N = transB ? B.nrow : B.ncol;
955  K = transA ? A.nrow : A.ncol;
956  if (K == (transB ? B.ncol : B.nrow))
957  return true;
958 
959  std::cerr <<"matrix::multiply: Incompatible matrices: A("
960  << A.nrow <<','<< A.ncol <<"), B("
961  << B.nrow <<','<< B.ncol <<")\n"
962  <<" when computing C = "
963  << (transA ? "A^t":"A") <<" * "
964  << (transB ? "B^t":"B") << std::endl;
965  ABORT_ON_INDEX_CHECK;
966  return false;
967  }
968 
971  bool compatible(const matrix<T>& A, const std::vector<T>& B,
972  bool transA, size_t& M, size_t& N, size_t& K)
973  {
974  M = transA ? A.ncol : A.nrow;
975  K = transA ? A.nrow : A.ncol;
976  N = K > 0 ? B.size()/K : 0;
977  if (N*K == B.size() && !B.empty())
978  return true;
979 
980  std::cerr <<"matrix::multiply: Incompatible matrices: A("
981  << A.nrow <<','<< A.ncol <<"), B(r*c="<< B.size() <<")\n"
982  <<" when computing C = "
983  << (transA ? "A^t":"A") <<" * B"<< std::endl;
984  ABORT_ON_INDEX_CHECK;
985  return false;
986  }
987 
990  bool compatible(const std::vector<T>& A, const matrix<T>& B,
991  bool transB, size_t& M, size_t& N, size_t& K)
992  {
993  N = transB ? B.nrow : B.ncol;
994  K = transB ? B.ncol : B.nrow;
995  M = K > 0 ? A.size() / K : 0;
996  if (M*K == A.size() && !A.empty())
997  return true;
998 
999  std::cerr <<"matrix::multiply: Incompatible matrices: A(r*c="<< A.size()
1000  <<"), B("<< B.nrow <<","<< B.ncol <<")\n"
1001  <<" when computing C = A * "
1002  << (transB ? "B^t":"B") << std::endl;
1003  ABORT_ON_INDEX_CHECK;
1004  return false;
1005  }
1006 
1008  bool compatible(const std::vector<T>& X, const std::vector<T>& Y)
1009  {
1010  if (X.size() == nrow && Y.size() == ncol)
1011  return true;
1012 
1013  std::cerr <<"matrix::outer_product: Incompatible matrix and vectors: A("
1014  << nrow <<','<< ncol <<"), X("
1015  << X.size() <<"), Y("<< Y.size() <<")\n"
1016  <<" when computing A += X*Y^t"
1017  << std::endl;
1018  ABORT_ON_INDEX_CHECK;
1019  return false;
1020  }
1021 
1022  protected:
1024  void clearIfNrowChanged(size_t n1, size_t, size_t) override
1025  {
1026  if (n1 != nrow) this->elem.clear();
1027  }
1028 
1029  private:
1030  size_t& nrow;
1031  size_t& ncol;
1032  };
1033 
1034 
1035 #ifdef HAS_BLAS
1036  //============================================================================
1037  //=== BLAS-implementation of the matrix/vector multiplication methods ====
1038  //============================================================================
1039 
1040  template<> inline
1042  {
1043  cblas_sscal(myVec.size(),c,myVec.data(),1);
1044  return *this;
1045  }
1046 
1047  template<> inline
1048  vector<double>& vector<double>::operator*=(double c)
1049  {
1050  cblas_dscal(myVec.size(),c,myVec.data(),1);
1051  return *this;
1052  }
1053 
1054  template<> inline
1055  float vector<float>::dot(const float* v, size_t nv,
1056  size_t o1, int i1, size_t o2, int i2) const
1057  {
1058  int n1 = i1 > 1 || i1 < -1 ? myVec.size()/abs(i1) : myVec.size()-o1;
1059  int n2 = i2 > 1 || i2 < -1 ? nv/abs(i2) : nv-o2;
1060  int n = n1 < n2 ? n1 : n2;
1061  return cblas_sdot(n,myVec.data()+o1,i1,v+o2,i2);
1062  }
1063 
1064  template<> inline
1065  double vector<double>::dot(const double* v, size_t nv,
1066  size_t o1, int i1, size_t o2, int i2) const
1067  {
1068  int n1 = i1 > 1 || i1 < -1 ? myVec.size()/abs(i1) : myVec.size()-o1;
1069  int n2 = i2 > 1 || i2 < -1 ? nv/abs(i2) : nv-o2;
1070  int n = n1 < n2 ? n1 : n2;
1071  return cblas_ddot(n,myVec.data()+o1,i1,v+o2,i2);
1072  }
1073 
1074  template<> inline
1075  float vector<float>::norm2(size_t off, int inc) const
1076  {
1077  int n = inc > 1 || inc < -1 ? myVec.size()/abs(inc) : myVec.size()-off;
1078  return cblas_snrm2(n,myVec.data()+off,inc);
1079  }
1080 
1081  template<> inline
1082  double vector<double>::norm2(size_t off, int inc) const
1083  {
1084  int n = inc > 1 || inc < -1 ? myVec.size()/abs(inc) : myVec.size()-off;
1085  return cblas_dnrm2(n,myVec.data()+off,inc);
1086  }
1087 
1088  template<> inline
1089  float vector<float>::normInf(size_t& off, int inc, bool sign) const
1090  {
1091  if (inc < 1 || myVec.empty())
1092  return 0.0f;
1093 
1094  const float* v = myVec.data() + off;
1095  off = 1 + cblas_isamax(myVec.size()/inc,v,inc);
1096  return sign ? v[(off-1)*inc] : fabsf(v[(off-1)*inc]);
1097  }
1098 
1099  template<> inline
1100  double vector<double>::normInf(size_t& off, int inc, bool sign) const
1101  {
1102  if (inc < 1 || myVec.empty())
1103  return 0.0;
1104 
1105  const double* v = myVec.data() + off;
1106  off = 1 + cblas_idamax(myVec.size()/inc,v,inc);
1107  return sign ? v[(off-1)*inc] : fabs(v[(off-1)*inc]);
1108  }
1109 
1110  template<> inline
1111  float vector<float>::asum(size_t off, int inc) const
1112  {
1113  int n = inc > 1 || inc < -1 ? myVec.size()/abs(inc) : myVec.size()-off;
1114  return cblas_sasum(n,myVec.data()+off,inc);
1115  }
1116 
1117  template<> inline
1118  double vector<double>::asum(size_t off, int inc) const
1119  {
1120  int n = inc > 1 || inc < -1 ? myVec.size()/abs(inc) : myVec.size()-off;
1121  return cblas_dasum(n,myVec.data()+off,inc);
1122  }
1123 
1124  template<> inline
1125  vector<float>& vector<float>::add(const std::vector<float>& X,
1126  const float& alfa,
1127  unsigned int ofsx, int stridex,
1128  unsigned int ofsy, int stridey)
1129  {
1130  if (myVec.empty() && stridex > 0 && stridey > 0)
1131  myVec.resize(ofsy+stridey*(X.size()-ofsx)/stridex);
1132 
1133  int nx = stridex == 0 ? 1 : 1 + (X.size()-ofsx-1)/abs(stridex);
1134  int ny = stridey == 0 ? 1 : 1 + (myVec.size()-ofsy-1)/abs(stridey);
1135  int n = nx < ny ? (stridex == 0 ? ny : nx) : (stridey == 0 ? nx : ny);
1136  if (n > 0)
1137  cblas_saxpy(n,alfa,X.data()+ofsx,stridex,myVec.data()+ofsy,stridey);
1138  return *this;
1139  }
1140 
1141  template<> inline
1142  vector<double>& vector<double>::add(const std::vector<double>& X,
1143  const double& alfa,
1144  unsigned int ofsx, int stridex,
1145  unsigned int ofsy, int stridey)
1146  {
1147  if (myVec.empty() && stridex > 0 && stridey > 0)
1148  myVec.resize(ofsy+stridey*(X.size()-ofsx)/stridex);
1149 
1150  int nx = stridex == 0 ? 1 : 1 + (X.size()-ofsx-1)/abs(stridex);
1151  int ny = stridey == 0 ? 1 : 1 + (myVec.size()-ofsy-1)/abs(stridey);
1152  int n = nx < ny ? (stridex == 0 ? ny : nx) : (stridey == 0 ? nx : ny);
1153  if (n > 0)
1154  cblas_daxpy(n,alfa,X.data()+ofsx,stridex,myVec.data()+ofsy,stridey);
1155  return *this;
1156  }
1157 
1158  template<> inline
1159  matrixBase<float>& matrixBase<float>::add(const matrixBase<float>& A,
1160  const float& alfa)
1161  {
1162  int n = this->size() < A.size() ? this->size() : A.size();
1163  if (n > 0)
1164  cblas_saxpy(n,alfa,A.ptr(),1,this->ptr(),1);
1165  return *this;
1166  }
1167 
1168  template<> inline
1169  matrixBase<double>& matrixBase<double>::add(const matrixBase<double>& A,
1170  const double& alfa)
1171  {
1172  int n = this->size() < A.size() ? this->size() : A.size();
1173  if (n > 0)
1174  cblas_daxpy(n,alfa,A.ptr(),1,this->ptr(),1);
1175  return *this;
1176  }
1177 
1178  template<> inline
1179  matrixBase<float>& matrixBase<float>::multiply(const float& c)
1180  {
1181  cblas_sscal(this->size(),c,this->ptr(),1);
1182  return *this;
1183  }
1184 
1185  template<> inline
1186  matrixBase<double>& matrixBase<double>::multiply(const double& c)
1187  {
1188  cblas_dscal(this->size(),c,this->ptr(),1);
1189  return *this;
1190  }
1191 
1192  template<> inline
1193  bool matrix<float>::multiply(const std::vector<float>& X,
1194  std::vector<float>& Y,
1195  bool transA, char addTo) const
1196  {
1197  if (!this->compatible(X,transA))
1198  return false;
1199  else if (!addTo || Y.empty())
1200  {
1201  Y.resize(transA ? ncol : nrow);
1202  if (addTo) std::fill(Y.begin(),Y.end(),0.0f);
1203  }
1204 
1205  cblas_sgemv(CblasColMajor,
1206  transA ? CblasTrans : CblasNoTrans,
1207  nrow, ncol, addTo < 0 ? -1.0f : 1.0f,
1208  this->ptr(), nrow,
1209  X.data(), 1, addTo ? 1.0f : 0.0f,
1210  Y.data(), 1);
1211 
1212  return true;
1213  }
1214 
1215  template<> inline
1216  bool matrix<double>::multiply(const std::vector<double>& X,
1217  std::vector<double>& Y,
1218  bool transA, char addTo) const
1219  {
1220  if (!this->compatible(X,transA))
1221  return false;
1222  else if (!addTo || Y.empty())
1223  {
1224  Y.resize(transA ? ncol : nrow);
1225  if (addTo) std::fill(Y.begin(),Y.end(),0.0);
1226  }
1227 
1228  cblas_dgemv(CblasColMajor,
1229  transA ? CblasTrans : CblasNoTrans,
1230  nrow, ncol, addTo < 0 ? -1.0 : 1.0,
1231  this->ptr(), nrow,
1232  X.data(), 1, addTo ? 1.0 : 0.0,
1233  Y.data(), 1);
1234 
1235  return true;
1236  }
1237 
1238  template<> inline
1239  bool matrix<float>::multiply(const std::vector<float>& X,
1240  std::vector<float>& Y,
1241  const float& alpha, const float& beta,
1242  bool transA, int stridex, int stridey,
1243  unsigned int ofsx, unsigned int ofsy) const
1244  {
1245  if (stridex == 0 || stridey == 0)
1246  {
1247  std::cerr <<"matrix::multiply: Stride must be non-zero ("
1248  << stridex <<", "<< stridey <<")"<< std::endl;
1249  ABORT_ON_INDEX_CHECK;
1250  return false;
1251  }
1252 
1253  if (ofsx == 0 && stridex == 1 && !this->compatible(X,transA))
1254  return false;
1255  else if (beta == 0.0f || Y.empty())
1256  {
1257  Y.resize(ofsy + 1 + ((transA ? ncol : nrow)-1)*abs(stridey));
1258  if (beta != 0.0f) std::fill(Y.begin(),Y.end(),0.0f);
1259  }
1260 
1261  cblas_sgemv(CblasColMajor,
1262  transA ? CblasTrans : CblasNoTrans,
1263  nrow, ncol, alpha,
1264  this->ptr(), nrow,
1265  X.data()+ofsx, stridex, beta,
1266  Y.data()+ofsy, stridey);
1267 
1268  return true;
1269  }
1270 
1271  template<> inline
1272  bool matrix<double>::multiply(const std::vector<double>& X,
1273  std::vector<double>& Y,
1274  const double& alpha, const double& beta,
1275  bool transA, int stridex, int stridey,
1276  unsigned int ofsx, unsigned int ofsy) const
1277  {
1278  if (stridex == 0 || stridey == 0)
1279  {
1280  std::cerr <<"matrix::multiply: Stride must be non-zero ("
1281  << stridex <<", "<< stridey <<")"<< std::endl;
1282  ABORT_ON_INDEX_CHECK;
1283  return false;
1284  }
1285 
1286  if (ofsx == 0 && stridex == 1 && !this->compatible(X,transA))
1287  return false;
1288  else if (beta == 0.0 || Y.empty())
1289  {
1290  Y.resize(ofsy + 1 + ((transA ? ncol : nrow)-1)*abs(stridey));
1291  if (beta != 0.0) std::fill(Y.begin(),Y.end(),0.0);
1292  }
1293 
1294  cblas_dgemv(CblasColMajor,
1295  transA ? CblasTrans : CblasNoTrans,
1296  nrow, ncol, alpha,
1297  this->ptr(), nrow,
1298  X.data()+ofsx, stridex, beta,
1299  Y.data()+ofsy, stridey);
1300 
1301  return true;
1302  }
1303 
1304  template<> inline
1305  matrix<float>& matrix<float>::multiply(const matrix<float>& A,
1306  const matrix<float>& B,
1307  bool transA, bool transB,
1308  bool addTo, const float& alpha)
1309  {
1310  size_t M, N, K;
1311  if (!this->compatible(A,B,transA,transB,M,N,K))
1312  {
1313  this->clear();
1314  return *this;
1315  }
1316  else if (!addTo || this->empty())
1317  this->resize(M,N);
1318 
1319  cblas_sgemm(CblasColMajor,
1320  transA ? CblasTrans : CblasNoTrans,
1321  transB ? CblasTrans : CblasNoTrans,
1322  M, N, K, alpha,
1323  A.ptr(), A.nrow,
1324  B.ptr(), B.nrow,
1325  addTo ? 1.0f : 0.0f,
1326  this->ptr(), nrow);
1327 
1328  return *this;
1329  }
1330 
1331  template<> inline
1332  matrix<double>& matrix<double>::multiply(const matrix<double>& A,
1333  const matrix<double>& B,
1334  bool transA, bool transB,
1335  bool addTo, const double& alpha)
1336  {
1337  size_t M, N, K;
1338  if (!this->compatible(A,B,transA,transB,M,N,K))
1339  {
1340  this->clear();
1341  return *this;
1342  }
1343  else if (!addTo || this->empty())
1344  this->resize(M,N);
1345 
1346  cblas_dgemm(CblasColMajor,
1347  transA ? CblasTrans : CblasNoTrans,
1348  transB ? CblasTrans : CblasNoTrans,
1349  M, N, K, alpha,
1350  A.ptr(), A.nrow,
1351  B.ptr(), B.nrow,
1352  addTo ? 1.0 : 0.0,
1353  this->ptr(), nrow);
1354 
1355  return *this;
1356  }
1357 
1358  template<> inline
1359  bool matrix<float>::multiplyMat(const matrix<float>& A,
1360  const std::vector<float>& B,
1361  bool transA, bool addTo)
1362  {
1363  size_t M, N, K;
1364  if (!this->compatible(A,B,transA,M,N,K))
1365  return false;
1366  else if (!addTo || this->empty())
1367  this->resize(M,N);
1368 
1369  cblas_sgemm(CblasColMajor,
1370  transA ? CblasTrans : CblasNoTrans, CblasNoTrans,
1371  M, N, K, 1.0f,
1372  A.ptr(), A.nrow,
1373  B.data(), K,
1374  addTo ? 1.0f : 0.0f,
1375  this->ptr(), nrow);
1376 
1377  return true;
1378  }
1379 
1380  template<> inline
1381  bool matrix<double>::multiplyMat(const matrix<double>& A,
1382  const std::vector<double>& B,
1383  bool transA, bool addTo)
1384  {
1385  size_t M, N, K;
1386  if (!this->compatible(A,B,transA,M,N,K))
1387  return false;
1388  else if (!addTo || this->empty())
1389  this->resize(M,N);
1390 
1391  cblas_dgemm(CblasColMajor,
1392  transA ? CblasTrans : CblasNoTrans, CblasNoTrans,
1393  M, N, K, 1.0,
1394  A.ptr(), A.nrow,
1395  B.data(), K,
1396  addTo ? 1.0 : 0.0,
1397  this->ptr(), nrow);
1398 
1399  return true;
1400  }
1401 
1402  template<> inline
1403  bool matrix<float>::multiplyMat(const std::vector<float>& A,
1404  const matrix<float>& B,
1405  bool transB, bool addTo)
1406  {
1407  size_t M, N, K;
1408  if (!this->compatible(A,B,transB,M,N,K))
1409  return false;
1410  else if (!addTo || this->empty())
1411  this->resize(M,N);
1412 
1413  cblas_sgemm(CblasColMajor,
1414  CblasNoTrans, transB ? CblasTrans : CblasNoTrans,
1415  M, N, K, 1.0f,
1416  A.data(), M,
1417  B.ptr(), B.nrow,
1418  addTo ? 1.0f : 0.0f,
1419  this->ptr(), nrow);
1420 
1421  return true;
1422  }
1423 
1424  template<> inline
1425  bool matrix<double>::multiplyMat(const std::vector<double>& A,
1426  const matrix<double>& B,
1427  bool transB, bool addTo)
1428  {
1429  size_t M, N, K;
1430  if (!this->compatible(A,B,transB,M,N,K))
1431  return false;
1432  else if (!addTo || this->empty())
1433  this->resize(M,N);
1434 
1435  cblas_dgemm(CblasColMajor,
1436  CblasNoTrans, transB ? CblasTrans : CblasNoTrans,
1437  M, N, K, 1.0,
1438  A.data(), M,
1439  B.ptr(), B.nrow,
1440  addTo ? 1.0 : 0.0,
1441  this->ptr(), nrow);
1442 
1443  return true;
1444  }
1445 
1446  template<> inline
1447  bool matrix<float>::outer_product(const std::vector<float>& X,
1448  const std::vector<float>& Y,
1449  bool addTo, float alpha)
1450  {
1451  if (!addTo)
1452  this->resize(X.size(),Y.size());
1453  else if (!this->compatible(X,Y))
1454  return false;
1455 
1456  cblas_sgemm(CblasColMajor,
1457  CblasNoTrans, CblasTrans,
1458  nrow, ncol, 1, alpha,
1459  X.data(), nrow,
1460  Y.data(), ncol,
1461  addTo ? 1.0f : 0.0f,
1462  this->ptr(), nrow);
1463 
1464  return true;
1465  }
1466 
1467  template<> inline
1468  bool matrix<double>::outer_product(const std::vector<double>& X,
1469  const std::vector<double>& Y,
1470  bool addTo, double alpha)
1471  {
1472  if (!addTo)
1473  this->resize(X.size(),Y.size());
1474  else if (!this->compatible(X,Y))
1475  return false;
1476 
1477  cblas_dgemm(CblasColMajor,
1478  CblasNoTrans, CblasTrans,
1479  nrow, ncol, 1, alpha,
1480  X.data(), nrow,
1481  Y.data(), ncol,
1482  addTo ? 1.0 : 0.0,
1483  this->ptr(), nrow);
1484 
1485  return true;
1486  }
1487 
1488 #else
1489  //============================================================================
1490  //=== Non-BLAS inlined implementations (slow...) =========================
1491  //============================================================================
1492 
1493  template<class T> inline
1495  {
1496  for (T& x : myVec)
1497  x *= c;
1498  return *this;
1499  }
1500 
1501  template<class T> inline
1502  T vector<T>::dot(const T* v, size_t nv,
1503  size_t o1, int i1, size_t o2, int i2) const
1504  {
1505  size_t i, j;
1506  T dotprod = T(0);
1507  for (i = o1, j = o2; i < myVec.size() && j < nv; i += i1, j += i2)
1508  dotprod += myVec[i] * v[j];
1509  return dotprod;
1510  }
1511 
1512  template<class T> inline
1513  T vector<T>::norm2(size_t off, int inc) const
1514  {
1515  double xsum = 0.0;
1516  if (inc < 1 || myVec.size() <= off)
1517  return xsum;
1518 
1519  // Warning: This might overflow or underflow for large/small values
1520  for (size_t i = off; i < myVec.size(); i += inc)
1521  xsum += myVec[i]*myVec[i];
1522  return sqrt(xsum);
1523  }
1524 
1525  template<class T> inline
1526  T vector<T>::normInf(size_t& off, int inc, bool sign) const
1527  {
1528  T xmax = T(0);
1529  if (inc < 1 || myVec.size() <= off)
1530  return xmax;
1531 
1532  T amax = T(0);
1533  for (size_t i = off; i < myVec.size(); i += inc)
1534  if (myVec[i] > amax)
1535  {
1536  off = 1+i/inc;
1537  xmax = amax = myVec[i];
1538  }
1539  else if (myVec[i] < -amax)
1540  {
1541  off = 1+i/inc;
1542  xmax = myVec[i];
1543  amax = -xmax;
1544  }
1545 
1546  return sign ? xmax : amax;
1547  }
1548 
1549  template<class T> inline
1550  T vector<T>::asum(size_t off, int inc) const
1551  {
1552  T xsum = T(0);
1553  if (inc < 1 || myVec.size() <= off)
1554  return xsum;
1555 
1556  for (size_t i = off; i < myVec.size(); i += inc)
1557  xsum += myVec[i] < T(0) ? -myVec[i] : myVec[i];
1558  return xsum;
1559  }
1560 
1561  template<class T> inline
1562  vector<T>& vector<T>::add(const std::vector<T>& X, const T& alfa,
1563  unsigned int ofsx, int stridex,
1564  unsigned int ofsy, int stridey)
1565  {
1566  if (stridex < 0 || stridey < 0 || stridex+stridey == 0)
1567  {
1568  std::cerr <<"vector::add: Negative stride not supported ("
1569  << stridex <<", "<< stridey <<")"<< std::endl;
1570  ABORT_ON_INDEX_CHECK;
1571  return *this;
1572  }
1573 
1574  std::vector<T>& Y = myVec;
1575  if (Y.empty() && stridex > 0)
1576  Y.resize(ofsy+stridey*(X.size()-ofsx)/stridex);
1577 
1578  for (; ofsx < X.size() && ofsy < Y.size(); ofsx += stridex, ofsy += stridey)
1579  Y[ofsy] += alfa*X[ofsx];
1580  return *this;
1581  }
1582 
1583  template<class T> inline
1585  {
1586  const vector<T>& X = A.elem;
1587  vector<T>& Y = this->elem;
1588 
1589  for (size_t i = 0; i < X.size() && i < Y.size(); i++)
1590  Y[i] += alfa*X[i];
1591  return *this;
1592  }
1593 
1594  template<class T> inline
1596  {
1597  for (T& x : this->elem)
1598  x *= c;
1599  return *this;
1600  }
1601 
1602  template<class T> inline
1603  bool matrix<T>::multiply(const std::vector<T>& X, std::vector<T>& Y,
1604  bool transA, char addTo) const
1605  {
1606  if (!this->compatible(X,transA))
1607  return false;
1608  else if (!addTo || Y.empty())
1609  {
1610  Y.resize(transA ? ncol : nrow);
1611  std::fill(Y.begin(),Y.end(),T(0));
1612  }
1613 
1614  for (size_t i = 0; i < Y.size(); i++)
1615  for (size_t j = 0; j < X.size(); j++)
1616  if (transA)
1617  Y[i] += THIS(j+1,i+1) * (addTo < 0 ? -X[j] : X[j]);
1618  else
1619  Y[i] += THIS(i+1,j+1) * (addTo < 0 ? -X[j] : X[j]);
1620 
1621  return true;
1622  }
1623 
1624  template<class T> inline
1625  bool matrix<T>::multiply(const std::vector<T>& X, std::vector<T>& Y,
1626  const T& alpha, const T& beta,
1627  bool transA, int stridex, int stridey,
1628  unsigned int ofsx, unsigned int ofsy) const
1629  {
1630  if (stridex <= 0 || stridey <= 0)
1631  {
1632  std::cerr <<"matrix::multiply: Non-positive stride not supported ("
1633  << stridex <<", "<< stridey <<")"<< std::endl;
1634  ABORT_ON_INDEX_CHECK;
1635  return false;
1636  }
1637 
1638  if (ofsx == 0 && stridex == 1 && !this->compatible(X,transA))
1639  return false;
1640  else if (beta == T(0) || Y.empty())
1641  {
1642  Y.resize(ofsy + 1 + ((transA ? ncol : nrow)-1)*stridey);
1643  std::fill(Y.begin(),Y.end(),T(0));
1644  }
1645  else if (beta != T(1))
1646  for (size_t i = ofsy; i < Y.size(); i += stridey)
1647  Y[i] *= beta;
1648 
1649  size_t a, b, i, j;
1650  for (a = 1, i = ofsy; i < Y.size(); a++, i += stridey)
1651  for (b = 1, j = ofsx; j < X.size(); b++, j += stridex)
1652  if (transA)
1653  Y[i] += alpha * THIS(b,a) * X[j];
1654  else
1655  Y[i] += alpha * THIS(a,b) * X[j];
1656 
1657  return true;
1658  }
1659 
1660  template<class T> inline
1662  const matrix<T>& B,
1663  bool transA, bool transB, bool addTo,
1664  const T& alpha)
1665  {
1666  size_t M, N, K;
1667  if (!this->compatible(A,B,transA,transB,M,N,K))
1668  {
1669  this->clear();
1670  return *this;
1671  }
1672  else if (!addTo || this->empty())
1673  this->resize(M,N,true);
1674 
1675  for (size_t i = 1; i <= M; i++)
1676  for (size_t j = 1; j <= N; j++)
1677  for (size_t k = 1; k <= K; k++)
1678  if (transA && transB)
1679  THIS(i,j) += alpha*A(k,i)*B(j,k);
1680  else if (transA)
1681  THIS(i,j) += alpha*A(k,i)*B(k,j);
1682  else if (transB)
1683  THIS(i,j) += alpha*A(i,k)*B(j,k);
1684  else
1685  THIS(i,j) += alpha*A(i,k)*B(k,j);
1686 
1687  return *this;
1688  }
1689 
1690  template<class T> inline
1691  bool matrix<T>::multiplyMat(const matrix<T>& A, const std::vector<T>& B,
1692  bool transA, bool addTo)
1693  {
1694  size_t M, N, K;
1695  if (!this->compatible(A,B,transA,M,N,K))
1696  return false;
1697  else if (!addTo || this->empty())
1698  this->resize(M,N,true);
1699 
1700  for (size_t i = 1; i <= M; i++)
1701  for (size_t j = 1; j <= N; j++)
1702  for (size_t k = 1; k <= K; k++)
1703  if (transA)
1704  THIS(i,j) += A(k,i)*B[k-1+K*(j-1)];
1705  else
1706  THIS(i,j) += A(i,k)*B[k-1+K*(j-1)];
1707 
1708  return true;
1709  }
1710 
1711  template<class T> inline
1712  bool matrix<T>::multiplyMat(const std::vector<T>& A, const matrix<T>& B,
1713  bool transB, bool addTo)
1714  {
1715  size_t M, N, K;
1716  if (!this->compatible(A,B,transB,M,N,K))
1717  return false;
1718  else if (!addTo || this->empty())
1719  this->resize(M,N,true);
1720 
1721  for (size_t i = 1; i <= M; i++)
1722  for (size_t j = 1; j <= N; j++)
1723  for (size_t k = 1; k <= K; k++)
1724  if (transB)
1725  THIS(i,j) += A[i-1+M*(k-1)]*B(j,k);
1726  else
1727  THIS(i,j) += A[i-1+M*(k-1)]*B(k,j);
1728 
1729  return true;
1730  }
1731 
1732  template<class T> inline
1733  bool matrix<T>::outer_product(const std::vector<T>& X,
1734  const std::vector<T>& Y,
1735  bool addTo, T alpha)
1736  {
1737  if (!addTo)
1738  this->resize(X.size(),Y.size());
1739  else if (!this->compatible(X,Y))
1740  return false;
1741 
1742  if (addTo)
1743  for (size_t j = 0; j < ncol; j++)
1744  for (size_t i = 0; i < nrow; i++)
1745  this->elem[i+nrow*j] += alpha*X[i]*Y[j];
1746  else
1747  for (size_t j = 0; j < ncol; j++)
1748  for (size_t i = 0; i < nrow; i++)
1749  this->elem[i+nrow*j] = alpha*X[i]*Y[j];
1750 
1751  return true;
1752  }
1753 
1754 #endif
1755 
1756  //============================================================================
1757  //=== Global operators ===================================================
1758  //============================================================================
1759 
1765  template<class T> inline T trunc(T v)
1766  {
1767  return v > T(zero_print_tol) || v < T(-zero_print_tol) || std::isnan(v) ?
1768  v : T(0);
1769  }
1770 
1772  template<class T> std::istream& operator>>(std::istream& s, vector<T>& X)
1773  {
1774  size_t n = 0;
1775  s >> n;
1776  X.resize(n,true);
1777  for (T& val : X)
1778  s >> val;
1779  return s;
1780  }
1781 
1783  template<class T> std::ostream& operator<<(std::ostream& s,
1784  const vector<T>& X)
1785  {
1786  if (X.size() < 1)
1787  s <<" (empty)";
1788  else for (size_t i = 0; i < X.size(); i++)
1789  s << ((i%nval_per_line) ? ' ':'\n') << trunc(X[i]);
1790 
1791  return s << std::endl;
1792  }
1793 
1795  template<class T> std::istream& operator>>(std::istream& s, matrix<T>& A)
1796  {
1797  size_t m = 0, n = 0;
1798  char c = 0;
1799  while (s.get(c) && isspace(c));
1800  bool symmetric = (c == 'S' || c == 's');
1801  bool columnori = (c == 'C' || c == 'c');
1802  if (symmetric)
1803  {
1804  s.ignore(10,':');
1805  s >> m;
1806  n = m;
1807  }
1808  else if (isalpha(c))
1809  {
1810  s.ignore(15,' ');
1811  s >> m >> n;
1812  }
1813  else
1814  {
1815  s.putback(c);
1816  s >> m >> n;
1817  }
1818  A.resize(m,n);
1819  for (size_t i = 1; i <= m; i++)
1820  {
1821  while (s.get(c) && isspace(c));
1822  if (c == 'R')
1823  s.ignore(10,':');
1824  else
1825  s.putback(c);
1826  for (size_t j = (symmetric ? i : 1); j <= n; j++)
1827  {
1828  s >> (columnori ? A(j,i) : A(i,j));
1829  if (symmetric && j > i)
1830  A(j,i) = A(i,j);
1831  }
1832  }
1833  return s;
1834  }
1835 
1841  template<class T> std::ostream& operator<<(std::ostream& s,
1842  const matrix<T>& A)
1843  {
1844  if (A.rows() < 1 || A.cols() < 1)
1845  return s <<" (empty)"<< std::endl;
1846 
1847  bool symm = A.isSymmetric(zero_print_tol);
1848  for (size_t i = 1; i <= A.rows(); i++)
1849  {
1850  size_t c1 = symm ? i : 1;
1851  s <<"\nRow "<< i <<": "<< trunc(A(i,c1));
1852  for (size_t j = c1+1; j <= A.cols(); j++)
1853  s <<' '<< trunc(A(i,j));
1854  }
1855 
1856  return s << std::endl;
1857  }
1858 
1860  template<class T> void writeMatlab(const char* label, const vector<T>& X,
1861  std::ostream& s = std::cout)
1862  {
1863  if (label)
1864  s << label <<" = [";
1865  else
1866  s <<"[";
1867 
1868  for (size_t i = 1; i <= X.size(); i++)
1869  s <<' '<< trunc(X(i));
1870  s <<" ];"<< std::endl;
1871  }
1872 
1874  template<class T> void writeMatlab(const char* label, const matrix<T>& A,
1875  std::ostream& s = std::cout)
1876  {
1877  if (label)
1878  s << label <<" = [";
1879  else
1880  s <<"[";
1881 
1882  size_t nsp = label ? 4 + strlen(label) : 1;
1883  for (size_t i = 1; i <= A.rows(); i++)
1884  {
1885  if (i > 1)
1886  {
1887  s <<";\n";
1888  for (size_t k = 0; k < nsp; k++) s <<' ';
1889  }
1890  for (size_t j = 1; j <= A.cols(); j++)
1891  s <<' '<< trunc(A(i,j));
1892  }
1893  s <<" ];"<< std::endl;
1894  }
1895 }
1896 
1897 #undef THIS
1898 #endif
BLAS support for various platforms.
static SystemMatrix * K
Pointer to coefficient matrix A.
Definition: EigSolver.C:91
static SystemMatrix * M
Pointer to coefficient matrix B.
Definition: EigSolver.C:92
virtual bool empty() const
Checks if the matrix is empty.
Definition: SystemMatrix.h:249
Common base class for multi-dimensional (2D and 3D) matrices.
Definition: matrix.h:308
size_t n[4]
Dimension of the matrix.
Definition: matrix.h:440
size_t dim(short int d=1) const
Query dimensions.
Definition: matrix.h:376
const T * ptr(size_t c=0) const
Reference through pointer.
Definition: matrix.h:397
matrixBase< T > & add(const matrixBase< T > &A, const T &alfa)
Add the given matrix A scaled by alfa to *this.
Definition: matrix.h:1584
matrixBase(vector< T > &vec)
Constructor using an external vector for matrix element storage.
Definition: matrix.h:313
void fill(const T *values, size_t n=0)
Fill the matrix with data from an array.
Definition: matrix.h:413
matrixBase()
The constructor is protected to allow sub-class instances only.
Definition: matrix.h:311
const vector< T > & toVec() const
Type casting to a one-dimensional utl::vector, for access.
Definition: matrix.h:385
vector< T > & elem
Actual matrix elements, stored column by column.
Definition: matrix.h:441
void fill(T s)
Fill the matrix with a scalar value.
Definition: matrix.h:411
matrixBase(const matrixBase< T > &mat, bool copyContent=true)
Copy constructor.
Definition: matrix.h:322
size_t size() const
Query total matrix size.
Definition: matrix.h:378
T norm2(int inc=1) const
Return the Euclidean norm of the matrix.
Definition: matrix.h:422
void clear()
Clears the matrix and sets its dimension to zero.
Definition: matrix.h:408
T * ptr(size_t c=0)
Access through pointer.
Definition: matrix.h:392
T asum(int inc=1) const
Return the sum of the absolute value of the matrix elements.
Definition: matrix.h:425
T sum(int inc=1) const
Return the sum of the matrix elements.
Definition: matrix.h:429
bool empty() const
Check if the matrix is empty.
Definition: matrix.h:380
bool zero(T tol=T(0)) const
Check if the matrix elements are all zero.
Definition: matrix.h:382
virtual void clearIfNrowChanged(size_t n1, size_t n2, size_t n3)=0
Clears the matrix content if the first dimension(s) changed.
void redim(size_t n_1, size_t n_2, size_t n_3, size_t n_4, bool forceClear)
Resize the matrix to dimension .
Definition: matrix.h:338
matrixBase(size_t n_1, size_t n_2, size_t n_3=1, size_t n_4=1)
Constructor creating a matrix of dimension .
Definition: matrix.h:316
vector< T > myElem
Internal matrix storage.
Definition: matrix.h:444
matrixBase< T > & multiply(const T &c)
Multiplication of this matrix by a scalar c.
Definition: matrix.h:1595
std::vector< T >::iterator begin()
Iterator to the start of the matrix elements.
Definition: matrix.h:403
std::vector< T >::iterator end()
Iterator to the end of the matrix elements.
Definition: matrix.h:405
Two-dimensional rectangular matrix with some algebraic operations.
Definition: matrix.h:456
matrix< T > & operator/=(T d)
Division by a scalar.
Definition: matrix.h:842
size_t cols() const
Query number of matrix columns.
Definition: matrix.h:560
bool multiply(const std::vector< T > &X, std::vector< T > &Y, const T &alpha, const T &beta=T(0), bool transA=false, int stridex=1, int stridey=1, unsigned int ofsx=0, unsigned int ofsy=0) const
Matrix-vector multiplication.
Definition: matrix.h:1625
bool multiply(const std::vector< T > &X, std::vector< T > &Y, bool transA=false, char addTo=0) const
Matrix-vector multiplication.
Definition: matrix.h:1603
matrix< T > & add(const matrix< T > &A, T alfa=T(1))
Add the given matrix A scaled by alfa to *this.
Definition: matrix.h:834
T trace() const
Return the trace of the matrix (sum of its diagonal elements).
Definition: matrix.h:741
matrix< T > & operator*=(T c)
Multiplication with a scalar.
Definition: matrix.h:840
std::vector< T > getColumn(size_t c) const
Extract a column from the matrix.
Definition: matrix.h:615
matrix(size_t r, size_t c)
Constructor creating a matrix of dimension .
Definition: matrix.h:464
void extractBlock(matrix< T > &block, size_t r, size_t c, bool addTo=false, bool transposed=false) const
Extract a block of the matrix to another matrix.
Definition: matrix.h:699
void fillRow(size_t r, const T *data)
Fill a row of the matrix.
Definition: matrix.h:661
void clearIfNrowChanged(size_t n1, size_t, size_t) override
Clears the content if the number of rows changed.
Definition: matrix.h:1024
T normInf() const
Return the infinite norm of the matrix.
Definition: matrix.h:921
size_t & ncol
Number of matrix columns.
Definition: matrix.h:1031
bool augmentRows(const matrix< T > &B)
Increase the number of rows by augmenting the given matrix.
Definition: matrix.h:527
bool augmentCols(const matrix< T > &B)
Increase the number of columns by augmenting the given matrix.
Definition: matrix.h:547
bool compatible(const std::vector< T > &X, const std::vector< T > &Y)
Check dimension compatibility for outer product multiplication.
Definition: matrix.h:1008
matrix< T > & multiply(const matrix< T > &A, const matrix< T > &B, bool transA=false, bool transB=false, bool addTo=false, const T &alpha=T(1))
Matrix-matrix multiplication.
Definition: matrix.h:1661
matrix< T > & multiply(T c)
Multiplication of this matrix by a scalar c.
Definition: matrix.h:844
T & operator()(size_t r, size_t c)
Index-1 based element access.
Definition: matrix.h:585
bool isSymmetric(T tol=T(0)) const
Check for symmetry.
Definition: matrix.h:813
vector< T > getRow(size_t r) const
Extract a row from the matrix.
Definition: matrix.h:602
matrix< T > & transpose()
Replace the current matrix by its transpose.
Definition: matrix.h:728
T det() const
Compute the determinant of a square matrix.
Definition: matrix.h:750
matrix< T > & operator=(const matrix< T > &A)
Assignment operator.
Definition: matrix.h:563
matrix< T > & operator+=(const matrix< T > &A)
Add the given matrix A to *this.
Definition: matrix.h:830
T inverse(T tol=T(0))
Compute the inverse of a square matrix.
Definition: matrix.h:773
void fillColumn(size_t c, const T *data)
Fill a column of the matrix.
Definition: matrix.h:654
matrix(const matrix< T > &mat, bool transposed=false)
Copy constructor, optionally creates the transpose of mat.
Definition: matrix.h:467
void fillBlock(const matrix< T > &block, size_t r, size_t c, bool transposed=false)
Fill a block of the matrix with another matrix.
Definition: matrix.h:671
bool compatible(const matrix< T > &A, const matrix< T > &B, bool transA, bool transB, size_t &M, size_t &N, size_t &K)
Check dimension compatibility for matrix-matrix multiplication.
Definition: matrix.h:950
size_t & nrow
Number of matrix rows.
Definition: matrix.h:1030
T colsum(size_t c) const
Return the sum of a matrix column.
Definition: matrix.h:745
const T & operator()(size_t r, size_t c) const
Index-1 based element reference.
Definition: matrix.h:594
void resize(size_t r, size_t c, bool forceClear=false)
Resize the matrix to dimension .
Definition: matrix.h:488
T rowsum(size_t r) const
Return the sum of a matrix row.
Definition: matrix.h:743
size_t rows() const
Query number of matrix rows.
Definition: matrix.h:558
matrix< T > & operator-=(const matrix< T > &A)
Subtract the given matrix A from *this.
Definition: matrix.h:832
virtual ~matrix()
Empty destructor.
Definition: matrix.h:480
matrix()
Constructor creating an empty matrix.
Definition: matrix.h:459
void addBlock(const matrix< T > &block, T s, size_t r, size_t c, bool transposed=false)
Add a scalar multiple of another matrix to a block of the matrix.
Definition: matrix.h:685
matrix< T > & diag(T d, size_t dim=0)
Create a diagonal matrix.
Definition: matrix.h:716
bool multiplyMat(const std::vector< T > &A, const matrix< T > &B, bool transB=false, bool addTo=false)
Matrix-matrix multiplication.
Definition: matrix.h:1712
bool outer_product(const std::vector< T > &X, const std::vector< T > &Y, bool addTo=false, T alpha=T(1))
Outer product between two vectors.
Definition: matrix.h:1733
void fillColumn(size_t c, const std::vector< T > &data)
Fill a column of the matrix.
Definition: matrix.h:646
matrix< T > & operator=(const std::vector< T > &X)
Overloaded assignment operator.
Definition: matrix.h:574
bool multiplyMat(const matrix< T > &A, const std::vector< T > &B, bool transA=false, bool addTo=false)
Matrix-matrix multiplication.
Definition: matrix.h:1691
matrix(vector< T > &vec)
Constructor using an external vector for matrix element storage.
Definition: matrix.h:461
bool compatible(const std::vector< T > &X, bool transA) const
Check dimension compatibility for matrix-vector multiplication.
Definition: matrix.h:935
matrix< T > & expandRows(int incRows)
Increase or decrease the number of rows in the matrix.
Definition: matrix.h:494
bool compatible(const matrix< T > &A, const std::vector< T > &B, bool transA, size_t &M, size_t &N, size_t &K)
Check dimension compatibility for matrix-matrix multiplication, when the matrix B is represented by a...
Definition: matrix.h:971
void fill(const std::vector< T > &v, size_t n, size_t m=0)
Fill the matrix with vector data.
Definition: matrix.h:628
bool compatible(const std::vector< T > &A, const matrix< T > &B, bool transB, size_t &M, size_t &N, size_t &K)
Check dimension compatibility for matrix-matrix multiplication, when the matrix A is represented by a...
Definition: matrix.h:990
A vector class with some added algebraic operations.
Definition: matrix.h:64
VecIter begin()
Start of the vector container, for update.
Definition: matrix.h:108
VecIter end()
End of the vector container, for update.
Definition: matrix.h:110
vector()
Constructor creating an empty vector.
Definition: matrix.h:67
vector(const T *values, size_t n)
Constructor creating a vector from a C-array.
Definition: matrix.h:71
const T & operator()(size_t i) const
Index-1 based element reference.
Definition: matrix.h:130
vector< T > & operator+=(const vector< T > &X)
Add the given vector X to *this.
Definition: matrix.h:182
void fill(T s)
Fill the vector with a scalar value.
Definition: matrix.h:137
void fill(const T *values, size_t n=0)
Fill the vector with data from an array.
Definition: matrix.h:139
T asum(size_t off=0, int inc=1) const
Return the sum of the absolute value of the vector elements.
Definition: matrix.h:1550
vector< T > & operator*=(const std::vector< T > &X)
Component-wise multiplication with a vector.
Definition: matrix.h:167
T max() const
Return the largest element of the vector.
Definition: matrix.h:247
typename std::vector< T >::iterator VecIter
Convenience alias for non-const iterators.
Definition: matrix.h:101
void push_back(T c)
Append a scalar value to the vector, increasing its size by one.
Definition: matrix.h:147
bool empty() const
Is the vector empty (zero size)?
Definition: matrix.h:90
T normInf(int inc=1) const
Return the infinite norm of the vector (no index offset).
Definition: matrix.h:244
void reserve(size_t n)
Pre-allocation of vector length to n.
Definition: matrix.h:293
void push_back(const T *p, const T *q)
Append a range of values increasing the size by q-p.
Definition: matrix.h:156
bool zero(T tol=T(0)) const
Is the vector elements all zero?
Definition: matrix.h:92
T & operator()(size_t i)
Index-1 based element access.
Definition: matrix.h:123
vector< T > & add(const std::vector< T > &X, const T &alfa=T(1), unsigned int ofsx=0, int stridex=1, unsigned int ofsy=0, int stridey=1)
Add the given vector X scaled by alfa to *this.
Definition: matrix.h:1562
vector< T > & operator/=(const std::vector< T > &X)
Component-wise division with a vector.
Definition: matrix.h:174
T norm2(size_t off=0, int inc=1) const
Return the Euclidean norm of the vector.
Definition: matrix.h:1513
ConstVecIter begin() const
Start of the vector container, for access.
Definition: matrix.h:104
void swap(vector< T > &vec)
Swap the content with another vector.
Definition: matrix.h:159
std::vector< T > myVec
Internal vector storage.
Definition: matrix.h:298
typename std::vector< T >::const_iterator ConstVecIter
Convenience alias for const iterators.
Definition: matrix.h:99
vector< T > & relax(T alfa, const std::vector< T > &X, const std::vector< T > &Y)
Perform where Z = *this.
Definition: matrix.h:203
vector< T > & operator-=(const vector< T > &X)
Subtract the given vector X from *this.
Definition: matrix.h:184
T & operator[](size_t i)
Index-0 based element access.
Definition: matrix.h:118
ConstVecIter end() const
End of the vector container, for access.
Definition: matrix.h:106
vector< T > & operator=(const std::vector< T > &X)
Overloaded assignment operator.
Definition: matrix.h:76
T min() const
Return the smallest element of the vector.
Definition: matrix.h:249
void push_back(ConstVecIter i1, ConstVecIter i2)
Append a range of values increasing the size by i2-i1.
Definition: matrix.h:150
vector(const std::vector< T > &X)
Overloaded copy constructor.
Definition: matrix.h:73
bool resize(size_t n, char forceClear=0)
Resize the vector to length n.
Definition: matrix.h:277
T sum(size_t off=0, int inc=1, size_t max=0) const
Return the sum of the vector elements.
Definition: matrix.h:260
void clear()
Clear the vector, setting its size to zero.
Definition: matrix.h:295
T dot(const T *v, size_t nv, size_t off1=0, int inc1=1, size_t off2=0, int inc2=1) const
Dot product between *this and another vector.
Definition: matrix.h:1502
vector< T > & operator*=(T c)
Multiplication with a scalar.
Definition: matrix.h:1494
vector< T > & relax(T alfa, const std::vector< T > &X)
Perform where Y = *this.
Definition: matrix.h:192
size_t size() const
Size of the vector.
Definition: matrix.h:88
vector< T > & operator/=(T d)
Division by a scalar.
Definition: matrix.h:164
const T & operator[](size_t i) const
Index-0 based element reference.
Definition: matrix.h:120
T normInf(size_t &off, int inc=1, bool sign=false) const
Return the infinite norm of the vector, or signed max value.
Definition: matrix.h:1526
vector(size_t n)
Constructor creating a vector of length n.
Definition: matrix.h:69
const T * ptr() const
Reference through pointer.
Definition: matrix.h:85
T * ptr()
Access through pointer.
Definition: matrix.h:83
T dot(const std::vector< T > &v, size_t off1=0, int inc1=1, size_t off2=0, int inc2=1) const
Dot product between *this and another vector.
Definition: matrix.h:225
General utility classes and functions.
Definition: SIMoptions.h:22
T trunc(T v)
Truncate a value to zero when it is less than a given threshold.
Definition: matrix.h:1765
double zero_print_tol
Zero tolerance for printing numbers.
Definition: MatVec.C:24
std::ostream & operator<<(std::ostream &s, const vector< T > &X)
Print the vector X to the stream s.
Definition: matrix.h:1783
int nval_per_line
Number of values to print per line.
Definition: MatVec.C:23
const char RETAIN
Flag for vector::resize() method telling it to retain its content.
Definition: matrix.h:55
std::istream & operator>>(std::istream &s, vector< T > &X)
Read the vector X from the stream s.
Definition: matrix.h:1772
void writeMatlab(const char *label, const vector< T > &X, std::ostream &s=std::cout)
Print the vector X to the stream s in matlab format.
Definition: matrix.h:1860
Global parameters for controlling the print of vectors and matrices.