Visual Servoing Platform  version 3.4.0
vpMatrix.h
1 /****************************************************************************
2  *
3  * ViSP, open source Visual Servoing Platform software.
4  * Copyright (C) 2005 - 2019 by Inria. All rights reserved.
5  *
6  * This software is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2 of the License, or
9  * (at your option) any later version.
10  * See the file LICENSE.txt at the root directory of this source
11  * distribution for additional information about the GNU GPL.
12  *
13  * For using ViSP with software that can not be combined with the GNU
14  * GPL, please contact Inria about acquiring a ViSP Professional
15  * Edition License.
16  *
17  * See http://visp.inria.fr for more information.
18  *
19  * This software was developed at:
20  * Inria Rennes - Bretagne Atlantique
21  * Campus Universitaire de Beaulieu
22  * 35042 Rennes Cedex
23  * France
24  *
25  * If you have questions regarding the use of this file, please contact
26  * Inria at visp@inria.fr
27  *
28  * This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
29  * WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
30  *
31  * Description:
32  * Matrix manipulation.
33  *
34  *****************************************************************************/
35 
36 #ifndef vpMatrix_H
37 #define vpMatrix_H
38 
39 #include <visp3/core/vpArray2D.h>
40 #include <visp3/core/vpConfig.h>
41 #include <visp3/core/vpException.h>
42 #include <visp3/core/vpForceTwistMatrix.h>
43 #include <visp3/core/vpHomogeneousMatrix.h>
44 #include <visp3/core/vpRotationMatrix.h>
45 #include <visp3/core/vpTime.h>
46 #include <visp3/core/vpVelocityTwistMatrix.h>
47 
48 #include <iostream>
49 #include <math.h>
50 
51 class vpRowVector;
52 class vpColVector;
56 class vpForceTwistMatrix;
57 
153 class VISP_EXPORT vpMatrix : public vpArray2D<double>
154 {
155 public:
160  typedef enum {
161  LU_DECOMPOSITION
162  } vpDetMethod;
163 
164 public:
169  vpMatrix() : vpArray2D<double>(0, 0) {}
170 
177  vpMatrix(unsigned int r, unsigned int c) : vpArray2D<double>(r, c) {}
178 
186  vpMatrix(unsigned int r, unsigned int c, double val) : vpArray2D<double>(r, c, val) {}
187  vpMatrix(const vpMatrix &M, unsigned int r, unsigned int c, unsigned int nrows, unsigned int ncols);
188 
201  vpMatrix(const vpArray2D<double> &A) : vpArray2D<double>(A) {}
202 
203  vpMatrix(const vpMatrix &A) : vpArray2D<double>(A) {}
204 
205 #if (VISP_CXX_STANDARD >= VISP_CXX_STANDARD_11)
206  vpMatrix(vpMatrix &&A);
207  explicit vpMatrix(const std::initializer_list<double> &list);
208  explicit vpMatrix(unsigned int nrows, unsigned int ncols, const std::initializer_list<double> &list);
209  explicit vpMatrix(const std::initializer_list<std::initializer_list<double> > &lists);
210 #endif
211 
213  virtual ~vpMatrix() {}
214 
219  void clear()
220  {
221  if (data != NULL) {
222  free(data);
223  data = NULL;
224  }
225 
226  if (rowPtrs != NULL) {
227  free(rowPtrs);
228  rowPtrs = NULL;
229  }
230  rowNum = colNum = dsize = 0;
231  }
232 
233  //-------------------------------------------------
234  // Setting a diagonal matrix
235  //-------------------------------------------------
246  static unsigned int getLapackMatrixMinSize() {
247  return m_lapack_min_size;
248  }
249 
262  static void setLapackMatrixMinSize(unsigned int min_size) {
263  m_lapack_min_size = min_size;
264  }
266 
267  //-------------------------------------------------
268  // Setting a diagonal matrix
269  //-------------------------------------------------
272  void diag(const double &val = 1.0);
273  void diag(const vpColVector &A);
274  // Initialize an identity matrix n-by-n
275  void eye();
276  void eye(unsigned int n);
277  // Initialize an identity matrix m-by-n
278  void eye(unsigned int m, unsigned int n);
280 
281  //---------------------------------
282  // Assignment
283  //---------------------------------
286  vpMatrix &operator<<(double *);
287  vpMatrix& operator<<(double val);
288  vpMatrix& operator,(double val);
290 #if (VISP_CXX_STANDARD >= VISP_CXX_STANDARD_11)
291  vpMatrix &operator=(const vpMatrix &A);
293 
294  vpMatrix& operator=(const std::initializer_list<double> &list);
295  vpMatrix& operator=(const std::initializer_list<std::initializer_list<double> > &lists);
296 #endif
297  vpMatrix &operator=(double x);
299 
300  //-------------------------------------------------
301  // Stacking
302  //-------------------------------------------------
305  // Stack the matrix A below the current one, copy if not initialized this =
306  // [ this A ]^T
307  void stack(const vpMatrix &A);
308  void stack(const vpRowVector &r);
309  void stack(const vpColVector &c);
310  // Stacks columns of a matrix in a vector
311  void stackColumns(vpColVector &out);
312 
313  // Stacks columns of a matrix in a vector
314  vpColVector stackColumns();
315 
316  // Stacks columns of a matrix in a vector
317  void stackRows(vpRowVector &out);
318 
319  // Stacks columns of a matrix in a vector
320  vpRowVector stackRows();
322 
323  //---------------------------------
324  // Matrix insertion
325  //---------------------------------
328  // Insert matrix A in the current matrix at the given position (r, c).
329  void insert(const vpMatrix &A, unsigned int r, unsigned int c);
331 
332  //-------------------------------------------------
333  // Columns, Rows, Diag extraction, SubMatrix
334  //-------------------------------------------------
337  vpMatrix extract(unsigned int r, unsigned int c, unsigned int nrows, unsigned int ncols) const;
338  vpColVector getCol(unsigned int j) const;
339  vpColVector getCol(unsigned int j, unsigned int i_begin, unsigned int size) const;
340  vpRowVector getRow(unsigned int i) const;
341  vpRowVector getRow(unsigned int i, unsigned int j_begin, unsigned int size) const;
342  vpColVector getDiag() const;
343  void init(const vpMatrix &M, unsigned int r, unsigned int c, unsigned int nrows, unsigned int ncols);
345 
346  //---------------------------------
347  // Matrix operations.
348  //---------------------------------
351  // return the determinant of the matrix.
352  double det(vpDetMethod method = LU_DECOMPOSITION) const;
353  double detByLU() const;
354 #ifdef VISP_HAVE_EIGEN3
355  double detByLUEigen3() const;
356 #endif
357 #if defined(VISP_HAVE_LAPACK)
358  double detByLULapack() const;
359 #endif
360 #if (VISP_HAVE_OPENCV_VERSION >= 0x020101)
361  double detByLUOpenCV() const;
362 #endif
363 
364  // Compute the exponential matrix of a square matrix
365  vpMatrix expm() const;
366 
367  // operation A = A + B
368  vpMatrix &operator+=(const vpMatrix &B);
369  // operation A = A - B
370  vpMatrix &operator-=(const vpMatrix &B);
371  vpMatrix operator*(const vpMatrix &B) const;
372  vpMatrix operator*(const vpRotationMatrix &R) const;
373  vpMatrix operator*(const vpHomogeneousMatrix &R) const;
374  vpMatrix operator*(const vpVelocityTwistMatrix &V) const;
375  vpMatrix operator*(const vpForceTwistMatrix &V) const;
376  // operation t_out = A * t (A is unchanged, t and t_out are translation
377  // vectors)
378  vpTranslationVector operator*(const vpTranslationVector &tv) const;
379  vpColVector operator*(const vpColVector &v) const;
380  vpMatrix operator+(const vpMatrix &B) const;
381  vpMatrix operator-(const vpMatrix &B) const;
382  vpMatrix operator-() const;
383 
385  vpMatrix &operator+=(double x);
387  vpMatrix &operator-=(double x);
389  vpMatrix &operator*=(double x);
391  vpMatrix &operator/=(double x);
392 
393  // Cij = Aij * x (A is unchanged)
394  vpMatrix operator*(double x) const;
395  // Cij = Aij / x (A is unchanged)
396  vpMatrix operator/(double x) const;
397 
403  double sum() const;
404  double sumSquare() const;
405 
406  //-------------------------------------------------
407  // Hadamard product
408  //-------------------------------------------------
410  vpMatrix hadamard(const vpMatrix &m) const;
411 
412  //-------------------------------------------------
413  // Kronecker product
414  //-------------------------------------------------
417  // Compute Kronecker product matrix
418  void kron(const vpMatrix &m1, vpMatrix &out) const;
419 
420  // Compute Kronecker product matrix
421  vpMatrix kron(const vpMatrix &m1) const;
423 
424  //-------------------------------------------------
425  // Transpose
426  //-------------------------------------------------
429  // Compute the transpose C = A^T
430  vpMatrix t() const;
431 
432  // Compute the transpose C = A^T
433  vpMatrix transpose() const;
434  void transpose(vpMatrix &At) const;
435 
436  vpMatrix AAt() const;
437  void AAt(vpMatrix &B) const;
438 
439  vpMatrix AtA() const;
440  void AtA(vpMatrix &B) const;
442 
443  //-------------------------------------------------
444  // Matrix inversion
445  //-------------------------------------------------
448  // inverse matrix A using the LU decomposition
449  vpMatrix inverseByLU() const;
450 
451 #if defined(VISP_HAVE_EIGEN3)
452  vpMatrix inverseByLUEigen3() const;
453 #endif
454 #if defined(VISP_HAVE_LAPACK)
455  vpMatrix inverseByLULapack() const;
456 #endif
457 #if (VISP_HAVE_OPENCV_VERSION >= 0x020101)
458  vpMatrix inverseByLUOpenCV() const;
459 #endif
460 
461  // inverse matrix A using the Cholesky decomposition (only for real
462  // symmetric matrices)
463  vpMatrix inverseByCholesky() const;
464 
465 #if defined(VISP_HAVE_LAPACK)
466  vpMatrix inverseByCholeskyLapack() const;
467 #endif
468 #if (VISP_HAVE_OPENCV_VERSION >= 0x020101)
469  vpMatrix inverseByCholeskyOpenCV() const;
470 #endif
471 
472  // inverse matrix A using the QR decomposition
473  vpMatrix inverseByQR() const;
474 #if defined(VISP_HAVE_LAPACK)
475  vpMatrix inverseByQRLapack() const;
476 #endif
477 
478  // inverse triangular matrix
479  vpMatrix inverseTriangular(bool upper = true) const;
480 
481  vpMatrix pseudoInverse(double svThreshold = 1e-6) const;
482  unsigned int pseudoInverse(vpMatrix &Ap, double svThreshold = 1e-6) const;
483  unsigned int pseudoInverse(vpMatrix &Ap, vpColVector &sv, double svThreshold = 1e-6) const;
484  unsigned int pseudoInverse(vpMatrix &Ap, vpColVector &sv, double svThreshold, vpMatrix &imA, vpMatrix &imAt) const;
485  unsigned int pseudoInverse(vpMatrix &Ap, vpColVector &sv, double svThreshold, vpMatrix &imA, vpMatrix &imAt, vpMatrix &kerAt) const;
486  vpMatrix pseudoInverse(int rank_in) const;
487  int pseudoInverse(vpMatrix &Ap, int rank_in) const;
488  int pseudoInverse(vpMatrix &Ap, vpColVector &sv, int rank_in) const;
489  int pseudoInverse(vpMatrix &Ap, vpColVector &sv, int rank_in, vpMatrix &imA, vpMatrix &imAt) const;
490  int pseudoInverse(vpMatrix &Ap, vpColVector &sv, int rank_in, vpMatrix &imA, vpMatrix &imAt, vpMatrix &kerAt) const;
491 
492 #if defined(VISP_HAVE_LAPACK)
493  vpMatrix pseudoInverseLapack(double svThreshold = 1e-6) const;
494  unsigned int pseudoInverseLapack(vpMatrix &Ap, double svThreshold = 1e-6) const;
495  unsigned int pseudoInverseLapack(vpMatrix &Ap, vpColVector &sv, double svThreshold = 1e-6) const;
496  unsigned int pseudoInverseLapack(vpMatrix &Ap, vpColVector &sv, double svThreshold, vpMatrix &imA, vpMatrix &imAt, vpMatrix &kerAt) const;
497  vpMatrix pseudoInverseLapack(int rank_in) const;
498  int pseudoInverseLapack(vpMatrix &Ap, int rank_in) const;
499  int pseudoInverseLapack(vpMatrix &Ap, vpColVector &sv, int rank_in) const;
500  int pseudoInverseLapack(vpMatrix &Ap, vpColVector &sv, int rank_in, vpMatrix &imA, vpMatrix &imAt, vpMatrix &kerAt) const;
501 #endif
502 #if defined(VISP_HAVE_EIGEN3)
503  vpMatrix pseudoInverseEigen3(double svThreshold = 1e-6) const;
504  unsigned int pseudoInverseEigen3(vpMatrix &Ap, double svThreshold = 1e-6) const;
505  unsigned int pseudoInverseEigen3(vpMatrix &Ap, vpColVector &sv, double svThreshold = 1e-6) const;
506  unsigned int pseudoInverseEigen3(vpMatrix &Ap, vpColVector &sv, double svThreshold, vpMatrix &imA, vpMatrix &imAt, vpMatrix &kerAt) const;
507  vpMatrix pseudoInverseEigen3(int rank_in) const;
508  int pseudoInverseEigen3(vpMatrix &Ap, int rank_in) const;
509  int pseudoInverseEigen3(vpMatrix &Ap, vpColVector &sv, int rank_in) const;
510  int pseudoInverseEigen3(vpMatrix &Ap, vpColVector &sv, int rank_in, vpMatrix &imA, vpMatrix &imAt, vpMatrix &kerAt) const;
511 #endif
512 #if (VISP_HAVE_OPENCV_VERSION >= 0x020101)
513  vpMatrix pseudoInverseOpenCV(double svThreshold = 1e-6) const;
514  unsigned int pseudoInverseOpenCV(vpMatrix &Ap, double svThreshold = 1e-6) const;
515  unsigned int pseudoInverseOpenCV(vpMatrix &Ap, vpColVector &sv, double svThreshold = 1e-6) const;
516  unsigned int pseudoInverseOpenCV(vpMatrix &Ap, vpColVector &sv, double svThreshold, vpMatrix &imA, vpMatrix &imAt, vpMatrix &kerAt) const;
517  vpMatrix pseudoInverseOpenCV(int rank_in) const;
518  int pseudoInverseOpenCV(vpMatrix &Ap, int rank_in) const;
519  int pseudoInverseOpenCV(vpMatrix &Ap, vpColVector &sv, int rank_in) const;
520  int pseudoInverseOpenCV(vpMatrix &Ap, vpColVector &sv, int rank_in, vpMatrix &imA, vpMatrix &imAt, vpMatrix &kerAt) const;
521 #endif
522 
523 
524  //-------------------------------------------------
525  // SVD decomposition
526  //-------------------------------------------------
527 
530  double cond(double svThreshold = 1e-6) const;
531  unsigned int kernel(vpMatrix &kerAt, double svThreshold = 1e-6) const;
532  unsigned int nullSpace(vpMatrix &kerA, double svThreshold = 1e-6) const;
533  unsigned int nullSpace(vpMatrix &kerA, int dim) const;
534 
535  // solve Ax=B using the SVD decomposition (usage A = solveBySVD(B,x) )
536  void solveBySVD(const vpColVector &B, vpColVector &x) const;
537  // solve Ax=B using the SVD decomposition (usage x=A.solveBySVD(B))
538  vpColVector solveBySVD(const vpColVector &B) const;
539 
540  // singular value decomposition SVD
541  void svd(vpColVector &w, vpMatrix &V);
542 #ifdef VISP_HAVE_EIGEN3
543  void svdEigen3(vpColVector &w, vpMatrix &V);
544 #endif
545 #if defined(VISP_HAVE_LAPACK)
546  void svdLapack(vpColVector &w, vpMatrix &V);
547 #endif
548 #if (VISP_HAVE_OPENCV_VERSION >= 0x020101) // Require opencv >= 2.1.1
549  void svdOpenCV(vpColVector &w, vpMatrix &V);
550 #endif
551 
552 
553  //-------------------------------------------------
554  // QR decomposition
555  //-------------------------------------------------
556 
559  unsigned int qr(vpMatrix &Q, vpMatrix &R, bool full = false, bool squareR = false, double tol = 1e-6) const;
560  unsigned int qrPivot(vpMatrix &Q, vpMatrix &R, vpMatrix &P, bool full = false, bool squareR = false, double tol = 1e-6) const;
561  void solveByQR(const vpColVector &b, vpColVector &x) const;
562  vpColVector solveByQR(const vpColVector &b) const;
564 
565  //-------------------------------------------------
566  // Eigen values and vectors
567  //-------------------------------------------------
571  // compute the eigen values using Lapack
572  vpColVector eigenValues() const;
573  void eigenValues(vpColVector &evalue, vpMatrix &evector) const;
575 
576  //-------------------------------------------------
577  // Norms
578  //-------------------------------------------------
581  double euclideanNorm() const;
582  double frobeniusNorm() const;
583  double inducedL2Norm() const;
584  double infinityNorm() const;
586 
587  //---------------------------------
588  // Printing
589  //---------------------------------
592  std::ostream &cppPrint(std::ostream &os, const std::string &matrixName = "A", bool octet = false) const;
593  std::ostream &csvPrint(std::ostream &os) const;
594  std::ostream &maplePrint(std::ostream &os) const;
595  std::ostream &matlabPrint(std::ostream &os) const;
596  int print(std::ostream &s, unsigned int length, const std::string &intro = "") const;
597  void printSize() const { std::cout << getRows() << " x " << getCols() << " "; }
599 
600  //------------------------------------------------------------------
601  // Static functionalities
602  //------------------------------------------------------------------
603 
604  //---------------------------------
605  // Setting a diagonal matrix with Static Public Member Functions
606  //---------------------------------
609  // Create a diagonal matrix with the element of a vector DAii = Ai
610  static void createDiagonalMatrix(const vpColVector &A, vpMatrix &DA);
612 
613  //---------------------------------
614  // Matrix insertion with Static Public Member Functions
615  //---------------------------------
618  // Insert matrix B in matrix A at the given position (r, c).
619  static vpMatrix insert(const vpMatrix &A, const vpMatrix &B, unsigned int r, unsigned int c);
620  // Insert matrix B in matrix A (not modified) at the given position (r, c),
621  // the result is given in matrix C.
622  static void insert(const vpMatrix &A, const vpMatrix &B, vpMatrix &C, unsigned int r, unsigned int c);
623 
624  //---------------------------------
625  // Stacking with Static Public Member Functions
626  //---------------------------------
629  // Juxtapose to matrices C = [ A B ]
630  static vpMatrix juxtaposeMatrices(const vpMatrix &A, const vpMatrix &B);
631  // Juxtapose to matrices C = [ A B ]
632  static void juxtaposeMatrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C);
633  // Stack two matrices C = [ A B ]^T
634  static vpMatrix stack(const vpMatrix &A, const vpMatrix &B);
635  static vpMatrix stack(const vpMatrix &A, const vpRowVector &r);
636  static vpMatrix stack(const vpMatrix &A, const vpColVector &c);
637 
638  // Stack two matrices C = [ A B ]^T
639  static void stack(const vpMatrix &A, const vpMatrix &B, vpMatrix &C);
640  static void stack(const vpMatrix &A, const vpRowVector &r, vpMatrix &C);
641  static void stack(const vpMatrix &A, const vpColVector &c, vpMatrix &C);
643 
644  //---------------------------------
645  // Matrix operations Static Public Member Functions
646  //---------------------------------
649  static void add2Matrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C);
650  static void add2Matrices(const vpColVector &A, const vpColVector &B, vpColVector &C);
651  static void add2WeightedMatrices(const vpMatrix &A, const double &wA, const vpMatrix &B, const double &wB,
652  vpMatrix &C);
653  static void computeHLM(const vpMatrix &H, const double &alpha, vpMatrix &HLM);
654  static void mult2Matrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C);
655  static void mult2Matrices(const vpMatrix &A, const vpMatrix &B, vpRotationMatrix &C);
656  static void mult2Matrices(const vpMatrix &A, const vpMatrix &B, vpHomogeneousMatrix &C);
657  static void mult2Matrices(const vpMatrix &A, const vpColVector &B, vpColVector &C);
658  static void multMatrixVector(const vpMatrix &A, const vpColVector &v, vpColVector &w);
659  static void negateMatrix(const vpMatrix &A, vpMatrix &C);
660  static void sub2Matrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C);
661  static void sub2Matrices(const vpColVector &A, const vpColVector &B, vpColVector &C);
663 
664  //---------------------------------
665  // Kronecker product Static Public Member Functions
666  //---------------------------------
669  // Compute Kronecker product matrix
670  static void kron(const vpMatrix &m1, const vpMatrix &m2, vpMatrix &out);
671 
672  // Compute Kronecker product matrix
673  static vpMatrix kron(const vpMatrix &m1, const vpMatrix &m2);
675 
676  //-------------------------------------------------
677  // 2D Convolution Static Public Member Functions
678  //-------------------------------------------------
680  static vpMatrix conv2(const vpMatrix &M, const vpMatrix &kernel, const std::string &mode="full");
681  static void conv2(const vpMatrix &M, const vpMatrix &kernel, vpMatrix &res, const std::string &mode="full");
682 
683  //---------------------------------
684  // Covariance computation Static Public Member Functions
685  //---------------------------------
688  static vpMatrix computeCovarianceMatrix(const vpMatrix &A, const vpColVector &x, const vpColVector &b);
689  static vpMatrix computeCovarianceMatrix(const vpMatrix &A, const vpColVector &x, const vpColVector &b,
690  const vpMatrix &w);
691  static vpMatrix computeCovarianceMatrixVVS(const vpHomogeneousMatrix &cMo, const vpColVector &deltaS,
692  const vpMatrix &Ls, const vpMatrix &W);
693  static vpMatrix computeCovarianceMatrixVVS(const vpHomogeneousMatrix &cMo, const vpColVector &deltaS,
694  const vpMatrix &Ls);
696 
697  //---------------------------------
698  // Matrix I/O Static Public Member Functions
699  //---------------------------------
713  static inline bool loadMatrix(const std::string &filename, vpArray2D<double> &M, bool binary = false,
714  char *header = NULL)
715  {
716  return vpArray2D<double>::load(filename, M, binary, header);
717  }
718 
729  static inline bool loadMatrixYAML(const std::string &filename, vpArray2D<double> &M, char *header = NULL)
730  {
731  return vpArray2D<double>::loadYAML(filename, M, header);
732  }
733 
748  static inline bool saveMatrix(const std::string &filename, const vpArray2D<double> &M, bool binary = false,
749  const char *header = "")
750  {
751  return vpArray2D<double>::save(filename, M, binary, header);
752  }
753 
766  static inline bool saveMatrixYAML(const std::string &filename, const vpArray2D<double> &M, const char *header = "")
767  {
768  return vpArray2D<double>::saveYAML(filename, M, header);
769  }
771 
772 #if defined(VISP_BUILD_DEPRECATED_FUNCTIONS)
773 
781  vp_deprecated void init() {}
782 
786  vp_deprecated void stackMatrices(const vpMatrix &A) { stack(A); }
791  vp_deprecated static vpMatrix stackMatrices(const vpMatrix &A, const vpMatrix &B) { return stack(A, B); }
796  vp_deprecated static void stackMatrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C) { stack(A, B, C); }
801  vp_deprecated static vpMatrix stackMatrices(const vpMatrix &A, const vpRowVector &B);
806  vp_deprecated static void stackMatrices(const vpMatrix &A, const vpRowVector &B, vpMatrix &C);
811  vp_deprecated static vpMatrix stackMatrices(const vpColVector &A, const vpColVector &B);
816  vp_deprecated static void stackMatrices(const vpColVector &A, const vpColVector &B, vpColVector &C);
817 
821  vp_deprecated void setIdentity(const double &val = 1.0);
822 
823  vp_deprecated vpRowVector row(unsigned int i);
824  vp_deprecated vpColVector column(unsigned int j);
825 
826  // Deprecated functions using GSL
827 #ifndef DOXYGEN_SHOULD_SKIP_THIS
828 
831  vp_deprecated double detByLUGsl() const {
832 #if defined(VISP_HAVE_LAPACK)
833  return detByLULapack();
834 #else
836  "Undefined detByLULapack(). Install Lapack 3rd party"));
837 #endif
838  }
839 
843  vp_deprecated vpMatrix inverseByLUGsl() const {
844 #if defined(VISP_HAVE_LAPACK)
845  return inverseByLULapack();
846 #else
848  "Undefined inverseByLULapack(). Install Lapack 3rd party"));
849 #endif
850  }
851 
855  vpMatrix inverseByCholeskyGsl() const {
856 #if defined(VISP_HAVE_LAPACK)
857  return inverseByCholeskyLapack();
858 #else
860  "Undefined inverseByCholeskyLapack(). Install Lapack 3rd party"));
861 #endif
862  }
863 
867  vpMatrix inverseByQRGsl() const {
868 #if defined(VISP_HAVE_LAPACK)
869  return inverseByQRLapack();
870 #else
872  "Undefined inverseByQRLapack(). Install Lapack 3rd party"));
873 #endif
874  }
875 
879  vpMatrix pseudoInverseGsl(double svThreshold = 1e-6) const {
880 #if defined(VISP_HAVE_LAPACK)
881  return pseudoInverseLapack(svThreshold);
882 #else
883  (void)svThreshold;
885  "Undefined pseudoInverseLapack(). Install Lapack 3rd party"));
886 #endif
887  }
888 
892  unsigned int pseudoInverseGsl(vpMatrix &Ap, double svThreshold = 1e-6) const {
893 #if defined(VISP_HAVE_LAPACK)
894  return pseudoInverseLapack(Ap, svThreshold);
895 #else
896  (void)Ap;
897  (void)svThreshold;
899  "Undefined pseudoInverseLapack(). Install Lapack 3rd party"));
900 #endif
901  }
902 
906  unsigned int pseudoInverseGsl(vpMatrix &Ap, vpColVector &sv, double svThreshold = 1e-6) const {
907 #if defined(VISP_HAVE_LAPACK)
908  return pseudoInverseLapack(Ap, sv, svThreshold);
909 #else
910  (void)Ap;
911  (void)sv;
912  (void)svThreshold;
914  "Undefined pseudoInverseLapack(). Install Lapack 3rd party"));
915 #endif
916  }
917 
921  unsigned int pseudoInverseGsl(vpMatrix &Ap, vpColVector &sv, double svThreshold, vpMatrix &imA, vpMatrix &imAt,
922  vpMatrix &kerAt) const {
923 #if defined(VISP_HAVE_LAPACK)
924  return pseudoInverseLapack(Ap, sv, svThreshold, imA, imAt, kerAt);
925 #else
926  (void)Ap;
927  (void)sv;
928  (void)svThreshold;
929  (void)imA;
930  (void)imAt;
931  (void)kerAt;
933  "Undefined pseudoInverseLapack(). Install Lapack 3rd party"));
934 #endif
935  }
936 
940  void svdGsl(vpColVector &w, vpMatrix &V) {
941 #if defined(VISP_HAVE_LAPACK)
942  svdLapack(w, V);
943 #else
944  (void)w;
945  (void)V;
947  "Undefined svdLapack(). Install Lapack 3rd party"));
948 #endif
949  }
950 
951 #endif // ifndef DOXYGEN_SHOULD_SKIP_THIS
952 
953 #endif
954 
955 private:
956  static unsigned int m_lapack_min_size;
957  static const unsigned int m_lapack_min_size_default;
958 
959 #if defined(VISP_HAVE_LAPACK)
960  static void blas_dgemm(char trans_a, char trans_b, unsigned int M_, unsigned int N_, unsigned int K_, double alpha,
961  double *a_data, unsigned int lda_, double *b_data, unsigned int ldb_, double beta, double *c_data,
962  unsigned int ldc_);
963  static void blas_dgemv(char trans, unsigned int M_, unsigned int N_, double alpha, double *a_data, unsigned int lda_,
964  double *x_data, int incx_, double beta, double *y_data, int incy_);
965  static void blas_dsyev(char jobz, char uplo, unsigned int n_, double *a_data, unsigned int lda_,
966  double *w_data, double *work_data, int lwork_, int &info_);
967 #endif
968 
969  static void computeCovarianceMatrixVVS(const vpHomogeneousMatrix &cMo, const vpColVector &deltaS, const vpMatrix &Ls,
970  vpMatrix &Js, vpColVector &deltaP);
971 };
972 
974 #if defined(VISP_USE_MSVC) && defined(visp_EXPORTS)
975 const __declspec(selectany) unsigned int vpMatrix::m_lapack_min_size_default = 0;
976 __declspec(selectany) unsigned int vpMatrix::m_lapack_min_size = vpMatrix::m_lapack_min_size_default;
977 #endif
978 
979 #ifndef DOXYGEN_SHOULD_SKIP_THIS
980 VISP_EXPORT
981 #endif
982 vpMatrix operator*(const double &x, const vpMatrix &A);
983 #endif
Implementation of a matrix and operations on matrices.
Definition: vpMatrix.h:153
static bool saveYAML(const std::string &filename, const vpArray2D< Type > &A, const char *header="")
Definition: vpArray2D.h:830
vpArray2D< Type > & operator=(Type x)
Set all the elements of the array to x.
Definition: vpArray2D.h:413
static bool loadYAML(const std::string &filename, vpArray2D< Type > &A, char *header=NULL)
Definition: vpArray2D.h:652
Implementation of an homogeneous matrix and operations on such kind of matrices.
Implementation of row vector and the associated operations.
Definition: vpRowVector.h:115
static bool save(const std::string &filename, const vpArray2D< Type > &A, bool binary=false, const char *header="")
Definition: vpArray2D.h:737
error that can be emited by ViSP classes.
Definition: vpException.h:71
vp_deprecated void stackMatrices(const vpMatrix &A)
Definition: vpMatrix.h:786
Implementation of a generic 2D array used as base class for matrices and vectors. ...
Definition: vpArray2D.h:131
unsigned int getCols() const
Definition: vpArray2D.h:279
void printSize() const
Definition: vpMatrix.h:597
static bool saveMatrix(const std::string &filename, const vpArray2D< double > &M, bool binary=false, const char *header="")
Definition: vpMatrix.h:748
vpMatrix(unsigned int r, unsigned int c, double val)
Definition: vpMatrix.h:186
vpMatrix()
Definition: vpMatrix.h:169
Implementation of a rotation matrix and operations on such kind of matrices.
static unsigned int getLapackMatrixMinSize()
Definition: vpMatrix.h:246
vpColVector operator*(const double &x, const vpColVector &v)
static vp_deprecated void stackMatrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C)
Definition: vpMatrix.h:796
vpMatrix(const vpMatrix &A)
Definition: vpMatrix.h:203
vp_deprecated void init()
Definition: vpMatrix.h:781
void clear()
Definition: vpMatrix.h:219
friend std::ostream & operator<<(std::ostream &s, const vpArray2D< Type > &A)
Definition: vpArray2D.h:493
unsigned int getRows() const
Definition: vpArray2D.h:289
vpArray2D< Type > hadamard(const vpArray2D< Type > &m) const
Definition: vpArray2D.h:932
static void setLapackMatrixMinSize(unsigned int min_size)
Definition: vpMatrix.h:262
static vp_deprecated vpMatrix stackMatrices(const vpMatrix &A, const vpMatrix &B)
Definition: vpMatrix.h:791
vpMatrix(unsigned int r, unsigned int c)
Definition: vpMatrix.h:177
Implementation of column vector and the associated operations.
Definition: vpColVector.h:130
virtual ~vpMatrix()
Destructor (Memory de-allocation)
Definition: vpMatrix.h:213
static bool saveMatrixYAML(const std::string &filename, const vpArray2D< double > &M, const char *header="")
Definition: vpMatrix.h:766
vpMatrix(const vpArray2D< double > &A)
Definition: vpMatrix.h:201
Class that consider the case of a translation vector.
static bool loadMatrix(const std::string &filename, vpArray2D< double > &M, bool binary=false, char *header=NULL)
Definition: vpMatrix.h:713
static bool load(const std::string &filename, vpArray2D< Type > &A, bool binary=false, char *header=NULL)
Definition: vpArray2D.h:540
static bool loadMatrixYAML(const std::string &filename, vpArray2D< double > &M, char *header=NULL)
Definition: vpMatrix.h:729