Visual Servoing Platform  version 3.6.1 under development (2024-07-27)
vpMatrix_svd.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 SVD decomposition.
32  */
33 
34 #include <visp3/core/vpColVector.h>
35 #include <visp3/core/vpConfig.h>
36 #include <visp3/core/vpException.h>
37 #include <visp3/core/vpMath.h>
38 #include <visp3/core/vpMatrix.h>
39 #include <visp3/core/vpMatrixException.h>
40 
41 #include <cmath> // std::fabs
42 #include <iostream>
43 
44 #ifdef VISP_HAVE_EIGEN3
45 #include <Eigen/SVD>
46 #endif
47 
48 #if defined(VISP_HAVE_OPENCV) // Require opencv >= 2.1.1
49 #include <opencv2/core/core.hpp>
50 #endif
51 
52 #ifdef VISP_HAVE_LAPACK
53 #ifdef VISP_HAVE_GSL
54 #include <gsl/gsl_linalg.h>
55 #endif
56 #ifdef VISP_HAVE_MKL
57 #include <mkl.h>
58 typedef MKL_INT integer;
59 #else
60 #if defined(VISP_HAVE_LAPACK_BUILT_IN)
61 typedef long int integer;
62 #else
63 typedef int integer;
64 #endif
65 extern "C" int dgesdd_(char *jobz, integer *m, integer *n, double *a, integer *lda, double *s, double *u, integer *ldu,
66  double *vt, integer *ldvt, double *work, integer *lwork, integer *iwork, integer *info);
67 
68 #include <stdio.h>
69 #include <string.h>
70 #endif
71 #endif
72 
74 
129 void vpMatrix::solveBySVD(const vpColVector &b, vpColVector &x) const { x = pseudoInverse(1e-6) * b; }
130 
185 {
186  vpColVector X(colNum);
187 
188  solveBySVD(B, X);
189  return X;
190 }
191 
260 {
261 #if defined(VISP_HAVE_LAPACK)
262  svdLapack(w, V);
263 #elif defined(VISP_HAVE_EIGEN3)
264  svdEigen3(w, V);
265 #elif defined(VISP_HAVE_OPENCV) // Require opencv >= 2.1.1
266  svdOpenCV(w, V);
267 #else
268  (void)w;
269  (void)V;
270  throw(vpException(vpException::fatalError, "Cannot compute SVD. Install Lapack, Eigen3 or OpenCV 3rd party"));
271 #endif
272 }
273 
274 #if defined(VISP_HAVE_OPENCV) // Require opencv >= 2.1.1
275 
346 {
347  int rows = (int)this->getRows();
348  int cols = (int)this->getCols();
349  cv::Mat m(rows, cols, CV_64F, this->data);
350  cv::SVD opencvSVD(m);
351  cv::Mat opencvV = opencvSVD.vt;
352  cv::Mat opencvW = opencvSVD.w;
353  V.resize((unsigned int)opencvV.rows, (unsigned int)opencvV.cols);
354  w.resize((unsigned int)(opencvW.rows * opencvW.cols));
355 
356  memcpy(V.data, opencvV.data, (size_t)(8 * opencvV.rows * opencvV.cols));
357  V = V.transpose();
358  memcpy(w.data, opencvW.data, (size_t)(8 * opencvW.rows * opencvW.cols));
359  this->resize((unsigned int)opencvSVD.u.rows, (unsigned int)opencvSVD.u.cols);
360  memcpy(this->data, opencvSVD.u.data, (size_t)(8 * opencvSVD.u.rows * opencvSVD.u.cols));
361 }
362 
363 #endif
364 
365 #if defined(VISP_HAVE_LAPACK)
437 {
438 #ifdef VISP_HAVE_GSL
439  {
440  // GSL cannot consider M < N. In that case we transpose input matrix
441  vpMatrix U;
442  unsigned int nc = getCols();
443  unsigned int nr = getRows();
444 
445  if (rowNum < colNum) {
446  U = this->transpose();
447  nc = getRows();
448  nr = getCols();
449  }
450  else {
451  nc = getCols();
452  nr = getRows();
453  }
454 
455  w.resize(nc);
456  V.resize(nc, nc);
457 
458  gsl_vector *work = gsl_vector_alloc(nc);
459 
460  gsl_matrix A;
461  A.size1 = nr;
462  A.size2 = nc;
463  A.tda = A.size2;
464  if (rowNum < colNum) {
465  A.data = U.data;
466  }
467  else {
468  A.data = this->data;
469  }
470  A.owner = 0;
471  A.block = 0;
472 
473  gsl_matrix V_;
474  V_.size1 = nc;
475  V_.size2 = nc;
476  V_.tda = V_.size2;
477  V_.data = V.data;
478  V_.owner = 0;
479  V_.block = 0;
480 
481  gsl_vector S;
482  S.size = nc;
483  S.stride = 1;
484  S.data = w.data;
485  S.owner = 0;
486  S.block = 0;
487 
488  gsl_linalg_SV_decomp(&A, &V_, &S, work);
489 
490  if (rowNum < colNum) {
491  (*this) = V;
492  V = U;
493  }
494 
495  gsl_vector_free(work);
496  }
497 #else
498  {
499  vpMatrix U;
500  unsigned int nc = getCols();
501  unsigned int nr = getRows();
502 
503  if (rowNum < colNum) {
504  U = this->transpose();
505  nc = getRows();
506  nr = getCols();
507  }
508  else {
509  nc = getCols();
510  nr = getRows();
511  }
512 
513  w.resize(nc);
514  V.resize(nc, nc);
515 
516  double *a = new double[static_cast<unsigned int>(nr * nc)];
517  if (rowNum < colNum) {
518  memcpy(a, U.data, this->getRows() * this->getCols() * sizeof(double));
519  }
520  else {
521  memcpy(a, this->data, this->getRows() * this->getCols() * sizeof(double));
522  }
523 
524  integer m = (integer)(nc);
525  integer n = (integer)(nr);
526  integer lda = nc;
527  integer ldu = nc;
528  integer ldvt = std::min<integer>(nr, nc);
529  integer info, lwork;
530 
531  double wkopt;
532  double *work;
533 
534  integer *iwork = new integer[8 * static_cast<integer>(std::min<integer>(nr, nc))];
535 
536  double *s = w.data;
537  double *u = V.data;
538  double *vt;
539  if (rowNum < colNum) {
540  vt = U.data;
541  }
542  else {
543  vt = this->data;
544  }
545 
546  lwork = -1;
547  dgesdd_((char *)"S", &m, &n, a, &lda, s, u, &ldu, vt, &ldvt, &wkopt, &lwork, iwork, &info);
548  lwork = (int)wkopt;
549  work = new double[static_cast<unsigned int>(lwork)];
550 
551  dgesdd_((char *)"S", &m, &n, a, &lda, s, u, &ldu, vt, &ldvt, work, &lwork, iwork, &info);
552 
553  if (info > 0) {
554  throw(vpMatrixException(vpMatrixException::fatalError, "The algorithm computing SVD failed to converge."));
555  }
556 
557  if (rowNum < colNum) {
558  (*this) = V.transpose();
559  V = U;
560  }
561  else {
562  V = V.transpose();
563  }
564  delete[] work;
565  delete[] iwork;
566  delete[] a;
567  }
568 #endif
569 }
570 #endif
571 
572 #ifdef VISP_HAVE_EIGEN3
644 {
645  if (rowNum < colNum) {
646  w.resize(rowNum);
647  V.resize(colNum, rowNum);
648  }
649  else {
650  w.resize(colNum);
651  V.resize(colNum, colNum);
652  }
653 
654  Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> > M(this->data, this->getRows(),
655  this->getCols());
656  Eigen::JacobiSVD<Eigen::MatrixXd> svd(M, Eigen::ComputeThinU | Eigen::ComputeThinV);
657 
658  Eigen::Map<Eigen::VectorXd> w_(w.data, w.size());
659 
660  if (rowNum < colNum) {
661  this->resize(rowNum, rowNum);
662  }
663  else {
664  this->resize(rowNum, colNum);
665  }
666  Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> > V_(V.data, V.getRows(), V.getCols());
667  Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> > U_(this->data, rowNum, colNum);
668  w_ = svd.singularValues();
669  V_ = svd.matrixV();
670  U_ = svd.matrixU();
671 }
672 #endif
673 END_VISP_NAMESPACE
unsigned int getCols() const
Definition: vpArray2D.h:337
double * data
Address of the first element of the data array.
Definition: vpArray2D.h:148
void resize(unsigned int nrows, unsigned int ncols, bool flagNullify=true, bool recopy_=true)
Definition: vpArray2D.h:362
unsigned int rowNum
Number of rows in the array.
Definition: vpArray2D.h:1098
unsigned int size() const
Return the number of elements of the 2D array.
Definition: vpArray2D.h:349
unsigned int getRows() const
Definition: vpArray2D.h:347
unsigned int colNum
Number of columns in the array.
Definition: vpArray2D.h:1100
Implementation of column vector and the associated operations.
Definition: vpColVector.h:191
void resize(unsigned int i, bool flagNullify=true)
Definition: vpColVector.h:1143
error that can be emitted by ViSP classes.
Definition: vpException.h:60
@ fatalError
Fatal error.
Definition: vpException.h:72
error that can be emitted by the vpMatrix class and its derivatives
Implementation of a matrix and operations on matrices.
Definition: vpMatrix.h:169
void svdLapack(vpColVector &w, vpMatrix &V)
void solveBySVD(const vpColVector &B, vpColVector &x) const
void svd(vpColVector &w, vpMatrix &V)
void svdOpenCV(vpColVector &w, vpMatrix &V)
vpMatrix pseudoInverse(double svThreshold=1e-6) const
void svdEigen3(vpColVector &w, vpMatrix &V)
vpMatrix transpose() const