Visual Servoing Platform  version 3.2.0 under development (2019-01-22)
testMatrixInverse.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  * Test various inversions.
33  *
34  * Authors:
35  * Fabien Spindler
36  *
37  *****************************************************************************/
38 
44 #include <cmath>
45 #include <fstream>
46 #include <iostream>
47 #include <stdio.h>
48 #include <stdlib.h>
49 #include <vector>
50 
51 #include <visp3/core/vpColVector.h>
52 #include <visp3/core/vpMatrix.h>
53 #include <visp3/core/vpTime.h>
54 #include <visp3/io/vpParseArgv.h>
55 
56 // List of allowed command line options
57 #define GETOPTARGS "cdn:i:pf:R:C:vh"
58 
67 void usage(const char *name, const char *badparam)
68 {
69  fprintf(stdout, "\n\
70 Test matrix inversions\n\
71 using LU, QR and Cholesky methods as well as Pseudo-inverse.\n\
72 Outputs a comparison of these methods.\n\
73 \n\
74 SYNOPSIS\n\
75  %s [-n <number of matrices>] [-f <plot filename>]\n\
76  [-R <number of rows>] [-C <number of columns>]\n\
77  [-i <number of iterations>] [-p] [-h]\n", name);
78 
79  fprintf(stdout, "\n\
80 OPTIONS: Default\n\
81  -n <number of matrices> \n\
82  Number of matrices inverted during each test loop.\n\
83 \n\
84  -i <number of iterations> \n\
85  Number of iterations of the test.\n\
86 \n\
87  -f <plot filename> \n\
88  Set output path for plot output.\n\
89  The plot logs the times of \n\
90  the different inversion methods: \n\
91  QR,LU,Cholesky and Pseudo-inverse.\n\
92 \n\
93  -R <number of rows>\n\
94  Number of rows of the automatically generated matrices \n\
95  we test on.\n\
96 \n\
97  -C <number of columns>\n\
98  Number of colums of the automatically generated matrices \n\
99  we test on.\n\
100 \n\
101  -p \n\
102  Plot into filename in the gnuplot format. \n\
103  If this option is used, tests results will be logged \n\
104  into a filename specified with -f.\n\
105 \n\
106  -h\n\
107  Print the help.\n\n");
108 
109  if (badparam) {
110  fprintf(stderr, "ERROR: \n");
111  fprintf(stderr, "\nBad parameter [%s]\n", badparam);
112  }
113 }
114 
122 bool getOptions(int argc, const char **argv, unsigned int &nb_matrices, unsigned int &nb_iterations,
123  bool &use_plot_file, std::string &plotfile, unsigned int &nbrows, unsigned int &nbcols, bool &verbose)
124 {
125  const char *optarg_;
126  int c;
127  while ((c = vpParseArgv::parse(argc, argv, GETOPTARGS, &optarg_)) > 1) {
128 
129  switch (c) {
130  case 'h':
131  usage(argv[0], NULL);
132  return false;
133  break;
134  case 'n':
135  nb_matrices = (unsigned int)atoi(optarg_);
136  break;
137  case 'i':
138  nb_iterations = (unsigned int)atoi(optarg_);
139  break;
140  case 'f':
141  plotfile = optarg_;
142  use_plot_file = true;
143  break;
144  case 'p':
145  use_plot_file = true;
146  break;
147  case 'R':
148  nbrows = (unsigned int)atoi(optarg_);
149  break;
150  case 'C':
151  nbcols = (unsigned int)atoi(optarg_);
152  break;
153  case 'v':
154  verbose = true;
155  break;
156  // add default options -c -d
157  case 'c':
158  break;
159  case 'd':
160  break;
161  default:
162  usage(argv[0], optarg_);
163  return false;
164  break;
165  }
166  }
167 
168  if ((c == 1) || (c == -1)) {
169  // standalone param or error
170  usage(argv[0], NULL);
171  std::cerr << "ERROR: " << std::endl;
172  std::cerr << " Bad argument " << optarg_ << std::endl << std::endl;
173  return false;
174  }
175 
176  return true;
177 }
178 
179 vpMatrix make_random_matrix(unsigned int nbrows, unsigned int nbcols)
180 {
181  vpMatrix A;
182  A.resize(nbrows, nbcols);
183 
184  for (unsigned int i = 0; i < A.getRows(); i++)
185  for (unsigned int j = 0; j < A.getCols(); j++)
186  A[i][j] = (double)rand() / (double)RAND_MAX;
187  return A;
188 }
189 
190 vpMatrix make_random_symmetric_positive_matrix(unsigned int n)
191 {
192  vpMatrix A;
193  A.resize(n, n);
194  vpMatrix I;
195  I.eye(n);
196 
197  for (unsigned int i = 0; i < A.getRows(); i++)
198  for (unsigned int j = 0; j < A.getCols(); j++)
199  A[i][j] = (double)rand() / (double)RAND_MAX;
200 
201  A = 0.5 * (A + A.t());
202  A = A + n * I;
203  return A;
204 }
205 
206 void create_bench_random_matrix(unsigned int nb_matrices, unsigned int nb_rows, unsigned int nb_cols, bool verbose,
207  std::vector<vpMatrix> &bench)
208 {
209  if (verbose)
210  std::cout << "Create a bench of " << nb_matrices << " " << nb_rows << " by " << nb_cols << " matrices" << std::endl;
211  bench.clear();
212  for (unsigned int i = 0; i < nb_matrices; i++) {
213  vpMatrix M;
214 #if defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_LAPACK) || (VISP_HAVE_OPENCV_VERSION >= 0x020101) || \
215  defined(VISP_HAVE_GSL)
216  double det = 0.;
217  // don't put singular matrices in the benchmark
218  for (M = make_random_matrix(nb_rows, nb_cols); std::fabs(det = M.AtA().det()) < .01;
219  M = make_random_matrix(nb_rows, nb_cols)) {
220  if (verbose) {
221  std::cout << " Generated random matrix AtA=" << std::endl << M.AtA() << std::endl;
222  std::cout << " Generated random matrix not invertible: det=" << det << ". Retrying..." << std::endl;
223  }
224  }
225 #else
226  M = make_random_matrix(nb_rows, nb_cols);
227 #endif
228  bench.push_back(M);
229  }
230 }
231 
232 void create_bench_symmetric_positive_matrix(unsigned int nb_matrices, unsigned int n, bool verbose,
233  std::vector<vpMatrix> &bench)
234 {
235  if (verbose)
236  std::cout << "Create a bench of " << nb_matrices << " " << n << " by " << n << " symmetric positive matrices"
237  << std::endl;
238  bench.clear();
239  for (unsigned int i = 0; i < nb_matrices; i++) {
240  vpMatrix M;
241 #if defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_LAPACK) || (VISP_HAVE_OPENCV_VERSION >= 0x020101) || \
242  defined(VISP_HAVE_GSL)
243  double det = 0.;
244  // don't put singular matrices in the benchmark
245  for (M = make_random_symmetric_positive_matrix(n); std::fabs(det = M.det()) < .01;
246  M = make_random_symmetric_positive_matrix(n)) {
247  if (verbose) {
248  std::cout << " Generated random symmetric positive matrix A=" << std::endl << M << std::endl;
249  std::cout << " Generated random symmetric positive matrix not "
250  "invertibleL: det="
251  << det << ". Retrying..." << std::endl;
252  }
253  }
254 #else
255  M = make_random_symmetric_positive_matrix(n);
256 #endif
257  bench.push_back(M);
258  }
259 }
260 
261 int test_inverse(const std::vector<vpMatrix> &bench, const std::vector<vpMatrix> &result)
262 {
263  for (unsigned int i = 0; i < bench.size(); i++) {
264  vpMatrix I = bench[i] * result[i];
265  if (std::fabs(I.euclideanNorm() - sqrt((double)bench[0].AtA().getRows())) > 1e-10) {
266  std::cout << "Bad inverse[" << i << "]: " << I.euclideanNorm() << " " << sqrt((double)bench[0].AtA().getRows())
267  << std::endl;
268  return EXIT_FAILURE;
269  }
270  }
271  return EXIT_SUCCESS;
272 }
273 
274 #if defined(VISP_HAVE_EIGEN3)
275 int test_inverse_lu_eigen3(bool verbose, const std::vector<vpMatrix> &bench, double &time)
276 {
277  if (verbose)
278  std::cout << "Test inverse by LU using Eigen3 3rd party" << std::endl;
279  // Compute inverse
280  if (verbose)
281  std::cout << " Inverting " << bench[0].AtA().getRows() << "x" << bench[0].AtA().getCols()
282  << " matrix using LU decomposition (Eigen3)." << std::endl;
283  std::vector<vpMatrix> result(bench.size());
284  double t = vpTime::measureTimeMs();
285  for (unsigned int i = 0; i < bench.size(); i++) {
286  result[i] = bench[i].AtA().inverseByLUEigen3() * bench[i].transpose();
287  }
288  time = vpTime::measureTimeMs() - t;
289 
290  // Test inverse
291  return test_inverse(bench, result);
292 }
293 #endif
294 
295 #if defined(VISP_HAVE_LAPACK)
296 int test_inverse_lu_lapack(bool verbose, const std::vector<vpMatrix> &bench, double &time)
297 {
298  if (verbose)
299  std::cout << "Test inverse by LU using Lapack 3rd party" << std::endl;
300  // Compute inverse
301  if (verbose)
302  std::cout << " Inverting " << bench[0].AtA().getRows() << "x" << bench[0].AtA().getCols()
303  << " matrix using LU decomposition (Lapack)." << std::endl;
304  std::vector<vpMatrix> result(bench.size());
305  double t = vpTime::measureTimeMs();
306  for (unsigned int i = 0; i < bench.size(); i++) {
307  result[i] = bench[i].AtA().inverseByLULapack() * bench[i].transpose();
308  }
309  time = vpTime::measureTimeMs() - t;
310 
311  // Test inverse
312  return test_inverse(bench, result);
313 }
314 
315 int test_inverse_cholesky_lapack(bool verbose, const std::vector<vpMatrix> &bench, double &time)
316 {
317  if (verbose)
318  std::cout << "Test inverse by Cholesky using Lapack 3rd party" << std::endl;
319  // Compute inverse
320  if (verbose)
321  std::cout << " Inverting " << bench[0].AtA().getRows() << "x" << bench[0].AtA().getCols()
322  << " matrix using cholesky decomposition (Lapack)." << std::endl;
323  std::vector<vpMatrix> result(bench.size());
324  double t = vpTime::measureTimeMs();
325  for (unsigned int i = 0; i < bench.size(); i++) {
326  result[i] = bench[i].AtA().inverseByCholeskyLapack() * bench[i].transpose();
327  }
328  time = vpTime::measureTimeMs() - t;
329 
330  // Test inverse
331  return test_inverse(bench, result);
332 }
333 
334 int test_inverse_qr_lapack(bool verbose, const std::vector<vpMatrix> &bench, double &time)
335 {
336  if (verbose)
337  std::cout << "Test inverse by QR using Lapack 3rd party" << std::endl;
338  // Compute inverse
339  if (verbose)
340  std::cout << " Inverting " << bench[0].AtA().getRows() << "x" << bench[0].AtA().getCols()
341  << " matrix using QR decomposition (Lapack)" << std::endl;
342  std::vector<vpMatrix> result(bench.size());
343  double t = vpTime::measureTimeMs();
344  for (unsigned int i = 0; i < bench.size(); i++) {
345  result[i] = bench[i].AtA().inverseByQRLapack() * bench[i].transpose();
346  }
347  time = vpTime::measureTimeMs() - t;
348 
349  // Test inverse
350  return test_inverse(bench, result);
351 }
352 #endif
353 
354 #if defined(VISP_HAVE_GSL)
355 int test_inverse_lu_gsl(bool verbose, const std::vector<vpMatrix> &bench, double &time)
356 {
357  if (verbose)
358  std::cout << "Test inverse by LU using GSL 3rd party" << std::endl;
359  // Compute inverse
360  if (verbose)
361  std::cout << " Inverting " << bench[0].AtA().getRows() << "x" << bench[0].AtA().getCols()
362  << " matrix using LU decomposition (GSL)" << std::endl;
363  std::vector<vpMatrix> result(bench.size());
364  double t = vpTime::measureTimeMs();
365  for (unsigned int i = 0; i < bench.size(); i++) {
366  result[i] = bench[i].AtA().inverseByLUGsl() * bench[i].transpose();
367  }
368  time = vpTime::measureTimeMs() - t;
369 
370  // Test inverse
371  return test_inverse(bench, result);
372 }
373 #endif
374 
375 #if (VISP_HAVE_OPENCV_VERSION >= 0x020101)
376 int test_inverse_lu_opencv(bool verbose, const std::vector<vpMatrix> &bench, double &time)
377 {
378  if (verbose)
379  std::cout << "Test inverse by LU using OpenCV 3rd party" << std::endl;
380  // Compute inverse
381  if (verbose)
382  std::cout << " Inverting " << bench[0].AtA().getRows() << "x" << bench[0].AtA().getCols()
383  << " matrix using LU decomposition (OpenCV)" << std::endl;
384  std::vector<vpMatrix> result(bench.size());
385  double t = vpTime::measureTimeMs();
386  for (unsigned int i = 0; i < bench.size(); i++) {
387  result[i] = bench[i].AtA().inverseByLUOpenCV() * bench[i].transpose();
388  }
389  time = vpTime::measureTimeMs() - t;
390 
391  // Test inverse
392  return test_inverse(bench, result);
393 }
394 
395 int test_inverse_cholesky_opencv(bool verbose, const std::vector<vpMatrix> &bench, double &time)
396 {
397  if (verbose)
398  std::cout << "Test inverse by Cholesky using OpenCV 3rd party" << std::endl;
399  // Compute inverse
400  if (verbose)
401  std::cout << " Inverting " << bench[0].AtA().getRows() << "x" << bench[0].AtA().getCols()
402  << " matrix using Cholesky decomposition (OpenCV)" << std::endl;
403  std::vector<vpMatrix> result(bench.size());
404  double t = vpTime::measureTimeMs();
405  for (unsigned int i = 0; i < bench.size(); i++) {
406  result[i] = bench[i].AtA().inverseByCholeskyOpenCV() * bench[i].transpose();
407  }
408  time = vpTime::measureTimeMs() - t;
409 
410  // Test inverse
411  return test_inverse(bench, result);
412 }
413 #endif
414 
415 #if defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_LAPACK) || (VISP_HAVE_OPENCV_VERSION >= 0x020101) || \
416  defined(VISP_HAVE_GSL)
417 // SVD is only available for these 3rd parties
418 int test_pseudo_inverse(bool verbose, const std::vector<vpMatrix> &bench, double &time)
419 {
420  if (verbose)
421  std::cout << "Test pseudo inverse using either Eigen3, Lapack, OpenCV or "
422  "GSL 3rd party"
423  << std::endl;
424  // Compute inverse
425  if (verbose)
426  std::cout << " Pseudo inverting " << bench[0].AtA().getRows() << "x" << bench[0].AtA().getCols() << " matrix"
427  << std::endl;
428  std::vector<vpMatrix> result(bench.size());
429  double t = vpTime::measureTimeMs();
430  for (unsigned int i = 0; i < bench.size(); i++) {
431  result[i] = bench[i].AtA().pseudoInverse() * bench[i].transpose();
432  }
433  time = vpTime::measureTimeMs() - t;
434 
435  // Test inverse
436  return test_inverse(bench, result);
437 }
438 #endif
439 
440 void save_time(const std::string &method, bool verbose, bool use_plot_file, std::ofstream &of, double time)
441 {
442  if (use_plot_file)
443  of << time << "\t";
444  if (verbose || !use_plot_file) {
445  std::cout << method << time << std::endl;
446  }
447 }
448 
449 int main(int argc, const char *argv[])
450 {
451  try {
452 #if defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_LAPACK) || (VISP_HAVE_OPENCV_VERSION >= 0x020101) || \
453  defined(VISP_HAVE_GSL)
454  unsigned int nb_matrices = 1000;
455  unsigned int nb_iterations = 10;
456  unsigned int nb_rows = 6;
457  unsigned int nb_cols = 6;
458  bool verbose = false;
459  std::string plotfile("plot-inv.csv");
460  bool use_plot_file = false;
461  std::ofstream of;
462 
463  // Read the command line options
464  if (getOptions(argc, argv, nb_matrices, nb_iterations, use_plot_file, plotfile, nb_rows, nb_cols, verbose) ==
465  false) {
466  exit(-1);
467  }
468 
469  if (use_plot_file) {
470  of.open(plotfile.c_str());
471  of << "iter"
472  << "\t";
473 
474 #if defined(VISP_HAVE_LAPACK)
475  of << "\"LU Lapack\""
476  << "\t";
477 #endif
478 #if defined(VISP_HAVE_EIGEN3)
479  of << "\"LU Eigen3\""
480  << "\t";
481 #endif
482 #if (VISP_HAVE_OPENCV_VERSION >= 0x020101)
483  of << "\"LU OpenCV\""
484  << "\t";
485 #endif
486 #if defined(VISP_HAVE_GSL)
487  of << "\"LU GSL\""
488  << "\t";
489 #endif
490 
491 #if defined(VISP_HAVE_LAPACK)
492  of << "\"Cholesky Lapack\""
493  << "\t";
494 #endif
495 
496 #if (VISP_HAVE_OPENCV_VERSION >= 0x020101)
497  of << "\"Cholesky OpenCV\""
498  << "\t";
499 #endif
500 
501 #if defined(VISP_HAVE_LAPACK)
502  of << "\"QR Lapack\""
503  << "\t";
504 #endif
505 
506 #if defined(VISP_HAVE_LAPACK) || (VISP_HAVE_OPENCV_VERSION >= 0x020101) || defined(VISP_HAVE_GSL)
507  of << "\"Pseudo inverse (Lapack, OpenCV, GSL)\""
508  << "\t";
509 #endif
510  of << std::endl;
511  }
512 
513  int ret = EXIT_SUCCESS;
514  for (unsigned int iter = 0; iter < nb_iterations; iter++) {
515  std::vector<vpMatrix> bench_random_matrices;
516  create_bench_random_matrix(nb_matrices, nb_rows, nb_cols, verbose, bench_random_matrices);
517  std::vector<vpMatrix> bench_symmetric_positive_matrices;
518  create_bench_symmetric_positive_matrix(nb_matrices, nb_rows, verbose, bench_symmetric_positive_matrices);
519 
520  if (use_plot_file)
521  of << iter << "\t";
522 
523  double time;
524 // LU decomposition
525 #if defined(VISP_HAVE_LAPACK)
526  ret += test_inverse_lu_lapack(verbose, bench_random_matrices, time);
527  save_time("Inverse by LU (Lapack): ", verbose, use_plot_file, of, time);
528 #endif
529 
530 #if defined(VISP_HAVE_EIGEN3)
531  ret += test_inverse_lu_eigen3(verbose, bench_random_matrices, time);
532  save_time("Inverse by LU (Eigen3): ", verbose, use_plot_file, of, time);
533 #endif
534 
535 #if (VISP_HAVE_OPENCV_VERSION >= 0x020101)
536  ret += test_inverse_lu_opencv(verbose, bench_random_matrices, time);
537  save_time("Inverse by LU (OpenCV): ", verbose, use_plot_file, of, time);
538 #endif
539 
540 #if defined(VISP_HAVE_GSL)
541  ret += test_inverse_lu_gsl(verbose, bench_random_matrices, time);
542  save_time("Inverse by LU (GSL): ", verbose, use_plot_file, of, time);
543 #endif
544 
545 // Cholesky for symmetric positive matrices
546 #if defined(VISP_HAVE_LAPACK)
547  ret += test_inverse_cholesky_lapack(verbose, bench_symmetric_positive_matrices, time);
548  save_time("Inverse by Cholesly (Lapack): ", verbose, use_plot_file, of, time);
549 #endif
550 
551 #if (VISP_HAVE_OPENCV_VERSION >= 0x020101)
552  ret += test_inverse_cholesky_opencv(verbose, bench_symmetric_positive_matrices, time);
553  save_time("Inverse by Cholesky (OpenCV): ", verbose, use_plot_file, of, time);
554 #endif
555 
556 // QR decomposition
557 #if defined(VISP_HAVE_LAPACK)
558  ret += test_inverse_qr_lapack(verbose, bench_random_matrices, time);
559  save_time("Inverse by QR (Lapack): ", verbose, use_plot_file, of, time);
560 #endif
561 
562 // Pseudo-inverse with SVD
563 #if defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_LAPACK) || (VISP_HAVE_OPENCV_VERSION >= 0x020101) || \
564  defined(VISP_HAVE_GSL)
565  ret += test_pseudo_inverse(verbose, bench_random_matrices, time);
566  save_time("Pseudo inverse (Lapack, Eigen3, OpenCV or GSL): ", verbose, use_plot_file, of, time);
567 #endif
568 
569  if (use_plot_file)
570  of << std::endl;
571  }
572  if (use_plot_file) {
573  of.close();
574  std::cout << "Result saved in " << plotfile << std::endl;
575  }
576 
577  if (ret == EXIT_SUCCESS) {
578  std::cout << "Test succeed" << std::endl;
579  } else {
580  std::cout << "Test failed" << std::endl;
581  }
582 
583  return ret;
584 #else
585  (void)argc;
586  (void)argv;
587  std::cout << "Test does nothing since you dont't have Eigen3, Lapack, "
588  "OpenCV or GSL 3rd party"
589  << std::endl;
590  return EXIT_SUCCESS;
591 #endif
592  } catch (const vpException &e) {
593  std::cout << "Catch an exception: " << e.getStringMessage() << std::endl;
594  return EXIT_FAILURE;
595  }
596 }
Implementation of a matrix and operations on matrices.
Definition: vpMatrix.h:104
void resize(const unsigned int nrows, const unsigned int ncols, const bool flagNullify=true, const bool recopy_=true)
Definition: vpArray2D.h:171
error that can be emited by ViSP classes.
Definition: vpException.h:71
unsigned int getCols() const
Definition: vpArray2D.h:146
VISP_EXPORT double measureTimeMs()
Definition: vpTime.cpp:88
static bool parse(int *argcPtr, const char **argv, vpArgvInfo *argTable, int flags)
Definition: vpParseArgv.cpp:69
vpMatrix AtA() const
Definition: vpMatrix.cpp:524
unsigned int getRows() const
Definition: vpArray2D.h:156
double det(vpDetMethod method=LU_DECOMPOSITION) const
Definition: vpMatrix.cpp:4906
vpMatrix t() const
Definition: vpMatrix.cpp:375
double euclideanNorm() const
Definition: vpMatrix.cpp:5096
const std::string & getStringMessage(void) const
Send a reference (constant) related the error message (can be empty).
Definition: vpException.cpp:92
void eye()
Definition: vpMatrix.cpp:360