Visual Servoing Platform  version 3.1.0
vpMatrix_lu.cpp
1 /****************************************************************************
2  *
3  * This file is part of the ViSP software.
4  * Copyright (C) 2005 - 2017 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 LU decomposition.
33  *
34  * Authors:
35  * Eric Marchand
36  *
37  *****************************************************************************/
38 
39 #include <visp3/core/vpConfig.h>
40 
41 #include <visp3/core/vpColVector.h>
42 #include <visp3/core/vpMath.h>
43 #include <visp3/core/vpMatrix.h>
44 
45 #ifdef VISP_HAVE_EIGEN3
46 #include <Eigen/LU>
47 #endif
48 
49 #ifdef VISP_HAVE_GSL
50 #include <gsl/gsl_linalg.h>
51 #include <gsl/gsl_permutation.h>
52 #endif
53 
54 #ifdef VISP_HAVE_LAPACK
55 #ifdef VISP_HAVE_LAPACK_BUILT_IN
56 typedef long int integer;
57 #else
58 typedef int integer;
59 #endif
60 
61 extern "C" int dgetrf_(integer *m, integer *n, double *a, integer *lda, integer *ipiv, integer *info);
62 extern "C" void dgetri_(integer *n, double *a, integer *lda, integer *ipiv, double *work, integer *lwork,
63  integer *info);
64 #endif
65 
66 #if (VISP_HAVE_OPENCV_VERSION >= 0x020101) // Require opencv >= 2.1.1
67 #include <opencv2/core/core.hpp>
68 #endif
69 
70 // Exception
71 #include <visp3/core/vpException.h>
72 #include <visp3/core/vpMatrixException.h>
73 
74 #include <cmath> // std::fabs
75 #include <limits> // numeric_limits
76 
77 /*--------------------------------------------------------------------
78  LU Decomposition related functions
79 -------------------------------------------------------------------- */
80 
131 {
132 #if defined(VISP_HAVE_LAPACK)
133  return inverseByLULapack();
134 #elif defined(VISP_HAVE_EIGEN3)
135  return inverseByLUEigen3();
136 #elif (VISP_HAVE_OPENCV_VERSION >= 0x020101)
137  return inverseByLUOpenCV();
138 #elif defined(VISP_HAVE_GSL)
139  return inverseByLUGsl();
140 #else
141  throw(vpException(vpException::fatalError, "Cannot compute matrix determinant. Install Eigen3, "
142  "Lapack, OpenCV or GSL 3rd party"));
143 #endif
144 }
145 
180 double vpMatrix::detByLU() const
181 {
182  if (rowNum == 2 && colNum == 2) {
183  return ((*this)[0][0] * (*this)[1][1] - (*this)[0][1] * (*this)[1][0]);
184  } else if (rowNum == 3 && colNum == 3) {
185  return ((*this)[0][0] * ((*this)[1][1] * (*this)[2][2] - (*this)[1][2] * (*this)[2][1]) -
186  (*this)[0][1] * ((*this)[1][0] * (*this)[2][2] - (*this)[1][2] * (*this)[2][0]) +
187  (*this)[0][2] * ((*this)[1][0] * (*this)[2][1] - (*this)[1][1] * (*this)[2][0]));
188  } else {
189 #if defined(VISP_HAVE_LAPACK)
190  return detByLULapack();
191 #elif defined(VISP_HAVE_EIGEN3)
192  return detByLUEigen3();
193 #elif (VISP_HAVE_OPENCV_VERSION >= 0x020101)
194  return detByLUOpenCV();
195 #elif defined(VISP_HAVE_GSL)
196  return detByLUGsl();
197 #else
198  throw(vpException(vpException::fatalError, "Cannot compute matrix determinant. Install Lapack, "
199  "Eigen3, OpenCV or GSL 3rd party"));
200 #endif
201  }
202 }
203 
204 #ifndef DOXYGEN_SHOULD_SKIP_THIS
205 
206 #if defined(VISP_HAVE_GSL)
207 
238 vpMatrix vpMatrix::inverseByLUGsl() const
239 {
240  if (rowNum != colNum) {
241  throw(vpException(vpException::fatalError, "Cannot inverse a non square matrix (%ux%u) by LU", rowNum, colNum));
242  }
243 
244  gsl_matrix *A = gsl_matrix_alloc(rowNum, colNum);
245 
246  // copy the input matrix to ensure the function doesn't modify its content
247  unsigned int tda = (unsigned int)A->tda;
248  for (unsigned int i = 0; i < rowNum; i++) {
249  unsigned int k = i * tda;
250  for (unsigned int j = 0; j < colNum; j++)
251  A->data[k + j] = (*this)[i][j];
252  }
253 
254  vpMatrix Ainv(rowNum, colNum);
255 
256  gsl_matrix inverse;
257  inverse.size1 = rowNum;
258  inverse.size2 = colNum;
259  inverse.tda = inverse.size2;
260  inverse.data = Ainv.data;
261  inverse.owner = 0;
262  inverse.block = 0;
263 
264  gsl_permutation *p = gsl_permutation_alloc(rowNum);
265  int s;
266 
267  // Do the LU decomposition on A and use it to solve the system
268  gsl_linalg_LU_decomp(A, p, &s);
269  gsl_linalg_LU_invert(A, p, &inverse);
270 
271  gsl_permutation_free(p);
272  gsl_matrix_free(A);
273 
274  return Ainv;
275 }
276 
302 double vpMatrix::detByLUGsl() const
303 {
304  double det = 0.;
305 
306  if (rowNum != colNum) {
307  throw(vpException(vpException::fatalError, "Cannot compute matrix determinant of a non square matrix (%ux%u)",
308  rowNum, colNum));
309  }
310 
311  gsl_matrix *A = gsl_matrix_alloc(rowNum, colNum);
312 
313  // copy the input matrix to ensure the function doesn't modify its content
314  unsigned int tda = (unsigned int)A->tda;
315  for (unsigned int i = 0; i < rowNum; i++) {
316  unsigned int k = i * tda;
317  for (unsigned int j = 0; j < colNum; j++)
318  A->data[k + j] = (*this)[i][j];
319  }
320 
321  gsl_permutation *p = gsl_permutation_alloc(rowNum);
322  int s;
323 
324  // Do the LU decomposition on A and use it to solve the system
325  gsl_linalg_LU_decomp(A, p, &s);
326  det = gsl_linalg_LU_det(A, s);
327 
328  gsl_permutation_free(p);
329  gsl_matrix_free(A);
330 
331  return det;
332 }
333 #endif
334 
335 #ifdef VISP_HAVE_LAPACK
336 
367 vpMatrix vpMatrix::inverseByLULapack() const
368 {
369  if (rowNum != colNum) {
370  throw(vpException(vpException::fatalError, "Cannot inverse a non square matrix (%ux%u) by LU", rowNum, colNum));
371  }
372 
373  integer dim = (integer)rowNum;
374  integer lda = dim;
375  integer info;
376  integer lwork = dim * dim;
377  integer *ipiv = new integer[dim + 1];
378  double *work = new double[lwork];
379 
380  vpMatrix A = *this;
381 
382  dgetrf_(&dim, &dim, A.data, &lda, &ipiv[1], &info);
383  if (info) {
384  delete[] ipiv;
385  delete[] work;
386  throw(vpException(vpException::fatalError, "Lapack LU decomposition failed; info=%d", info));
387  }
388 
389  dgetri_(&dim, A.data, &dim, &ipiv[1], work, &lwork, &info);
390 
391  delete[] ipiv;
392  delete[] work;
393 
394  return A;
395 }
396 
422 double vpMatrix::detByLULapack() const
423 {
424  if (rowNum != colNum) {
425  throw(vpException(vpException::fatalError, "Cannot compute matrix determinant of a non square matrix (%ux%u)",
426  rowNum, colNum));
427  }
428 
429  integer dim = (integer)rowNum;
430  integer lda = dim;
431  integer info;
432  integer *ipiv = new integer[dim + 1];
433 
434  vpMatrix A = *this;
435 
436  dgetrf_(&dim, &dim, A.data, &lda, &ipiv[1], &info);
437  if (info) {
438  delete[] ipiv;
439  throw(vpException(vpException::fatalError, "Lapack LU decomposition failed; info=%d", info));
440  }
441 
442  double det = A[0][0];
443  for (unsigned int i = 1; i < rowNum; i++) {
444  det *= A[i][i];
445  }
446 
447  double sign = 1.;
448  for (int i = 1; i <= dim; i++) {
449  if (ipiv[i] != i)
450  sign = -sign;
451  }
452 
453  det *= sign;
454 
455  delete[] ipiv;
456 
457  return det;
458 }
459 #endif
460 
461 #if (VISP_HAVE_OPENCV_VERSION >= 0x020101)
462 
493 vpMatrix vpMatrix::inverseByLUOpenCV() const
494 {
495  if (rowNum != colNum) {
496  throw(vpException(vpException::fatalError, "Cannot inverse a non square matrix (%ux%u) by LU", rowNum, colNum));
497  }
498 
499  cv::Mat M(rowNum, colNum, CV_64F, this->data);
500 
501  cv::Mat Minv = M.inv(cv::DECOMP_LU);
502 
503  vpMatrix A(rowNum, colNum);
504  memcpy(A.data, Minv.data, (size_t)(8 * Minv.rows * Minv.cols));
505 
506  return A;
507 }
508 
534 double vpMatrix::detByLUOpenCV() const
535 {
536  double det = 0.;
537 
538  if (rowNum != colNum) {
539  throw(vpException(vpException::fatalError, "Cannot compute matrix determinant of a non square matrix (%ux%u)",
540  rowNum, colNum));
541  }
542 
543  cv::Mat M(rowNum, colNum, CV_64F, this->data);
544  det = cv::determinant(M);
545 
546  return (det);
547 }
548 #endif
549 
550 #if defined(VISP_HAVE_EIGEN3)
551 
583 vpMatrix vpMatrix::inverseByLUEigen3() const
584 {
585  if (rowNum != colNum) {
586  throw(vpException(vpException::fatalError, "Cannot inverse a non square matrix (%ux%u) by LU", rowNum, colNum));
587  }
588  vpMatrix A(this->getRows(), this->getCols());
589 
590  Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> > M(this->data, this->getRows(),
591  this->getCols());
592  Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> > A_(A.data, this->getRows(),
593  this->getCols());
594 
595  A_ = M.inverse();
596 
597  return A;
598 }
599 
625 double vpMatrix::detByLUEigen3() const
626 {
627  if (rowNum != colNum) {
628  throw(vpException(vpException::fatalError, "Cannot compute matrix determinant of a non square matrix (%ux%u)",
629  rowNum, colNum));
630  }
631 
632  Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> > M(this->data, this->getRows(),
633  this->getCols());
634 
635  return M.determinant();
636 }
637 #endif
638 
639 #endif // #ifndef DOXYGEN_SHOULD_SKIP_THIS
Implementation of a matrix and operations on matrices.
Definition: vpMatrix.h:104
double det(vpDetMethod method=LU_DECOMPOSITION) const
Definition: vpMatrix.cpp:4848
vpMatrix inverseByLU() const
error that can be emited by ViSP classes.
Definition: vpException.h:71
unsigned int getRows() const
Definition: vpArray2D.h:156
Type * data
Address of the first element of the data array.
Definition: vpArray2D.h:84
unsigned int getCols() const
Definition: vpArray2D.h:146
unsigned int rowNum
Number of rows in the array.
Definition: vpArray2D.h:74
unsigned int colNum
Number of columns in the array.
Definition: vpArray2D.h:76
double detByLU() const