Visual Servoing Platform  version 3.5.1 under development (2023-05-30)
vpMatrix_svd.cpp
1 /****************************************************************************
2  *
3  * ViSP, open source Visual Servoing Platform software.
4  * Copyright (C) 2005 - 2019 by Inria. All rights reserved.
5  *
6  * This software is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2 of the License, or
9  * (at your option) any later version.
10  * See the file LICENSE.txt at the root directory of this source
11  * distribution for additional information about the GNU GPL.
12  *
13  * For using ViSP with software that can not be combined with the GNU
14  * GPL, please contact Inria about acquiring a ViSP Professional
15  * Edition License.
16  *
17  * See http://visp.inria.fr for more information.
18  *
19  * This software was developed at:
20  * Inria Rennes - Bretagne Atlantique
21  * Campus Universitaire de Beaulieu
22  * 35042 Rennes Cedex
23  * France
24  *
25  * If you have questions regarding the use of this file, please contact
26  * Inria at visp@inria.fr
27  *
28  * This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
29  * WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
30  *
31  * Description:
32  * Matrix SVD decomposition.
33  *
34  * Authors:
35  * Eric Marchand
36  *
37  *****************************************************************************/
38 
39 #include <visp3/core/vpColVector.h>
40 #include <visp3/core/vpConfig.h>
41 #include <visp3/core/vpDebug.h>
42 #include <visp3/core/vpException.h>
43 #include <visp3/core/vpMath.h>
44 #include <visp3/core/vpMatrix.h>
45 #include <visp3/core/vpMatrixException.h>
46 
47 #include <cmath> // std::fabs
48 #include <iostream>
49 #include <limits> // numeric_limits
50 
51 #ifdef VISP_HAVE_EIGEN3
52 #include <Eigen/SVD>
53 #endif
54 
55 #if (VISP_HAVE_OPENCV_VERSION >= 0x020101) // Require opencv >= 2.1.1
56 #include <opencv2/core/core.hpp>
57 #endif
58 
59 #ifdef VISP_HAVE_LAPACK
60 #ifdef VISP_HAVE_GSL
61 #include <gsl/gsl_linalg.h>
62 #endif
63 #ifdef VISP_HAVE_MKL
64 #include <mkl.h>
65 typedef MKL_INT integer;
66 #else
67 #if defined(VISP_HAVE_LAPACK_BUILT_IN)
68 typedef long int integer;
69 #else
70 typedef int integer;
71 #endif
72 extern "C" int dgesdd_(char *jobz, integer *m, integer *n, double *a, integer *lda, double *s, double *u, integer *ldu,
73  double *vt, integer *ldvt, double *work, integer *lwork, integer *iwork, integer *info);
74 
75 #include <stdio.h>
76 #include <string.h>
77 #endif
78 #endif
79 /*---------------------------------------------------------------------
80 
81 SVD related functions
82 
83 ---------------------------------------------------------------------*/
84 
85 #if (VISP_HAVE_OPENCV_VERSION >= 0x020101) // Require opencv >= 2.1.1
86 
154 {
155  int rows = (int)this->getRows();
156  int cols = (int)this->getCols();
157  cv::Mat m(rows, cols, CV_64F, this->data);
158  cv::SVD opencvSVD(m);
159  cv::Mat opencvV = opencvSVD.vt;
160  cv::Mat opencvW = opencvSVD.w;
161  V.resize((unsigned int)opencvV.rows, (unsigned int)opencvV.cols);
162  w.resize((unsigned int)(opencvW.rows * opencvW.cols));
163 
164  memcpy(V.data, opencvV.data, (size_t)(8 * opencvV.rows * opencvV.cols));
165  V = V.transpose();
166  memcpy(w.data, opencvW.data, (size_t)(8 * opencvW.rows * opencvW.cols));
167  this->resize((unsigned int)opencvSVD.u.rows, (unsigned int)opencvSVD.u.cols);
168  memcpy(this->data, opencvSVD.u.data, (size_t)(8 * opencvSVD.u.rows * opencvSVD.u.cols));
169 }
170 
171 #endif
172 
173 #if defined(VISP_HAVE_LAPACK)
241 {
242 #ifdef VISP_HAVE_GSL
243  {
244  // GSL cannot consider M < N. In that case we transpose input matrix
245  vpMatrix U;
246  unsigned int nc = getCols();
247  unsigned int nr = getRows();
248 
249  if (rowNum < colNum) {
250  U = this->transpose();
251  nc = getRows();
252  nr = getCols();
253  } else {
254  nc = getCols();
255  nr = getRows();
256  }
257 
258  w.resize(nc);
259  V.resize(nc, nc);
260 
261  gsl_vector *work = gsl_vector_alloc(nc);
262 
263  gsl_matrix A;
264  A.size1 = nr;
265  A.size2 = nc;
266  A.tda = A.size2;
267  if (rowNum < colNum) {
268  A.data = U.data;
269  } else {
270  A.data = this->data;
271  }
272  A.owner = 0;
273  A.block = 0;
274 
275  gsl_matrix V_;
276  V_.size1 = nc;
277  V_.size2 = nc;
278  V_.tda = V_.size2;
279  V_.data = V.data;
280  V_.owner = 0;
281  V_.block = 0;
282 
283  gsl_vector S;
284  S.size = nc;
285  S.stride = 1;
286  S.data = w.data;
287  S.owner = 0;
288  S.block = 0;
289 
290  gsl_linalg_SV_decomp(&A, &V_, &S, work);
291 
292  if (rowNum < colNum) {
293  (*this) = V.transpose();
294  V = U;
295  }
296 
297  gsl_vector_free(work);
298  }
299 #else
300  {
301  w.resize(this->getCols());
302  V.resize(this->getCols(), this->getCols());
303 
304  integer m = (integer)(this->getCols());
305  integer n = (integer)(this->getRows());
306  integer lda = m;
307  integer ldu = m;
308  integer ldvt = (std::min)(m, n);
309  integer info, lwork;
310 
311  double wkopt;
312  double *work;
313 
314  integer *iwork = new integer[8 * static_cast<integer>((std::min)(n, m))];
315 
316  double *s = w.data;
317  double *a = new double[static_cast<unsigned int>(lda * n)];
318  memcpy(a, this->data, this->getRows() * this->getCols() * sizeof(double));
319  double *u = V.data;
320  double *vt = this->data;
321 
322  lwork = -1;
323  dgesdd_((char *)"S", &m, &n, a, &lda, s, u, &ldu, vt, &ldvt, &wkopt, &lwork, iwork, &info);
324  lwork = (int)wkopt;
325  work = new double[static_cast<unsigned int>(lwork)];
326 
327  dgesdd_((char *)"S", &m, &n, a, &lda, s, u, &ldu, vt, &ldvt, work, &lwork, iwork, &info);
328 
329  if (info > 0) {
330  throw(vpMatrixException(vpMatrixException::fatalError, "The algorithm computing SVD failed to converge."));
331  }
332 
333  V = V.transpose();
334  delete[] work;
335  delete[] iwork;
336  delete[] a;
337  }
338 #endif
339 }
340 #endif
341 
342 #ifdef VISP_HAVE_EIGEN3
410 {
411  w.resize(this->getCols());
412  V.resize(this->getCols(), this->getCols());
413 
414  Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> > M(this->data, this->getRows(),
415  this->getCols());
416 
417  Eigen::JacobiSVD<Eigen::MatrixXd> svd(M, Eigen::ComputeThinU | Eigen::ComputeThinV);
418 
419  Eigen::Map<Eigen::VectorXd> w_(w.data, w.size());
420  Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> > V_(V.data, V.getRows(),
421  V.getCols());
422  Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> > U_(this->data, this->getRows(),
423  this->getCols());
424  w_ = svd.singularValues();
425  V_ = svd.matrixV();
426  U_ = svd.matrixU();
427 }
428 #endif
unsigned int getCols() const
Definition: vpArray2D.h:282
double * data
Address of the first element of the data array.
Definition: vpArray2D.h:146
void resize(unsigned int nrows, unsigned int ncols, bool flagNullify=true, bool recopy_=true)
Definition: vpArray2D.h:307
unsigned int rowNum
Number of rows in the array.
Definition: vpArray2D.h:136
unsigned int size() const
Return the number of elements of the 2D array.
Definition: vpArray2D.h:294
unsigned int getRows() const
Definition: vpArray2D.h:292
unsigned int colNum
Number of columns in the array.
Definition: vpArray2D.h:138
Implementation of column vector and the associated operations.
Definition: vpColVector.h:172
void resize(unsigned int i, bool flagNullify=true)
Definition: vpColVector.h:357
@ fatalError
Fatal error.
Definition: vpException.h:96
error that can be emited by the vpMatrix class and its derivates
Implementation of a matrix and operations on matrices.
Definition: vpMatrix.h:154
void svdLapack(vpColVector &w, vpMatrix &V)
void svd(vpColVector &w, vpMatrix &V)
Definition: vpMatrix.cpp:2021
void svdOpenCV(vpColVector &w, vpMatrix &V)
void svdEigen3(vpColVector &w, vpMatrix &V)
vpMatrix transpose() const
Definition: vpMatrix.cpp:469