Imagine++
Matrix.h
1  // ===========================================================================
2 // Imagine++ Libraries
3 // Copyright (C) Imagine
4 // For detailed information: http://imagine.enpc.fr/software
5 // ===========================================================================
6 
7 
8 namespace Imagine {
11 
12  template<typename T> class SymMatrix;
13 
22  template <typename T>
23  class Matrix : public MultiArray<T,2> {
24  typedef MultiArray<T,2> Base;
25  public:
27  typedef typename Base::iterator iterator;
30 
36  Matrix() : Base() {}
44  Matrix(int M,int N) : Base(M,N) {}
57  Matrix(T* t,int M,int N,bool handleDelete=false) : Base(t,M,N,handleDelete) {}
64  Matrix(const Base& A) : Base(A) {}
71  Matrix(const SymMatrix<T>& A);
74  ~Matrix() {}
82  Matrix& operator=(const Matrix& A) { Base::operator=(A); return *this; }
89  Matrix clone() const
90  {
91  Matrix A(nrow(),ncol());
92  equalize((int)this->totalSize(),this->data(),1,A.data(),1);
93  return A;
94  }
101  int nrow() const { return this->width(); }
108  int ncol() const { return this->height(); }
116  Matrix& fill(T x) { Base::fill(x); return *this; }
124  static Matrix Zero(int M,int N){
125  Matrix Z(M,N);
126  memset(Z.data(),0,Z.totalSize()*sizeof(T));
127  return Z;
128  }
136  static Matrix Identity(int N){
137  Matrix I=Matrix::Zero(N,N);
138  for (int i=0;i<N;i++)
139  I(i,i)=T(1);
140  return I;
141  }
149  Vector<T> getCol(int j) const {
150  assert(j>=0 && j<ncol());
151  Vector<T> v(nrow());
152  equalize(nrow(),this->data()+nrow()*j,1,v.data(),1);
153  return v;
154  }
162  Vector<T> getRow(int i) const {
163  assert(i>=0 && i<nrow());
164  Vector<T> v(ncol());
165  equalize(ncol(),this->data()+i,nrow(),v.data(),1);
166  return v;
167  }
175  int n=std::min(nrow(),ncol());
176  Vector<T> v(n);
177  for (int i=0;i<n;i++)
178  v[i]=(*this)(i,i);
179  return v;
180  }
189  Matrix& setCol(int j, const Vector<T>& v) {
190  assert(v.size()==(size_t)nrow() && j>=0 && j<ncol());
191  equalize(nrow(),v.data(),1,this->data()+nrow()*j,1);
192  return *this;
193  }
202  Matrix& setRow(int i, const Vector<T>& v) {
203  assert(v.size()==(size_t)ncol() && i>=0 && i<nrow());
204  equalize(ncol(),v.data(),1,this->data()+i,nrow());
205  return *this;
206  }
215  int n=int(v.size());
216  assert(n==std::min(nrow(),ncol()));
217  for (int i=0;i<n;i++)
218  (*this)(i,i)=v[i];
219  return *this;
220  }
230  assert (j>=0 && j<ncol());
231  return Vector<T>(this->data()+(size_t)nrow()*j,nrow());
232  }
243  Vector<T> getSubColRef(int j, int i0, int m) {
244  assert (j>=0 && j<ncol() && i0>=0 && m>=0 && i0+m<=nrow());
245  return Vector<T>(this->data()+i0+(size_t)+nrow()*j,m);
246  }
256  Matrix getSubMat(int i0, int m, int j0, int n) const {
257  assert (i0>=0 && i0+m<=nrow() && j0>=0 && j0+n<=ncol());
258  Matrix retMat(m,n);
259  for(int j=0;j<n;j++)
260  equalize(m,this->data()+(j+j0)*nrow()+i0, 1, retMat.data()+j*m, 1);
261  return retMat;
262  }
272  Matrix getColsRef(int j0,int n) {
273  assert (j0>=0 && n>=0 && j0+n<=ncol());
274  return Matrix(this->data()+(size_t)nrow()*j0, nrow(), n);
275  }
283  Vector<T> operator *(const Vector<T>& v) const {
284  assert((size_t)ncol()==v.size());
285  Vector<T> y(nrow());
286  y.fill(T(0));
287  affinity(nrow(), ncol(), T(1), this->data(), v.data(), T(0), y.data(), '0');
288  return y;
289  }
298  friend Vector<T> tmult(const Matrix& A,const Vector<T>& v) {
299  assert((size_t)A.nrow()==v.size());
300  Vector<T> y(A.ncol());
301  y.fill(T(0));
302  affinity(A.nrow(), A.ncol(), T(1), A.data(), v.data(), T(0), y.data(), 'T');
303  return y;
304  }
312  friend Matrix transpose(const Matrix& A) {
313  Matrix r(A.ncol(),A.nrow());
314  for(int i=0;i<A.nrow();i++) for(int j=0;j<A.ncol();j++) r(j,i)=A(i,j);
315  return r;
316  }
324  friend Matrix inverse(const Matrix& A) {
325  assert(A.nrow()==A.ncol());
326  Matrix invA=A.clone();
327  int *pivots=new int[A.nrow()];
328  int info,n=A.ncol();
329  // LU Inverse
330  matrixInverse(&n,invA.data(),&n,pivots,&info);
331  if (info) {
332  std::cerr<<"Cannot invert matrix"<<std::endl;
333  invA.fill(T(0));
334  delete [] pivots;
335  return invA;
336  }
337  delete[] pivots;
338  return invA;
339  }
347  Matrix operator *(const Matrix& B) const {
348  assert(ncol()==B.nrow());
349  int p=ncol();
350  Matrix C(nrow(),B.ncol());
351  C.fill(T(0));
352  matrixProduct(C.nrow(),C.ncol(),p,T(1),this->data(),nrow(),B.data(),B.nrow(),T(0),C.data(), 'N', 'N');
353  return C;
354  }
363  friend Matrix tmult(const Matrix& A,const Matrix& B) {
364  assert(A.nrow()==B.nrow());
365  int p=A.nrow();
366  Matrix C(A.ncol(),B.ncol());
367  C.fill(T(0));
368  matrixProduct(C.nrow(),C.ncol(),p,T(1),A.data(),A.nrow(),B.data(),B.nrow(),T(0),C.data(),'T','N');
369  return C;
370  }
379  friend Matrix multt(const Matrix& A,const Matrix& B) {
380  assert(A.ncol()==B.ncol());
381  int p=A.ncol();
382  Matrix C(A.nrow(),B.nrow());
383  C.fill(T(0));
384  matrixProduct(C.nrow(),C.ncol(),p,T(1),A.data(),A.nrow(),B.data(),B.nrow(),T(0),C.data(),'N','T');
385  return C;
386  }
395  friend Matrix tmultt(const Matrix& A,const Matrix& B) {
396  assert(A.nrow()==B.ncol());
397  int p=A.nrow();
398  Matrix C(A.ncol(),B.nrow());
399  C.fill(T(0));
400  matrixProduct(C.nrow(),C.ncol(),p,T(1),A.data(),A.nrow(),B.data(),B.nrow(),T(0),C.data(),'T','T');
401  return C;
402  }
410  Matrix operator *(const SymMatrix<T>& B) const;
418  Matrix operator +(const Matrix& B) const {
419  Matrix C=clone();
420  return (C+=B);
421  }
429  Matrix operator -(const Matrix& B) const {
430  Matrix C=clone();
431  return (C-=B);
432  }
440  Matrix& operator +=(const Matrix& B) {
441  assert(ncol()==B.ncol());
442  assert(nrow()==B.nrow());
443  combine((int)this->totalSize(), 1, B.data(), this->data());
444  return *this;
445  }
453  Matrix& operator -=(const Matrix& B) {
454  assert(ncol()==B.ncol());
455  assert(nrow()==B.nrow());
456  combine((int)this->totalSize(), -1, B.data(), this->data());
457  return *this;
458  }
465  Matrix operator - () const {
466  return (*this)*T(-1);
467  }
475  Matrix operator *(T x) const {
476  Matrix C=clone();
477  return (C*=x);
478  }
486  Matrix operator /(T x) const {
487  Matrix C=clone();
488  return (C/=x);
489  }
498  multiply((int)this->totalSize(), x, this->data());
499  return *this;
500  }
509  multiply((int)this->totalSize(), 1 / x, this->data());
510  return *this;
511  }
520  friend inline Matrix operator *(T x, const Matrix& B) {return B*x;}
529  friend T dot(const Matrix& a,const Matrix& b) {
530  assert(a.nrow()==b.nrow() && a.ncol()==b.ncol());
531  return scalarProduct((int)a.totalSize(), a.data(), b.data());
532  }
546  friend T norm(const Matrix& A,char type='F') {
547  int nr=A.nrow(),nc=A.ncol();
548  T n = matrixNorm(&type,&nr,&nc,(T*)A.data(),&nr);
549  return n;
550  }
558  friend T rcond(const Matrix& A)
559  {
560  if (A.nrow()<A.ncol())
561  return rcond(transpose(A));
562  // QR decomposition
563  Matrix tmpA=A.clone();
564  int info,m=A.nrow(),n=A.ncol();
565  // Estimation of condition number
566  T rcond;
567  conditionNumber(&n,tmpA.data(),&m,&rcond,&info);
568  if (info!=0) {
569  std::cerr << "Bad cond number " << std::endl;
570  rcond = T(-1);
571  }
572  return rcond;
573  }
584  friend void svd(const Matrix& A,Matrix &U,Vector<T> &S, Matrix &Vt,bool all=false) {
585  Matrix cpy=A.clone();
586  U=Matrix(A.nrow(),A.nrow());
587  Vt=Matrix(A.ncol(),A.ncol());
588  S=Vector<T>(std::min(A.nrow(),A.ncol()));
589  int nr=A.nrow(),nc=A.ncol();
590  int info;
591  singularValuesDecomposition(all?"A":"S",&nr,&nc,cpy.data(),S.data(),U.data(),Vt.data(),&info);
592  if (info) {
593  std::cerr<<"Cannot compute SVD"<<std::endl;
594  return;
595  }
596  }
610  friend Vector<T> linSolve(const Matrix& A,const Vector<T>& b)
611  {
612  assert(b.size() == (size_t)A.nrow());
613  if(A.nrow()==A.ncol())
614  {
615  // Linear system is not over- or under- determined
616  Matrix tmpA=A.clone();
617  Vector<T> tmpb=b.clone();
618  // LU factorization
619  int *pivots=new int[A.nrow()];
620  int info,n=A.nrow();
621  LUSystemSolve(&n,tmpA.data(),pivots,tmpb.data(),&info);
622  if (info) {
623  std::cerr<<"Cannot solve linear system"<<std::endl;
624  tmpb.fill(T(0));
625  delete [] pivots;
626  return tmpb;
627  tmpb.fill(T(0));
628  std::cerr<<"Cannot solve linear system"<<std::endl;
629  }
630  delete[] pivots;
631  return tmpb;
632  } else {
633  // Least-squares solution in the over-determined case, minimum norm solution in the under-determined one
634  int info,m=A.nrow(),n=A.ncol(),o=std::max(m,n);
635  Matrix tmpA=A.clone();
636  Vector<T> tmpb(o);
637  for (int i=0;i<m;i++) tmpb[i] = b[i];
638  QRSystemSolve(&m,&n,tmpA.data(),tmpb.data(),&info);
639  if (info) {
640  tmpb.fill(T(0));
641  std::cerr<<"Cannot solve linear system"<<std::endl;
642  }
643  return tmpb.getSubVect(0,A.ncol());
644  }
645  }
655  friend Matrix pseudoInverse(const Matrix& A,T tolrel=0) {
656  int m=A.nrow(),n=A.ncol();
657  if (n > m)
658  return transpose(pseudoInverse(transpose(A)));
659  Matrix U,Vt;
660  Vector<T> S;
661  svd(A,U,S,Vt);
662  if (tolrel==0)
663  tolrel=std::numeric_limits<T>::epsilon();
664  T tol = S[0] * tolrel;
665  for (int i=0;i<n;i++) {
666  T s=(S[i]>tol)?T(1)/S[i]:T(0);
667  for (int j=0;j<n;j++)
668  Vt(i,j)*=s;
669  }
670  return transpose(U.getColsRef(0,n)*Vt);
671  }
683  friend bool QR(const Matrix& A,Matrix& Q, Matrix& R, bool all=false){
684  assert(A.nrow()>=A.ncol());
685  int info,m=A.nrow(),n=A.ncol(),qn=(all?m:n);
686  if(all) {
687  // Copy A in beginning of Q
688  Q = Matrix(A.nrow(),A.nrow());
689  equalize((int)A.totalSize(),A.data(),1,Q.data(),1);
690  R=Matrix(A.nrow(),A.ncol());
691  QRFactorization(&m,&qn,&n,Q.data(),R.data(),&info);
692  } else
693  {
694  Q = A.clone();
695  R=Matrix(A.ncol(),A.ncol());
696  QRFactorization(&m,&qn,&n,Q.data(),R.data(),&info);
697  }
698  if(info!=1) {
699  std::cerr << "Failed QR" << std::endl;
700  return false;
701  }
702  return true;
703  }
715  friend Matrix cholesky(const Matrix& A,bool low=true) {
716  assert(A.nrow()==A.ncol());
717  Matrix b=A.clone();
718  int n=A.nrow();
719  int info;
720  choleskyDecomposition(low?"L":"U",&n,b.data(),&info);
721  if (info) {
722  std::cerr << "Bad cholesky" << std::endl;
723  b.fill(T(0));
724  return b;
725  }
726  if (low) {
727  for (int j=1; j<n; j++)
728  for (int i=0; i<j; i++)
729  b(i,j)=0;
730  } else {
731  for (int j=1; j<n; j++)
732  for (int i=0; i<j; i++)
733  b(j,i)=0;
734  }
735  return b;
736  }
737 
738 
747  friend T det(const Matrix& A) {
748  assert(A.nrow()==A.ncol());
749  Matrix b=A.clone();
750  // LU decomposition
751  int info,m=A.nrow(),n=A.ncol();
752  T d = determinant(&m,&n,b.data(),&info);
753  if (info!=0) {
754  std::cerr<<"Cannot compute LU decomposition"<<std::endl;
755  return 0;
756  }
757  return d;
758  }
768  friend inline std::ostream& operator<<(std::ostream& out,const Matrix& A) {
769  out<<A.nrow()<<" "<<A.ncol()<< std::endl;
770  for (int i=0;i<A.nrow();i++) {
771  for (int j=0;j<A.ncol();j++) {
772  out<<A(i,j);
773  if (j<A.ncol()-1) out<<" ";
774  }
775  out<<std::endl;
776  }
777  return out;
778  }
788  friend inline std::istream& operator>>(std::istream& in,Matrix& A) {
789  int m,n;
790  in>>m>>n;
791  A.setSize(m,n);
792  for (int i=0;i<A.nrow();i++)
793  for (int j=0;j<A.ncol();j++)
794  in>>A(i,j);
795  return in;
796  }
797  };
798 
808  template <typename T>
809  inline Matrix<T> Diagonal(const Vector<T>& d){
810  int n=int(d.size());
811  Matrix<T> D=Matrix<T>::Zero(n,n);
812  for (int i=0;i<n;i++)
813  D(i,i)=d[i];
814  return D;
815  }
816 
824  template <typename T>
825  Matrix<T> outerProduct(const Vector<T>& v1,const Vector<T>& v2) {
826  assert(v1.size()==v2.size());
827  int n=(int)v1.size();
828  Matrix<T> A(n,n);
829  A.fill(T(0));
830  matrixProduct(n,n,1,T(1),v1.data(),n,v2.data(),1,T(0),A.data(),'N','N');
831  return A;
832  }
833 
847  template <typename T,int N>
849  if (N<=3)
850  return inverse(A); // Direct inversion
851  FMatrix<T,N,N> invA(A);
852  // LU
853  int *pivots=new int[N];
854  int info,n=N;
855  // Inverse
856  matrixInverse(&n,invA.data(),&n,pivots,&info);
857  if (info) {
858  std::cerr<<"Cannot invert matrix"<<std::endl;
859  invA.fill(T(0));
860  }
861  delete[] pivots;
862  return invA;
863  }
864 
875  template <typename T,int M,int N>
876  void svd(const FMatrix<T,M,N>& A, FMatrix<T,M,M> &U, FVector<T,M> &S, FMatrix<T,N,N> &Vt, bool all=false) {
877  FMatrix<T,M,N> cpy(A);
878  int nr=M,nc=N;
879  int info;
880  singularValuesDecomposition(all?"A":"S",&nr,&nc,cpy.data(),S.data(),U.data(),Vt.data(),&info);
881  if (info) {
882  std::cerr<<"Cannot compute SVD"<<std::endl;
883  return;
884  }
885  }
886 
892  template <typename T>
893  void eigenvalues(const Matrix<T>& A, Vector<T>& wr, Vector<T>& wi)
894  {
895  assert(A.ncol()==A.nrow());
896  Matrix<T> cpy = A.clone();
897  wr = Vector<T>(A.nrow());
898  wi = Vector<T>(A.nrow());
899 
900  int n = A.nrow();
901  Vector<T> vl(A.nrow()), vr(A.nrow());
902  int info;
903 
904  eigenValues(&n, cpy.data(), wr.data(), wi.data(), vr.data(), &info);
905  if (info)
906  std::cerr << "Cannot compute eigenvalues" << std::endl;
907  }
908 
922  template <typename T,int M,int N>
924  {
925  if(M==N)
926  {
927  // Linear system is not over- or under- determined
928  FMatrix<T,M,N> tmpA(A);
929  FVector<T,M> tmpb(b);
930  // LU factorization
931  int *pivots=new int[M];
932  int info,n=M;
933  // Solve linear system using factorization
934  int nrhs=1;
935  LUSystemSolve(&n,tmpA.data(),pivots,tmpb.data(),&info);
936  if (info) {
937  std::cerr<<"Cannot solve linear system"<<std::endl;
938  delete [] pivots;
939  return FVector<T,N>(T(0));
940  tmpb.fill(T(0));
941  std::cerr<<"Cannot solve linear system"<<std::endl;
942  }
943  delete[] pivots;
944  return FVector<T,N>(tmpb.data()); //Safe since M==N but necessary for instantiation if M!=N
945  }
946  // Least-squares solution in the over-determined case, minimum norm solution in the under-determined one
947  int info,m=M,n=N,o=std::max(m,n);
948  int nrhs = 1;
949  FMatrix<T,M,N> tmpA(A);
950  if(M > N) { // Over-determined system: least squares solution
951  FVector<T,M> tmpb;
952  for (int i=0;i<m;i++) tmpb[i] = b[i];
953  QRSystemSolve(&m,&n,tmpA.data(),tmpb.data(),&info);
954  if (info) {
955  tmpb.fill(T(0));
956  std::cerr<<"Cannot solve linear system"<<std::endl;
957  }
958  return FVector<T,N>( tmpb.data() );
959  }
960  // Under-determined system (M<N): minimum norm solution
961  FVector<T,N> tmpb;
962  for (int i=0;i<m;i++) tmpb[i] = b[i];
963  QRSystemSolve(&m,&n,tmpA.data(),tmpb.data(),&info);
964  if (info) {
965  tmpb.fill(T(0));
966  std::cerr<<"Cannot solve linear system"<<std::endl;
967  }
968  return tmpb;
969  }
984  template <typename T,int M,int N>
986  if (N > M)
987  return transpose(pseudoInverse(transpose(A)));
989  FVector<T,M> S;
990  svd(A,U,S,Vt);
991  if (tolrel==0)
992  tolrel=std::numeric_limits<T>::epsilon();
993  T tol = S[0] * tolrel;
994  for (int i=0;i<N;i++) {
995  T s=(S[i]>tol)?T(1)/S[i]:T(0);
996  for (int j=0;j<N;j++)
997  Vt(i,j)*=s;
998  }
999  return transpose(FMatrix<T,M,N>(U.data())*Vt);
1000  }
1011  template <typename T,int M,int N>
1013  assert(M>=N);
1014  int info,m=M,n=N;
1015  Q = A;
1016  QRFactorization(&m,&n,&n,Q.data(),R.data(),&info);
1017  if(info!=1) {
1018  std::cerr << "Failed QR" << std::endl;
1019  return false;
1020  }
1021  return true;
1022  }
1033  template <typename T,int M,int N>
1035  assert(M>=N);
1036  int info,m=M,n=N;
1037  equalize(M*N,A.data(),1,Q.data(),1);
1038  QRFactorization(&m,&m,&n,Q.data(),R.data(),&info);
1039  if(info!=1) {
1040  std::cerr << "Failed QR" << std::endl;
1041  return false;
1042  }
1043  return true;
1044  }
1056  template <typename T,int N>
1057  FMatrix<T,N,N> cholesky(const FMatrix<T,N,N>& A,bool low=true) {
1058  FMatrix<T,N,N> b(A);
1059  int n=N;
1060  int info;
1061  choleskyDecomposition(low?"L":"U",&n,b.data(),&info);
1062  if (info) {
1063  std::cerr << "Bad cholesky" << std::endl;
1064  b.fill(T(0));
1065  return b;
1066  }
1067  if (low) {
1068  for (int j=1; j<n; j++)
1069  for (int i=0; i<j; i++)
1070  b(i,j)=0;
1071  } else {
1072  for (int j=1; j<n; j++)
1073  for (int i=0; i<j; i++)
1074  b(j,i)=0;
1075  }
1076  return b;
1077  }
1078 
1088  template <typename T,int N>
1090  int n=N;
1091  FMatrix<T,N,N> b(A);
1092  int info;
1093  // Try cholesky factorization -> only for positive matrices!
1094  choleskyDecomposition((char*)"L",&n,b.data(),&info);
1095  if (!info) {
1096  T d=T(1);
1097  for (int i=0; i<n; i++) d*=b(i,i);
1098  return d*d;
1099  }
1100  // Try SVD
1101  FMatrix<T,N,N> U,V;
1102  FVector<T,N> s;
1103  svd(A,U,s,V);
1104  T d=T(1);
1105  for (int i=0; i<N; i++) d*=s[i];
1106  return d;
1107  }
1108 
1110 }
FMatrix< T, N, N > inverseFMatrix(const FMatrix< T, N, N > &A)
Inverse.
Definition: Matrix.h:848
Vector< T > getCol(int j) const
Get column.
Definition: Matrix.h:149
Vector & fill(T x)
Filling.
Definition: Vector.h:93
friend Matrix multt(const Matrix &A, const Matrix &B)
Product.
Definition: Matrix.h:379
const T * data() const
Data pointer (read).
Definition: FArray.h:188
Vector of fixed size.
Definition: FVector.h:17
Symmetric Matrix of variable size.
Definition: Matrix.h:12
Vector of variable size.
Definition: Vector.h:24
friend Matrix tmult(const Matrix &A, const Matrix &B)
Product.
Definition: Matrix.h:363
Matrix & setCol(int j, const Vector< T > &v)
Set column.
Definition: Matrix.h:189
friend void svd(const Matrix &A, Matrix &U, Vector< T > &S, Matrix &Vt, bool all=false)
SVD.
Definition: Matrix.h:584
Matrix operator/(T x) const
Scalar division.
Definition: Matrix.h:486
size_t size() const
Size.
Definition: Array.h:194
FVector & fill(const T &v)
Filling.
Definition: FVector.h:75
friend T rcond(const Matrix &A)
Reciprocal condition number.
Definition: Matrix.h:558
friend Matrix pseudoInverse(const Matrix &A, T tolrel=0)
Pseudo inverse.
Definition: Matrix.h:655
void setSize(const Coords< dim > &sz)
Change sizes.
Definition: MultiArray.h:146
friend Matrix tmultt(const Matrix &A, const Matrix &B)
Product.
Definition: Matrix.h:395
Base::iterator iterator
Iterator type.
Definition: Matrix.h:27
Matrix< T > outerProduct(const Vector< T > &v1, const Vector< T > &v2)
Outer product matrix.
Definition: Matrix.h:825
Matrix & setDiagonal(const Vector< T > &v)
Set diagonal.
Definition: Matrix.h:214
Vector< T > getRow(int i) const
Get row.
Definition: Matrix.h:162
~Matrix()
Destructor.
Definition: Matrix.h:74
size_t totalSize() const
Total Size.
Definition: MultiArray.h:218
int height() const
Size alias 1.
Definition: MultiArray.h:239
Matrix & operator/=(T x)
Scalar in place division.
Definition: Matrix.h:508
Matrix operator-() const
Opposite.
Definition: Matrix.h:465
Matrix getSubMat(int i0, int m, int j0, int n) const
Get sub matrix.
Definition: Matrix.h:256
FMatrix< T, M, M > Diagonal(const FVector< T, M > &d)
Diagonal.
Definition: FMatrix.h:743
MultiArray & operator=(const MultiArray &A)
Assignment.
Definition: MultiArray.h:175
friend Matrix inverse(const Matrix &A)
Inverse.
Definition: Matrix.h:324
bool QR(const FMatrix< T, M, N > &A, FMatrix< T, M, N > &Q, FMatrix< T, N, N > &R)
QR decomposition.
Definition: Matrix.h:1012
friend Matrix cholesky(const Matrix &A, bool low=true)
Cholesky decomposition.
Definition: Matrix.h:715
FMatrix< T, N, M > pseudoInverse(const FMatrix< T, M, N > &A, T tolrel=0)
Pseudo inverse.
Definition: Matrix.h:985
const T * data() const
Data pointer (read).
Definition: FMatrix.h:182
Vector< T > getSubColRef(int j, int i0, int m)
Get sub column by reference.
Definition: Matrix.h:243
Vector< T > getDiagonal() const
Get diagonal.
Definition: Matrix.h:174
Vector getSubVect(size_t offset, size_t size) const
Sub vector.
Definition: Vector.h:114
Matrix & operator*=(T x)
Scalar in place multiplication.
Definition: Matrix.h:497
friend T det(const Matrix &A)
Determinant.
Definition: Matrix.h:747
int nrow() const
Number of rows.
Definition: Matrix.h:101
T detFMatrix(const FMatrix< T, N, N > &A)
Determinant.
Definition: Matrix.h:1089
FMatrix & fill(const T &v)
Filling.
Definition: FMatrix.h:78
FMatrix< T, N, N > cholesky(const FMatrix< T, N, N > &A, bool low=true)
Cholesky decomposition.
Definition: Matrix.h:1057
bool QRAll(const FMatrix< T, M, N > &A, FMatrix< T, M, M > &Q, FMatrix< T, M, N > &R)
QR decomposition.
Definition: Matrix.h:1034
Matrix(const Base &A)
Copy constructor.
Definition: Matrix.h:64
Matrix getColsRef(int j0, int n)
Get part of columns by reference.
Definition: Matrix.h:272
Matrix()
Empty constructor.
Definition: Matrix.h:36
T * data()
Data pointer (read/write).
Definition: Array.h:217
Matrix & operator=(const Matrix &A)
Assignment.
Definition: Matrix.h:82
Matrix operator+(const Matrix &B) const
Addition.
Definition: Matrix.h:418
nD array of variable size.
Definition: MultiArray.h:21
friend std::ostream & operator<<(std::ostream &out, const Matrix &A)
ASCII write.
Definition: Matrix.h:768
FMatrix< T, M, M > inverse(const FMatrix< T, M, M > &A)
Inverse.
Definition: FMatrix.h:678
friend std::istream & operator>>(std::istream &in, Matrix &A)
ASCII read.
Definition: Matrix.h:788
Matrix & setRow(int i, const Vector< T > &v)
Set row.
Definition: Matrix.h:202
FMatrix< T, N, M > transpose(const FMatrix< T, M, N > &A)
Transpose.
Definition: FMatrix.h:564
Matrix clone() const
Cloning.
Definition: Matrix.h:89
MultiArray & fill(T x)
Filling.
Definition: MultiArray.h:204
friend bool QR(const Matrix &A, Matrix &Q, Matrix &R, bool all=false)
QR decomposition.
Definition: Matrix.h:683
Base::const_iterator const_iterator
Const iterator type.
Definition: Matrix.h:29
friend Vector< T > tmult(const Matrix &A, const Vector< T > &v)
Product with vector.
Definition: Matrix.h:298
T * iterator
Iterator type.
Definition: MultiArray.h:45
friend T dot(const Matrix &a, const Matrix &b)
Scalar product.
Definition: Matrix.h:529
Matrix(T *t, int M, int N, bool handleDelete=false)
Constructor (pre-allocated).
Definition: Matrix.h:57
int width() const
Size alias 0.
Definition: MultiArray.h:232
static Matrix Zero(int M, int N)
Zero matrix.
Definition: Matrix.h:124
void eigenvalues(const Matrix< T > &A, Vector< T > &wr, Vector< T > &wi)
Eigenvalues calculation Given square matrix, only eigenvalues are obtained.
Definition: Matrix.h:893
int ncol() const
Number of columns.
Definition: Matrix.h:108
Matrix & operator+=(const Matrix &B)
In place Addition.
Definition: Matrix.h:440
Matrix(int M, int N)
Constructor (known size).
Definition: Matrix.h:44
friend Matrix transpose(const Matrix &A)
Transpose.
Definition: Matrix.h:312
static Matrix Identity(int N)
Identity.
Definition: Matrix.h:136
friend Vector< T > linSolve(const Matrix &A, const Vector< T > &b)
Linear system.
Definition: Matrix.h:610
void svd(const FMatrix< T, M, N > &A, FMatrix< T, M, M > &U, FVector< T, M > &S, FMatrix< T, N, N > &Vt, bool all=false)
SVD.
Definition: Matrix.h:876
Imagine++ namespace.
Definition: Array.h:7
Matrix of variable size.
Definition: Matrix.h:23
Matrix & operator-=(const Matrix &B)
In place Substraction.
Definition: Matrix.h:453
Vector clone() const
Cloning.
Definition: Vector.h:80
Matrix & fill(T x)
Filling.
Definition: Matrix.h:116
FVector< T, N > linSolve(const FMatrix< T, M, N > &A, const FVector< T, M > &b)
Linear system.
Definition: Matrix.h:923
Vector< T > getColRef(int j)
Get column by reference.
Definition: Matrix.h:229
const T * const_iterator
Const iterator type.
Definition: MultiArray.h:47
friend T norm(const Matrix &A, char type='F')
Norm.
Definition: Matrix.h:546
Matrix of fixed dimension.
Definition: FMatrix.h:17
Vector< T > operator*(const Vector< T > &v) const
Product with vector.
Definition: Matrix.h:283