Visual Servoing Platform  version 3.6.1 under development (2025-02-05)
vpMatrix.cpp
1 /*
2  * ViSP, open source Visual Servoing Platform software.
3  * Copyright (C) 2005 - 2024 by Inria. All rights reserved.
4  *
5  * This software is free software; you can redistribute it and/or modify
6  * it under the terms of the GNU General Public License as published by
7  * the Free Software Foundation; either version 2 of the License, or
8  * (at your option) any later version.
9  * See the file LICENSE.txt at the root directory of this source
10  * distribution for additional information about the GNU GPL.
11  *
12  * For using ViSP with software that can not be combined with the GNU
13  * GPL, please contact Inria about acquiring a ViSP Professional
14  * Edition License.
15  *
16  * See https://visp.inria.fr for more information.
17  *
18  * This software was developed at:
19  * Inria Rennes - Bretagne Atlantique
20  * Campus Universitaire de Beaulieu
21  * 35042 Rennes Cedex
22  * France
23  *
24  * If you have questions regarding the use of this file, please contact
25  * Inria at visp@inria.fr
26  *
27  * This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
28  * WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
29  *
30  * Description:
31  * Matrix manipulation.
32  */
33 
39 #include <algorithm>
40 #include <assert.h>
41 #include <cmath> // std::fabs
42 #include <fstream>
43 #include <limits> // numeric_limits
44 #include <sstream>
45 #include <stdio.h>
46 #include <stdlib.h>
47 #include <string.h>
48 #include <string>
49 #include <vector>
50 
51 #include <visp3/core/vpCPUFeatures.h>
52 #include <visp3/core/vpColVector.h>
53 #include <visp3/core/vpConfig.h>
54 #include <visp3/core/vpException.h>
55 #include <visp3/core/vpMath.h>
56 #include <visp3/core/vpMatrix.h>
57 #include <visp3/core/vpTranslationVector.h>
58 
59 #if defined(VISP_HAVE_SIMDLIB)
60 #include <Simd/SimdLib.h>
61 #endif
62 
63 #ifdef VISP_HAVE_LAPACK
64 #ifdef VISP_HAVE_GSL
65 #include <gsl/gsl_eigen.h>
66 #include <gsl/gsl_linalg.h>
67 #include <gsl/gsl_math.h>
68 #elif defined(VISP_HAVE_MKL)
69 #include <mkl.h>
70 #endif
71 #endif
72 
73 BEGIN_VISP_NAMESPACE
74 
75 #ifdef VISP_HAVE_LAPACK
76 #ifdef VISP_HAVE_GSL
77 #elif defined(VISP_HAVE_MKL)
78 typedef MKL_INT integer;
79 
80 void vpMatrix::blas_dsyev(char jobz, char uplo, unsigned int n_, double *a_data, unsigned int lda_, double *w_data,
81  double *work_data, int lwork_, int &info_)
82 {
83  MKL_INT n = static_cast<MKL_INT>(n_);
84  MKL_INT lda = static_cast<MKL_INT>(lda_);
85  MKL_INT lwork = static_cast<MKL_INT>(lwork_);
86  MKL_INT info = static_cast<MKL_INT>(info_);
87 
88  dsyev(&jobz, &uplo, &n, a_data, &lda, w_data, work_data, &lwork, &info);
89 }
90 
91 #else
92 #if defined(VISP_HAVE_LAPACK_BUILT_IN)
93 typedef long int integer;
94 #else
95 typedef int integer;
96 #endif
97 extern "C" integer dsyev_(char *jobz, char *uplo, integer *n, double *a, integer *lda, double *w, double *WORK,
98  integer *lwork, integer *info);
99 
100 void vpMatrix::blas_dsyev(char jobz, char uplo, unsigned int n_, double *a_data, unsigned int lda_, double *w_data,
101  double *work_data, int lwork_, int &info_)
102 {
103  integer n = static_cast<integer>(n_);
104  integer lda = static_cast<integer>(lda_);
105  integer lwork = static_cast<integer>(lwork_);
106  integer info = static_cast<integer>(info_);
107 
108  dsyev_(&jobz, &uplo, &n, a_data, &lda, w_data, work_data, &lwork, &info);
109 
110  lwork_ = static_cast<int>(lwork);
111  info_ = static_cast<int>(info);
112 }
113 #endif
114 #endif
115 
116 #if !defined(VISP_USE_MSVC) || (defined(VISP_USE_MSVC) && !defined(VISP_BUILD_SHARED_LIBS))
117 const unsigned int vpMatrix::m_lapack_min_size_default = 0;
118 unsigned int vpMatrix::m_lapack_min_size = vpMatrix::m_lapack_min_size_default;
119 #endif
120 
121 // Prototypes of specific functions
122 vpMatrix subblock(const vpMatrix &M, unsigned int col, unsigned int row);
123 
128 vpMatrix::vpMatrix(const vpMatrix &M, unsigned int r, unsigned int c, unsigned int nrows, unsigned int ncols)
129  : vpArray2D<double>()
130 {
131  if (((r + nrows) > M.rowNum) || ((c + ncols) > M.colNum)) {
133  "Cannot construct a sub matrix (%dx%d) starting at "
134  "position (%d,%d) that is not contained in the "
135  "original matrix (%dx%d)",
136  nrows, ncols, r, c, M.rowNum, M.colNum));
137  }
138 
139  init(M, r, c, nrows, ncols);
140 }
141 
147 {
148  *this = M;
149 }
150 
156 {
157  *this = V;
158 }
159 
165 {
166  *this = F;
167 }
168 
174 {
175  *this = R;
176 }
177 
183 {
184  *this = v;
185 }
186 
192 {
193  *this = v;
194 }
195 
201 {
202  *this = t;
203 }
204 
205 #if (VISP_CXX_STANDARD >= VISP_CXX_STANDARD_11)
206 vpMatrix::vpMatrix(vpMatrix &&A) : vpArray2D<double>(std::move(A))
207 { }
208 
236 vpMatrix::vpMatrix(const std::initializer_list<double> &list) : vpArray2D<double>(list) { }
237 
264 vpMatrix::vpMatrix(unsigned int nrows, unsigned int ncols, const std::initializer_list<double> &list)
265  : vpArray2D<double>(nrows, ncols, list)
266 { }
267 
292 vpMatrix::vpMatrix(const std::initializer_list<std::initializer_list<double> > &lists) : vpArray2D<double>(lists) { }
293 #endif
294 
306 vpMatrix vpMatrix::view(double *data, unsigned int rows, unsigned int cols)
307 {
308  vpMatrix M;
309  vpArray2D<double>::view(M, data, rows, cols);
310  return M;
311 }
312 
363 void vpMatrix::init(const vpMatrix &M, unsigned int r, unsigned int c, unsigned int nrows, unsigned int ncols)
364 {
365  unsigned int rnrows = r + nrows;
366  unsigned int cncols = c + ncols;
367 
368  if (rnrows > M.getRows()) {
369  throw(vpException(vpException::dimensionError, "Bad row dimension (%d > %d) used to initialize vpMatrix", rnrows,
370  M.getRows()));
371  }
372  if (cncols > M.getCols()) {
373  throw(vpException(vpException::dimensionError, "Bad column dimension (%d > %d) used to initialize vpMatrix", cncols,
374  M.getCols()));
375  }
376  resize(nrows, ncols, false, false);
377 
378  if (this->rowPtrs == nullptr) { // Fix coverity scan: explicit null dereferenced
379  return; // Noting to do
380  }
381  for (unsigned int i = 0; i < nrows; ++i) {
382  memcpy((*this)[i], &M[i + r][c], ncols * sizeof(double));
383  }
384 }
385 
430 vpMatrix vpMatrix::extract(unsigned int r, unsigned int c, unsigned int nrows, unsigned int ncols) const
431 {
432  unsigned int rnrows = r + nrows;
433  unsigned int cncols = c + ncols;
434 
435  if (rnrows > getRows()) {
436  throw(vpException(vpException::dimensionError, "Bad row dimension (%d > %d) used to initialize vpMatrix", rnrows,
437  getRows()));
438  }
439  if (cncols > getCols()) {
440  throw(vpException(vpException::dimensionError, "Bad column dimension (%d > %d) used to initialize vpMatrix", cncols,
441  getCols()));
442  }
443 
444  vpMatrix M;
445  M.resize(nrows, ncols, false, false);
446  for (unsigned int i = 0; i < nrows; ++i) {
447  memcpy(M[i], &(*this)[i + r][c], ncols * sizeof(double));
448  }
449 
450  return M;
451 }
452 
498 vpColVector vpMatrix::getCol(unsigned int j, unsigned int i_begin, unsigned int column_size) const
499 {
500  if (((i_begin + column_size) > getRows()) || (j >= getCols())) {
501  throw(vpException(vpException::dimensionError, "Unable to extract column %u from the %ux%u matrix", j, getRows(),
502  getCols()));
503  }
504  vpColVector c(column_size);
505  for (unsigned int i = 0; i < column_size; ++i) {
506  c[i] = (*this)[i_begin + i][j];
507  }
508  return c;
509 }
510 
554 vpColVector vpMatrix::getCol(unsigned int j) const { return getCol(j, 0, rowNum); }
555 
596 vpRowVector vpMatrix::getRow(unsigned int i) const { return getRow(i, 0, colNum); }
597 
641 vpRowVector vpMatrix::getRow(unsigned int i, unsigned int j_begin, unsigned int row_size) const
642 {
643  if (((j_begin + row_size) > colNum) || (i >= rowNum)) {
644  throw(vpException(vpException::dimensionError, "Unable to extract a row vector from the matrix"));
645  }
646  vpRowVector r(row_size);
647  if ((r.data != nullptr) && (data != nullptr)) {
648  memcpy(r.data, (*this)[i] + j_begin, row_size * sizeof(double));
649  }
650 
651  return r;
652 }
653 
694 {
695  unsigned int min_size = std::min<unsigned int>(rowNum, colNum);
697 
698  if (min_size > 0) {
699  diag.resize(min_size, false);
700 
701  for (unsigned int i = 0; i < min_size; ++i) {
702  diag[i] = (*this)[i][i];
703  }
704  }
705 
706  return diag;
707 }
708 
721 vpMatrix vpMatrix::insert(const vpMatrix &A, const vpMatrix &B, unsigned int r, unsigned int c)
722 {
724 
725  vpArray2D<double>::insert(A, B, C, r, c);
726 
727  return vpMatrix(C);
728 }
729 
743 void vpMatrix::insert(const vpMatrix &A, const vpMatrix &B, vpMatrix &C, unsigned int r, unsigned int c)
744 {
745  vpArray2D<double> C_array;
746 
747  vpArray2D<double>::insert(A, B, C_array, r, c);
748 
749  C = C_array;
750 }
751 
764 {
765  vpMatrix C;
766 
767  juxtaposeMatrices(A, B, C);
768 
769  return C;
770 }
771 
785 {
786  unsigned int nca = A.getCols();
787  unsigned int ncb = B.getCols();
788 
789  if (nca != 0) {
790  if (A.getRows() != B.getRows()) {
791  throw(vpException(vpException::dimensionError, "Cannot juxtapose (%dx%d) matrix with (%dx%d) matrix", A.getRows(),
792  A.getCols(), B.getRows(), B.getCols()));
793  }
794  }
795 
796  if ((B.getRows() == 0) || ((nca + ncb) == 0)) {
797  std::cerr << "B.getRows() == 0 || nca+ncb == 0" << std::endl;
798  return;
799  }
800 
801  C.resize(B.getRows(), nca + ncb, false, false);
802 
803  C.insert(A, 0, 0);
804  C.insert(B, 0, nca);
805 }
806 
807 //--------------------------------------------------------------------
808 // Output
809 //--------------------------------------------------------------------
810 
830 int vpMatrix::print(std::ostream &s, unsigned int length, const std::string &intro) const
831 {
832  typedef std::string::size_type size_type;
833 
834  unsigned int m = getRows();
835  unsigned int n = getCols();
836 
837  std::vector<std::string> values(m * n);
838  std::ostringstream oss;
839  std::ostringstream ossFixed;
840  std::ios_base::fmtflags original_flags = oss.flags();
841 
842  ossFixed.setf(std::ios::fixed, std::ios::floatfield);
843 
844  size_type maxBefore = 0; // the length of the integral part
845  size_type maxAfter = 0; // number of decimals plus
846  // one place for the decimal point
847  for (unsigned int i = 0; i < m; ++i) {
848  for (unsigned int j = 0; j < n; ++j) {
849  oss.str("");
850  oss << (*this)[i][j];
851  if (oss.str().find("e") != std::string::npos) {
852  ossFixed.str("");
853  ossFixed << (*this)[i][j];
854  oss.str(ossFixed.str());
855  }
856 
857  values[(i * n) + j] = oss.str();
858  size_type thislen = values[(i * n) + j].size();
859  size_type p = values[(i * n) + j].find('.');
860 
861  if (p == std::string::npos) {
862  maxBefore = vpMath::maximum(maxBefore, thislen);
863  // maxAfter remains the same
864  }
865  else {
866  maxBefore = vpMath::maximum(maxBefore, p);
867  maxAfter = vpMath::maximum(maxAfter, thislen - p);
868  }
869  }
870  }
871 
872  size_type totalLength = length;
873  // increase totalLength according to maxBefore
874  totalLength = vpMath::maximum(totalLength, maxBefore);
875  // decrease maxAfter according to totalLength
876  maxAfter = std::min<size_type>(maxAfter, totalLength - maxBefore);
877 
878  if (!intro.empty()) {
879  s << intro;
880  }
881  s << "[" << m << "," << n << "]=\n";
882 
883  for (unsigned int i = 0; i < m; ++i) {
884  s << " ";
885  for (unsigned int j = 0; j < n; ++j) {
886  size_type p = values[(i * n) + j].find('.');
887  s.setf(std::ios::right, std::ios::adjustfield);
888  s.width(static_cast<std::streamsize>(maxBefore));
889  s << values[(i * n) + j].substr(0, p).c_str();
890 
891  if (maxAfter > 0) {
892  s.setf(std::ios::left, std::ios::adjustfield);
893  if (p != std::string::npos) {
894  s.width(static_cast<std::streamsize>(maxAfter));
895  s << values[(i * n) + j].substr(p, maxAfter).c_str();
896  }
897  else {
898  s.width(static_cast<std::streamsize>(maxAfter));
899  s << ".0";
900  }
901  }
902 
903  s << ' ';
904  }
905  s << std::endl;
906  }
907 
908  s.flags(original_flags); // restore s to standard state
909 
910  return static_cast<int>(maxBefore + maxAfter);
911 }
912 
953 std::ostream &vpMatrix::matlabPrint(std::ostream &os) const
954 {
955  unsigned int this_rows = this->getRows();
956  unsigned int this_col = this->getCols();
957  os << "[ ";
958  for (unsigned int i = 0; i < this_rows; ++i) {
959  for (unsigned int j = 0; j < this_col; ++j) {
960  os << (*this)[i][j] << ", ";
961  }
962  if (this->getRows() != (i + 1)) {
963  os << ";" << std::endl;
964  }
965  else {
966  os << "]" << std::endl;
967  }
968  }
969  return os;
970 }
971 
1004 std::ostream &vpMatrix::maplePrint(std::ostream &os) const
1005 {
1006  unsigned int this_rows = this->getRows();
1007  unsigned int this_col = this->getCols();
1008  os << "([ " << std::endl;
1009  for (unsigned int i = 0; i < this_rows; ++i) {
1010  os << "[";
1011  for (unsigned int j = 0; j < this_col; ++j) {
1012  os << (*this)[i][j] << ", ";
1013  }
1014  os << "]," << std::endl;
1015  }
1016  os << "])" << std::endl;
1017  return os;
1018 }
1019 
1051 std::ostream &vpMatrix::csvPrint(std::ostream &os) const
1052 {
1053  unsigned int this_rows = this->getRows();
1054  unsigned int this_col = this->getCols();
1055  for (unsigned int i = 0; i < this_rows; ++i) {
1056  for (unsigned int j = 0; j < this_col; ++j) {
1057  os << (*this)[i][j];
1058  if (!(j == (this->getCols() - 1))) {
1059  os << ", ";
1060  }
1061  }
1062  os << std::endl;
1063  }
1064  return os;
1065 }
1066 
1106 std::ostream &vpMatrix::cppPrint(std::ostream &os, const std::string &matrixName, bool octet) const
1107 {
1108  os << "vpMatrix " << matrixName << " (" << this->getRows() << ", " << this->getCols() << "); " << std::endl;
1109 
1110  unsigned int this_rows = this->getRows();
1111  unsigned int this_col = this->getCols();
1112  for (unsigned int i = 0; i < this_rows; ++i) {
1113  for (unsigned int j = 0; j < this_col; ++j) {
1114  if (!octet) {
1115  os << matrixName << "[" << i << "][" << j << "] = " << (*this)[i][j] << "; " << std::endl;
1116  }
1117  else {
1118  for (unsigned int k = 0; k < sizeof(double); ++k) {
1119  os << "((unsigned char*)&(" << matrixName << "[" << i << "][" << j << "]) )[" << k << "] = 0x" << std::hex
1120  << static_cast<unsigned int>(((unsigned char *)&((*this)[i][j]))[k]) << "; " << std::endl;
1121  }
1122  }
1123  }
1124  os << std::endl;
1125  }
1126  return os;
1127 }
1128 
1139 void vpMatrix::insert(const vpMatrix &A, unsigned int r, unsigned int c)
1140 {
1141  if (((r + A.getRows()) <= rowNum) && ((c + A.getCols()) <= colNum)) {
1142  if ((A.colNum == colNum) && (data != nullptr) && (A.data != nullptr) && (A.data != data)) {
1143  memcpy(data + (r * colNum), A.data, sizeof(double) * A.size());
1144  }
1145  else if ((data != nullptr) && (A.data != nullptr) && (A.data != data)) {
1146  unsigned int a_rows = A.getRows();
1147  for (unsigned int i = r; i < (r + a_rows); ++i) {
1148  memcpy(data + (i * colNum) + c, A.data + ((i - r) * A.colNum), sizeof(double) * A.colNum);
1149  }
1150  }
1151  }
1152  else {
1153  throw vpException(vpException::dimensionError, "Cannot insert (%dx%d) matrix in (%dx%d) matrix at position (%d,%d)",
1154  A.getRows(), A.getCols(), rowNum, colNum, r, c);
1155  }
1156 }
1157 
1199 {
1200  vpColVector evalue(rowNum); // Eigen values
1201 
1202  if (rowNum != colNum) {
1203  throw(vpException(vpException::dimensionError, "Cannot compute eigen values on a non square matrix (%dx%d)", rowNum,
1204  colNum));
1205  }
1206 
1207  // Check if the matrix is symmetric: At - A = 0
1208  vpMatrix At_A = (*this).t() - (*this);
1209  for (unsigned int i = 0; i < rowNum; ++i) {
1210  for (unsigned int j = 0; j < rowNum; ++j) {
1211  // in comment: if (At_A[i][j] != 0) {
1212  if (std::fabs(At_A[i][j]) > std::numeric_limits<double>::epsilon()) {
1213  throw(vpException(vpException::fatalError, "Cannot compute eigen values on a non symmetric matrix"));
1214  }
1215  }
1216  }
1217 
1218 #if defined(VISP_HAVE_LAPACK)
1219 #if defined(VISP_HAVE_GSL) /* be careful of the copy below */
1220  {
1221  gsl_vector *eval = gsl_vector_alloc(rowNum);
1222  gsl_matrix *evec = gsl_matrix_alloc(rowNum, colNum);
1223 
1224  gsl_eigen_symmv_workspace *w = gsl_eigen_symmv_alloc(rowNum);
1225  gsl_matrix *m = gsl_matrix_alloc(rowNum, colNum);
1226 
1227  unsigned int Atda = (unsigned int)m->tda;
1228  for (unsigned int i = 0; i < rowNum; i++) {
1229  unsigned int k = i * Atda;
1230  for (unsigned int j = 0; j < colNum; j++)
1231  m->data[k + j] = (*this)[i][j];
1232  }
1233  gsl_eigen_symmv(m, eval, evec, w);
1234 
1235  gsl_eigen_symmv_sort(eval, evec, GSL_EIGEN_SORT_ABS_ASC);
1236 
1237  for (unsigned int i = 0; i < rowNum; i++) {
1238  evalue[i] = gsl_vector_get(eval, i);
1239  }
1240 
1241  gsl_eigen_symmv_free(w);
1242  gsl_vector_free(eval);
1243  gsl_matrix_free(m);
1244  gsl_matrix_free(evec);
1245  }
1246 #else
1247  {
1248  const char jobz = 'N';
1249  const char uplo = 'U';
1250  vpMatrix A = (*this);
1251  vpColVector WORK;
1252  int lwork = -1;
1253  int info = 0;
1254  double wkopt;
1255  vpMatrix::blas_dsyev(jobz, uplo, rowNum, A.data, colNum, evalue.data, &wkopt, lwork, info);
1256  lwork = static_cast<int>(wkopt);
1257  WORK.resize(lwork);
1258  vpMatrix::blas_dsyev(jobz, uplo, rowNum, A.data, colNum, evalue.data, WORK.data, lwork, info);
1259  }
1260 #endif
1261 #else
1262  throw(vpException(vpException::functionNotImplementedError, "Eigen values computation is not implemented. "
1263  "You should install Lapack 3rd party"));
1264 #endif
1265  return evalue;
1266 }
1267 
1321 void vpMatrix::eigenValues(vpColVector &evalue, vpMatrix &evector) const
1322 {
1323  if (rowNum != colNum) {
1324  throw(vpException(vpException::dimensionError, "Cannot compute eigen values on a non square matrix (%dx%d)", rowNum,
1325  colNum));
1326  }
1327 
1328  // Check if the matrix is symmetric: At - A = 0
1329  vpMatrix At_A = (*this).t() - (*this);
1330  for (unsigned int i = 0; i < rowNum; ++i) {
1331  for (unsigned int j = 0; j < rowNum; ++j) {
1332  // -- in comment: if (At_A[i][j] != 0) {
1333  if (std::fabs(At_A[i][j]) > std::numeric_limits<double>::epsilon()) {
1334  throw(vpException(vpException::fatalError, "Cannot compute eigen values on a non symmetric matrix"));
1335  }
1336  }
1337  }
1338 
1339  // Resize the output matrices
1340  evalue.resize(rowNum);
1341  evector.resize(rowNum, colNum);
1342 
1343 #if defined(VISP_HAVE_LAPACK)
1344 #if defined(VISP_HAVE_GSL) /* be careful of the copy below */
1345  {
1346  gsl_vector *eval = gsl_vector_alloc(rowNum);
1347  gsl_matrix *evec = gsl_matrix_alloc(rowNum, colNum);
1348 
1349  gsl_eigen_symmv_workspace *w = gsl_eigen_symmv_alloc(rowNum);
1350  gsl_matrix *m = gsl_matrix_alloc(rowNum, colNum);
1351 
1352  unsigned int Atda = (unsigned int)m->tda;
1353  for (unsigned int i = 0; i < rowNum; i++) {
1354  unsigned int k = i * Atda;
1355  for (unsigned int j = 0; j < colNum; j++)
1356  m->data[k + j] = (*this)[i][j];
1357  }
1358  gsl_eigen_symmv(m, eval, evec, w);
1359 
1360  gsl_eigen_symmv_sort(eval, evec, GSL_EIGEN_SORT_ABS_ASC);
1361 
1362  for (unsigned int i = 0; i < rowNum; i++) {
1363  evalue[i] = gsl_vector_get(eval, i);
1364  }
1365  Atda = (unsigned int)evec->tda;
1366  for (unsigned int i = 0; i < rowNum; i++) {
1367  unsigned int k = i * Atda;
1368  for (unsigned int j = 0; j < rowNum; j++) {
1369  evector[i][j] = evec->data[k + j];
1370  }
1371  }
1372 
1373  gsl_eigen_symmv_free(w);
1374  gsl_vector_free(eval);
1375  gsl_matrix_free(m);
1376  gsl_matrix_free(evec);
1377  }
1378 #else // defined(VISP_HAVE_GSL)
1379  {
1380  const char jobz = 'V';
1381  const char uplo = 'U';
1382  vpMatrix A = (*this);
1383  vpColVector WORK;
1384  int lwork = -1;
1385  int info = 0;
1386  double wkopt;
1387  vpMatrix::blas_dsyev(jobz, uplo, rowNum, A.data, colNum, evalue.data, &wkopt, lwork, info);
1388  lwork = static_cast<int>(wkopt);
1389  WORK.resize(lwork);
1390  vpMatrix::blas_dsyev(jobz, uplo, rowNum, A.data, colNum, evalue.data, WORK.data, lwork, info);
1391  evector = A.t();
1392  }
1393 #endif // defined(VISP_HAVE_GSL)
1394 #else
1395  throw(vpException(vpException::functionNotImplementedError, "Eigen values computation is not implemented. "
1396  "You should install Lapack 3rd party"));
1397 #endif
1398 }
1399 
1418 unsigned int vpMatrix::kernel(vpMatrix &kerAt, double svThreshold) const
1419 {
1420  unsigned int nbline = getRows();
1421  unsigned int nbcol = getCols();
1422 
1423  vpMatrix U, V; // Copy of the matrix, SVD function is destructive
1424  vpColVector sv;
1425  sv.resize(nbcol, false); // singular values
1426  V.resize(nbcol, nbcol, false); // V matrix of singular value decomposition
1427 
1428  // Copy and resize matrix to have at least as many rows as columns
1429  // kernel is computed in svd method only if the matrix has more rows than
1430  // columns
1431 
1432  if (nbline < nbcol) {
1433  U.resize(nbcol, nbcol, true);
1434  }
1435  else {
1436  U.resize(nbline, nbcol, false);
1437  }
1438 
1439  U.insert(*this, 0, 0);
1440 
1441  U.svd(sv, V);
1442 
1443  // Compute the highest singular value and rank of the matrix
1444  double maxsv = 0;
1445  for (unsigned int i = 0; i < nbcol; ++i) {
1446  if (sv[i] > maxsv) {
1447  maxsv = sv[i];
1448  }
1449  }
1450 
1451  unsigned int rank = 0;
1452  for (unsigned int i = 0; i < nbcol; ++i) {
1453  if (sv[i] >(maxsv * svThreshold)) {
1454  ++rank;
1455  }
1456  }
1457 
1458  kerAt.resize(nbcol - rank, nbcol);
1459  if (rank != nbcol) {
1460  unsigned int k = 0;
1461  for (unsigned int j = 0; j < nbcol; ++j) {
1462  // if( v.col(j) in kernel and non zero )
1463  if ((sv[j] <= (maxsv * svThreshold)) &&
1464  (std::fabs(V.getCol(j).sumSquare()) > std::numeric_limits<double>::epsilon())) {
1465  unsigned int v_rows = V.getRows();
1466  for (unsigned int i = 0; i < v_rows; ++i) {
1467  kerAt[k][i] = V[i][j];
1468  }
1469  ++k;
1470  }
1471  }
1472  }
1473 
1474  return rank;
1475 }
1476 
1493 unsigned int vpMatrix::nullSpace(vpMatrix &kerA, double svThreshold) const
1494 {
1495  unsigned int nbrow = getRows();
1496  unsigned int nbcol = getCols();
1497 
1498  vpMatrix U, V; // Copy of the matrix, SVD function is destructive
1499  vpColVector sv;
1500  sv.resize(nbcol, false); // singular values
1501  V.resize(nbcol, nbcol, false); // V matrix of singular value decomposition
1502 
1503  // Copy and resize matrix to have at least as many rows as columns
1504  // kernel is computed in svd method only if the matrix has more rows than
1505  // columns
1506 
1507  if (nbrow < nbcol) {
1508  U.resize(nbcol, nbcol, true);
1509  }
1510  else {
1511  U.resize(nbrow, nbcol, false);
1512  }
1513 
1514  U.insert(*this, 0, 0);
1515 
1516  U.svd(sv, V);
1517 
1518  // Compute the highest singular value and rank of the matrix
1519  double maxsv = sv[0];
1520 
1521  unsigned int rank = 0;
1522  for (unsigned int i = 0; i < nbcol; ++i) {
1523  if (sv[i] >(maxsv * svThreshold)) {
1524  ++rank;
1525  }
1526  }
1527 
1528  kerA.resize(nbcol, nbcol - rank);
1529  if (rank != nbcol) {
1530  unsigned int k = 0;
1531  for (unsigned int j = 0; j < nbcol; ++j) {
1532  // if( v.col(j) in kernel and non zero )
1533  if (sv[j] <= (maxsv * svThreshold)) {
1534  for (unsigned int i = 0; i < nbcol; ++i) {
1535  kerA[i][k] = V[i][j];
1536  }
1537  ++k;
1538  }
1539  }
1540  }
1541 
1542  return (nbcol - rank);
1543 }
1544 
1561 unsigned int vpMatrix::nullSpace(vpMatrix &kerA, int dim) const
1562 {
1563  unsigned int nbrow = getRows();
1564  unsigned int nbcol = getCols();
1565  unsigned int dim_ = static_cast<unsigned int>(dim);
1566 
1567  vpMatrix U, V; // Copy of the matrix, SVD function is destructive
1568  vpColVector sv;
1569  sv.resize(nbcol, false); // singular values
1570  V.resize(nbcol, nbcol, false); // V matrix of singular value decomposition
1571 
1572  // Copy and resize matrix to have at least as many rows as columns
1573  // kernel is computed in svd method only if the matrix has more rows than
1574  // columns
1575 
1576  if (nbrow < nbcol) {
1577  U.resize(nbcol, nbcol, true);
1578  }
1579  else {
1580  U.resize(nbrow, nbcol, false);
1581  }
1582 
1583  U.insert(*this, 0, 0);
1584 
1585  U.svd(sv, V);
1586 
1587  kerA.resize(nbcol, dim_);
1588  if (dim_ != 0) {
1589  unsigned int rank = nbcol - dim_;
1590  for (unsigned int k = 0; k < dim_; ++k) {
1591  unsigned int j = k + rank;
1592  for (unsigned int i = 0; i < nbcol; ++i) {
1593  kerA[i][k] = V[i][j];
1594  }
1595  }
1596  }
1597 
1598  double maxsv = sv[0];
1599  unsigned int rank = 0;
1600  for (unsigned int i = 0; i < nbcol; ++i) {
1601  if (sv[i] >(maxsv * 1e-6)) {
1602  ++rank;
1603  }
1604  }
1605  return (nbcol - rank);
1606 }
1607 
1643 double vpMatrix::det(vpDetMethod method) const
1644 {
1645  double det = 0.;
1646 
1647  if (method == LU_DECOMPOSITION) {
1648  det = this->detByLU();
1649  }
1650 
1651  return det;
1652 }
1653 
1655 {
1656 #if defined(VISP_HAVE_EIGEN3)
1657  return choleskyByEigen3();
1658 #elif defined(VISP_HAVE_LAPACK)
1659  return choleskyByLapack();
1660 #elif defined(VISP_HAVE_OPENCV)
1661  return choleskyByOpenCV();
1662 #else
1663  throw(vpException(vpException::fatalError, "Cannot compute matrix Chloesky decomposition. "
1664  "Install Lapack, Eigen3 or OpenCV 3rd party"));
1665 #endif
1666 }
1667 
1676 {
1677  if (colNum != rowNum) {
1678  throw(vpException(vpException::dimensionError, "Cannot compute the exponential of a non square (%dx%d) matrix",
1679  rowNum, colNum));
1680  }
1681  else {
1682 #ifdef VISP_HAVE_GSL
1683  size_t size_ = rowNum * colNum;
1684  double *b = new double[size_];
1685  for (size_t i = 0; i < size_; i++)
1686  b[i] = 0.;
1687  gsl_matrix_view m = gsl_matrix_view_array(this->data, rowNum, colNum);
1688  gsl_matrix_view em = gsl_matrix_view_array(b, rowNum, colNum);
1689  gsl_linalg_exponential_ss(&m.matrix, &em.matrix, 0);
1690  // gsl_matrix_fprintf(stdout, &em.matrix, "%g");
1691  vpMatrix expA;
1692  expA.resize(rowNum, colNum, false);
1693  memcpy(expA.data, b, size_ * sizeof(double));
1694 
1695  delete[] b;
1696  return expA;
1697 #else
1698  vpMatrix v_expE(rowNum, colNum, false);
1699  vpMatrix v_expD(rowNum, colNum, false);
1700  vpMatrix v_expX(rowNum, colNum, false);
1701  vpMatrix v_expcX(rowNum, colNum, false);
1702  vpMatrix v_eye(rowNum, colNum, false);
1703 
1704  v_eye.eye();
1705  vpMatrix exp(*this);
1706 
1707  // -- in comment: double f;
1708  int e;
1709  double c = 0.5;
1710  int q = 6;
1711  int p = 1;
1712 
1713  double nA = 0;
1714  for (unsigned int i = 0; i < rowNum; ++i) {
1715  double sum = 0;
1716  for (unsigned int j = 0; j < colNum; ++j) {
1717  sum += fabs((*this)[i][j]);
1718  }
1719  if ((sum > nA) || (i == 0)) {
1720  nA = sum;
1721  }
1722  }
1723 
1724  /* f = */ frexp(nA, &e);
1725  // -- in comment: double s = (0 > e+1)?0:e+1;
1726  double s = e + 1;
1727 
1728  double sca = 1.0 / pow(2.0, s);
1729  exp = sca * exp;
1730  v_expX = *this;
1731  v_expE = (c * exp) + v_eye;
1732  v_expD = (-c * exp) + v_eye;
1733  for (int k = 2; k <= q; ++k) {
1734  c = (c * (static_cast<double>((q - k) + 1))) / (static_cast<double>(k * (((2.0 * q) - k) + 1)));
1735  v_expcX = exp * v_expX;
1736  v_expX = v_expcX;
1737  v_expcX = c * v_expX;
1738  v_expE = v_expE + v_expcX;
1739  if (p) {
1740  v_expD = v_expD + v_expcX;
1741  }
1742  else {
1743  v_expD = v_expD - v_expcX;
1744  }
1745  p = !p;
1746  }
1747  v_expX = v_expD.inverseByLU();
1748  exp = v_expX * v_expE;
1749  for (int k = 1; k <= s; ++k) {
1750  v_expE = exp * exp;
1751  exp = v_expE;
1752  }
1753  return exp;
1754 #endif
1755  }
1756 }
1757 
1758 /**************************************************************************************************************/
1759 /**************************************************************************************************************/
1760 
1761 // Specific functions
1762 
1763 /*
1764  input:: matrix M(nCols,nRows), nCols > 3, nRows > 3 , nCols == nRows.
1765 
1766  output:: the complement matrix of the element (rowNo,colNo).
1767  This is the matrix obtained from M after elimenating the row rowNo and column
1768  colNo
1769 
1770  example:
1771  1 2 3
1772  M = 4 5 6
1773  7 8 9
1774  1 3
1775  subblock(M, 1, 1) give the matrix 7 9
1776 */
1777 vpMatrix subblock(const vpMatrix &M, unsigned int col, unsigned int row)
1778 {
1779  vpMatrix M_comp;
1780  M_comp.resize(M.getRows() - 1, M.getCols() - 1, false);
1781 
1782  unsigned int m_rows = M.getRows();
1783  unsigned int m_col = M.getCols();
1784  for (unsigned int i = 0; i < col; ++i) {
1785  for (unsigned int j = 0; j < row; ++j) {
1786  M_comp[i][j] = M[i][j];
1787  }
1788  for (unsigned int j = row + 1; j < m_rows; ++j) {
1789  M_comp[i][j - 1] = M[i][j];
1790  }
1791  }
1792  for (unsigned int i = col + 1; i < m_col; ++i) {
1793  for (unsigned int j = 0; j < row; ++j) {
1794  M_comp[i - 1][j] = M[i][j];
1795  }
1796  for (unsigned int j = row + 1; j < m_rows; ++j) {
1797  M_comp[i - 1][j - 1] = M[i][j];
1798  }
1799  }
1800  return M_comp;
1801 }
1802 
1812 double vpMatrix::cond(double svThreshold) const
1813 {
1814  unsigned int nbline = getRows();
1815  unsigned int nbcol = getCols();
1816 
1817  vpMatrix U, V; // Copy of the matrix, SVD function is destructive
1818  vpColVector sv;
1819  sv.resize(nbcol); // singular values
1820  V.resize(nbcol, nbcol, false); // V matrix of singular value decomposition
1821 
1822  // Copy and resize matrix to have at least as many rows as columns
1823  // kernel is computed in svd method only if the matrix has more rows than
1824  // columns
1825 
1826  if (nbline < nbcol) {
1827  U.resize(nbcol, nbcol, true);
1828  }
1829  else {
1830  U.resize(nbline, nbcol, false);
1831  }
1832 
1833  U.insert(*this, 0, 0);
1834 
1835  U.svd(sv, V);
1836 
1837  // Compute the highest singular value
1838  double maxsv = 0;
1839  for (unsigned int i = 0; i < nbcol; ++i) {
1840  if (sv[i] > maxsv) {
1841  maxsv = sv[i];
1842  }
1843  }
1844 
1845  // Compute the rank of the matrix
1846  unsigned int rank = 0;
1847  for (unsigned int i = 0; i < nbcol; ++i) {
1848  if (sv[i] >(maxsv * svThreshold)) {
1849  ++rank;
1850  }
1851  }
1852 
1853  // Compute the lowest singular value
1854  double minsv = maxsv;
1855  for (unsigned int i = 0; i < rank; ++i) {
1856  if (sv[i] < minsv) {
1857  minsv = sv[i];
1858  }
1859  }
1860 
1861  if (std::fabs(minsv) > std::numeric_limits<double>::epsilon()) {
1862  return maxsv / minsv;
1863  }
1864  else {
1865  return std::numeric_limits<double>::infinity();
1866  }
1867 }
1868 
1875 void vpMatrix::computeHLM(const vpMatrix &H, const double &alpha, vpMatrix &HLM)
1876 {
1877  if (H.getCols() != H.getRows()) {
1878  throw(vpException(vpException::dimensionError, "Cannot compute HLM on a non square matrix (%dx%d)", H.getRows(),
1879  H.getCols()));
1880  }
1881 
1882  HLM = H;
1883  unsigned int h_col = H.getCols();
1884  for (unsigned int i = 0; i < h_col; ++i) {
1885  HLM[i][i] += alpha * H[i][i];
1886  }
1887 }
1888 
1897 {
1898  double norm = 0.0;
1899  for (unsigned int i = 0; i < dsize; ++i) {
1900  double x = *(data + i);
1901  norm += x * x;
1902  }
1903 
1904  return sqrt(norm);
1905 }
1906 
1916 {
1917  if (this->dsize != 0) {
1918  vpMatrix v;
1919  vpColVector w;
1920 
1921  vpMatrix M = *this;
1922 
1923  M.svd(w, v);
1924 
1925  double max = w[0];
1926  unsigned int maxRank = std::min<unsigned int>(this->getCols(), this->getRows());
1927  // The maximum reachable rank is either the number of columns or the number of rows
1928  // of the matrix.
1929  unsigned int boundary = std::min<unsigned int>(maxRank, w.size());
1930  // boundary is here to ensure that the number of singular values used for the com-
1931  // putation of the euclidean norm of the matrix is not greater than the maximum
1932  // reachable rank. Indeed, some svd library pad the singular values vector with 0s
1933  // if the input matrix is non-square.
1934  for (unsigned int i = 0; i < boundary; ++i) {
1935  if (max < w[i]) {
1936  max = w[i];
1937  }
1938  }
1939  return max;
1940  }
1941  else {
1942  return 0.;
1943  }
1944 }
1945 
1957 {
1958  double norm = 0.0;
1959  for (unsigned int i = 0; i < rowNum; ++i) {
1960  double x = 0;
1961  for (unsigned int j = 0; j < colNum; ++j) {
1962  x += fabs(*(*(rowPtrs + i) + j));
1963  }
1964  if (x > norm) {
1965  norm = x;
1966  }
1967  }
1968  return norm;
1969 }
1970 
1977 double vpMatrix::sumSquare() const
1978 {
1979  double sum_square = 0.0;
1980  double x;
1981 
1982  for (unsigned int i = 0; i < rowNum; ++i) {
1983  for (unsigned int j = 0; j < colNum; ++j) {
1984  x = rowPtrs[i][j];
1985  sum_square += x * x;
1986  }
1987  }
1988 
1989  return sum_square;
1990 }
1991 #if defined(VISP_BUILD_DEPRECATED_FUNCTIONS)
2001 double vpMatrix::euclideanNorm() const { return frobeniusNorm(); }
2002 
2021 vpRowVector vpMatrix::row(unsigned int i)
2022 {
2023  vpRowVector c(getCols());
2024 
2025  for (unsigned int j = 0; j < getCols(); ++j) {
2026  c[j] = (*this)[i - 1][j];
2027  }
2028  return c;
2029 }
2030 
2048 vpColVector vpMatrix::column(unsigned int j)
2049 {
2050  vpColVector c(getRows());
2051 
2052  for (unsigned int i = 0; i < getRows(); ++i) {
2053  c[i] = (*this)[i][j - 1];
2054  }
2055  return c;
2056 }
2057 
2064 void vpMatrix::setIdentity(const double &val)
2065 {
2066  for (unsigned int i = 0; i < rowNum; ++i) {
2067  for (unsigned int j = 0; j < colNum; ++j) {
2068  if (i == j) {
2069  (*this)[i][j] = val;
2070  }
2071  else {
2072  (*this)[i][j] = 0;
2073  }
2074  }
2075  }
2076 }
2077 
2078 #endif //#if defined(VISP_BUILD_DEPRECATED_FUNCTIONS)
2079 
2080 END_VISP_NAMESPACE
Implementation of a generic 2D array used as base class for matrices and vectors.
Definition: vpArray2D.h:145
unsigned int getCols() const
Definition: vpArray2D.h:417
double * data
Address of the first element of the data array.
Definition: vpArray2D.h:148
double ** rowPtrs
Address of the first element of each rows.
Definition: vpArray2D.h:1190
void insert(const vpArray2D< Type > &A, unsigned int r, unsigned int c)
Definition: vpArray2D.h:580
void resize(unsigned int nrows, unsigned int ncols, bool flagNullify=true, bool recopy_=true)
Definition: vpArray2D.h:442
unsigned int rowNum
Number of rows in the array.
Definition: vpArray2D.h:1186
static vpArray2D< Type > view(const vpArray2D< Type > &A)
Creates a view of the Matrix A. A view shares the same underlying memory as the original array....
Definition: vpArray2D.h:323
unsigned int dsize
Current array size (rowNum * colNum)
Definition: vpArray2D.h:1192
unsigned int size() const
Return the number of elements of the 2D array.
Definition: vpArray2D.h:429
unsigned int getRows() const
Definition: vpArray2D.h:427
unsigned int colNum
Number of columns in the array.
Definition: vpArray2D.h:1188
Implementation of column vector and the associated operations.
Definition: vpColVector.h:191
double sumSquare() const
void resize(unsigned int i, bool flagNullify=true)
Definition: vpColVector.h:1148
error that can be emitted by ViSP classes.
Definition: vpException.h:60
@ functionNotImplementedError
Function not implemented.
Definition: vpException.h:66
@ dimensionError
Bad dimension.
Definition: vpException.h:71
@ fatalError
Fatal error.
Definition: vpException.h:72
Implementation of an homogeneous matrix and operations on such kind of matrices.
static Type maximum(const Type &a, const Type &b)
Definition: vpMath.h:254
Implementation of a matrix and operations on matrices.
Definition: vpMatrix.h:169
vpColVector eigenValues() const
Definition: vpMatrix.cpp:1198
std::ostream & cppPrint(std::ostream &os, const std::string &matrixName="A", bool octet=false) const
Definition: vpMatrix.cpp:1106
double cond(double svThreshold=1e-6) const
Definition: vpMatrix.cpp:1812
vpMatrix expm() const
Definition: vpMatrix.cpp:1675
int print(std::ostream &s, unsigned int length, const std::string &intro="") const
Definition: vpMatrix.cpp:830
std::ostream & maplePrint(std::ostream &os) const
Definition: vpMatrix.cpp:1004
void init(const vpMatrix &M, unsigned int r, unsigned int c, unsigned int nrows, unsigned int ncols)
Definition: vpMatrix.cpp:363
vpMatrix cholesky() const
Definition: vpMatrix.cpp:1654
unsigned int kernel(vpMatrix &kerAt, double svThreshold=1e-6) const
Definition: vpMatrix.cpp:1418
vpMatrix inverseByLU() const
static vpMatrix juxtaposeMatrices(const vpMatrix &A, const vpMatrix &B)
Definition: vpMatrix.cpp:763
void svd(vpColVector &w, vpMatrix &V)
double sum() const
void diag(const double &val=1.0)
vpRowVector getRow(unsigned int i) const
Definition: vpMatrix.cpp:596
vpMatrix()
Definition: vpMatrix.h:185
vpMatrix choleskyByOpenCV() const
Compute the Cholesky decomposition of a Hermitian positive-definite matrix using OpenCV library.
vpMatrix choleskyByEigen3() const
Compute the Cholesky decomposition of a Hermitian positive-definite matrix using Eigen3 library.
double inducedL2Norm() const
Definition: vpMatrix.cpp:1915
double infinityNorm() const
Definition: vpMatrix.cpp:1956
double frobeniusNorm() const
Definition: vpMatrix.cpp:1896
vpColVector getDiag() const
Definition: vpMatrix.cpp:693
static void computeHLM(const vpMatrix &H, const double &alpha, vpMatrix &HLM)
Definition: vpMatrix.cpp:1875
vpColVector getCol(unsigned int j) const
Definition: vpMatrix.cpp:554
double detByLU() const
vpMatrix choleskyByLapack() const
Compute the Cholesky decomposition of a Hermitian positive-definite matrix using Lapack library.
double det(vpDetMethod method=LU_DECOMPOSITION) const
Definition: vpMatrix.cpp:1643
double sumSquare() const
Definition: vpMatrix.cpp:1977
void insert(const vpMatrix &A, unsigned int r, unsigned int c)
Definition: vpMatrix.cpp:1139
@ LU_DECOMPOSITION
Definition: vpMatrix.h:177
vpMatrix t() const
std::ostream & csvPrint(std::ostream &os) const
Definition: vpMatrix.cpp:1051
static vpMatrix view(double *data, unsigned int rows, unsigned int cols)
Create a matrix view of a raw data array. The view can modify the contents of the raw data array,...
Definition: vpMatrix.cpp:306
unsigned int nullSpace(vpMatrix &kerA, double svThreshold=1e-6) const
Definition: vpMatrix.cpp:1493
std::ostream & matlabPrint(std::ostream &os) const
Definition: vpMatrix.cpp:953
vpMatrix extract(unsigned int r, unsigned int c, unsigned int nrows, unsigned int ncols) const
Definition: vpMatrix.cpp:430
Implementation of a rotation matrix and operations on such kind of matrices.
Implementation of row vector and the associated operations.
Definition: vpRowVector.h:124
Class that consider the case of a translation vector.