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
66 void usage(
const char *name,
const char *badparam)
69 Test matrix inversions\n\
70 using LU, QR and Cholesky methods as well as Pseudo-inverse.\n\
71 Outputs a comparison of these methods.\n\
74 %s [-n <number of matrices>] [-f <plot filename>]\n\
75 [-R <number of rows>] [-C <number of columns>]\n\
76 [-i <number of iterations>] [-p] [-h]\n",
81 -n <number of matrices> \n\
82 Number of matrices inverted during each test loop.\n\
84 -i <number of iterations> \n\
85 Number of iterations of the test.\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\
93 -R <number of rows>\n\
94 Number of rows of the automatically generated matrices \n\
97 -C <number of columns>\n\
98 Number of colums of the automatically generated matrices \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\
107 Print the help.\n\n");
110 fprintf(stderr,
"ERROR: \n");
111 fprintf(stderr,
"\nBad parameter [%s]\n", badparam);
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)
131 usage(argv[0],
nullptr);
135 nb_matrices = (
unsigned int)atoi(optarg_);
138 nb_iterations = (
unsigned int)atoi(optarg_);
142 use_plot_file =
true;
145 use_plot_file =
true;
148 nbrows = (
unsigned int)atoi(optarg_);
151 nbcols = (
unsigned int)atoi(optarg_);
162 usage(argv[0], optarg_);
168 if ((c == 1) || (c == -1)) {
170 usage(argv[0],
nullptr);
171 std::cerr <<
"ERROR: " << std::endl;
172 std::cerr <<
" Bad argument " << optarg_ << std::endl << std::endl;
179 vpMatrix make_random_matrix(
unsigned int nbrows,
unsigned int nbcols)
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;
190 vpMatrix make_random_symmetric_positive_matrix(
unsigned int n)
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;
201 A = 0.5 * (A + A.
t());
206 vpMatrix make_random_triangular_matrix(
unsigned int nbrows)
211 for (
unsigned int i = 0; i < A.
getRows(); i++) {
212 for (
unsigned int j = i; j < A.
getCols(); j++) {
213 A[i][j] =
static_cast<double>(rand()) /
static_cast<double>(RAND_MAX);
223 void create_bench_random_matrix(
unsigned int nb_matrices,
unsigned int nb_rows,
unsigned int nb_cols,
bool verbose,
224 std::vector<vpMatrix> &bench)
227 std::cout <<
"Create a bench of " << nb_matrices <<
" " << nb_rows <<
" by " << nb_cols <<
" matrices" << std::endl;
229 for (
unsigned int i = 0; i < nb_matrices; i++) {
231 #if defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_LAPACK) || defined(VISP_HAVE_OPENCV)
234 for (M = make_random_matrix(nb_rows, nb_cols); std::fabs(det = M.
AtA().
det()) < .01;
235 M = make_random_matrix(nb_rows, nb_cols)) {
237 std::cout <<
" Generated random matrix AtA=" << std::endl << M.
AtA() << std::endl;
238 std::cout <<
" Generated random matrix not invertible: det=" << det <<
". Retrying..." << std::endl;
242 M = make_random_matrix(nb_rows, nb_cols);
248 void create_bench_symmetric_positive_matrix(
unsigned int nb_matrices,
unsigned int n,
bool verbose,
249 std::vector<vpMatrix> &bench)
252 std::cout <<
"Create a bench of " << nb_matrices <<
" " << n <<
" by " << n <<
" symmetric positive matrices"
255 for (
unsigned int i = 0; i < nb_matrices; i++) {
257 #if defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_LAPACK) || defined(VISP_HAVE_OPENCV)
260 for (M = make_random_symmetric_positive_matrix(n); std::fabs(det = M.
det()) < .01;
261 M = make_random_symmetric_positive_matrix(n)) {
263 std::cout <<
" Generated random symmetric positive matrix A=" << std::endl << M << std::endl;
264 std::cout <<
" Generated random symmetric positive matrix not "
266 << det <<
". Retrying..." << std::endl;
270 M = make_random_symmetric_positive_matrix(n);
276 void create_bench_random_triangular_matrix(
unsigned int nb_matrices,
unsigned int n,
bool verbose,
277 std::vector<vpMatrix> &bench)
280 std::cout <<
"Create a bench of " << nb_matrices <<
" " << n <<
" by " << n <<
" triangular matrices" << std::endl;
282 for (
unsigned int i = 0; i < nb_matrices; i++) {
284 #if defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_LAPACK) || defined(VISP_HAVE_OPENCV)
287 for (M = make_random_triangular_matrix(n); std::fabs(det = M.
det()) < .01; M = make_random_triangular_matrix(n)) {
289 std::cout <<
" Generated random symmetric positive matrix A=" << std::endl << M << std::endl;
290 std::cout <<
" Generated random symmetric positive matrix not "
292 << det <<
". Retrying..." << std::endl;
296 M = make_random_triangular_matrix(n);
302 int test_inverse(
const std::vector<vpMatrix> &bench,
const std::vector<vpMatrix> &result)
304 double epsilon = 1e-10;
305 for (
unsigned int i = 0; i < bench.size(); i++) {
307 if (std::fabs(I.frobeniusNorm() - sqrt(
static_cast<double>(bench[0].AtA().getRows()))) > epsilon) {
308 std::cout <<
"Bad inverse[" << i <<
"]: " << I.frobeniusNorm() <<
" " << sqrt((
double)bench[0].AtA().getRows())
316 int test_inverse_lu_small(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
319 std::cout <<
"Test inverse by LU on small matrices" << std::endl;
322 std::cout <<
" Inverting " << bench[0].getRows() <<
"x" << bench[0].getCols() <<
" small matrix." << std::endl;
323 std::vector<vpMatrix> result(bench.size());
325 for (
unsigned int i = 0; i < bench.size(); i++) {
326 result[i] = bench[i].inverseByLU();
331 return test_inverse(bench, result);
334 #if defined(VISP_HAVE_EIGEN3)
335 int test_inverse_lu_eigen3(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
338 std::cout <<
"Test inverse by LU using Eigen3 3rd party" << std::endl;
341 std::cout <<
" Inverting " << bench[0].AtA().getRows() <<
"x" << bench[0].AtA().getCols()
342 <<
" matrix using LU decomposition (Eigen3)." << std::endl;
343 std::vector<vpMatrix> result(bench.size());
345 for (
unsigned int i = 0; i < bench.size(); i++) {
346 result[i] = bench[i].AtA().inverseByLUEigen3() * bench[i].transpose();
351 return test_inverse(bench, result);
355 #if defined(VISP_HAVE_LAPACK)
356 int test_inverse_lu_lapack(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
359 std::cout <<
"Test inverse by LU using Lapack 3rd party" << std::endl;
362 std::cout <<
" Inverting " << bench[0].AtA().getRows() <<
"x" << bench[0].AtA().getCols()
363 <<
" matrix using LU decomposition (Lapack)." << std::endl;
364 std::vector<vpMatrix> result(bench.size());
366 for (
unsigned int i = 0; i < bench.size(); i++) {
367 result[i] = bench[i].AtA().inverseByLULapack() * bench[i].transpose();
372 return test_inverse(bench, result);
375 int test_inverse_cholesky_lapack(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
378 std::cout <<
"Test inverse by Cholesky using Lapack 3rd party" << std::endl;
381 std::cout <<
" Inverting " << bench[0].AtA().getRows() <<
"x" << bench[0].AtA().getCols()
382 <<
" matrix using cholesky decomposition (Lapack)." << std::endl;
383 std::vector<vpMatrix> result(bench.size());
385 for (
unsigned int i = 0; i < bench.size(); i++) {
386 result[i] = bench[i].AtA().inverseByCholeskyLapack() * bench[i].transpose();
391 return test_inverse(bench, result);
394 int test_inverse_qr_lapack(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
397 std::cout <<
"Test inverse by QR using Lapack 3rd party" << std::endl;
400 std::cout <<
" Inverting " << bench[0].AtA().getRows() <<
"x" << bench[0].AtA().getCols()
401 <<
" matrix using QR decomposition (Lapack)" << std::endl;
402 std::vector<vpMatrix> result(bench.size());
404 for (
unsigned int i = 0; i < bench.size(); i++) {
405 result[i] = bench[i].AtA().inverseByQRLapack() * bench[i].transpose();
410 return test_inverse(bench, result);
414 #if defined(VISP_HAVE_OPENCV)
415 int test_inverse_lu_opencv(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
418 std::cout <<
"Test inverse by LU using OpenCV 3rd party" << std::endl;
421 std::cout <<
" Inverting " << bench[0].AtA().getRows() <<
"x" << bench[0].AtA().getCols()
422 <<
" matrix using LU decomposition (OpenCV)" << std::endl;
423 std::vector<vpMatrix> result(bench.size());
425 for (
unsigned int i = 0; i < bench.size(); i++) {
426 result[i] = bench[i].AtA().inverseByLUOpenCV() * bench[i].transpose();
431 return test_inverse(bench, result);
434 int test_inverse_cholesky_opencv(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
437 std::cout <<
"Test inverse by Cholesky using OpenCV 3rd party" << std::endl;
440 std::cout <<
" Inverting " << bench[0].AtA().getRows() <<
"x" << bench[0].AtA().getCols()
441 <<
" matrix using Cholesky decomposition (OpenCV)" << std::endl;
442 std::vector<vpMatrix> result(bench.size());
444 for (
unsigned int i = 0; i < bench.size(); i++) {
445 result[i] = bench[i].AtA().inverseByCholeskyOpenCV() * bench[i].transpose();
450 return test_inverse(bench, result);
454 #if defined(VISP_HAVE_LAPACK) || defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_OPENCV)
456 int test_pseudo_inverse(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
459 std::cout <<
"Test pseudo inverse using either Lapack, Eigen3 or OpenCV 3rd party" << std::endl;
462 std::cout <<
" Pseudo inverting " << bench[0].AtA().getRows() <<
"x" << bench[0].AtA().getCols() <<
" matrix"
464 std::vector<vpMatrix> result(bench.size());
466 for (
unsigned int i = 0; i < bench.size(); i++) {
467 result[i] = bench[i].AtA().pseudoInverse() * bench[i].transpose();
472 return test_inverse(bench, result);
475 int test_inverse_triangular(
bool verbose,
const std::vector<vpMatrix> &bench,
double &time)
478 std::cout <<
"Test inverse triangular using Lapack" << std::endl;
481 std::cout <<
" Triangular inverse " << bench[0].getRows() <<
"x" << bench[0].getCols() <<
" matrix" << std::endl;
482 std::vector<vpMatrix> result(bench.size());
484 for (
unsigned int i = 0; i < bench.size(); i++) {
485 result[i] = bench[i].inverseTriangular(
true);
490 return test_inverse(bench, result);
494 void save_time(
const std::string &method,
bool verbose,
bool use_plot_file, std::ofstream &of,
double time)
498 if (verbose || !use_plot_file) {
499 std::cout << method << time << std::endl;
503 int main(
int argc,
const char *argv[])
506 unsigned int nb_matrices = 1000;
507 unsigned int nb_iterations = 10;
508 unsigned int nb_rows = 6;
509 unsigned int nb_cols = 6;
510 bool verbose =
false;
511 std::string plotfile(
"plot-inv.csv");
512 bool use_plot_file =
false;
516 if (getOptions(argc, argv, nb_matrices, nb_iterations, use_plot_file, plotfile, nb_rows, nb_cols, verbose) ==
522 of.open(plotfile.c_str());
526 #if defined(VISP_HAVE_LAPACK)
527 of <<
"\"LU Lapack\""
530 #if defined(VISP_HAVE_EIGEN3)
531 of <<
"\"LU Eigen3\""
534 #if defined(VISP_HAVE_OPENCV)
535 of <<
"\"LU OpenCV\""
539 #if defined(VISP_HAVE_LAPACK)
540 of <<
"\"Cholesky Lapack\""
544 #if defined(VISP_HAVE_OPENCV)
545 of <<
"\"Cholesky OpenCV\""
549 #if defined(VISP_HAVE_LAPACK)
550 of <<
"\"QR Lapack\""
554 #if defined(VISP_HAVE_LAPACK) || defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_OPENCV)
555 of <<
"\"Pseudo inverse (Lapack, Eigen3 or OpenCV)\""
561 int ret = EXIT_SUCCESS;
562 for (
unsigned int iter = 0; iter < nb_iterations; iter++) {
563 std::vector<vpMatrix> bench_random_matrices_11;
564 create_bench_random_matrix(nb_matrices, 1, 1, verbose, bench_random_matrices_11);
565 std::vector<vpMatrix> bench_random_matrices_22;
566 create_bench_random_matrix(nb_matrices, 2, 2, verbose, bench_random_matrices_22);
567 std::vector<vpMatrix> bench_random_matrices_33;
568 create_bench_random_matrix(nb_matrices, 3, 3, verbose, bench_random_matrices_33);
569 #if defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_LAPACK) || defined(VISP_HAVE_OPENCV)
570 std::vector<vpMatrix> bench_random_matrices;
571 create_bench_random_matrix(nb_matrices, nb_rows, nb_cols, verbose, bench_random_matrices);
572 std::vector<vpMatrix> bench_symmetric_positive_matrices;
573 create_bench_symmetric_positive_matrix(nb_matrices, nb_rows, verbose, bench_symmetric_positive_matrices);
574 std::vector<vpMatrix> bench_triangular_matrices;
575 create_bench_random_triangular_matrix(nb_matrices, nb_rows, verbose, bench_triangular_matrices);
584 ret += test_inverse_lu_small(verbose, bench_random_matrices_11, time);
585 save_time(
"Inverse by LU 1x1: ", verbose, use_plot_file, of, time);
587 ret += test_inverse_lu_small(verbose, bench_random_matrices_22, time);
588 save_time(
"Inverse by LU 2x2: ", verbose, use_plot_file, of, time);
590 ret += test_inverse_lu_small(verbose, bench_random_matrices_33, time);
591 save_time(
"Inverse by LU 3x3: ", verbose, use_plot_file, of, time);
594 #if defined(VISP_HAVE_LAPACK)
595 ret += test_inverse_lu_lapack(verbose, bench_random_matrices, time);
596 save_time(
"Inverse by LU (Lapack): ", verbose, use_plot_file, of, time);
599 #if defined(VISP_HAVE_EIGEN3)
600 ret += test_inverse_lu_eigen3(verbose, bench_random_matrices, time);
601 save_time(
"Inverse by LU (Eigen3): ", verbose, use_plot_file, of, time);
604 #if defined(VISP_HAVE_OPENCV)
605 ret += test_inverse_lu_opencv(verbose, bench_random_matrices, time);
606 save_time(
"Inverse by LU (OpenCV): ", verbose, use_plot_file, of, time);
610 #if defined(VISP_HAVE_LAPACK)
611 ret += test_inverse_cholesky_lapack(verbose, bench_symmetric_positive_matrices, time);
612 save_time(
"Inverse by Cholesly (Lapack): ", verbose, use_plot_file, of, time);
615 #if defined(VISP_HAVE_OPENCV)
616 ret += test_inverse_cholesky_opencv(verbose, bench_symmetric_positive_matrices, time);
617 save_time(
"Inverse by Cholesky (OpenCV): ", verbose, use_plot_file, of, time);
621 #if defined(VISP_HAVE_LAPACK)
622 ret += test_inverse_qr_lapack(verbose, bench_random_matrices, time);
623 save_time(
"Inverse by QR (Lapack): ", verbose, use_plot_file, of, time);
627 #if defined(VISP_HAVE_LAPACK) || defined(VISP_HAVE_EIGEN3) || defined(VISP_HAVE_OPENCV)
628 ret += test_pseudo_inverse(verbose, bench_random_matrices, time);
629 save_time(
"Pseudo inverse (Lapack, Eigen3, OpenCV): ", verbose, use_plot_file, of, time);
633 #if defined(VISP_HAVE_LAPACK)
634 ret += test_inverse_triangular(verbose, bench_triangular_matrices, time);
635 save_time(
"Triangular inverse (Lapack): ", verbose, use_plot_file, of, time);
643 std::cout <<
"Result saved in " << plotfile << std::endl;
646 if (ret == EXIT_SUCCESS) {
647 std::cout <<
"Test succeed" << std::endl;
650 std::cout <<
"Test failed" << 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.
const std::string & getStringMessage() const
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()