Visual Servoing Platform version 3.7.0
Loading...
Searching...
No Matches
catchMatrixCholesky.cpp
1/*
2 * ViSP, open source Visual Servoing Platform software.
3 * Copyright (C) 2005 - 2024 by Inria. All rights reserved.
4 *
5 * This software is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or
8 * (at your option) any later version.
9 * See the file LICENSE.txt at the root directory of this source
10 * distribution for additional information about the GNU GPL.
11 *
12 * For using ViSP with software that can not be combined with the GNU
13 * GPL, please contact Inria about acquiring a ViSP Professional
14 * Edition License.
15 *
16 * See https://visp.inria.fr for more information.
17 *
18 * This software was developed at:
19 * Inria Rennes - Bretagne Atlantique
20 * Campus Universitaire de Beaulieu
21 * 35042 Rennes Cedex
22 * France
23 *
24 * If you have questions regarding the use of this file, please contact
25 * Inria at visp@inria.fr
26 *
27 * This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
28 * WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
29 *
30 * Description:
31 * Test matrix condition number.
32 */
33
38
39#include <visp3/core/vpConfig.h>
40
41#if defined(VISP_HAVE_CATCH2)
42#include <iostream>
43#include <stdlib.h>
44#include <visp3/core/vpMatrix.h>
45
46#include <catch_amalgamated.hpp>
47
48#ifdef ENABLE_VISP_NAMESPACE
49using namespace VISP_NAMESPACE_NAME;
50#endif
51
52bool equal(const vpMatrix &A, const vpMatrix &B)
53{
54 unsigned int rowsA = A.getRows(), rowsB = B.getRows();
55 unsigned int colsA = A.getCols(), colsB = B.getCols();
56 bool areEqual = (rowsA == rowsB) && (colsA == colsB);
57
58 if (!areEqual) {
59 return false;
60 }
61
62 for (unsigned int r = 0; r < rowsA; ++r) {
63 for (unsigned int c = 0; c < colsA; ++c) {
64 areEqual = areEqual && vpMath::equal(A[r][c], B[r][c]);
65 if (!areEqual) {
66 return false;
67 }
68 }
69 }
70 return areEqual;
71}
72
73bool testCholeskyDecomposition(const vpMatrix &M, const vpMatrix &gtResult, const vpMatrix &result,
74 const std::string &title, const bool &verbose)
75{
76 if (verbose) {
77 std::cout << "-------------------------" << std::endl;
78 std::cout << "Test: " << title << std::endl;
79 std::cout << std::endl;
80 }
81 unsigned int gtRows = gtResult.getRows(), resRows = result.getRows();
82 unsigned int gtCols = gtResult.getCols(), resCols = result.getCols();
83 bool areEqual = (gtRows == resRows) && (gtCols == resCols);
84 if (!areEqual) {
85 if (verbose) {
86 std::cout << "Failed: dimensions mismatch (" << gtRows << " x " << gtCols << ")";
87 std::cout << " vs (" << resRows << " x " << resCols << ")" << std::endl;
88 }
89 return false;
90 }
91
92 areEqual = equal(gtResult, result);
93
94 if ((!areEqual) && verbose) {
95 std::cout << "Failed: L matrices differ." << std::endl;
96 std::cout << "Result =\n" << result << std::endl;
97 std::cout << "GT =\n" << gtResult << std::endl;
98 return false;
99 }
100
101 vpMatrix LLt = result * result.transpose();
102 areEqual = equal(M, LLt);
103
104 if ((!areEqual) && verbose) {
105 std::cout << "Failed: LL^T differ from M." << std::endl;
106 std::cout << "LL^T =\n" << LLt << std::endl;
107 std::cout << "GT M =\n" << M << std::endl;
108 }
109 if (areEqual && verbose) {
110 std::cout << "Test " << title << " succeeded" << std::endl;
111 }
112
113 return areEqual;
114}
115
116
117TEST_CASE("3 x 3 input", "[cholesky]")
118{
119 vpMatrix M(3, 3);
120 M[0][0] = 4; M[0][1] = 12; M[0][2] = -16;
121 M[1][0] = 12; M[1][1] = 37; M[1][2] = -43;
122 M[2][0] = -16; M[2][1] = -43; M[2][2] = 98;
123
124 vpMatrix gtL(3, 3, 0.);
125 gtL[0][0] = 2; // M[0][1] = 0; M[0][2] = 0;
126 gtL[1][0] = 6; gtL[1][1] = 1; // M[1][2] = 0;
127 gtL[2][0] = -8; gtL[2][1] = 5; gtL[2][2] = 3;
128
129#if defined(VISP_HAVE_LAPACK)
130 SECTION("LAPACK")
131 {
133 CHECK(testCholeskyDecomposition(M, gtL, L, "Test Cholesky's decomposition using Lapack", true));
134 }
135#endif
136
137#if defined(VISP_HAVE_EIGEN3)
138 SECTION("EIGEN3")
139 {
141 CHECK(testCholeskyDecomposition(M, gtL, L, "Test Cholesky's decomposition using Eigen3", true));
142 }
143#endif
144
145// There is a bug with OpenCV 2.4.9 and OpenCV 3.1.0
146// See https://answers.opencv.org/question/99704/cholesky-function-in-opencv/
147#if defined(VISP_HAVE_OPENCV) && (VISP_HAVE_OPENCV_VERSION != 0x020409) && (VISP_HAVE_OPENCV_VERSION != 0x030100)
148 SECTION("OPENCV")
149 {
151 CHECK(testCholeskyDecomposition(M, gtL, L, "Test Cholesky's decomposition using OpenCV", true));
152 }
153#endif
154}
155
156TEST_CASE("4 x 4 input", "[cholesky]")
157{
158 vpMatrix M(4, 4);
159 M[0][0] = 4.16; M[0][1] = -3.12; M[0][2] = 0.56; M[0][3] = -0.10;
160 M[1][0] = -3.12; M[1][1] = 5.03; M[1][2] = -0.83; M[1][3] = 1.18;
161 M[2][0] = 0.56; M[2][1] = -0.83; M[2][2] = 0.76; M[2][3] = 0.34;
162 M[3][0] = -0.10; M[3][1] = 1.18; M[3][2] = 0.34; M[3][3] = 1.18;
163
164 vpMatrix gtL(4, 4, 0.);
165 gtL[0][0] = 2.0396;
166 gtL[1][0] = -1.5297; gtL[1][1] = 1.6401;
167 gtL[2][0] = 0.2746; gtL[2][1] = -0.2500; gtL[2][2] = 0.7887;
168 gtL[3][0] = -0.0490; gtL[3][1] = 0.6737; gtL[3][2] = 0.6617; gtL[3][3] = 0.5347;
169
170#if defined(VISP_HAVE_LAPACK)
171 SECTION("LAPACK")
172 {
174 CHECK(testCholeskyDecomposition(M, gtL, L, "Test Cholesky's decomposition using Lapack", true));
175 }
176#endif
177
178#if defined(VISP_HAVE_EIGEN3)
179 SECTION("EIGEN3")
180 {
182 CHECK(testCholeskyDecomposition(M, gtL, L, "Test Cholesky's decomposition using Eigen3", true));
183 }
184#endif
185
186// There is a bug with OpenCV 2.4.9 and OpenCV 3.1.0
187// See https://answers.opencv.org/question/99704/cholesky-function-in-opencv/
188#if defined(VISP_HAVE_OPENCV) && (VISP_HAVE_OPENCV_VERSION != 0x020409) && (VISP_HAVE_OPENCV_VERSION != 0x030100)
189 SECTION("OPENCV")
190 {
192 CHECK(testCholeskyDecomposition(M, gtL, L, "Test Cholesky's decomposition using OpenCV", true));
193 }
194#endif
195}
196
197int main(int argc, char *argv[])
198{
199 Catch::Session session;
200 session.applyCommandLine(argc, argv);
201 int numFailed = session.run();
202 return numFailed;
203}
204#else
205#include <iostream>
206
207int main()
208{
209 std::cout << "Test ignored: Catch2 is not available" << std::endl;
210 return EXIT_SUCCESS;
211}
212#endif
unsigned int getCols() const
Definition vpArray2D.h:423
unsigned int getRows() const
Definition vpArray2D.h:433
static bool equal(double x, double y, double threshold=0.001)
Definition vpMath.h:470
Implementation of a matrix and operations on matrices.
Definition vpMatrix.h:175
vpMatrix choleskyByOpenCV() const
Compute the Cholesky decomposition of a Hermitian positive-definite matrix using OpenCV library.
vpMatrix choleskyByEigen3() const
Compute the Cholesky decomposition of a Hermitian positive-definite matrix using Eigen3 library.
vpMatrix choleskyByLapack() const
Compute the Cholesky decomposition of a Hermitian positive-definite matrix using Lapack library.
vpMatrix transpose() const