Visual Servoing Platform  version 3.6.1 under development (2024-04-16)
vpMatrix.cpp
1 /****************************************************************************
2  *
3  * ViSP, open source Visual Servoing Platform software.
4  * Copyright (C) 2005 - 2023 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 https://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 
40 #include <algorithm>
41 #include <assert.h>
42 #include <cmath> // std::fabs
43 #include <fstream>
44 #include <limits> // numeric_limits
45 #include <sstream>
46 #include <stdio.h>
47 #include <stdlib.h>
48 #include <string.h>
49 #include <string>
50 #include <vector>
51 
52 #include <visp3/core/vpCPUFeatures.h>
53 #include <visp3/core/vpColVector.h>
54 #include <visp3/core/vpConfig.h>
55 #include <visp3/core/vpDebug.h>
56 #include <visp3/core/vpException.h>
57 #include <visp3/core/vpMath.h>
58 #include <visp3/core/vpMatrix.h>
59 #include <visp3/core/vpTranslationVector.h>
60 
61 #if defined(VISP_HAVE_SIMDLIB)
62 #include <Simd/SimdLib.h>
63 #endif
64 
65 #ifdef VISP_HAVE_LAPACK
66 #ifdef VISP_HAVE_GSL
67 #include <gsl/gsl_eigen.h>
68 #include <gsl/gsl_linalg.h>
69 #include <gsl/gsl_math.h>
70 #elif defined(VISP_HAVE_MKL)
71 #include <mkl.h>
72 typedef MKL_INT integer;
73 
74 void vpMatrix::blas_dsyev(char jobz, char uplo, unsigned int n_, double *a_data, unsigned int lda_, double *w_data,
75  double *work_data, int lwork_, int &info_)
76 {
77  MKL_INT n = static_cast<MKL_INT>(n_);
78  MKL_INT lda = static_cast<MKL_INT>(lda_);
79  MKL_INT lwork = static_cast<MKL_INT>(lwork_);
80  MKL_INT info = static_cast<MKL_INT>(info_);
81 
82  dsyev(&jobz, &uplo, &n, a_data, &lda, w_data, work_data, &lwork, &info);
83 }
84 
85 #else
86 #if defined(VISP_HAVE_LAPACK_BUILT_IN)
87 typedef long int integer;
88 #else
89 typedef int integer;
90 #endif
91 extern "C" integer dsyev_(char *jobz, char *uplo, integer *n, double *a, integer *lda, double *w, double *WORK,
92  integer *lwork, integer *info);
93 
94 void vpMatrix::blas_dsyev(char jobz, char uplo, unsigned int n_, double *a_data, unsigned int lda_, double *w_data,
95  double *work_data, int lwork_, int &info_)
96 {
97  integer n = static_cast<integer>(n_);
98  integer lda = static_cast<integer>(lda_);
99  integer lwork = static_cast<integer>(lwork_);
100  integer info = static_cast<integer>(info_);
101 
102  dsyev_(&jobz, &uplo, &n, a_data, &lda, w_data, work_data, &lwork, &info);
103 
104  lwork_ = static_cast<int>(lwork);
105  info_ = static_cast<int>(info);
106 }
107 #endif
108 #endif
109 
110 #if !defined(VISP_USE_MSVC) || (defined(VISP_USE_MSVC) && !defined(VISP_BUILD_SHARED_LIBS))
111 const unsigned int vpMatrix::m_lapack_min_size_default = 0;
112 unsigned int vpMatrix::m_lapack_min_size = vpMatrix::m_lapack_min_size_default;
113 #endif
114 
115 // Prototypes of specific functions
116 vpMatrix subblock(const vpMatrix &M, unsigned int col, unsigned int row);
117 
118 void compute_pseudo_inverse(const vpMatrix &U, const vpColVector &sv, const vpMatrix &V, unsigned int nrows,
119  unsigned int ncols, double svThreshold, vpMatrix &Ap, int &rank_out, int *rank_in,
120  vpMatrix *imA, vpMatrix *imAt, vpMatrix *kerAt)
121 {
122  Ap.resize(ncols, nrows, true, false);
123 
124  // compute the highest singular value and the rank of h
125  double maxsv = sv[0];
126 
127  rank_out = 0;
128 
129  unsigned int sv_size = sv.size();
130  for (unsigned int i = 0; i < sv_size; ++i) {
131  if (sv[i] > (maxsv * svThreshold)) {
132  rank_out++;
133  }
134  }
135 
136  unsigned int rank = static_cast<unsigned int>(rank_out);
137  if (rank_in) {
138  rank = static_cast<unsigned int>(*rank_in);
139  }
140 
141  for (unsigned int i = 0; i < ncols; ++i) {
142  for (unsigned int j = 0; j < nrows; ++j) {
143  for (unsigned int k = 0; k < rank; ++k) {
144  Ap[i][j] += (V[i][k] * U[j][k]) / sv[k];
145  }
146  }
147  }
148 
149  // Compute im(A)
150  if (imA) {
151  imA->resize(nrows, rank);
152 
153  for (unsigned int i = 0; i < nrows; ++i) {
154  for (unsigned int j = 0; j < rank; ++j) {
155  (*imA)[i][j] = U[i][j];
156  }
157  }
158  }
159 
160  // Compute im(At)
161  if (imAt) {
162  imAt->resize(ncols, rank);
163  for (unsigned int i = 0; i < ncols; ++i) {
164  for (unsigned int j = 0; j < rank; ++j) {
165  (*imAt)[i][j] = V[i][j];
166  }
167  }
168  }
169 
170  // Compute ker(At)
171  if (kerAt) {
172  kerAt->resize(ncols - rank, ncols);
173  if (rank != ncols) {
174  unsigned int v_rows = V.getRows();
175  for (unsigned int k = 0; k < (ncols - rank); ++k) {
176  unsigned j = k + rank;
177  for (unsigned int i = 0; i < v_rows; ++i) {
178  (*kerAt)[k][i] = V[i][j];
179  }
180  }
181  }
182  }
183 }
184 
190 vpMatrix::vpMatrix(const vpMatrix &M, unsigned int r, unsigned int c, unsigned int nrows, unsigned int ncols)
191  : vpArray2D<double>()
192 {
193  if (((r + nrows) > M.rowNum) || ((c + ncols) > M.colNum)) {
195  "Cannot construct a sub matrix (%dx%d) starting at "
196  "position (%d,%d) that is not contained in the "
197  "original matrix (%dx%d)",
198  nrows, ncols, r, c, M.rowNum, M.colNum));
199  }
200 
201  init(M, r, c, nrows, ncols);
202 }
203 
204 #if (VISP_CXX_STANDARD >= VISP_CXX_STANDARD_11)
206 {
207  rowNum = A.rowNum;
208  colNum = A.colNum;
209  rowPtrs = A.rowPtrs;
210  dsize = A.dsize;
211  data = A.data;
212 
213  A.rowNum = 0;
214  A.colNum = 0;
215  A.rowPtrs = nullptr;
216  A.dsize = 0;
217  A.data = nullptr;
218 }
219 
243 vpMatrix::vpMatrix(const std::initializer_list<double> &list) : vpArray2D<double>(list) { }
244 
267 vpMatrix::vpMatrix(unsigned int nrows, unsigned int ncols, const std::initializer_list<double> &list)
268  : vpArray2D<double>(nrows, ncols, list)
269 { }
270 
291 vpMatrix::vpMatrix(const std::initializer_list<std::initializer_list<double> > &lists) : vpArray2D<double>(lists) { }
292 #endif
293 
340 void vpMatrix::init(const vpMatrix &M, unsigned int r, unsigned int c, unsigned int nrows, unsigned int ncols)
341 {
342  unsigned int rnrows = r + nrows;
343  unsigned int cncols = c + ncols;
344 
345  if (rnrows > M.getRows()) {
346  throw(vpException(vpException::dimensionError, "Bad row dimension (%d > %d) used to initialize vpMatrix", rnrows,
347  M.getRows()));
348  }
349  if (cncols > M.getCols()) {
350  throw(vpException(vpException::dimensionError, "Bad column dimension (%d > %d) used to initialize vpMatrix", cncols,
351  M.getCols()));
352  }
353  resize(nrows, ncols, false, false);
354 
355  if (this->rowPtrs == nullptr) { // Fix coverity scan: explicit null dereferenced
356  return; // Noting to do
357  }
358  for (unsigned int i = 0; i < nrows; ++i) {
359  memcpy((*this)[i], &M[i + r][c], ncols * sizeof(double));
360  }
361 }
362 
404 vpMatrix vpMatrix::extract(unsigned int r, unsigned int c, unsigned int nrows, unsigned int ncols) const
405 {
406  unsigned int rnrows = r + nrows;
407  unsigned int cncols = c + ncols;
408 
409  if (rnrows > getRows()) {
410  throw(vpException(vpException::dimensionError, "Bad row dimension (%d > %d) used to initialize vpMatrix", rnrows,
411  getRows()));
412  }
413  if (cncols > getCols()) {
414  throw(vpException(vpException::dimensionError, "Bad column dimension (%d > %d) used to initialize vpMatrix", cncols,
415  getCols()));
416  }
417 
418  vpMatrix M;
419  M.resize(nrows, ncols, false, false);
420  for (unsigned int i = 0; i < nrows; ++i) {
421  memcpy(M[i], &(*this)[i + r][c], ncols * sizeof(double));
422  }
423 
424  return M;
425 }
426 
431 void vpMatrix::eye(unsigned int n) { eye(n, n); }
432 
437 void vpMatrix::eye(unsigned int m, unsigned int n)
438 {
439  resize(m, n);
440 
441  eye();
442 }
443 
449 {
450  for (unsigned int i = 0; i < rowNum; ++i) {
451  for (unsigned int j = 0; j < colNum; ++j) {
452  if (i == j) {
453  (*this)[i][j] = 1.0;
454  }
455  else {
456  (*this)[i][j] = 0;
457  }
458  }
459  }
460 }
461 
465 vpMatrix vpMatrix::t() const { return transpose(); }
466 
473 {
474  vpMatrix At;
475  transpose(At);
476  return At;
477 }
478 
485 {
486  At.resize(colNum, rowNum, false, false);
487 
488  if ((rowNum <= 16) || (colNum <= 16)) {
489  for (unsigned int i = 0; i < rowNum; ++i) {
490  for (unsigned int j = 0; j < colNum; ++j) {
491  At[j][i] = (*this)[i][j];
492  }
493  }
494  }
495  else {
496 #if defined(VISP_HAVE_SIMDLIB)
497  SimdMatTranspose(data, rowNum, colNum, At.data);
498 #else
499  // https://stackoverflow.com/a/21548079
500  const int tileSize = 32;
501  for (unsigned int i = 0; i < rowNum; i += tileSize) {
502  for (unsigned int j = 0; j < colNum; j++) {
503  for (unsigned int b = 0; b < static_cast<unsigned int>(tileSize) && i + b < rowNum; b++) {
504  At[j][i + b] = (*this)[i + b][j];
505  }
506  }
507  }
508 #endif
509  }
510 }
511 
518 {
519  vpMatrix B;
520 
521  AAt(B);
522 
523  return B;
524 }
525 
537 void vpMatrix::AAt(vpMatrix &B) const
538 {
539  if ((B.rowNum != rowNum) || (B.colNum != rowNum)) {
540  B.resize(rowNum, rowNum, false, false);
541  }
542 
543  // If available use Lapack only for large matrices
544  bool useLapack = ((rowNum > vpMatrix::m_lapack_min_size) || (colNum > vpMatrix::m_lapack_min_size));
545 #if !(defined(VISP_HAVE_LAPACK) && !defined(VISP_HAVE_LAPACK_BUILT_IN) && !defined(VISP_HAVE_GSL))
546  useLapack = false;
547 #endif
548 
549  if (useLapack) {
550 #if defined(VISP_HAVE_LAPACK) && !defined(VISP_HAVE_LAPACK_BUILT_IN) && !defined(VISP_HAVE_GSL)
551  const double alpha = 1.0;
552  const double beta = 0.0;
553  const char transa = 't';
554  const char transb = 'n';
555 
556  vpMatrix::blas_dgemm(transa, transb, rowNum, rowNum, colNum, alpha, data, colNum, data, colNum, beta, B.data,
557  rowNum);
558 #endif
559  }
560  else {
561  // compute A*A^T
562  for (unsigned int i = 0; i < rowNum; ++i) {
563  for (unsigned int j = i; j < rowNum; ++j) {
564  double *pi = rowPtrs[i]; // row i
565  double *pj = rowPtrs[j]; // row j
566 
567  // sum (row i .* row j)
568  double ssum = 0;
569  for (unsigned int k = 0; k < colNum; ++k) {
570  ssum += *(pi++) * *(pj++);
571  }
572 
573  B[i][j] = ssum; // upper triangle
574  if (i != j) {
575  B[j][i] = ssum; // lower triangle
576  }
577  }
578  }
579  }
580 }
581 
593 void vpMatrix::AtA(vpMatrix &B) const
594 {
595  if ((B.rowNum != colNum) || (B.colNum != colNum)) {
596  B.resize(colNum, colNum, false, false);
597  }
598 
599  // If available use Lapack only for large matrices
600  bool useLapack = ((rowNum > vpMatrix::m_lapack_min_size) || (colNum > vpMatrix::m_lapack_min_size));
601 #if !(defined(VISP_HAVE_LAPACK) && !defined(VISP_HAVE_LAPACK_BUILT_IN) && !defined(VISP_HAVE_GSL))
602  useLapack = false;
603 #endif
604 
605  if (useLapack) {
606 #if defined(VISP_HAVE_LAPACK) && !defined(VISP_HAVE_LAPACK_BUILT_IN) && !defined(VISP_HAVE_GSL)
607  const double alpha = 1.0;
608  const double beta = 0.0;
609  const char transa = 'n';
610  const char transb = 't';
611 
612  vpMatrix::blas_dgemm(transa, transb, colNum, colNum, rowNum, alpha, data, colNum, data, colNum, beta, B.data,
613  colNum);
614 #endif
615  }
616  else {
617  for (unsigned int i = 0; i < colNum; ++i) {
618  double *Bi = B[i];
619  for (unsigned int j = 0; j < i; ++j) {
620  double *ptr = data;
621  double s = 0;
622  for (unsigned int k = 0; k < rowNum; ++k) {
623  s += (*(ptr + i)) * (*(ptr + j));
624  ptr += colNum;
625  }
626  *Bi++ = s;
627  B[j][i] = s;
628  }
629  double *ptr = data;
630  double s = 0;
631  for (unsigned int k = 0; k < rowNum; ++k) {
632  s += (*(ptr + i)) * (*(ptr + i));
633  ptr += colNum;
634  }
635  *Bi = s;
636  }
637  }
638 }
639 
646 {
647  vpMatrix B;
648 
649  AtA(B);
650 
651  return B;
652 }
653 
671 {
672  resize(A.getRows(), A.getCols(), false, false);
673 
674  if ((data != nullptr) && (A.data != nullptr) && (data != A.data)) {
675  memcpy(data, A.data, dsize * sizeof(double));
676  }
677 
678  return *this;
679 }
680 
682 {
683  resize(A.getRows(), A.getCols(), false, false);
684 
685  if ((data != nullptr) && (A.data != nullptr) && (data != A.data)) {
686  memcpy(data, A.data, dsize * sizeof(double));
687  }
688 
689  return *this;
690 }
691 
692 #if (VISP_CXX_STANDARD >= VISP_CXX_STANDARD_11)
694 {
695  if (this != &other) {
696  if (data) {
697  free(data);
698  }
699  if (rowPtrs) {
700  free(rowPtrs);
701  }
702 
703  rowNum = other.rowNum;
704  colNum = other.colNum;
705  rowPtrs = other.rowPtrs;
706  dsize = other.dsize;
707  data = other.data;
708 
709  other.rowNum = 0;
710  other.colNum = 0;
711  other.rowPtrs = nullptr;
712  other.dsize = 0;
713  other.data = nullptr;
714  }
715 
716  return *this;
717 }
718 
745 vpMatrix &vpMatrix::operator=(const std::initializer_list<double> &list)
746 {
747  if (dsize != static_cast<unsigned int>(list.size())) {
748  resize(1, static_cast<unsigned int>(list.size()), false, false);
749  }
750 
751  std::copy(list.begin(), list.end(), data);
752 
753  return *this;
754 }
755 
779 vpMatrix &vpMatrix::operator=(const std::initializer_list<std::initializer_list<double> > &lists)
780 {
781  unsigned int nrows = static_cast<unsigned int>(lists.size()), ncols = 0;
782  for (auto &l : lists) {
783  if (static_cast<unsigned int>(l.size()) > ncols) {
784  ncols = static_cast<unsigned int>(l.size());
785  }
786  }
787 
788  resize(nrows, ncols, false, false);
789  auto it = lists.begin();
790  for (unsigned int i = 0; i < rowNum; ++i, ++it) {
791  std::copy(it->begin(), it->end(), rowPtrs[i]);
792  }
793 
794  return *this;
795 }
796 #endif
797 
800 {
801  std::fill(data, data + (rowNum * colNum), x);
802  return *this;
803 }
804 
811 {
812  for (unsigned int i = 0; i < rowNum; ++i) {
813  for (unsigned int j = 0; j < colNum; ++j) {
814  rowPtrs[i][j] = *x++;
815  }
816  }
817  return *this;
818 }
819 
821 {
822  resize(1, 1, false, false);
823  rowPtrs[0][0] = val;
824  return *this;
825 }
826 
828 {
829  resize(1, colNum + 1, false, false);
830  rowPtrs[0][colNum - 1] = val;
831  return *this;
832 }
833 
870 {
871  unsigned int rows = A.getRows();
872  this->resize(rows, rows);
873 
874  (*this) = 0;
875  for (unsigned int i = 0; i < rows; ++i) {
876  (*this)[i][i] = A[i];
877  }
878 }
879 
910 void vpMatrix::diag(const double &val)
911 {
912  (*this) = 0;
913  unsigned int min_ = (rowNum < colNum) ? rowNum : colNum;
914  for (unsigned int i = 0; i < min_; ++i) {
915  (*this)[i][i] = val;
916  }
917 }
918 
930 {
931  unsigned int rows = A.getRows();
932  DA.resize(rows, rows, true);
933 
934  for (unsigned int i = 0; i < rows; ++i) {
935  DA[i][i] = A[i];
936  }
937 }
938 
944 {
945  vpTranslationVector t_out;
946 
947  if ((rowNum != 3) || (colNum != 3)) {
948  throw(vpException(vpException::dimensionError, "Cannot multiply a (%dx%d) matrix by a (%dx%d) translation vector",
949  rowNum, colNum, tv.getRows(), tv.getCols()));
950  }
951 
952  for (unsigned int j = 0; j < 3; ++j) {
953  t_out[j] = 0;
954  }
955 
956  for (unsigned int j = 0; j < 3; ++j) {
957  double tj = tv[j]; // optimization em 5/12/2006
958  for (unsigned int i = 0; i < 3; ++i) {
959  t_out[i] += rowPtrs[i][j] * tj;
960  }
961  }
962  return t_out;
963 }
964 
970 {
971  vpColVector v_out;
972  vpMatrix::multMatrixVector(*this, v, v_out);
973  return v_out;
974 }
975 
985 {
986  if (A.colNum != v.getRows()) {
987  throw(vpException(vpException::dimensionError, "Cannot multiply a (%dx%d) matrix by a (%d) column vector",
988  A.getRows(), A.getCols(), v.getRows()));
989  }
990 
991  if (A.rowNum != w.rowNum) {
992  w.resize(A.rowNum, false);
993  }
994 
995  // If available use Lapack only for large matrices
996  bool useLapack = ((A.rowNum > vpMatrix::m_lapack_min_size) || (A.colNum > vpMatrix::m_lapack_min_size));
997 #if !(defined(VISP_HAVE_LAPACK) && !defined(VISP_HAVE_LAPACK_BUILT_IN) && !defined(VISP_HAVE_GSL))
998  useLapack = false;
999 #endif
1000 
1001  if (useLapack) {
1002 #if defined(VISP_HAVE_LAPACK) && !defined(VISP_HAVE_LAPACK_BUILT_IN) && !defined(VISP_HAVE_GSL)
1003  double alpha = 1.0;
1004  double beta = 0.0;
1005  char trans = 't';
1006  int incr = 1;
1007 
1008  vpMatrix::blas_dgemv(trans, A.colNum, A.rowNum, alpha, A.data, A.colNum, v.data, incr, beta, w.data, incr);
1009 #endif
1010  }
1011  else {
1012  w = 0.0;
1013  for (unsigned int j = 0; j < A.colNum; ++j) {
1014  double vj = v[j]; // optimization em 5/12/2006
1015  for (unsigned int i = 0; i < A.rowNum; ++i) {
1016  w[i] += A.rowPtrs[i][j] * vj;
1017  }
1018  }
1019  }
1020 }
1021 
1022 //---------------------------------
1023 // Matrix operations.
1024 //---------------------------------
1025 
1036 {
1037  if ((A.rowNum != C.rowNum) || (B.colNum != C.colNum)) {
1038  C.resize(A.rowNum, B.colNum, false, false);
1039  }
1040 
1041  if (A.colNum != B.rowNum) {
1042  throw(vpException(vpException::dimensionError, "Cannot multiply (%dx%d) matrix by (%dx%d) matrix", A.getRows(),
1043  A.getCols(), B.getRows(), B.getCols()));
1044  }
1045 
1046  // If available use Lapack only for large matrices
1047  bool useLapack = ((A.rowNum > vpMatrix::m_lapack_min_size) || (A.colNum > vpMatrix::m_lapack_min_size) ||
1048  (B.colNum > vpMatrix::m_lapack_min_size));
1049 #if !(defined(VISP_HAVE_LAPACK) && !defined(VISP_HAVE_LAPACK_BUILT_IN) && !defined(VISP_HAVE_GSL))
1050  useLapack = false;
1051 #endif
1052 
1053  if (useLapack) {
1054 #if defined(VISP_HAVE_LAPACK) && !defined(VISP_HAVE_LAPACK_BUILT_IN) && !defined(VISP_HAVE_GSL)
1055  const double alpha = 1.0;
1056  const double beta = 0.0;
1057  const char trans = 'n';
1058  vpMatrix::blas_dgemm(trans, trans, B.colNum, A.rowNum, A.colNum, alpha, B.data, B.colNum, A.data, A.colNum, beta,
1059  C.data, B.colNum);
1060 #endif
1061  }
1062  else {
1063  // 5/12/06 some "very" simple optimization to avoid indexation
1064  const unsigned int BcolNum = B.colNum;
1065  const unsigned int BrowNum = B.rowNum;
1066  double **BrowPtrs = B.rowPtrs;
1067  for (unsigned int i = 0; i < A.rowNum; ++i) {
1068  const double *rowptri = A.rowPtrs[i];
1069  double *ci = C[i];
1070  for (unsigned int j = 0; j < BcolNum; ++j) {
1071  double s = 0;
1072  for (unsigned int k = 0; k < BrowNum; ++k) {
1073  s += rowptri[k] * BrowPtrs[k][j];
1074  }
1075  ci[j] = s;
1076  }
1077  }
1078  }
1079 }
1080 
1094 {
1095  if ((A.colNum != 3) || (A.rowNum != 3) || (B.colNum != 3) || (B.rowNum != 3)) {
1097  "Cannot multiply (%dx%d) matrix by (%dx%d) matrix as a "
1098  "rotation matrix",
1099  A.getRows(), A.getCols(), B.getRows(), B.getCols()));
1100  }
1101  // 5/12/06 some "very" simple optimization to avoid indexation
1102  const unsigned int BcolNum = B.colNum;
1103  const unsigned int BrowNum = B.rowNum;
1104  double **BrowPtrs = B.rowPtrs;
1105  for (unsigned int i = 0; i < A.rowNum; ++i) {
1106  const double *rowptri = A.rowPtrs[i];
1107  double *ci = C[i];
1108  for (unsigned int j = 0; j < BcolNum; ++j) {
1109  double s = 0;
1110  for (unsigned int k = 0; k < BrowNum; ++k) {
1111  s += rowptri[k] * BrowPtrs[k][j];
1112  }
1113  ci[j] = s;
1114  }
1115  }
1116 }
1117 
1131 {
1132  if ((A.colNum != 4) || (A.rowNum != 4) || (B.colNum != 4) || (B.rowNum != 4)) {
1134  "Cannot multiply (%dx%d) matrix by (%dx%d) matrix as a "
1135  "rotation matrix",
1136  A.getRows(), A.getCols(), B.getRows(), B.getCols()));
1137  }
1138  // Considering perfMatrixMultiplication.cpp benchmark,
1139  // using either MKL, OpenBLAS, or Netlib can slow down this function with respect to the naive code.
1140  // Lapack usage needs to be validated again.
1141  // If available use Lapack only for large matrices.
1142  // Using SSE2 doesn't speed up.
1143  bool useLapack = ((A.rowNum > vpMatrix::m_lapack_min_size) || (A.colNum > vpMatrix::m_lapack_min_size) ||
1144  (B.colNum > vpMatrix::m_lapack_min_size));
1145 #if !(defined(VISP_HAVE_LAPACK) && !defined(VISP_HAVE_LAPACK_BUILT_IN) && !defined(VISP_HAVE_GSL))
1146  useLapack = false;
1147 #endif
1148 
1149  if (useLapack) {
1150 #if defined(VISP_HAVE_LAPACK) && !defined(VISP_HAVE_LAPACK_BUILT_IN) && !defined(VISP_HAVE_GSL)
1151  const double alpha = 1.0;
1152  const double beta = 0.0;
1153  const char trans = 'n';
1154  vpMatrix::blas_dgemm(trans, trans, B.colNum, A.rowNum, A.colNum, alpha, B.data, B.colNum, A.data, A.colNum, beta,
1155  C.data, B.colNum);
1156 #endif
1157  }
1158  else {
1159  // 5/12/06 some "very" simple optimization to avoid indexation
1160  const unsigned int BcolNum = B.colNum;
1161  const unsigned int BrowNum = B.rowNum;
1162  double **BrowPtrs = B.rowPtrs;
1163  for (unsigned int i = 0; i < A.rowNum; ++i) {
1164  const double *rowptri = A.rowPtrs[i];
1165  double *ci = C[i];
1166  for (unsigned int j = 0; j < BcolNum; ++j) {
1167  double s = 0;
1168  for (unsigned int k = 0; k < BrowNum; ++k) {
1169  s += rowptri[k] * BrowPtrs[k][j];
1170  }
1171  ci[j] = s;
1172  }
1173  }
1174  }
1175 }
1176 
1190 {
1191  vpMatrix::multMatrixVector(A, B, C);
1192 }
1193 
1199 {
1200  vpMatrix C;
1201 
1202  vpMatrix::mult2Matrices(*this, B, C);
1203 
1204  return C;
1205 }
1206 
1212 {
1213  if (colNum != R.getRows()) {
1214  throw(vpException(vpException::dimensionError, "Cannot multiply (%dx%d) matrix by (3x3) rotation matrix", rowNum,
1215  colNum));
1216  }
1217  vpMatrix C;
1218  C.resize(rowNum, 3, false, false);
1219 
1220  unsigned int RcolNum = R.getCols();
1221  unsigned int RrowNum = R.getRows();
1222  for (unsigned int i = 0; i < rowNum; ++i) {
1223  double *rowptri = rowPtrs[i];
1224  double *ci = C[i];
1225  for (unsigned int j = 0; j < RcolNum; ++j) {
1226  double s = 0;
1227  for (unsigned int k = 0; k < RrowNum; ++k) {
1228  s += rowptri[k] * R[k][j];
1229  }
1230  ci[j] = s;
1231  }
1232  }
1233 
1234  return C;
1235 }
1236 
1242 {
1243  if (colNum != M.getRows()) {
1244  throw(vpException(vpException::dimensionError, "Cannot multiply (%dx%d) matrix by (3x3) rotation matrix", rowNum,
1245  colNum));
1246  }
1247  vpMatrix C;
1248  C.resize(rowNum, 4, false, false);
1249 
1250  const unsigned int McolNum = M.getCols();
1251  const unsigned int MrowNum = M.getRows();
1252  for (unsigned int i = 0; i < rowNum; ++i) {
1253  const double *rowptri = rowPtrs[i];
1254  double *ci = C[i];
1255  for (unsigned int j = 0; j < McolNum; ++j) {
1256  double s = 0;
1257  for (unsigned int k = 0; k < MrowNum; ++k) {
1258  s += rowptri[k] * M[k][j];
1259  }
1260  ci[j] = s;
1261  }
1262  }
1263 
1264  return C;
1265 }
1266 
1272 {
1273  if (colNum != V.getRows()) {
1274  throw(vpException(vpException::dimensionError, "Cannot multiply (%dx%d) matrix by (6x6) velocity twist matrix",
1275  rowNum, colNum));
1276  }
1277  vpMatrix M;
1278  M.resize(rowNum, 6, false, false);
1279 
1280  // Considering perfMatrixMultiplication.cpp benchmark,
1281  // using either MKL, OpenBLAS, or Netlib can slow down this function with respect to the naive code.
1282  // Lapack usage needs to be validated again.
1283  // If available use Lapack only for large matrices.
1284  // Speed up obtained using SSE2.
1285  bool useLapack = ((rowNum > vpMatrix::m_lapack_min_size) || (colNum > vpMatrix::m_lapack_min_size) ||
1286  (V.colNum > vpMatrix::m_lapack_min_size));
1287 #if !(defined(VISP_HAVE_LAPACK) && !defined(VISP_HAVE_LAPACK_BUILT_IN) && !defined(VISP_HAVE_GSL))
1288  useLapack = false;
1289 #endif
1290 
1291  if (useLapack) {
1292 #if defined(VISP_HAVE_LAPACK) && !defined(VISP_HAVE_LAPACK_BUILT_IN) && !defined(VISP_HAVE_GSL)
1293  const double alpha = 1.0;
1294  const double beta = 0.0;
1295  const char trans = 'n';
1296  vpMatrix::blas_dgemm(trans, trans, V.colNum, rowNum, colNum, alpha, V.data, V.colNum, data, colNum, beta, M.data,
1297  M.colNum);
1298 #endif
1299  }
1300  else {
1301 #if defined(VISP_HAVE_SIMDLIB)
1302  SimdMatMulTwist(data, rowNum, V.data, M.data);
1303 #else
1304  unsigned int VcolNum = V.getCols();
1305  unsigned int VrowNum = V.getRows();
1306  for (unsigned int i = 0; i < rowNum; ++i) {
1307  double *rowptri = rowPtrs[i];
1308  double *ci = M[i];
1309  for (unsigned int j = 0; j < VcolNum; ++j) {
1310  double s = 0;
1311  for (unsigned int k = 0; k < VrowNum; ++k) {
1312  s += rowptri[k] * V[k][j];
1313  }
1314  ci[j] = s;
1315  }
1316  }
1317 #endif
1318  }
1319 
1320  return M;
1321 }
1322 
1328 {
1329  if (colNum != V.getRows()) {
1330  throw(vpException(vpException::dimensionError, "Cannot multiply (%dx%d) matrix by (6x6) force/torque twist matrix",
1331  rowNum, colNum));
1332  }
1333  vpMatrix M;
1334  M.resize(rowNum, 6, false, false);
1335 
1336  // Considering perfMatrixMultiplication.cpp benchmark,
1337  // using either MKL, OpenBLAS, or Netlib can slow down this function with respect to the naive code.
1338  // Lapack usage needs to be validated again.
1339  // If available use Lapack only for large matrices.
1340  // Speed up obtained using SSE2.
1341  bool useLapack = ((rowNum > vpMatrix::m_lapack_min_size) || (colNum > vpMatrix::m_lapack_min_size) ||
1342  (V.getCols() > vpMatrix::m_lapack_min_size));
1343 #if !(defined(VISP_HAVE_LAPACK) && !defined(VISP_HAVE_LAPACK_BUILT_IN) && !defined(VISP_HAVE_GSL))
1344  useLapack = false;
1345 #endif
1346 
1347  if (useLapack) {
1348 #if defined(VISP_HAVE_LAPACK) && !defined(VISP_HAVE_LAPACK_BUILT_IN) && !defined(VISP_HAVE_GSL)
1349  const double alpha = 1.0;
1350  const double beta = 0.0;
1351  const char trans = 'n';
1352  vpMatrix::blas_dgemm(trans, trans, V.getCols(), rowNum, colNum, alpha, V.data, V.getCols(), data, colNum, beta,
1353  M.data, M.colNum);
1354 #endif
1355  }
1356  else {
1357 #if defined(VISP_HAVE_SIMDLIB)
1358  SimdMatMulTwist(data, rowNum, V.data, M.data);
1359 #else
1360  unsigned int VcolNum = V.getCols();
1361  unsigned int VrowNum = V.getRows();
1362  for (unsigned int i = 0; i < rowNum; ++i) {
1363  double *rowptri = rowPtrs[i];
1364  double *ci = M[i];
1365  for (unsigned int j = 0; j < VcolNum; ++j) {
1366  double s = 0;
1367  for (unsigned int k = 0; k < VrowNum; ++k) {
1368  s += rowptri[k] * V[k][j];
1369  }
1370  ci[j] = s;
1371  }
1372  }
1373 #endif
1374  }
1375 
1376  return M;
1377 }
1378 
1389 void vpMatrix::add2WeightedMatrices(const vpMatrix &A, const double &wA, const vpMatrix &B, const double &wB,
1390  vpMatrix &C)
1391 {
1392  if ((A.rowNum != C.rowNum) || (B.colNum != C.colNum)) {
1393  C.resize(A.rowNum, B.colNum, false, false);
1394  }
1395 
1396  if ((A.colNum != B.getCols()) || (A.rowNum != B.getRows())) {
1397  throw(vpException(vpException::dimensionError, "Cannot add (%dx%d) matrix with (%dx%d) matrix", A.getRows(),
1398  A.getCols(), B.getRows(), B.getCols()));
1399  }
1400 
1401  double **ArowPtrs = A.rowPtrs;
1402  double **BrowPtrs = B.rowPtrs;
1403  double **CrowPtrs = C.rowPtrs;
1404 
1405  for (unsigned int i = 0; i < A.rowNum; ++i) {
1406  for (unsigned int j = 0; j < A.colNum; ++j) {
1407  CrowPtrs[i][j] = (wB * BrowPtrs[i][j]) + (wA * ArowPtrs[i][j]);
1408  }
1409  }
1410 }
1411 
1421 void vpMatrix::add2Matrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C)
1422 {
1423  if ((A.rowNum != C.rowNum) || (B.colNum != C.colNum)) {
1424  C.resize(A.rowNum, B.colNum, false, false);
1425  }
1426 
1427  if ((A.colNum != B.getCols()) || (A.rowNum != B.getRows())) {
1428  throw(vpException(vpException::dimensionError, "Cannot add (%dx%d) matrix with (%dx%d) matrix", A.getRows(),
1429  A.getCols(), B.getRows(), B.getCols()));
1430  }
1431 
1432  double **ArowPtrs = A.rowPtrs;
1433  double **BrowPtrs = B.rowPtrs;
1434  double **CrowPtrs = C.rowPtrs;
1435 
1436  for (unsigned int i = 0; i < A.rowNum; ++i) {
1437  for (unsigned int j = 0; j < A.colNum; ++j) {
1438  CrowPtrs[i][j] = BrowPtrs[i][j] + ArowPtrs[i][j];
1439  }
1440  }
1441 }
1442 
1456 {
1457  if ((A.rowNum != C.rowNum) || (B.colNum != C.colNum)) {
1458  C.resize(A.rowNum);
1459  }
1460 
1461  if ((A.colNum != B.getCols()) || (A.rowNum != B.getRows())) {
1462  throw(vpException(vpException::dimensionError, "Cannot add (%dx%d) matrix with (%dx%d) matrix", A.getRows(),
1463  A.getCols(), B.getRows(), B.getCols()));
1464  }
1465 
1466  double **ArowPtrs = A.rowPtrs;
1467  double **BrowPtrs = B.rowPtrs;
1468  double **CrowPtrs = C.rowPtrs;
1469 
1470  for (unsigned int i = 0; i < A.rowNum; ++i) {
1471  for (unsigned int j = 0; j < A.colNum; ++j) {
1472  CrowPtrs[i][j] = BrowPtrs[i][j] + ArowPtrs[i][j];
1473  }
1474  }
1475 }
1476 
1482 {
1483  vpMatrix C;
1484  vpMatrix::add2Matrices(*this, B, C);
1485  return C;
1486 }
1487 
1504 {
1505  if ((A.rowNum != C.rowNum) || (A.colNum != C.colNum)) {
1506  C.resize(A.rowNum);
1507  }
1508 
1509  if ((A.colNum != B.getCols()) || (A.rowNum != B.getRows())) {
1510  throw(vpException(vpException::dimensionError, "Cannot subtract (%dx%d) matrix to (%dx%d) matrix", A.getRows(),
1511  A.getCols(), B.getRows(), B.getCols()));
1512  }
1513 
1514  double **ArowPtrs = A.rowPtrs;
1515  double **BrowPtrs = B.rowPtrs;
1516  double **CrowPtrs = C.rowPtrs;
1517 
1518  for (unsigned int i = 0; i < A.rowNum; ++i) {
1519  for (unsigned int j = 0; j < A.colNum; ++j) {
1520  CrowPtrs[i][j] = ArowPtrs[i][j] - BrowPtrs[i][j];
1521  }
1522  }
1523 }
1524 
1537 void vpMatrix::sub2Matrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C)
1538 {
1539  if ((A.rowNum != C.rowNum) || (A.colNum != C.colNum)) {
1540  C.resize(A.rowNum, A.colNum, false, false);
1541  }
1542 
1543  if ((A.colNum != B.getCols()) || (A.rowNum != B.getRows())) {
1544  throw(vpException(vpException::dimensionError, "Cannot subtract (%dx%d) matrix to (%dx%d) matrix", A.getRows(),
1545  A.getCols(), B.getRows(), B.getCols()));
1546  }
1547 
1548  double **ArowPtrs = A.rowPtrs;
1549  double **BrowPtrs = B.rowPtrs;
1550  double **CrowPtrs = C.rowPtrs;
1551 
1552  for (unsigned int i = 0; i < A.rowNum; ++i) {
1553  for (unsigned int j = 0; j < A.colNum; ++j) {
1554  CrowPtrs[i][j] = ArowPtrs[i][j] - BrowPtrs[i][j];
1555  }
1556  }
1557 }
1558 
1564 {
1565  vpMatrix C;
1566  vpMatrix::sub2Matrices(*this, B, C);
1567  return C;
1568 }
1569 
1571 
1573 {
1574  if ((colNum != B.getCols()) || (rowNum != B.getRows())) {
1575  throw(vpException(vpException::dimensionError, "Cannot add (%dx%d) matrix to (%dx%d) matrix", rowNum, colNum,
1576  B.getRows(), B.getCols()));
1577  }
1578 
1579  double **BrowPtrs = B.rowPtrs;
1580 
1581  for (unsigned int i = 0; i < rowNum; ++i) {
1582  for (unsigned int j = 0; j < colNum; ++j) {
1583  rowPtrs[i][j] += BrowPtrs[i][j];
1584  }
1585  }
1586 
1587  return *this;
1588 }
1589 
1592 {
1593  if ((colNum != B.getCols()) || (rowNum != B.getRows())) {
1594  throw(vpException(vpException::dimensionError, "Cannot subtract (%dx%d) matrix to (%dx%d) matrix", rowNum, colNum,
1595  B.getRows(), B.getCols()));
1596  }
1597 
1598  double **BrowPtrs = B.rowPtrs;
1599  for (unsigned int i = 0; i < rowNum; ++i) {
1600  for (unsigned int j = 0; j < colNum; ++j) {
1601  rowPtrs[i][j] -= BrowPtrs[i][j];
1602  }
1603  }
1604 
1605  return *this;
1606 }
1607 
1618 {
1619  if ((A.rowNum != C.rowNum) || (A.colNum != C.colNum)) {
1620  C.resize(A.rowNum, A.colNum, false, false);
1621  }
1622 
1623  double **ArowPtrs = A.rowPtrs;
1624  double **CrowPtrs = C.rowPtrs;
1625 
1626  for (unsigned int i = 0; i < A.rowNum; ++i) {
1627  for (unsigned int j = 0; j < A.colNum; ++j) {
1628  CrowPtrs[i][j] = -ArowPtrs[i][j];
1629  }
1630  }
1631 }
1632 
1638 {
1639  vpMatrix C;
1640  vpMatrix::negateMatrix(*this, C);
1641  return C;
1642 }
1643 
1644 double vpMatrix::sum() const
1645 {
1646  double s = 0.0;
1647  for (unsigned int i = 0; i < rowNum; ++i) {
1648  for (unsigned int j = 0; j < colNum; ++j) {
1649  s += rowPtrs[i][j];
1650  }
1651  }
1652 
1653  return s;
1654 }
1655 
1656 //---------------------------------
1657 // Matrix/vector operations.
1658 //---------------------------------
1659 
1660 //---------------------------------
1661 // Matrix/real operations.
1662 //---------------------------------
1663 
1668 vpMatrix operator*(const double &x, const vpMatrix &B)
1669 {
1670  if (std::fabs(x - 1.) < std::numeric_limits<double>::epsilon()) {
1671  return B;
1672  }
1673 
1674  unsigned int Brow = B.getRows();
1675  unsigned int Bcol = B.getCols();
1676 
1677  vpMatrix C;
1678  C.resize(Brow, Bcol, false, false);
1679 
1680  for (unsigned int i = 0; i < Brow; ++i) {
1681  for (unsigned int j = 0; j < Bcol; ++j) {
1682  C[i][j] = B[i][j] * x;
1683  }
1684  }
1685 
1686  return C;
1687 }
1688 
1694 {
1695  if (std::fabs(x - 1.) < std::numeric_limits<double>::epsilon()) {
1696  return *this;
1697  }
1698 
1699  vpMatrix M;
1700  M.resize(rowNum, colNum, false, false);
1701 
1702  for (unsigned int i = 0; i < rowNum; ++i) {
1703  for (unsigned int j = 0; j < colNum; ++j) {
1704  M[i][j] = rowPtrs[i][j] * x;
1705  }
1706  }
1707 
1708  return M;
1709 }
1710 
1713 {
1714  if (std::fabs(x - 1.) < std::numeric_limits<double>::epsilon()) {
1715  return *this;
1716  }
1717 
1718  if (std::fabs(x) < std::numeric_limits<double>::epsilon()) {
1719  throw vpException(vpException::divideByZeroError, "Divide matrix by zero scalar");
1720  }
1721 
1722  vpMatrix C;
1723  C.resize(rowNum, colNum, false, false);
1724 
1725  double xinv = 1 / x;
1726 
1727  for (unsigned int i = 0; i < rowNum; ++i) {
1728  for (unsigned int j = 0; j < colNum; ++j) {
1729  C[i][j] = rowPtrs[i][j] * xinv;
1730  }
1731  }
1732 
1733  return C;
1734 }
1735 
1738 {
1739  for (unsigned int i = 0; i < rowNum; ++i) {
1740  for (unsigned int j = 0; j < colNum; ++j) {
1741  rowPtrs[i][j] += x;
1742  }
1743  }
1744 
1745  return *this;
1746 }
1747 
1750 {
1751  for (unsigned int i = 0; i < rowNum; ++i) {
1752  for (unsigned int j = 0; j < colNum; ++j) {
1753  rowPtrs[i][j] -= x;
1754  }
1755  }
1756 
1757  return *this;
1758 }
1759 
1765 {
1766  if (std::fabs(x - 1.) < std::numeric_limits<double>::epsilon()) {
1767  return *this;
1768  }
1769 
1770  for (unsigned int i = 0; i < rowNum; ++i) {
1771  for (unsigned int j = 0; j < colNum; ++j) {
1772  rowPtrs[i][j] *= x;
1773  }
1774  }
1775 
1776  return *this;
1777 }
1778 
1781 {
1782  if (std::fabs(x - 1.) < std::numeric_limits<double>::epsilon()) {
1783  return *this;
1784  }
1785 
1786  if (std::fabs(x) < std::numeric_limits<double>::epsilon()) {
1787  throw vpException(vpException::divideByZeroError, "Divide matrix by zero scalar");
1788  }
1789 
1790  double xinv = 1 / x;
1791 
1792  for (unsigned int i = 0; i < rowNum; ++i) {
1793  for (unsigned int j = 0; j < colNum; ++j) {
1794  rowPtrs[i][j] *= xinv;
1795  }
1796  }
1797 
1798  return *this;
1799 }
1800 
1801 //----------------------------------------------------------------
1802 // Matrix Operation
1803 //----------------------------------------------------------------
1804 
1810 {
1811  if ((out.rowNum != (colNum * rowNum)) || (out.colNum != 1)) {
1812  out.resize(colNum * rowNum, false, false);
1813  }
1814 
1815  double *optr = out.data;
1816  for (unsigned int j = 0; j < colNum; ++j) {
1817  for (unsigned int i = 0; i < rowNum; ++i) {
1818  *(optr++) = rowPtrs[i][j];
1819  }
1820  }
1821 }
1822 
1828 {
1829  vpColVector out(colNum * rowNum);
1830  stackColumns(out);
1831  return out;
1832 }
1833 
1839 {
1840  if ((out.getRows() != 1) || (out.getCols() != (colNum * rowNum))) {
1841  out.resize(colNum * rowNum, false, false);
1842  }
1843 
1844  memcpy(out.data, data, sizeof(double) * out.getCols());
1845 }
1851 {
1852  vpRowVector out(colNum * rowNum);
1853  stackRows(out);
1854  return out;
1855 }
1856 
1864 {
1865  if ((m.getRows() != rowNum) || (m.getCols() != colNum)) {
1866  throw(vpException(vpException::dimensionError, "In Hadamard product: bad dimension of input matrix"));
1867  }
1868 
1869  vpMatrix out;
1870  out.resize(rowNum, colNum, false, false);
1871 
1872 #if defined(VISP_HAVE_SIMDLIB)
1873  SimdVectorHadamard(data, m.data, dsize, out.data);
1874 #else
1875  for (unsigned int i = 0; i < dsize; ++i) {
1876  out.data[i] = data[i] * m.data[i];
1877  }
1878 #endif
1879 
1880  return out;
1881 }
1882 
1889 void vpMatrix::kron(const vpMatrix &m1, const vpMatrix &m2, vpMatrix &out)
1890 {
1891  unsigned int r1 = m1.getRows();
1892  unsigned int c1 = m1.getCols();
1893  unsigned int r2 = m2.getRows();
1894  unsigned int c2 = m2.getCols();
1895 
1896  out.resize(r1 * r2, c1 * c2, false, false);
1897 
1898  for (unsigned int r = 0; r < r1; ++r) {
1899  for (unsigned int c = 0; c < c1; ++c) {
1900  double alpha = m1[r][c];
1901  double *m2ptr = m2[0];
1902  unsigned int roffset = r * r2;
1903  unsigned int coffset = c * c2;
1904  for (unsigned int rr = 0; rr < r2; ++rr) {
1905  for (unsigned int cc = 0; cc < c2; ++cc) {
1906  out[roffset + rr][coffset + cc] = alpha * *(m2ptr++);
1907  }
1908  }
1909  }
1910  }
1911 }
1912 
1919 void vpMatrix::kron(const vpMatrix &m, vpMatrix &out) const { kron(*this, m, out); }
1920 
1928 {
1929  unsigned int r1 = m1.getRows();
1930  unsigned int c1 = m1.getCols();
1931  unsigned int r2 = m2.getRows();
1932  unsigned int c2 = m2.getCols();
1933 
1934  vpMatrix out;
1935  out.resize(r1 * r2, c1 * c2, false, false);
1936 
1937  for (unsigned int r = 0; r < r1; ++r) {
1938  for (unsigned int c = 0; c < c1; ++c) {
1939  double alpha = m1[r][c];
1940  double *m2ptr = m2[0];
1941  unsigned int roffset = r * r2;
1942  unsigned int coffset = c * c2;
1943  for (unsigned int rr = 0; rr < r2; ++rr) {
1944  for (unsigned int cc = 0; cc < c2; ++cc) {
1945  out[roffset + rr][coffset + cc] = alpha * *(m2ptr++);
1946  }
1947  }
1948  }
1949  }
1950  return out;
1951 }
1952 
1958 vpMatrix vpMatrix::kron(const vpMatrix &m) const { return kron(*this, m); }
1959 
2010 void vpMatrix::solveBySVD(const vpColVector &b, vpColVector &x) const { x = pseudoInverse(1e-6) * b; }
2011 
2062 {
2063  vpColVector X(colNum);
2064 
2065  solveBySVD(B, X);
2066  return X;
2067 }
2068 
2133 {
2134 #if defined(VISP_HAVE_LAPACK)
2135  svdLapack(w, V);
2136 #elif defined(VISP_HAVE_EIGEN3)
2137  svdEigen3(w, V);
2138 #elif defined(VISP_HAVE_OPENCV) // Require opencv >= 2.1.1
2139  svdOpenCV(w, V);
2140 #else
2141  (void)w;
2142  (void)V;
2143  throw(vpException(vpException::fatalError, "Cannot compute SVD. Install Lapack, Eigen3 or OpenCV 3rd party"));
2144 #endif
2145 }
2146 
2201 unsigned int vpMatrix::pseudoInverse(vpMatrix &Ap, double svThreshold) const
2202 {
2203 #if defined(VISP_HAVE_LAPACK)
2204  return pseudoInverseLapack(Ap, svThreshold);
2205 #elif defined(VISP_HAVE_EIGEN3)
2206  return pseudoInverseEigen3(Ap, svThreshold);
2207 #elif defined(VISP_HAVE_OPENCV) // Require opencv >= 2.1.1
2208  return pseudoInverseOpenCV(Ap, svThreshold);
2209 #else
2210  (void)Ap;
2211  (void)svThreshold;
2212  throw(vpException(vpException::fatalError, "Cannot compute pseudo-inverse. "
2213  "Install Lapack, Eigen3 or OpenCV 3rd party"));
2214 #endif
2215 }
2216 
2277 int vpMatrix::pseudoInverse(vpMatrix &Ap, int rank_in) const
2278 {
2279 #if defined(VISP_HAVE_LAPACK)
2280  return pseudoInverseLapack(Ap, rank_in);
2281 #elif defined(VISP_HAVE_EIGEN3)
2282  return pseudoInverseEigen3(Ap, rank_in);
2283 #elif defined(VISP_HAVE_OPENCV) // Require opencv >= 2.1.1
2284  return pseudoInverseOpenCV(Ap, rank_in);
2285 #else
2286  (void)Ap;
2287  (void)rank_in;
2288  throw(vpException(vpException::fatalError, "Cannot compute pseudo-inverse. "
2289  "Install Lapack, Eigen3 or OpenCV 3rd party"));
2290 #endif
2291 }
2292 
2343 vpMatrix vpMatrix::pseudoInverse(double svThreshold) const
2344 {
2345 #if defined(VISP_HAVE_LAPACK)
2346  return pseudoInverseLapack(svThreshold);
2347 #elif defined(VISP_HAVE_EIGEN3)
2348  return pseudoInverseEigen3(svThreshold);
2349 #elif defined(VISP_HAVE_OPENCV) // Require opencv >= 2.1.1
2350  return pseudoInverseOpenCV(svThreshold);
2351 #else
2352  (void)svThreshold;
2353  throw(vpException(vpException::fatalError, "Cannot compute pseudo-inverse. "
2354  "Install Lapack, Eigen3 or OpenCV 3rd party"));
2355 #endif
2356 }
2357 
2409 {
2410 #if defined(VISP_HAVE_LAPACK)
2411  return pseudoInverseLapack(rank_in);
2412 #elif defined(VISP_HAVE_EIGEN3)
2413  return pseudoInverseEigen3(rank_in);
2414 #elif defined(VISP_HAVE_OPENCV) // Require opencv >= 2.1.1
2415  return pseudoInverseOpenCV(rank_in);
2416 #else
2417  (void)rank_in;
2418  throw(vpException(vpException::fatalError, "Cannot compute pseudo-inverse. "
2419  "Install Lapack, Eigen3 or OpenCV 3rd party"));
2420 #endif
2421 }
2422 
2423 #ifndef DOXYGEN_SHOULD_SKIP_THIS
2424 #if defined(VISP_HAVE_LAPACK)
2461 vpMatrix vpMatrix::pseudoInverseLapack(double svThreshold) const
2462 {
2463  unsigned int nrows = getRows();
2464  unsigned int ncols = getCols();
2465  int rank_out;
2466 
2467  vpMatrix U, V, Ap;
2468  vpColVector sv;
2469 
2470  U = *this;
2471  U.svdLapack(sv, V);
2472 
2473  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, nullptr, nullptr, nullptr, nullptr);
2474 
2475  return Ap;
2476 }
2477 
2518 unsigned int vpMatrix::pseudoInverseLapack(vpMatrix &Ap, double svThreshold) const
2519 {
2520  unsigned int nrows = getRows();
2521  unsigned int ncols = getCols();
2522  int rank_out;
2523 
2524  vpMatrix U, V;
2525  vpColVector sv;
2526 
2527  U = *this;
2528  U.svdLapack(sv, V);
2529 
2530  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, nullptr, nullptr, nullptr, nullptr);
2531 
2532  return static_cast<unsigned int>(rank_out);
2533 }
2581 unsigned int vpMatrix::pseudoInverseLapack(vpMatrix &Ap, vpColVector &sv, double svThreshold) const
2582 {
2583  unsigned int nrows = getRows();
2584  unsigned int ncols = getCols();
2585  int rank_out;
2586 
2587  vpMatrix U, V;
2588 
2589  Ap.resize(ncols, nrows, false, false);
2590 
2591  U = *this;
2592  U.svdLapack(sv, V);
2593 
2594  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, nullptr, nullptr, nullptr, nullptr);
2595 
2596  return static_cast<unsigned int>(rank_out);
2597 }
2598 
2705 unsigned int vpMatrix::pseudoInverseLapack(vpMatrix &Ap, vpColVector &sv, double svThreshold, vpMatrix &imA,
2706  vpMatrix &imAt, vpMatrix &kerAt) const
2707 {
2708  unsigned int nrows = getRows();
2709  unsigned int ncols = getCols();
2710  int rank_out;
2711  vpMatrix U, V;
2712 
2713  if (nrows < ncols) {
2714  U.resize(ncols, ncols, true);
2715  U.insert(*this, 0, 0);
2716  }
2717  else {
2718  U = *this;
2719  }
2720 
2721  U.svdLapack(sv, V);
2722 
2723  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, nullptr, &imA, &imAt, &kerAt);
2724 
2725  return static_cast<unsigned int>(rank_out);
2726 }
2727 
2764 vpMatrix vpMatrix::pseudoInverseLapack(int rank_in) const
2765 {
2766  unsigned int nrows = getRows();
2767  unsigned int ncols = getCols();
2768  int rank_out;
2769  double svThreshold = 1e-26;
2770 
2771  vpMatrix U, V, Ap;
2772  vpColVector sv;
2773 
2774  U = *this;
2775  U.svdLapack(sv, V);
2776 
2777  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, &rank_in, nullptr, nullptr, nullptr);
2778 
2779  return Ap;
2780 }
2781 
2828 int vpMatrix::pseudoInverseLapack(vpMatrix &Ap, int rank_in) const
2829 {
2830  unsigned int nrows = getRows();
2831  unsigned int ncols = getCols();
2832  int rank_out;
2833  double svThreshold = 1e-26;
2834 
2835  vpMatrix U, V;
2836  vpColVector sv;
2837 
2838  U = *this;
2839  U.svdLapack(sv, V);
2840 
2841  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, &rank_in, nullptr, nullptr, nullptr);
2842 
2843  return rank_out;
2844 }
2845 
2899 int vpMatrix::pseudoInverseLapack(vpMatrix &Ap, vpColVector &sv, int rank_in) const
2900 {
2901  unsigned int nrows = getRows();
2902  unsigned int ncols = getCols();
2903  int rank_out;
2904  double svThreshold = 1e-26;
2905  vpMatrix U, V;
2906 
2907  Ap.resize(ncols, nrows, false, false);
2908 
2909  U = *this;
2910  U.svdLapack(sv, V);
2911 
2912  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, &rank_in, nullptr, nullptr, nullptr);
2913 
2914  return rank_out;
2915 }
2916 
3028 int vpMatrix::pseudoInverseLapack(vpMatrix &Ap, vpColVector &sv, int rank_in, vpMatrix &imA, vpMatrix &imAt,
3029  vpMatrix &kerAt) const
3030 {
3031  unsigned int nrows = getRows();
3032  unsigned int ncols = getCols();
3033  int rank_out;
3034  double svThreshold = 1e-26;
3035  vpMatrix U, V;
3036 
3037  if (nrows < ncols) {
3038  U.resize(ncols, ncols, true);
3039  U.insert(*this, 0, 0);
3040  }
3041  else {
3042  U = *this;
3043  }
3044 
3045  U.svdLapack(sv, V);
3046 
3047  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, &rank_in, &imA, &imAt, &kerAt);
3048 
3049  return rank_out;
3050 }
3051 
3052 #endif
3053 #if defined(VISP_HAVE_EIGEN3)
3090 vpMatrix vpMatrix::pseudoInverseEigen3(double svThreshold) const
3091 {
3092  unsigned int nrows = getRows();
3093  unsigned int ncols = getCols();
3094  int rank_out;
3095  vpMatrix U, V, Ap;
3096  vpColVector sv;
3097 
3098  U = *this;
3099  U.svdEigen3(sv, V);
3100 
3101  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, nullptr, nullptr, nullptr, nullptr);
3102 
3103  return Ap;
3104 }
3105 
3146 unsigned int vpMatrix::pseudoInverseEigen3(vpMatrix &Ap, double svThreshold) const
3147 {
3148  unsigned int nrows = getRows();
3149  unsigned int ncols = getCols();
3150  int rank_out;
3151  vpMatrix U, V;
3152  vpColVector sv;
3153 
3154  U = *this;
3155  U.svdEigen3(sv, V);
3156 
3157  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, nullptr, nullptr, nullptr, nullptr);
3158 
3159  return static_cast<unsigned int>(rank_out);
3160 }
3161 
3209 unsigned int vpMatrix::pseudoInverseEigen3(vpMatrix &Ap, vpColVector &sv, double svThreshold) const
3210 {
3211  unsigned int nrows = getRows();
3212  unsigned int ncols = getCols();
3213  int rank_out;
3214  vpMatrix U, V;
3215 
3216  Ap.resize(ncols, nrows, false, false);
3217 
3218  U = *this;
3219  U.svdEigen3(sv, V);
3220 
3221  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, nullptr, nullptr, nullptr, nullptr);
3222 
3223  return static_cast<unsigned int>(rank_out);
3224 }
3225 
3332 unsigned int vpMatrix::pseudoInverseEigen3(vpMatrix &Ap, vpColVector &sv, double svThreshold, vpMatrix &imA,
3333  vpMatrix &imAt, vpMatrix &kerAt) const
3334 {
3335  unsigned int nrows = getRows();
3336  unsigned int ncols = getCols();
3337  int rank_out;
3338  vpMatrix U, V;
3339 
3340  if (nrows < ncols) {
3341  U.resize(ncols, ncols, true);
3342  U.insert(*this, 0, 0);
3343  }
3344  else {
3345  U = *this;
3346  }
3347 
3348  U.svdEigen3(sv, V);
3349 
3350  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, nullptr, &imA, &imAt, &kerAt);
3351 
3352  return static_cast<unsigned int>(rank_out);
3353 }
3354 
3391 vpMatrix vpMatrix::pseudoInverseEigen3(int rank_in) const
3392 {
3393  unsigned int nrows = getRows();
3394  unsigned int ncols = getCols();
3395  int rank_out;
3396  double svThreshold = 1e-26;
3397  vpMatrix U, V, Ap;
3398  vpColVector sv;
3399 
3400  U = *this;
3401  U.svdEigen3(sv, V);
3402 
3403  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, &rank_in, nullptr, nullptr, nullptr);
3404 
3405  return Ap;
3406 }
3407 
3454 int vpMatrix::pseudoInverseEigen3(vpMatrix &Ap, int rank_in) const
3455 {
3456  unsigned int nrows = getRows();
3457  unsigned int ncols = getCols();
3458  int rank_out;
3459  double svThreshold = 1e-26;
3460  vpMatrix U, V;
3461  vpColVector sv;
3462 
3463  Ap.resize(ncols, nrows, false, false);
3464 
3465  U = *this;
3466  U.svdEigen3(sv, V);
3467 
3468  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, &rank_in, nullptr, nullptr, nullptr);
3469 
3470  return rank_out;
3471 }
3472 
3526 int vpMatrix::pseudoInverseEigen3(vpMatrix &Ap, vpColVector &sv, int rank_in) const
3527 {
3528  unsigned int nrows = getRows();
3529  unsigned int ncols = getCols();
3530  int rank_out;
3531  double svThreshold = 1e-26;
3532  vpMatrix U, V;
3533 
3534  Ap.resize(ncols, nrows, false, false);
3535 
3536  U = *this;
3537  U.svdEigen3(sv, V);
3538 
3539  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, &rank_in, nullptr, nullptr, nullptr);
3540 
3541  return rank_out;
3542 }
3543 
3655 int vpMatrix::pseudoInverseEigen3(vpMatrix &Ap, vpColVector &sv, int rank_in, vpMatrix &imA, vpMatrix &imAt,
3656  vpMatrix &kerAt) const
3657 {
3658  unsigned int nrows = getRows();
3659  unsigned int ncols = getCols();
3660  int rank_out;
3661  double svThreshold = 1e-26;
3662  vpMatrix U, V;
3663 
3664  if (nrows < ncols) {
3665  U.resize(ncols, ncols, true);
3666  U.insert(*this, 0, 0);
3667  }
3668  else {
3669  U = *this;
3670  }
3671  U.svdEigen3(sv, V);
3672 
3673  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, &rank_in, &imA, &imAt, &kerAt);
3674 
3675  return rank_out;
3676 }
3677 
3678 #endif
3679 #if defined(VISP_HAVE_OPENCV)
3716 vpMatrix vpMatrix::pseudoInverseOpenCV(double svThreshold) const
3717 {
3718  unsigned int nrows = getRows();
3719  unsigned int ncols = getCols();
3720  int rank_out;
3721  vpMatrix U, V, Ap;
3722  vpColVector sv;
3723 
3724  U = *this;
3725  U.svdOpenCV(sv, V);
3726 
3727  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, nullptr, nullptr, nullptr, nullptr);
3728 
3729  return Ap;
3730 }
3731 
3772 unsigned int vpMatrix::pseudoInverseOpenCV(vpMatrix &Ap, double svThreshold) const
3773 {
3774  unsigned int nrows = getRows();
3775  unsigned int ncols = getCols();
3776  int rank_out;
3777  vpMatrix U, V;
3778  vpColVector sv;
3779 
3780  Ap.resize(ncols, nrows, false, false);
3781 
3782  U = *this;
3783  U.svdOpenCV(sv, V);
3784 
3785  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, nullptr, nullptr, nullptr, nullptr);
3786 
3787  return static_cast<unsigned int>(rank_out);
3788 }
3836 unsigned int vpMatrix::pseudoInverseOpenCV(vpMatrix &Ap, vpColVector &sv, double svThreshold) const
3837 {
3838  unsigned int nrows = getRows();
3839  unsigned int ncols = getCols();
3840  int rank_out;
3841  vpMatrix U, V;
3842 
3843  Ap.resize(ncols, nrows, false, false);
3844 
3845  U = *this;
3846  U.svdOpenCV(sv, V);
3847 
3848  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, nullptr, nullptr, nullptr, nullptr);
3849 
3850  return static_cast<unsigned int>(rank_out);
3851 }
3852 
3959 unsigned int vpMatrix::pseudoInverseOpenCV(vpMatrix &Ap, vpColVector &sv, double svThreshold, vpMatrix &imA,
3960  vpMatrix &imAt, vpMatrix &kerAt) const
3961 {
3962  unsigned int nrows = getRows();
3963  unsigned int ncols = getCols();
3964  int rank_out;
3965  vpMatrix U, V;
3966 
3967  if (nrows < ncols) {
3968  U.resize(ncols, ncols, true);
3969  U.insert(*this, 0, 0);
3970  }
3971  else {
3972  U = *this;
3973  }
3974 
3975  U.svdOpenCV(sv, V);
3976 
3977  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, nullptr, &imA, &imAt, &kerAt);
3978 
3979  return static_cast<unsigned int>(rank_out);
3980 }
3981 
4018 vpMatrix vpMatrix::pseudoInverseOpenCV(int rank_in) const
4019 {
4020  unsigned int nrows = getRows();
4021  unsigned int ncols = getCols();
4022  int rank_out;
4023  double svThreshold = 1e-26;
4024  vpMatrix U, V, Ap;
4025  vpColVector sv;
4026 
4027  U = *this;
4028  U.svdOpenCV(sv, V);
4029 
4030  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, &rank_in, nullptr, nullptr, nullptr);
4031 
4032  return Ap;
4033 }
4034 
4081 int vpMatrix::pseudoInverseOpenCV(vpMatrix &Ap, int rank_in) const
4082 {
4083  unsigned int nrows = getRows();
4084  unsigned int ncols = getCols();
4085  int rank_out;
4086  double svThreshold = 1e-26;
4087  vpMatrix U, V;
4088  vpColVector sv;
4089 
4090  Ap.resize(ncols, nrows, false, false);
4091 
4092  U = *this;
4093  U.svdOpenCV(sv, V);
4094 
4095  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, &rank_in, nullptr, nullptr, nullptr);
4096 
4097  return rank_out;
4098 }
4099 
4153 int vpMatrix::pseudoInverseOpenCV(vpMatrix &Ap, vpColVector &sv, int rank_in) const
4154 {
4155  unsigned int nrows = getRows();
4156  unsigned int ncols = getCols();
4157  int rank_out;
4158  double svThreshold = 1e-26;
4159  vpMatrix U, V;
4160 
4161  Ap.resize(ncols, nrows, false, false);
4162 
4163  U = *this;
4164  U.svdOpenCV(sv, V);
4165 
4166  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, &rank_in, nullptr, nullptr, nullptr);
4167 
4168  return rank_out;
4169 }
4170 
4282 int vpMatrix::pseudoInverseOpenCV(vpMatrix &Ap, vpColVector &sv, int rank_in, vpMatrix &imA, vpMatrix &imAt,
4283  vpMatrix &kerAt) const
4284 {
4285  unsigned int nrows = getRows();
4286  unsigned int ncols = getCols();
4287  int rank_out;
4288  double svThreshold = 1e-26;
4289  vpMatrix U, V;
4290 
4291  if (nrows < ncols) {
4292  U.resize(ncols, ncols, true);
4293  U.insert(*this, 0, 0);
4294  }
4295  else {
4296  U = *this;
4297  }
4298 
4299  U.svdOpenCV(sv, V);
4300 
4301  compute_pseudo_inverse(U, sv, V, nrows, ncols, svThreshold, Ap, rank_out, &rank_in, &imA, &imAt, &kerAt);
4302 
4303  return rank_out;
4304 }
4305 
4306 #endif
4307 #endif // #ifndef DOXYGEN_SHOULD_SKIP_THIS
4308 
4370 unsigned int vpMatrix::pseudoInverse(vpMatrix &Ap, vpColVector &sv, double svThreshold) const
4371 {
4372 #if defined(VISP_HAVE_LAPACK)
4373  return pseudoInverseLapack(Ap, sv, svThreshold);
4374 #elif defined(VISP_HAVE_EIGEN3)
4375  return pseudoInverseEigen3(Ap, sv, svThreshold);
4376 #elif defined(VISP_HAVE_OPENCV) // Require opencv >= 2.1.1
4377  return pseudoInverseOpenCV(Ap, sv, svThreshold);
4378 #else
4379  (void)Ap;
4380  (void)sv;
4381  (void)svThreshold;
4382  throw(vpException(vpException::fatalError, "Cannot compute pseudo-inverse. "
4383  "Install Lapack, Eigen3 or OpenCV 3rd party"));
4384 #endif
4385 }
4386 
4453 int vpMatrix::pseudoInverse(vpMatrix &Ap, vpColVector &sv, int rank_in) const
4454 {
4455 #if defined(VISP_HAVE_LAPACK)
4456  return pseudoInverseLapack(Ap, sv, rank_in);
4457 #elif defined(VISP_HAVE_EIGEN3)
4458  return pseudoInverseEigen3(Ap, sv, rank_in);
4459 #elif defined(VISP_HAVE_OPENCV) // Require opencv >= 2.1.1
4460  return pseudoInverseOpenCV(Ap, sv, rank_in);
4461 #else
4462  (void)Ap;
4463  (void)sv;
4464  (void)rank_in;
4465  throw(vpException(vpException::fatalError, "Cannot compute pseudo-inverse. "
4466  "Install Lapack, Eigen3 or OpenCV 3rd party"));
4467 #endif
4468 }
4543 unsigned int vpMatrix::pseudoInverse(vpMatrix &Ap, vpColVector &sv, double svThreshold, vpMatrix &imA, vpMatrix &imAt) const
4544 {
4545  vpMatrix kerAt;
4546  return pseudoInverse(Ap, sv, svThreshold, imA, imAt, kerAt);
4547 }
4548 
4627 int vpMatrix::pseudoInverse(vpMatrix &Ap, vpColVector &sv, int rank_in, vpMatrix &imA, vpMatrix &imAt) const
4628 {
4629  vpMatrix kerAt;
4630  return pseudoInverse(Ap, sv, rank_in, imA, imAt, kerAt);
4631 }
4632 
4734 unsigned int vpMatrix::pseudoInverse(vpMatrix &Ap, vpColVector &sv, double svThreshold, vpMatrix &imA, vpMatrix &imAt,
4735  vpMatrix &kerAt) const
4736 {
4737 #if defined(VISP_HAVE_LAPACK)
4738  return pseudoInverseLapack(Ap, sv, svThreshold, imA, imAt, kerAt);
4739 #elif defined(VISP_HAVE_EIGEN3)
4740  return pseudoInverseEigen3(Ap, sv, svThreshold, imA, imAt, kerAt);
4741 #elif defined(VISP_HAVE_OPENCV) // Require opencv >= 2.1.1
4742  return pseudoInverseOpenCV(Ap, sv, svThreshold, imA, imAt, kerAt);
4743 #else
4744  (void)Ap;
4745  (void)sv;
4746  (void)svThreshold;
4747  (void)imA;
4748  (void)imAt;
4749  (void)kerAt;
4750  throw(vpException(vpException::fatalError, "Cannot compute pseudo-inverse. "
4751  "Install Lapack, Eigen3 or OpenCV 3rd party"));
4752 #endif
4753 }
4754 
4861 int vpMatrix::pseudoInverse(vpMatrix &Ap, vpColVector &sv, int rank_in, vpMatrix &imA, vpMatrix &imAt,
4862  vpMatrix &kerAt) const
4863 {
4864 #if defined(VISP_HAVE_LAPACK)
4865  return pseudoInverseLapack(Ap, sv, rank_in, imA, imAt, kerAt);
4866 #elif defined(VISP_HAVE_EIGEN3)
4867  return pseudoInverseEigen3(Ap, sv, rank_in, imA, imAt, kerAt);
4868 #elif defined(VISP_HAVE_OPENCV) // Require opencv >= 2.1.1
4869  return pseudoInverseOpenCV(Ap, sv, rank_in, imA, imAt, kerAt);
4870 #else
4871  (void)Ap;
4872  (void)sv;
4873  (void)rank_in;
4874  (void)imA;
4875  (void)imAt;
4876  (void)kerAt;
4877  throw(vpException(vpException::fatalError, "Cannot compute pseudo-inverse. "
4878  "Install Lapack, Eigen3 or OpenCV 3rd party"));
4879 #endif
4880 }
4881 
4923 vpColVector vpMatrix::getCol(unsigned int j, unsigned int i_begin, unsigned int column_size) const
4924 {
4925  if (((i_begin + column_size) > getRows()) || (j >= getCols())) {
4926  throw(vpException(vpException::dimensionError, "Unable to extract column %u from the %ux%u matrix", j, getRows(),
4927  getCols()));
4928  }
4929  vpColVector c(column_size);
4930  for (unsigned int i = 0; i < column_size; ++i) {
4931  c[i] = (*this)[i_begin + i][j];
4932  }
4933  return c;
4934 }
4935 
4975 vpColVector vpMatrix::getCol(unsigned int j) const { return getCol(j, 0, rowNum); }
4976 
5013 vpRowVector vpMatrix::getRow(unsigned int i) const { return getRow(i, 0, colNum); }
5014 
5054 vpRowVector vpMatrix::getRow(unsigned int i, unsigned int j_begin, unsigned int row_size) const
5055 {
5056  if (((j_begin + row_size) > colNum) || (i >= rowNum)) {
5057  throw(vpException(vpException::dimensionError, "Unable to extract a row vector from the matrix"));
5058  }
5059  vpRowVector r(row_size);
5060  if ((r.data != nullptr) && (data != nullptr)) {
5061  memcpy(r.data, (*this)[i] + j_begin, row_size * sizeof(double));
5062  }
5063 
5064  return r;
5065 }
5066 
5103 {
5104  unsigned int min_size = std::min<unsigned int>(rowNum, colNum);
5105  vpColVector diag;
5106 
5107  if (min_size > 0) {
5108  diag.resize(min_size, false);
5109 
5110  for (unsigned int i = 0; i < min_size; ++i) {
5111  diag[i] = (*this)[i][i];
5112  }
5113  }
5114 
5115  return diag;
5116 }
5117 
5129 {
5130  vpMatrix C;
5131 
5132  vpMatrix::stack(A, B, C);
5133 
5134  return C;
5135 }
5136 
5148 void vpMatrix::stack(const vpMatrix &A, const vpMatrix &B, vpMatrix &C)
5149 {
5150  unsigned int nra = A.getRows();
5151  unsigned int nrb = B.getRows();
5152 
5153  if (nra != 0) {
5154  if (A.getCols() != B.getCols()) {
5155  throw(vpException(vpException::dimensionError, "Cannot stack (%dx%d) matrix with (%dx%d) matrix", A.getRows(),
5156  A.getCols(), B.getRows(), B.getCols()));
5157  }
5158  }
5159 
5160  if ((A.data != nullptr) && (A.data == C.data)) {
5161  std::cerr << "A and C must be two different objects!" << std::endl;
5162  return;
5163  }
5164 
5165  if ((B.data != nullptr) && (B.data == C.data)) {
5166  std::cerr << "B and C must be two different objects!" << std::endl;
5167  return;
5168  }
5169 
5170  C.resize(nra + nrb, B.getCols(), false, false);
5171 
5172  if ((C.data != nullptr) && (A.data != nullptr) && (A.size() > 0)) {
5173  // Copy A in C
5174  memcpy(C.data, A.data, sizeof(double) * A.size());
5175  }
5176 
5177  if ((C.data != nullptr) && (B.data != nullptr) && (B.size() > 0)) {
5178  // Copy B in C
5179  memcpy(C.data + A.size(), B.data, sizeof(double) * B.size());
5180  }
5181 }
5182 
5193 {
5194  vpMatrix C;
5195  vpMatrix::stack(A, r, C);
5196 
5197  return C;
5198 }
5199 
5211 void vpMatrix::stack(const vpMatrix &A, const vpRowVector &r, vpMatrix &C)
5212 {
5213  if ((A.data != nullptr) && (A.data == C.data)) {
5214  std::cerr << "A and C must be two different objects!" << std::endl;
5215  return;
5216  }
5217 
5218  C = A;
5219  C.stack(r);
5220 }
5221 
5232 {
5233  vpMatrix C;
5234  vpMatrix::stack(A, c, C);
5235 
5236  return C;
5237 }
5238 
5250 void vpMatrix::stack(const vpMatrix &A, const vpColVector &c, vpMatrix &C)
5251 {
5252  if ((A.data != nullptr) && (A.data == C.data)) {
5253  std::cerr << "A and C must be two different objects!" << std::endl;
5254  return;
5255  }
5256 
5257  C = A;
5258  C.stack(c);
5259 }
5260 
5273 vpMatrix vpMatrix::insert(const vpMatrix &A, const vpMatrix &B, unsigned int r, unsigned int c)
5274 {
5276 
5277  vpArray2D<double>::insert(A, B, C, r, c);
5278 
5279  return vpMatrix(C);
5280 }
5281 
5295 void vpMatrix::insert(const vpMatrix &A, const vpMatrix &B, vpMatrix &C, unsigned int r, unsigned int c)
5296 {
5297  vpArray2D<double> C_array;
5298 
5299  vpArray2D<double>::insert(A, B, C_array, r, c);
5300 
5301  C = C_array;
5302 }
5303 
5316 {
5317  vpMatrix C;
5318 
5319  juxtaposeMatrices(A, B, C);
5320 
5321  return C;
5322 }
5323 
5337 {
5338  unsigned int nca = A.getCols();
5339  unsigned int ncb = B.getCols();
5340 
5341  if (nca != 0) {
5342  if (A.getRows() != B.getRows()) {
5343  throw(vpException(vpException::dimensionError, "Cannot juxtapose (%dx%d) matrix with (%dx%d) matrix", A.getRows(),
5344  A.getCols(), B.getRows(), B.getCols()));
5345  }
5346  }
5347 
5348  if ((B.getRows() == 0) || ((nca + ncb) == 0)) {
5349  std::cerr << "B.getRows() == 0 || nca+ncb == 0" << std::endl;
5350  return;
5351  }
5352 
5353  C.resize(B.getRows(), nca + ncb, false, false);
5354 
5355  C.insert(A, 0, 0);
5356  C.insert(B, 0, nca);
5357 }
5358 
5359 //--------------------------------------------------------------------
5360 // Output
5361 //--------------------------------------------------------------------
5362 
5382 int vpMatrix::print(std::ostream &s, unsigned int length, const std::string &intro) const
5383 {
5384  typedef std::string::size_type size_type;
5385 
5386  unsigned int m = getRows();
5387  unsigned int n = getCols();
5388 
5389  std::vector<std::string> values(m * n);
5390  std::ostringstream oss;
5391  std::ostringstream ossFixed;
5392  std::ios_base::fmtflags original_flags = oss.flags();
5393 
5394  ossFixed.setf(std::ios::fixed, std::ios::floatfield);
5395 
5396  size_type maxBefore = 0; // the length of the integral part
5397  size_type maxAfter = 0; // number of decimals plus
5398  // one place for the decimal point
5399  for (unsigned int i = 0; i < m; ++i) {
5400  for (unsigned int j = 0; j < n; ++j) {
5401  oss.str("");
5402  oss << (*this)[i][j];
5403  if (oss.str().find("e") != std::string::npos) {
5404  ossFixed.str("");
5405  ossFixed << (*this)[i][j];
5406  oss.str(ossFixed.str());
5407  }
5408 
5409  values[(i * n) + j] = oss.str();
5410  size_type thislen = values[(i * n) + j].size();
5411  size_type p = values[(i * n) + j].find('.');
5412 
5413  if (p == std::string::npos) {
5414  maxBefore = vpMath::maximum(maxBefore, thislen);
5415  // maxAfter remains the same
5416  }
5417  else {
5418  maxBefore = vpMath::maximum(maxBefore, p);
5419  maxAfter = vpMath::maximum(maxAfter, thislen - p);
5420  }
5421  }
5422  }
5423 
5424  size_type totalLength = length;
5425  // increase totalLength according to maxBefore
5426  totalLength = vpMath::maximum(totalLength, maxBefore);
5427  // decrease maxAfter according to totalLength
5428  maxAfter = std::min<size_type>(maxAfter, totalLength - maxBefore);
5429 
5430  if (!intro.empty()) {
5431  s << intro;
5432  }
5433  s << "[" << m << "," << n << "]=\n";
5434 
5435  for (unsigned int i = 0; i < m; ++i) {
5436  s << " ";
5437  for (unsigned int j = 0; j < n; ++j) {
5438  size_type p = values[(i * n) + j].find('.');
5439  s.setf(std::ios::right, std::ios::adjustfield);
5440  s.width(static_cast<std::streamsize>(maxBefore));
5441  s << values[(i * n) + j].substr(0, p).c_str();
5442 
5443  if (maxAfter > 0) {
5444  s.setf(std::ios::left, std::ios::adjustfield);
5445  if (p != std::string::npos) {
5446  s.width(static_cast<std::streamsize>(maxAfter));
5447  s << values[(i * n) + j].substr(p, maxAfter).c_str();
5448  }
5449  else {
5450  s.width(static_cast<std::streamsize>(maxAfter));
5451  s << ".0";
5452  }
5453  }
5454 
5455  s << ' ';
5456  }
5457  s << std::endl;
5458  }
5459 
5460  s.flags(original_flags); // restore s to standard state
5461 
5462  return static_cast<int>(maxBefore + maxAfter);
5463 }
5464 
5501 std::ostream &vpMatrix::matlabPrint(std::ostream &os) const
5502 {
5503  unsigned int this_rows = this->getRows();
5504  unsigned int this_col = this->getCols();
5505  os << "[ ";
5506  for (unsigned int i = 0; i < this_rows; ++i) {
5507  for (unsigned int j = 0; j < this_col; ++j) {
5508  os << (*this)[i][j] << ", ";
5509  }
5510  if (this->getRows() != (i + 1)) {
5511  os << ";" << std::endl;
5512  }
5513  else {
5514  os << "]" << std::endl;
5515  }
5516  }
5517  return os;
5518 }
5519 
5548 std::ostream &vpMatrix::maplePrint(std::ostream &os) const
5549 {
5550  unsigned int this_rows = this->getRows();
5551  unsigned int this_col = this->getCols();
5552  os << "([ " << std::endl;
5553  for (unsigned int i = 0; i < this_rows; ++i) {
5554  os << "[";
5555  for (unsigned int j = 0; j < this_col; ++j) {
5556  os << (*this)[i][j] << ", ";
5557  }
5558  os << "]," << std::endl;
5559  }
5560  os << "])" << std::endl;
5561  return os;
5562 }
5563 
5591 std::ostream &vpMatrix::csvPrint(std::ostream &os) const
5592 {
5593  unsigned int this_rows = this->getRows();
5594  unsigned int this_col = this->getCols();
5595  for (unsigned int i = 0; i < this_rows; ++i) {
5596  for (unsigned int j = 0; j < this_col; ++j) {
5597  os << (*this)[i][j];
5598  if (!(j == (this->getCols() - 1))) {
5599  os << ", ";
5600  }
5601  }
5602  os << std::endl;
5603  }
5604  return os;
5605 }
5606 
5642 std::ostream &vpMatrix::cppPrint(std::ostream &os, const std::string &matrixName, bool octet) const
5643 {
5644  os << "vpMatrix " << matrixName << " (" << this->getRows() << ", " << this->getCols() << "); " << std::endl;
5645 
5646  unsigned int this_rows = this->getRows();
5647  unsigned int this_col = this->getCols();
5648  for (unsigned int i = 0; i < this_rows; ++i) {
5649  for (unsigned int j = 0; j < this_col; ++j) {
5650  if (!octet) {
5651  os << matrixName << "[" << i << "][" << j << "] = " << (*this)[i][j] << "; " << std::endl;
5652  }
5653  else {
5654  for (unsigned int k = 0; k < sizeof(double); ++k) {
5655  os << "((unsigned char*)&(" << matrixName << "[" << i << "][" << j << "]) )[" << k << "] = 0x" << std::hex
5656  << static_cast<unsigned int>(((unsigned char *)&((*this)[i][j]))[k]) << "; " << std::endl;
5657  }
5658  }
5659  }
5660  os << std::endl;
5661  }
5662  return os;
5663 }
5664 
5670 {
5671  if (rowNum == 0) {
5672  *this = A;
5673  }
5674  else if (A.getRows() > 0) {
5675  if (colNum != A.getCols()) {
5676  throw(vpException(vpException::dimensionError, "Cannot stack (%dx%d) matrix with (%dx%d) matrix", rowNum, colNum,
5677  A.getRows(), A.getCols()));
5678  }
5679 
5680  unsigned int rowNumOld = rowNum;
5681  resize(rowNum + A.getRows(), colNum, false, false);
5682  insert(A, rowNumOld, 0);
5683  }
5684 }
5685 
5702 {
5703  if (rowNum == 0) {
5704  *this = r;
5705  }
5706  else {
5707  if (colNum != r.getCols()) {
5708  throw(vpException(vpException::dimensionError, "Cannot stack (%dx%d) matrix with (1x%d) row vector", rowNum,
5709  colNum, r.getCols()));
5710  }
5711 
5712  if (r.size() == 0) {
5713  return;
5714  }
5715 
5716  unsigned int oldSize = size();
5717  resize(rowNum + 1, colNum, false, false);
5718 
5719  if ((data != nullptr) && (r.data != nullptr) && (data != r.data)) {
5720  // Copy r in data
5721  memcpy(data + oldSize, r.data, sizeof(double) * r.size());
5722  }
5723  }
5724 }
5725 
5743 {
5744  if (colNum == 0) {
5745  *this = c;
5746  }
5747  else {
5748  if (rowNum != c.getRows()) {
5749  throw(vpException(vpException::dimensionError, "Cannot stack (%dx%d) matrix with (%dx1) column vector", rowNum,
5750  colNum, c.getRows()));
5751  }
5752 
5753  if (c.size() == 0) {
5754  return;
5755  }
5756 
5757  vpMatrix tmp = *this;
5758  unsigned int oldColNum = colNum;
5759  resize(rowNum, colNum + 1, false, false);
5760 
5761  if ((data != nullptr) && (tmp.data != nullptr) && (data != tmp.data)) {
5762  // Copy c in data
5763  for (unsigned int i = 0; i < rowNum; ++i) {
5764  memcpy(data + (i * colNum), tmp.data + (i * oldColNum), sizeof(double) * oldColNum);
5765  rowPtrs[i][oldColNum] = c[i];
5766  }
5767  }
5768  }
5769 }
5770 
5781 void vpMatrix::insert(const vpMatrix &A, unsigned int r, unsigned int c)
5782 {
5783  if (((r + A.getRows()) <= rowNum) && ((c + A.getCols()) <= colNum)) {
5784  if ((A.colNum == colNum) && (data != nullptr) && (A.data != nullptr) && (A.data != data)) {
5785  memcpy(data + (r * colNum), A.data, sizeof(double) * A.size());
5786  }
5787  else if ((data != nullptr) && (A.data != nullptr) && (A.data != data)) {
5788  unsigned int a_rows = A.getRows();
5789  for (unsigned int i = r; i < (r + a_rows); ++i) {
5790  memcpy(data + (i * colNum) + c, A.data + ((i - r) * A.colNum), sizeof(double) * A.colNum);
5791  }
5792  }
5793  }
5794  else {
5795  throw vpException(vpException::dimensionError, "Cannot insert (%dx%d) matrix in (%dx%d) matrix at position (%d,%d)",
5796  A.getRows(), A.getCols(), rowNum, colNum, r, c);
5797  }
5798 }
5799 
5837 {
5838  vpColVector evalue(rowNum); // Eigen values
5839 
5840  if (rowNum != colNum) {
5841  throw(vpException(vpException::dimensionError, "Cannot compute eigen values on a non square matrix (%dx%d)", rowNum,
5842  colNum));
5843  }
5844 
5845  // Check if the matrix is symmetric: At - A = 0
5846  vpMatrix At_A = (*this).t() - (*this);
5847  for (unsigned int i = 0; i < rowNum; ++i) {
5848  for (unsigned int j = 0; j < rowNum; ++j) {
5849  // in comment: if (At_A[i][j] != 0) {
5850  if (std::fabs(At_A[i][j]) > std::numeric_limits<double>::epsilon()) {
5851  throw(vpException(vpException::fatalError, "Cannot compute eigen values on a non symmetric matrix"));
5852  }
5853  }
5854  }
5855 
5856 #if defined(VISP_HAVE_LAPACK)
5857 #if defined(VISP_HAVE_GSL) /* be careful of the copy below */
5858  {
5859  gsl_vector *eval = gsl_vector_alloc(rowNum);
5860  gsl_matrix *evec = gsl_matrix_alloc(rowNum, colNum);
5861 
5862  gsl_eigen_symmv_workspace *w = gsl_eigen_symmv_alloc(rowNum);
5863  gsl_matrix *m = gsl_matrix_alloc(rowNum, colNum);
5864 
5865  unsigned int Atda = (unsigned int)m->tda;
5866  for (unsigned int i = 0; i < rowNum; i++) {
5867  unsigned int k = i * Atda;
5868  for (unsigned int j = 0; j < colNum; j++)
5869  m->data[k + j] = (*this)[i][j];
5870  }
5871  gsl_eigen_symmv(m, eval, evec, w);
5872 
5873  gsl_eigen_symmv_sort(eval, evec, GSL_EIGEN_SORT_ABS_ASC);
5874 
5875  for (unsigned int i = 0; i < rowNum; i++) {
5876  evalue[i] = gsl_vector_get(eval, i);
5877  }
5878 
5879  gsl_eigen_symmv_free(w);
5880  gsl_vector_free(eval);
5881  gsl_matrix_free(m);
5882  gsl_matrix_free(evec);
5883  }
5884 #else
5885  {
5886  const char jobz = 'N';
5887  const char uplo = 'U';
5888  vpMatrix A = (*this);
5889  vpColVector WORK;
5890  int lwork = -1;
5891  int info = 0;
5892  double wkopt;
5893  vpMatrix::blas_dsyev(jobz, uplo, rowNum, A.data, colNum, evalue.data, &wkopt, lwork, info);
5894  lwork = static_cast<int>(wkopt);
5895  WORK.resize(lwork);
5896  vpMatrix::blas_dsyev(jobz, uplo, rowNum, A.data, colNum, evalue.data, WORK.data, lwork, info);
5897  }
5898 #endif
5899 #else
5900  {
5901  throw(vpException(vpException::functionNotImplementedError, "Eigen values computation is not implemented. "
5902  "You should install Lapack 3rd party"));
5903  }
5904 #endif
5905  return evalue;
5906 }
5907 
5957 void vpMatrix::eigenValues(vpColVector &evalue, vpMatrix &evector) const
5958 {
5959  if (rowNum != colNum) {
5960  throw(vpException(vpException::dimensionError, "Cannot compute eigen values on a non square matrix (%dx%d)", rowNum,
5961  colNum));
5962  }
5963 
5964  // Check if the matrix is symmetric: At - A = 0
5965  vpMatrix At_A = (*this).t() - (*this);
5966  for (unsigned int i = 0; i < rowNum; ++i) {
5967  for (unsigned int j = 0; j < rowNum; ++j) {
5968  // -- in comment: if (At_A[i][j] != 0) {
5969  if (std::fabs(At_A[i][j]) > std::numeric_limits<double>::epsilon()) {
5970  throw(vpException(vpException::fatalError, "Cannot compute eigen values on a non symmetric matrix"));
5971  }
5972  }
5973  }
5974 
5975  // Resize the output matrices
5976  evalue.resize(rowNum);
5977  evector.resize(rowNum, colNum);
5978 
5979 #if defined(VISP_HAVE_LAPACK)
5980 #if defined(VISP_HAVE_GSL) /* be careful of the copy below */
5981  {
5982  gsl_vector *eval = gsl_vector_alloc(rowNum);
5983  gsl_matrix *evec = gsl_matrix_alloc(rowNum, colNum);
5984 
5985  gsl_eigen_symmv_workspace *w = gsl_eigen_symmv_alloc(rowNum);
5986  gsl_matrix *m = gsl_matrix_alloc(rowNum, colNum);
5987 
5988  unsigned int Atda = (unsigned int)m->tda;
5989  for (unsigned int i = 0; i < rowNum; i++) {
5990  unsigned int k = i * Atda;
5991  for (unsigned int j = 0; j < colNum; j++)
5992  m->data[k + j] = (*this)[i][j];
5993  }
5994  gsl_eigen_symmv(m, eval, evec, w);
5995 
5996  gsl_eigen_symmv_sort(eval, evec, GSL_EIGEN_SORT_ABS_ASC);
5997 
5998  for (unsigned int i = 0; i < rowNum; i++) {
5999  evalue[i] = gsl_vector_get(eval, i);
6000  }
6001  Atda = (unsigned int)evec->tda;
6002  for (unsigned int i = 0; i < rowNum; i++) {
6003  unsigned int k = i * Atda;
6004  for (unsigned int j = 0; j < rowNum; j++) {
6005  evector[i][j] = evec->data[k + j];
6006  }
6007  }
6008 
6009  gsl_eigen_symmv_free(w);
6010  gsl_vector_free(eval);
6011  gsl_matrix_free(m);
6012  gsl_matrix_free(evec);
6013  }
6014 #else // defined(VISP_HAVE_GSL)
6015  {
6016  const char jobz = 'V';
6017  const char uplo = 'U';
6018  vpMatrix A = (*this);
6019  vpColVector WORK;
6020  int lwork = -1;
6021  int info = 0;
6022  double wkopt;
6023  vpMatrix::blas_dsyev(jobz, uplo, rowNum, A.data, colNum, evalue.data, &wkopt, lwork, info);
6024  lwork = static_cast<int>(wkopt);
6025  WORK.resize(lwork);
6026  vpMatrix::blas_dsyev(jobz, uplo, rowNum, A.data, colNum, evalue.data, WORK.data, lwork, info);
6027  evector = A.t();
6028  }
6029 #endif // defined(VISP_HAVE_GSL)
6030 #else
6031  {
6032  throw(vpException(vpException::functionNotImplementedError, "Eigen values computation is not implemented. "
6033  "You should install Lapack 3rd party"));
6034  }
6035 #endif
6036 }
6037 
6056 unsigned int vpMatrix::kernel(vpMatrix &kerAt, double svThreshold) const
6057 {
6058  unsigned int nbline = getRows();
6059  unsigned int nbcol = getCols();
6060 
6061  vpMatrix U, V; // Copy of the matrix, SVD function is destructive
6062  vpColVector sv;
6063  sv.resize(nbcol, false); // singular values
6064  V.resize(nbcol, nbcol, false); // V matrix of singular value decomposition
6065 
6066  // Copy and resize matrix to have at least as many rows as columns
6067  // kernel is computed in svd method only if the matrix has more rows than
6068  // columns
6069 
6070  if (nbline < nbcol) {
6071  U.resize(nbcol, nbcol, true);
6072  }
6073  else {
6074  U.resize(nbline, nbcol, false);
6075  }
6076 
6077  U.insert(*this, 0, 0);
6078 
6079  U.svd(sv, V);
6080 
6081  // Compute the highest singular value and rank of the matrix
6082  double maxsv = 0;
6083  for (unsigned int i = 0; i < nbcol; ++i) {
6084  if (sv[i] > maxsv) {
6085  maxsv = sv[i];
6086  }
6087  }
6088 
6089  unsigned int rank = 0;
6090  for (unsigned int i = 0; i < nbcol; ++i) {
6091  if (sv[i] > (maxsv * svThreshold)) {
6092  rank++;
6093  }
6094  }
6095 
6096  kerAt.resize(nbcol - rank, nbcol);
6097  if (rank != nbcol) {
6098  for (unsigned int j = 0, k = 0; j < nbcol; ++j) {
6099  // if( v.col(j) in kernel and non zero )
6100  if ((sv[j] <= (maxsv * svThreshold)) &&
6101  (std::fabs(V.getCol(j).sumSquare()) > std::numeric_limits<double>::epsilon())) {
6102  unsigned int v_rows = V.getRows();
6103  for (unsigned int i = 0; i < v_rows; ++i) {
6104  kerAt[k][i] = V[i][j];
6105  }
6106  k++;
6107  }
6108  }
6109  }
6110 
6111  return rank;
6112 }
6113 
6130 unsigned int vpMatrix::nullSpace(vpMatrix &kerA, double svThreshold) const
6131 {
6132  unsigned int nbrow = getRows();
6133  unsigned int nbcol = getCols();
6134 
6135  vpMatrix U, V; // Copy of the matrix, SVD function is destructive
6136  vpColVector sv;
6137  sv.resize(nbcol, false); // singular values
6138  V.resize(nbcol, nbcol, false); // V matrix of singular value decomposition
6139 
6140  // Copy and resize matrix to have at least as many rows as columns
6141  // kernel is computed in svd method only if the matrix has more rows than
6142  // columns
6143 
6144  if (nbrow < nbcol) {
6145  U.resize(nbcol, nbcol, true);
6146  }
6147  else {
6148  U.resize(nbrow, nbcol, false);
6149  }
6150 
6151  U.insert(*this, 0, 0);
6152 
6153  U.svd(sv, V);
6154 
6155  // Compute the highest singular value and rank of the matrix
6156  double maxsv = sv[0];
6157 
6158  unsigned int rank = 0;
6159  for (unsigned int i = 0; i < nbcol; ++i) {
6160  if (sv[i] > (maxsv * svThreshold)) {
6161  rank++;
6162  }
6163  }
6164 
6165  kerA.resize(nbcol, nbcol - rank);
6166  if (rank != nbcol) {
6167  for (unsigned int j = 0, k = 0; j < nbcol; ++j) {
6168  // if( v.col(j) in kernel and non zero )
6169  if (sv[j] <= (maxsv * svThreshold)) {
6170  for (unsigned int i = 0; i < nbcol; ++i) {
6171  kerA[i][k] = V[i][j];
6172  }
6173  k++;
6174  }
6175  }
6176  }
6177 
6178  return (nbcol - rank);
6179 }
6180 
6197 unsigned int vpMatrix::nullSpace(vpMatrix &kerA, int dim) const
6198 {
6199  unsigned int nbrow = getRows();
6200  unsigned int nbcol = getCols();
6201  unsigned int dim_ = static_cast<unsigned int>(dim);
6202 
6203  vpMatrix U, V; // Copy of the matrix, SVD function is destructive
6204  vpColVector sv;
6205  sv.resize(nbcol, false); // singular values
6206  V.resize(nbcol, nbcol, false); // V matrix of singular value decomposition
6207 
6208  // Copy and resize matrix to have at least as many rows as columns
6209  // kernel is computed in svd method only if the matrix has more rows than
6210  // columns
6211 
6212  if (nbrow < nbcol) {
6213  U.resize(nbcol, nbcol, true);
6214  }
6215  else {
6216  U.resize(nbrow, nbcol, false);
6217  }
6218 
6219  U.insert(*this, 0, 0);
6220 
6221  U.svd(sv, V);
6222 
6223  kerA.resize(nbcol, dim_);
6224  if (dim_ != 0) {
6225  unsigned int rank = nbcol - dim_;
6226  for (unsigned int k = 0; k < dim_; ++k) {
6227  unsigned int j = k + rank;
6228  for (unsigned int i = 0; i < nbcol; ++i) {
6229  kerA[i][k] = V[i][j];
6230  }
6231  }
6232  }
6233 
6234  double maxsv = sv[0];
6235  unsigned int rank = 0;
6236  for (unsigned int i = 0; i < nbcol; ++i) {
6237  if (sv[i] > (maxsv * 1e-6)) {
6238  rank++;
6239  }
6240  }
6241  return (nbcol - rank);
6242 }
6243 
6275 double vpMatrix::det(vpDetMethod method) const
6276 {
6277  double det = 0.;
6278 
6279  if (method == LU_DECOMPOSITION) {
6280  det = this->detByLU();
6281  }
6282 
6283  return det;
6284 }
6285 
6294 {
6295  if (colNum != rowNum) {
6296  throw(vpException(vpException::dimensionError, "Cannot compute the exponential of a non square (%dx%d) matrix",
6297  rowNum, colNum));
6298  }
6299  else {
6300 #ifdef VISP_HAVE_GSL
6301  size_t size_ = rowNum * colNum;
6302  double *b = new double[size_];
6303  for (size_t i = 0; i < size_; i++)
6304  b[i] = 0.;
6305  gsl_matrix_view m = gsl_matrix_view_array(this->data, rowNum, colNum);
6306  gsl_matrix_view em = gsl_matrix_view_array(b, rowNum, colNum);
6307  gsl_linalg_exponential_ss(&m.matrix, &em.matrix, 0);
6308  // gsl_matrix_fprintf(stdout, &em.matrix, "%g");
6309  vpMatrix expA;
6310  expA.resize(rowNum, colNum, false);
6311  memcpy(expA.data, b, size_ * sizeof(double));
6312 
6313  delete[] b;
6314  return expA;
6315 #else
6316  vpMatrix v_expE(rowNum, colNum, false);
6317  vpMatrix v_expD(rowNum, colNum, false);
6318  vpMatrix v_expX(rowNum, colNum, false);
6319  vpMatrix v_expcX(rowNum, colNum, false);
6320  vpMatrix v_eye(rowNum, colNum, false);
6321 
6322  v_eye.eye();
6323  vpMatrix exp(*this);
6324 
6325  // -- in comment: double f;
6326  int e;
6327  double c = 0.5;
6328  int q = 6;
6329  int p = 1;
6330 
6331  double nA = 0;
6332  for (unsigned int i = 0; i < rowNum; ++i) {
6333  double sum = 0;
6334  for (unsigned int j = 0; j < colNum; ++j) {
6335  sum += fabs((*this)[i][j]);
6336  }
6337  if ((sum > nA) || (i == 0)) {
6338  nA = sum;
6339  }
6340  }
6341 
6342  /* f = */ frexp(nA, &e);
6343  // -- in comment: double s = (0 > e+1)?0:e+1;
6344  double s = e + 1;
6345 
6346  double sca = 1.0 / pow(2.0, s);
6347  exp = sca * exp;
6348  v_expX = *this;
6349  v_expE = (c * exp) + v_eye;
6350  v_expD = (-c * exp) + v_eye;
6351  for (int k = 2; k <= q; ++k) {
6352  c = (c * (static_cast<double>((q - k) + 1))) / (static_cast<double>(k * (((2 * q) - k) + 1)));
6353  v_expcX = exp * v_expX;
6354  v_expX = v_expcX;
6355  v_expcX = c * v_expX;
6356  v_expE = v_expE + v_expcX;
6357  if (p) {
6358  v_expD = v_expD + v_expcX;
6359  }
6360  else {
6361  v_expD = v_expD - v_expcX;
6362  }
6363  p = !p;
6364  }
6365  v_expX = v_expD.inverseByLU();
6366  exp = v_expX * v_expE;
6367  for (int k = 1; k <= s; ++k) {
6368  v_expE = exp * exp;
6369  exp = v_expE;
6370  }
6371  return exp;
6372 #endif
6373  }
6374 }
6375 
6376 /**************************************************************************************************************/
6377 /**************************************************************************************************************/
6378 
6379 // Specific functions
6380 
6381 /*
6382  input:: matrix M(nCols,nRows), nCols > 3, nRows > 3 , nCols == nRows.
6383 
6384  output:: the complement matrix of the element (rowNo,colNo).
6385  This is the matrix obtained from M after elimenating the row rowNo and column
6386  colNo
6387 
6388  example:
6389  1 2 3
6390  M = 4 5 6
6391  7 8 9
6392  1 3
6393  subblock(M, 1, 1) give the matrix 7 9
6394 */
6395 vpMatrix subblock(const vpMatrix &M, unsigned int col, unsigned int row)
6396 {
6397  vpMatrix M_comp;
6398  M_comp.resize(M.getRows() - 1, M.getCols() - 1, false);
6399 
6400  unsigned int m_rows = M.getRows();
6401  unsigned int m_col = M.getCols();
6402  for (unsigned int i = 0; i < col; ++i) {
6403  for (unsigned int j = 0; j < row; ++j) {
6404  M_comp[i][j] = M[i][j];
6405  }
6406  for (unsigned int j = row + 1; j < m_rows; ++j) {
6407  M_comp[i][j - 1] = M[i][j];
6408  }
6409  }
6410  for (unsigned int i = col + 1; i < m_col; ++i) {
6411  for (unsigned int j = 0; j < row; ++j) {
6412  M_comp[i - 1][j] = M[i][j];
6413  }
6414  for (unsigned int j = row + 1; j < m_rows; ++j) {
6415  M_comp[i - 1][j - 1] = M[i][j];
6416  }
6417  }
6418  return M_comp;
6419 }
6420 
6430 double vpMatrix::cond(double svThreshold) const
6431 {
6432  unsigned int nbline = getRows();
6433  unsigned int nbcol = getCols();
6434 
6435  vpMatrix U, V; // Copy of the matrix, SVD function is destructive
6436  vpColVector sv;
6437  sv.resize(nbcol); // singular values
6438  V.resize(nbcol, nbcol, false); // V matrix of singular value decomposition
6439 
6440  // Copy and resize matrix to have at least as many rows as columns
6441  // kernel is computed in svd method only if the matrix has more rows than
6442  // columns
6443 
6444  if (nbline < nbcol) {
6445  U.resize(nbcol, nbcol, true);
6446  }
6447  else {
6448  U.resize(nbline, nbcol, false);
6449  }
6450 
6451  U.insert(*this, 0, 0);
6452 
6453  U.svd(sv, V);
6454 
6455  // Compute the highest singular value
6456  double maxsv = 0;
6457  for (unsigned int i = 0; i < nbcol; ++i) {
6458  if (sv[i] > maxsv) {
6459  maxsv = sv[i];
6460  }
6461  }
6462 
6463  // Compute the rank of the matrix
6464  unsigned int rank = 0;
6465  for (unsigned int i = 0; i < nbcol; ++i) {
6466  if (sv[i] > (maxsv * svThreshold)) {
6467  rank++;
6468  }
6469  }
6470 
6471  // Compute the lowest singular value
6472  double minsv = maxsv;
6473  for (unsigned int i = 0; i < rank; ++i) {
6474  if (sv[i] < minsv) {
6475  minsv = sv[i];
6476  }
6477  }
6478 
6479  if (std::fabs(minsv) > std::numeric_limits<double>::epsilon()) {
6480  return maxsv / minsv;
6481  }
6482  else {
6483  return std::numeric_limits<double>::infinity();
6484  }
6485 }
6486 
6493 void vpMatrix::computeHLM(const vpMatrix &H, const double &alpha, vpMatrix &HLM)
6494 {
6495  if (H.getCols() != H.getRows()) {
6496  throw(vpException(vpException::dimensionError, "Cannot compute HLM on a non square matrix (%dx%d)", H.getRows(),
6497  H.getCols()));
6498  }
6499 
6500  HLM = H;
6501  unsigned int h_col = H.getCols();
6502  for (unsigned int i = 0; i < h_col; ++i) {
6503  HLM[i][i] += alpha * H[i][i];
6504  }
6505 }
6506 
6515 {
6516  double norm = 0.0;
6517  for (unsigned int i = 0; i < dsize; ++i) {
6518  double x = *(data + i);
6519  norm += x * x;
6520  }
6521 
6522  return sqrt(norm);
6523 }
6524 
6534 {
6535  if (this->dsize != 0) {
6536  vpMatrix v;
6537  vpColVector w;
6538 
6539  vpMatrix M = *this;
6540 
6541  M.svd(w, v);
6542 
6543  double max = w[0];
6544  unsigned int maxRank = std::min<unsigned int>(this->getCols(), this->getRows());
6545  // The maximum reachable rank is either the number of columns or the number of rows
6546  // of the matrix.
6547  unsigned int boundary = std::min<unsigned int>(maxRank, w.size());
6548  // boundary is here to ensure that the number of singular values used for the com-
6549  // putation of the euclidean norm of the matrix is not greater than the maximum
6550  // reachable rank. Indeed, some svd library pad the singular values vector with 0s
6551  // if the input matrix is non-square.
6552  for (unsigned int i = 0; i < boundary; ++i) {
6553  if (max < w[i]) {
6554  max = w[i];
6555  }
6556  }
6557  return max;
6558  }
6559  else {
6560  return 0.;
6561  }
6562 }
6563 
6575 {
6576  double norm = 0.0;
6577  for (unsigned int i = 0; i < rowNum; ++i) {
6578  double x = 0;
6579  for (unsigned int j = 0; j < colNum; ++j) {
6580  x += fabs(*(*(rowPtrs + i) + j));
6581  }
6582  if (x > norm) {
6583  norm = x;
6584  }
6585  }
6586  return norm;
6587 }
6588 
6595 double vpMatrix::sumSquare() const
6596 {
6597  double sum_square = 0.0;
6598  double x;
6599 
6600  for (unsigned int i = 0; i < rowNum; ++i) {
6601  for (unsigned int j = 0; j < colNum; ++j) {
6602  x = rowPtrs[i][j];
6603  sum_square += x * x;
6604  }
6605  }
6606 
6607  return sum_square;
6608 }
6609 #if defined(VISP_BUILD_DEPRECATED_FUNCTIONS)
6619 double vpMatrix::euclideanNorm() const { return frobeniusNorm(); }
6620 
6622 {
6623  return (vpMatrix)(vpColVector::stack(A, B));
6624 }
6625 
6627 {
6628  vpColVector::stack(A, B, C);
6629 }
6630 
6632 
6633 void vpMatrix::stackMatrices(const vpMatrix &A, const vpRowVector &B, vpMatrix &C) { vpMatrix::stack(A, B, C); }
6634 
6654 {
6655  vpRowVector c(getCols());
6656 
6657  for (unsigned int j = 0; j < getCols(); ++j) {
6658  c[j] = (*this)[i - 1][j];
6659  }
6660  return c;
6661 }
6662 
6681 {
6682  vpColVector c(getRows());
6683 
6684  for (unsigned int i = 0; i < getRows(); ++i) {
6685  c[i] = (*this)[i][j - 1];
6686  }
6687  return c;
6688 }
6689 
6696 void vpMatrix::setIdentity(const double &val)
6697 {
6698  for (unsigned int i = 0; i < rowNum; ++i) {
6699  for (unsigned int j = 0; j < colNum; ++j) {
6700  if (i == j) {
6701  (*this)[i][j] = val;
6702  }
6703  else {
6704  (*this)[i][j] = 0;
6705  }
6706  }
6707  }
6708 }
6709 
6710 #endif //#if defined(VISP_BUILD_DEPRECATED_FUNCTIONS)
Implementation of a generic 2D array used as base class for matrices and vectors.
Definition: vpArray2D.h:126
unsigned int getCols() const
Definition: vpArray2D.h:314
double * data
Address of the first element of the data array.
Definition: vpArray2D.h:139
double ** rowPtrs
Address of the first element of each rows.
Definition: vpArray2D.h:133
void insert(const vpArray2D< Type > &A, unsigned int r, unsigned int c)
Definition: vpArray2D.h:471
void resize(unsigned int nrows, unsigned int ncols, bool flagNullify=true, bool recopy_=true)
Definition: vpArray2D.h:339
unsigned int rowNum
Number of rows in the array.
Definition: vpArray2D.h:129
unsigned int dsize
Current array size (rowNum * colNum)
Definition: vpArray2D.h:135
unsigned int size() const
Return the number of elements of the 2D array.
Definition: vpArray2D.h:326
unsigned int getRows() const
Definition: vpArray2D.h:324
unsigned int colNum
Number of columns in the array.
Definition: vpArray2D.h:131
Implementation of column vector and the associated operations.
Definition: vpColVector.h:163
double sumSquare() const
void stack(double d)
void resize(unsigned int i, bool flagNullify=true)
Definition: vpColVector.h:1056
error that can be emitted by ViSP classes.
Definition: vpException.h:59
@ functionNotImplementedError
Function not implemented.
Definition: vpException.h:78
@ dimensionError
Bad dimension.
Definition: vpException.h:83
@ fatalError
Fatal error.
Definition: vpException.h:84
@ divideByZeroError
Division by zero.
Definition: vpException.h:82
Implementation of an homogeneous matrix and operations on such kind of matrices.
static Type maximum(const Type &a, const Type &b)
Definition: vpMath.h:252
Implementation of a matrix and operations on matrices.
Definition: vpMatrix.h:146
void svdLapack(vpColVector &w, vpMatrix &V)
vpColVector eigenValues() const
Definition: vpMatrix.cpp:5836
vpMatrix pseudoInverseOpenCV(double svThreshold=1e-6) const
std::ostream & cppPrint(std::ostream &os, const std::string &matrixName="A", bool octet=false) const
Definition: vpMatrix.cpp:5642
vpMatrix & operator<<(double *)
Definition: vpMatrix.cpp:810
double cond(double svThreshold=1e-6) const
Definition: vpMatrix.cpp:6430
vpMatrix expm() const
Definition: vpMatrix.cpp:6293
vpMatrix pseudoInverseEigen3(double svThreshold=1e-6) const
vpMatrix hadamard(const vpMatrix &m) const
Definition: vpMatrix.cpp:1863
vpMatrix operator-() const
Definition: vpMatrix.cpp:1637
vp_deprecated void setIdentity(const double &val=1.0)
Definition: vpMatrix.cpp:6696
vpMatrix & operator/=(double x)
Divide all the element of the matrix by x : Aij = Aij / x.
Definition: vpMatrix.cpp:1780
int print(std::ostream &s, unsigned int length, const std::string &intro="") const
Definition: vpMatrix.cpp:5382
std::ostream & maplePrint(std::ostream &os) const
Definition: vpMatrix.cpp:5548
void eye()
Definition: vpMatrix.cpp:448
unsigned int kernel(vpMatrix &kerAt, double svThreshold=1e-6) const
Definition: vpMatrix.cpp:6056
vpMatrix inverseByLU() const
vpMatrix operator*(const vpMatrix &B) const
Definition: vpMatrix.cpp:1198
vpMatrix operator/(double x) const
Cij = Aij / x (A is unchanged)
Definition: vpMatrix.cpp:1712
void stack(const vpMatrix &A)
Definition: vpMatrix.cpp:5669
vp_deprecated void stackMatrices(const vpMatrix &A)
Definition: vpMatrix.h:994
vpMatrix & operator-=(const vpMatrix &B)
Operation A = A - B.
Definition: vpMatrix.cpp:1591
vpMatrix & operator*=(double x)
Multiply all the element of the matrix by x : Aij = Aij * x.
Definition: vpMatrix.cpp:1764
vpMatrix operator*(const double &x, const vpMatrix &B)
Definition: vpMatrix.cpp:1668
static vpMatrix juxtaposeMatrices(const vpMatrix &A, const vpMatrix &B)
Definition: vpMatrix.cpp:5315
vp_deprecated double euclideanNorm() const
Definition: vpMatrix.cpp:6619
vpColVector stackColumns()
Definition: vpMatrix.cpp:1827
vp_deprecated vpColVector column(unsigned int j)
Definition: vpMatrix.cpp:6680
static void sub2Matrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C)
Definition: vpMatrix.cpp:1537
void svd(vpColVector &w, vpMatrix &V)
Definition: vpMatrix.cpp:2132
vpRowVector stackRows()
Definition: vpMatrix.cpp:1850
vpMatrix operator+(const vpMatrix &B) const
Definition: vpMatrix.cpp:1481
double sum() const
Definition: vpMatrix.cpp:1644
void diag(const double &val=1.0)
Definition: vpMatrix.cpp:910
vpRowVector getRow(unsigned int i) const
Definition: vpMatrix.cpp:5013
void svdOpenCV(vpColVector &w, vpMatrix &V)
vpMatrix t() const
Definition: vpMatrix.cpp:465
vpMatrix()
Definition: vpMatrix.h:162
vpMatrix & operator+=(const vpMatrix &B)
Operation A = A + B.
Definition: vpMatrix.cpp:1572
static void add2WeightedMatrices(const vpMatrix &A, const double &wA, const vpMatrix &B, const double &wB, vpMatrix &C)
Definition: vpMatrix.cpp:1389
double inducedL2Norm() const
Definition: vpMatrix.cpp:6533
double infinityNorm() const
Definition: vpMatrix.cpp:6574
vpMatrix & operator=(const vpArray2D< double > &A)
Definition: vpMatrix.cpp:670
double frobeniusNorm() const
Definition: vpMatrix.cpp:6514
vpColVector getDiag() const
Definition: vpMatrix.cpp:5102
vpMatrix AtA() const
Definition: vpMatrix.cpp:645
static void createDiagonalMatrix(const vpColVector &A, vpMatrix &DA)
Definition: vpMatrix.cpp:929
void solveBySVD(const vpColVector &B, vpColVector &x) const
Definition: vpMatrix.cpp:2010
static void computeHLM(const vpMatrix &H, const double &alpha, vpMatrix &HLM)
Definition: vpMatrix.cpp:6493
vpMatrix pseudoInverseLapack(double svThreshold=1e-6) const
vpMatrix pseudoInverse(double svThreshold=1e-6) const
Definition: vpMatrix.cpp:2343
vpMatrix & operator,(double val)
Definition: vpMatrix.cpp:827
static void multMatrixVector(const vpMatrix &A, const vpColVector &v, vpColVector &w)
Definition: vpMatrix.cpp:984
vp_deprecated vpRowVector row(unsigned int i)
Definition: vpMatrix.cpp:6653
vpColVector getCol(unsigned int j) const
Definition: vpMatrix.cpp:4975
double detByLU() const
double det(vpDetMethod method=LU_DECOMPOSITION) const
Definition: vpMatrix.cpp:6275
double sumSquare() const
Definition: vpMatrix.cpp:6595
vp_deprecated void init()
Definition: vpMatrix.h:989
static void add2Matrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C)
Definition: vpMatrix.cpp:1421
void insert(const vpMatrix &A, unsigned int r, unsigned int c)
Definition: vpMatrix.cpp:5781
void svdEigen3(vpColVector &w, vpMatrix &V)
void kron(const vpMatrix &m1, vpMatrix &out) const
Definition: vpMatrix.cpp:1919
vpMatrix AAt() const
Definition: vpMatrix.cpp:517
vpMatrix transpose() const
Definition: vpMatrix.cpp:472
static void mult2Matrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C)
Definition: vpMatrix.cpp:1035
static void negateMatrix(const vpMatrix &A, vpMatrix &C)
Definition: vpMatrix.cpp:1617
@ LU_DECOMPOSITION
Definition: vpMatrix.h:154
std::ostream & csvPrint(std::ostream &os) const
Definition: vpMatrix.cpp:5591
unsigned int nullSpace(vpMatrix &kerA, double svThreshold=1e-6) const
Definition: vpMatrix.cpp:6130
std::ostream & matlabPrint(std::ostream &os) const
Definition: vpMatrix.cpp:5501
vpMatrix extract(unsigned int r, unsigned int c, unsigned int nrows, unsigned int ncols) const
Definition: vpMatrix.cpp:404
Implementation of a rotation matrix and operations on such kind of matrices.
Implementation of row vector and the associated operations.
Definition: vpRowVector.h:107
void resize(unsigned int i, bool flagNullify=true)
Definition: vpRowVector.h:259
Class that consider the case of a translation vector.