46#include <visp3/core/vpColVector.h>
47#include <visp3/core/vpMatrix.h>
48#include <visp3/core/vpTime.h>
49#include <visp3/io/vpParseArgv.h>
52#define GETOPTARGS "cdn:i:pf:R:C:vh"
54#ifdef ENABLE_VISP_NAMESPACE
58void usage(
const char *name,
const char *badparam);
59bool getOptions(
int argc,
const char **argv,
unsigned int &nb_matrices,
unsigned int &nb_iterations,
60 bool &use_plot_file, std::string &plotfile,
unsigned int &nbrows,
unsigned int &nbcols,
bool &verbose);
62vpMatrix make_random_matrix(
unsigned int nbrows,
unsigned int nbcols);
63vpMatrix make_random_symmetric_positive_matrix(
unsigned int n);
64vpMatrix make_random_triangular_matrix(
unsigned int nbrows);
65void create_bench_random_matrix(
unsigned int nb_matrices,
unsigned int nb_rows,
unsigned int nb_cols,
bool verbose,
66 std::vector<vpMatrix> &bench);
67void create_bench_symmetric_positive_matrix(
unsigned int nb_matrices,
unsigned int n,
bool verbose,
68 std::vector<vpMatrix> &bench);
69void create_bench_random_triangular_matrix(
unsigned int nb_matrices,
unsigned int n,
bool verbose,
70 std::vector<vpMatrix> &bench);
71int test_inverse(
const std::vector<vpMatrix> &bench,
const std::vector<vpMatrix> &result);
72int test_inverse_lu_small(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time);
73void save_time(
const std::string &method,
bool verbose,
bool use_plot_file, std::ofstream &of,
double time);
74#if defined(VISP_HAVE_EIGEN3)
75int test_inverse_lu_eigen3(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time);
78#if defined(VISP_HAVE_LAPACK)
79int test_inverse_lu_lapack(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time);
80int test_inverse_cholesky_lapack(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time);
81int test_inverse_qr_lapack(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time);
83#if defined(VISP_HAVE_OPENCV)
84int test_inverse_lu_opencv(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time);
85int test_inverse_cholesky_opencv(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time);
87#if defined(VISP_HAVE_LAPACK) || defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_OPENCV)
88int test_pseudo_inverse(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time);
89int test_inverse_triangular(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time);
98void usage(
const char *name,
const char *badparam)
101Test matrix inversions\n\
102using LU, QR and Cholesky methods as well as Pseudo-inverse.\n\
103Outputs a comparison of these methods.\n\
106 %s [-n <number of matrices>] [-f <plot filename>]\n\
107 [-R <number of rows>] [-C <number of columns>]\n\
108 [-i <number of iterations>] [-p] [-h]\n",
113 -n <number of matrices> \n\
114 Number of matrices inverted during each test loop.\n\
116 -i <number of iterations> \n\
117 Number of iterations of the test.\n\
119 -f <plot filename> \n\
120 Set output path for plot output.\n\
121 The plot logs the times of \n\
122 the different inversion methods: \n\
123 QR,LU,Cholesky and Pseudo-inverse.\n\
125 -R <number of rows>\n\
126 Number of rows of the automatically generated matrices \n\
129 -C <number of columns>\n\
130 Number of colums of the automatically generated matrices \n\
134 Plot into filename in the gnuplot format. \n\
135 If this option is used, tests results will be logged \n\
136 into a filename specified with -f.\n\
139 Print the help.\n\n");
142 fprintf(stderr,
"ERROR: \n");
143 fprintf(stderr,
"\nBad parameter [%s]\n", badparam);
154bool getOptions(
int argc,
const char **argv,
unsigned int &nb_matrices,
unsigned int &nb_iterations,
155 bool &use_plot_file, std::string &plotfile,
unsigned int &nbrows,
unsigned int &nbcols,
bool &verbose)
163 usage(argv[0],
nullptr);
166 nb_matrices =
static_cast<unsigned int>(atoi(optarg_));
169 nb_iterations =
static_cast<unsigned int>(atoi(optarg_));
173 use_plot_file =
true;
176 use_plot_file =
true;
179 nbrows =
static_cast<unsigned int>(atoi(optarg_));
182 nbcols =
static_cast<unsigned int>(atoi(optarg_));
193 usage(argv[0], optarg_);
198 if ((c == 1) || (c == -1)) {
200 usage(argv[0],
nullptr);
201 std::cerr <<
"ERROR: " << std::endl;
202 std::cerr <<
" Bad argument " << optarg_ << std::endl << std::endl;
209vpMatrix make_random_matrix(
unsigned int nbrows,
unsigned int nbcols)
214 for (
unsigned int i = 0;
i < A.
getRows();
i++)
215 for (
unsigned int j = 0;
j < A.
getCols();
j++)
216 A[i][j] =
static_cast<double>(rand()) /
static_cast<double>(RAND_MAX);
220vpMatrix make_random_symmetric_positive_matrix(
unsigned int n)
227 for (
unsigned int i = 0;
i < A.
getRows();
i++)
228 for (
unsigned int j = 0;
j < A.
getCols();
j++)
229 A[i][j] =
static_cast<double>(rand()) /
static_cast<double>(RAND_MAX);
231 A = 0.5 * (A + A.
t());
236vpMatrix make_random_triangular_matrix(
unsigned int nbrows)
241 for (
unsigned int i = 0;
i < A.
getRows();
i++) {
242 for (
unsigned int j = i;
j < A.
getCols();
j++) {
243 A[
i][
j] =
static_cast<double>(rand()) /
static_cast<double>(RAND_MAX);
253void create_bench_random_matrix(
unsigned int nb_matrices,
unsigned int nb_rows,
unsigned int nb_cols,
bool verbose,
254 std::vector<vpMatrix> &bench)
257 std::cout <<
"Create a bench of " << nb_matrices <<
" " << nb_rows <<
" by " << nb_cols <<
" matrices" << std::endl;
259 for (
unsigned int i = 0;
i < nb_matrices;
i++) {
261#if defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_LAPACK) || defined(VISP_HAVE_OPENCV)
264 for (M = make_random_matrix(nb_rows, nb_cols); std::fabs(det = M.
AtA().
det()) < .01;
265 M = make_random_matrix(nb_rows, nb_cols)) {
267 std::cout <<
" Generated random matrix AtA=" << std::endl << M.
AtA() << std::endl;
268 std::cout <<
" Generated random matrix not invertible: det=" << det <<
". Retrying..." << std::endl;
272 M = make_random_matrix(nb_rows, nb_cols);
278void create_bench_symmetric_positive_matrix(
unsigned int nb_matrices,
unsigned int n,
bool verbose,
279 std::vector<vpMatrix> &bench)
282 std::cout <<
"Create a bench of " << nb_matrices <<
" " << n <<
" by " << n <<
" symmetric positive matrices"
285 for (
unsigned int i = 0;
i < nb_matrices;
i++) {
287#if defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_LAPACK) || defined(VISP_HAVE_OPENCV)
290 for (M = make_random_symmetric_positive_matrix(n); std::fabs(det = M.
det()) < .01;
291 M = make_random_symmetric_positive_matrix(n)) {
293 std::cout <<
" Generated random symmetric positive matrix A=" << std::endl << M << std::endl;
294 std::cout <<
" Generated random symmetric positive matrix not "
296 << det <<
". Retrying..." << std::endl;
300 M = make_random_symmetric_positive_matrix(n);
306void create_bench_random_triangular_matrix(
unsigned int nb_matrices,
unsigned int n,
bool verbose,
307 std::vector<vpMatrix> &bench)
310 std::cout <<
"Create a bench of " << nb_matrices <<
" " << n <<
" by " << n <<
" triangular matrices" << std::endl;
312 for (
unsigned int i = 0;
i < nb_matrices;
i++) {
314#if defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_LAPACK) || defined(VISP_HAVE_OPENCV)
317 for (M = make_random_triangular_matrix(n); std::fabs(det = M.
det()) < .01; M = make_random_triangular_matrix(n)) {
319 std::cout <<
" Generated random symmetric positive matrix A=" << std::endl << M << std::endl;
320 std::cout <<
" Generated random symmetric positive matrix not "
322 << det <<
". Retrying..." << std::endl;
326 M = make_random_triangular_matrix(n);
332int test_inverse(
const std::vector<vpMatrix> &bench,
const std::vector<vpMatrix> &result)
334 double epsilon = 1
e-10;
335 for (
unsigned int i = 0;
i < bench.size();
i++) {
337 if (std::fabs(I.frobeniusNorm() - sqrt(
static_cast<double>(bench[0].AtA().getRows()))) > epsilon) {
338 std::cout <<
"Bad inverse[" <<
i <<
"]: " << I.frobeniusNorm() <<
" " << sqrt(
static_cast<double>(bench[0].AtA().getRows()))
346int test_inverse_lu_small(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
349 std::cout <<
"Test inverse by LU on small matrices" << std::endl;
352 std::cout <<
" Inverting " << bench[0].getRows() <<
"x" << bench[0].getCols() <<
" small matrix." << std::endl;
353 std::vector<vpMatrix> result(bench.size());
355 for (
unsigned int i = 0;
i < bench.size();
i++) {
356 result[
i] = bench[
i].inverseByLU();
361 return test_inverse(bench, result);
364#if defined(VISP_HAVE_EIGEN3)
365int test_inverse_lu_eigen3(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
368 std::cout <<
"Test inverse by LU using Eigen3 3rd party" << std::endl;
371 std::cout <<
" Inverting " << bench[0].AtA().getRows() <<
"x" << bench[0].AtA().getCols()
372 <<
" matrix using LU decomposition (Eigen3)." << std::endl;
373 std::vector<vpMatrix> result(bench.size());
375 for (
unsigned int i = 0;
i < bench.size();
i++) {
376 result[
i] = bench[
i].AtA().inverseByLUEigen3() * bench[
i].transpose();
381 return test_inverse(bench, result);
385#if defined(VISP_HAVE_LAPACK)
386int test_inverse_lu_lapack(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
389 std::cout <<
"Test inverse by LU using Lapack 3rd party" << std::endl;
392 std::cout <<
" Inverting " << bench[0].AtA().getRows() <<
"x" << bench[0].AtA().getCols()
393 <<
" matrix using LU decomposition (Lapack)." << std::endl;
394 std::vector<vpMatrix> result(bench.size());
396 for (
unsigned int i = 0;
i < bench.size();
i++) {
397 result[
i] = bench[
i].AtA().inverseByLULapack() * bench[
i].transpose();
402 return test_inverse(bench, result);
405int test_inverse_cholesky_lapack(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
408 std::cout <<
"Test inverse by Cholesky using Lapack 3rd party" << std::endl;
411 std::cout <<
" Inverting " << bench[0].AtA().getRows() <<
"x" << bench[0].AtA().getCols()
412 <<
" matrix using cholesky decomposition (Lapack)." << std::endl;
413 std::vector<vpMatrix> result(bench.size());
415 for (
unsigned int i = 0;
i < bench.size();
i++) {
416 result[
i] = bench[
i].AtA().inverseByCholeskyLapack() * bench[
i].transpose();
421 return test_inverse(bench, result);
424int test_inverse_qr_lapack(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
427 std::cout <<
"Test inverse by QR using Lapack 3rd party" << std::endl;
430 std::cout <<
" Inverting " << bench[0].AtA().getRows() <<
"x" << bench[0].AtA().getCols()
431 <<
" matrix using QR decomposition (Lapack)" << std::endl;
432 std::vector<vpMatrix> result(bench.size());
434 for (
unsigned int i = 0;
i < bench.size();
i++) {
435 result[
i] = bench[
i].AtA().inverseByQRLapack() * bench[
i].transpose();
440 return test_inverse(bench, result);
444#if defined(VISP_HAVE_OPENCV)
445int test_inverse_lu_opencv(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
448 std::cout <<
"Test inverse by LU using OpenCV 3rd party" << std::endl;
451 std::cout <<
" Inverting " << bench[0].AtA().getRows() <<
"x" << bench[0].AtA().getCols()
452 <<
" matrix using LU decomposition (OpenCV)" << std::endl;
453 std::vector<vpMatrix> result(bench.size());
455 for (
unsigned int i = 0;
i < bench.size();
i++) {
456 result[
i] = bench[
i].AtA().inverseByLUOpenCV() * bench[
i].transpose();
461 return test_inverse(bench, result);
464int test_inverse_cholesky_opencv(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
467 std::cout <<
"Test inverse by Cholesky using OpenCV 3rd party" << std::endl;
470 std::cout <<
" Inverting " << bench[0].AtA().getRows() <<
"x" << bench[0].AtA().getCols()
471 <<
" matrix using Cholesky decomposition (OpenCV)" << std::endl;
472 std::vector<vpMatrix> result(bench.size());
474 for (
unsigned int i = 0;
i < bench.size();
i++) {
475 result[
i] = bench[
i].AtA().inverseByCholeskyOpenCV() * bench[
i].transpose();
480 return test_inverse(bench, result);
484#if defined(VISP_HAVE_LAPACK) || defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_OPENCV)
486int test_pseudo_inverse(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
489 std::cout <<
"Test pseudo inverse using either Lapack, Eigen3 or OpenCV 3rd party" << std::endl;
492 std::cout <<
" Pseudo inverting " << bench[0].AtA().getRows() <<
"x" << bench[0].AtA().getCols() <<
" matrix"
494 std::vector<vpMatrix> result(bench.size());
496 for (
unsigned int i = 0;
i < bench.size();
i++) {
497 result[
i] = bench[
i].AtA().pseudoInverse() * bench[
i].transpose();
502 return test_inverse(bench, result);
505int test_inverse_triangular(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
508 std::cout <<
"Test inverse triangular using Lapack" << std::endl;
511 std::cout <<
" Triangular inverse " << bench[0].getRows() <<
"x" << bench[0].getCols() <<
" matrix" << std::endl;
512 std::vector<vpMatrix> result(bench.size());
514 for (
unsigned int i = 0;
i < bench.size();
i++) {
515 result[
i] = bench[
i].inverseTriangular(
true);
520 return test_inverse(bench, result);
524void save_time(
const std::string &method,
bool verbose,
bool use_plot_file, std::ofstream &of,
double time)
528 if (verbose || !use_plot_file) {
529 std::cout << method << time << std::endl;
533int main(
int argc,
const char *argv[])
536 unsigned int nb_matrices = 1000;
537 unsigned int nb_iterations = 10;
538 unsigned int nb_rows = 6;
539 unsigned int nb_cols = 6;
540 bool verbose =
false;
541 std::string plotfile(
"plot-inv.csv");
542 bool use_plot_file =
false;
546 if (getOptions(argc, argv, nb_matrices, nb_iterations, use_plot_file, plotfile, nb_rows, nb_cols, verbose) ==
552 of.open(plotfile.c_str());
556#if defined(VISP_HAVE_LAPACK)
557 of <<
"\"LU Lapack\""
560#if defined(VISP_HAVE_EIGEN3)
561 of <<
"\"LU Eigen3\""
564#if defined(VISP_HAVE_OPENCV)
565 of <<
"\"LU OpenCV\""
569#if defined(VISP_HAVE_LAPACK)
570 of <<
"\"Cholesky Lapack\""
574#if defined(VISP_HAVE_OPENCV)
575 of <<
"\"Cholesky OpenCV\""
579#if defined(VISP_HAVE_LAPACK)
580 of <<
"\"QR Lapack\""
584#if defined(VISP_HAVE_LAPACK) || defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_OPENCV)
585 of <<
"\"Pseudo inverse (Lapack, Eigen3 or OpenCV)\""
591 int ret = EXIT_SUCCESS;
592 for (
unsigned int iter = 0;
iter < nb_iterations;
iter++) {
593 std::vector<vpMatrix> bench_random_matrices_11;
594 create_bench_random_matrix(nb_matrices, 1, 1, verbose, bench_random_matrices_11);
595 std::vector<vpMatrix> bench_random_matrices_22;
596 create_bench_random_matrix(nb_matrices, 2, 2, verbose, bench_random_matrices_22);
597 std::vector<vpMatrix> bench_random_matrices_33;
598 create_bench_random_matrix(nb_matrices, 3, 3, verbose, bench_random_matrices_33);
599#if defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_LAPACK) || defined(VISP_HAVE_OPENCV)
600 std::vector<vpMatrix> bench_random_matrices;
601 create_bench_random_matrix(nb_matrices, nb_rows, nb_cols, verbose, bench_random_matrices);
602 std::vector<vpMatrix> bench_symmetric_positive_matrices;
603 create_bench_symmetric_positive_matrix(nb_matrices, nb_rows, verbose, bench_symmetric_positive_matrices);
604 std::vector<vpMatrix> bench_triangular_matrices;
605 create_bench_random_triangular_matrix(nb_matrices, nb_rows, verbose, bench_triangular_matrices);
614 ret += test_inverse_lu_small(verbose, bench_random_matrices_11, time);
615 save_time(
"Inverse by LU 1x1: ", verbose, use_plot_file, of, time);
617 ret += test_inverse_lu_small(verbose, bench_random_matrices_22, time);
618 save_time(
"Inverse by LU 2x2: ", verbose, use_plot_file, of, time);
620 ret += test_inverse_lu_small(verbose, bench_random_matrices_33, time);
621 save_time(
"Inverse by LU 3x3: ", verbose, use_plot_file, of, time);
624#if defined(VISP_HAVE_LAPACK)
625 ret += test_inverse_lu_lapack(verbose, bench_random_matrices, time);
626 save_time(
"Inverse by LU (Lapack): ", verbose, use_plot_file, of, time);
629#if defined(VISP_HAVE_EIGEN3)
630 ret += test_inverse_lu_eigen3(verbose, bench_random_matrices, time);
631 save_time(
"Inverse by LU (Eigen3): ", verbose, use_plot_file, of, time);
634#if defined(VISP_HAVE_OPENCV)
635 ret += test_inverse_lu_opencv(verbose, bench_random_matrices, time);
636 save_time(
"Inverse by LU (OpenCV): ", verbose, use_plot_file, of, time);
640#if defined(VISP_HAVE_LAPACK)
641 ret += test_inverse_cholesky_lapack(verbose, bench_symmetric_positive_matrices, time);
642 save_time(
"Inverse by Cholesly (Lapack): ", verbose, use_plot_file, of, time);
645#if defined(VISP_HAVE_OPENCV)
646 ret += test_inverse_cholesky_opencv(verbose, bench_symmetric_positive_matrices, time);
647 save_time(
"Inverse by Cholesky (OpenCV): ", verbose, use_plot_file, of, time);
651#if defined(VISP_HAVE_LAPACK)
652 ret += test_inverse_qr_lapack(verbose, bench_random_matrices, time);
653 save_time(
"Inverse by QR (Lapack): ", verbose, use_plot_file, of, time);
657#if defined(VISP_HAVE_LAPACK) || defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_OPENCV)
658 ret += test_pseudo_inverse(verbose, bench_random_matrices, time);
659 save_time(
"Pseudo inverse (Lapack, Eigen3, OpenCV): ", verbose, use_plot_file, of, time);
663#if defined(VISP_HAVE_LAPACK)
664 ret += test_inverse_triangular(verbose, bench_triangular_matrices, time);
665 save_time(
"Triangular inverse (Lapack): ", verbose, use_plot_file, of, time);
673 std::cout <<
"Result saved in " << plotfile << std::endl;
676 if (ret == EXIT_SUCCESS) {
677 std::cout <<
"Test succeed" << std::endl;
680 std::cout <<
"Test failed" << std::endl;
686 std::cout <<
"Catch an exception: " <<
e.getStringMessage() << std::endl;
unsigned int getCols() const
void resize(unsigned int nrows, unsigned int ncols, bool flagNullify=true, bool recopy_=true)
unsigned int getRows() const
error that can be emitted by ViSP classes.
Implementation of a matrix and operations on matrices.
double det(vpDetMethod method=LU_DECOMPOSITION) const
static bool parse(int *argcPtr, const char **argv, vpArgvInfo *argTable, int flags)
VISP_EXPORT double measureTimeMs()