Visual Servoing Platform  version 3.6.1 under development (2024-10-02)
testMatrixCholesky.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 
39 #include <visp3/core/vpConfig.h>
40 
41 #ifdef VISP_HAVE_CATCH2
42 #include <iostream>
43 #include <stdlib.h>
44 #include <visp3/core/vpMatrix.h>
45 
46 #define CATCH_CONFIG_RUNNER
47 #include <catch.hpp>
48 
49 #ifdef ENABLE_VISP_NAMESPACE
50 using namespace VISP_NAMESPACE_NAME;
51 #endif
52 
53 bool equal(const vpMatrix &A, const vpMatrix &B)
54 {
55  unsigned int rowsA = A.getRows(), rowsB = B.getRows();
56  unsigned int colsA = A.getCols(), colsB = B.getCols();
57  bool areEqual = (rowsA == rowsB) && (colsA == colsB);
58 
59  if (!areEqual) {
60  return false;
61  }
62 
63  for (unsigned int r = 0; r < rowsA; ++r) {
64  for (unsigned int c = 0; c < colsA; ++c) {
65  areEqual = areEqual && vpMath::equal(A[r][c], B[r][c]);
66  if (!areEqual) {
67  return false;
68  }
69  }
70  }
71  return areEqual;
72 }
73 
74 bool testCholeskyDecomposition(const vpMatrix &M, const vpMatrix &gtResult, const vpMatrix &result,
75  const std::string &title, const bool &verbose)
76 {
77  if (verbose) {
78  std::cout << "-------------------------" << std::endl;
79  std::cout << "Test: " << title << std::endl;
80  std::cout << std::endl;
81  }
82  unsigned int gtRows = gtResult.getRows(), resRows = result.getRows();
83  unsigned int gtCols = gtResult.getCols(), resCols = result.getCols();
84  bool areEqual = (gtRows == resRows) && (gtCols == resCols);
85  if (!areEqual) {
86  if (verbose) {
87  std::cout << "Failed: dimensions mismatch (" << gtRows << " x " << gtCols << ")";
88  std::cout << " vs (" << resRows << " x " << resCols << ")" << std::endl;
89  }
90  return false;
91  }
92 
93  areEqual = equal(gtResult, result);
94 
95  if ((!areEqual) && verbose) {
96  std::cout << "Failed: L matrices differ." << std::endl;
97  std::cout << "Result =\n" << result << std::endl;
98  std::cout << "GT =\n" << gtResult << std::endl;
99  return false;
100  }
101 
102  vpMatrix LLt = result * result.transpose();
103  areEqual = equal(M, LLt);
104 
105  if ((!areEqual) && verbose) {
106  std::cout << "Failed: LL^T differ from M." << std::endl;
107  std::cout << "LL^T =\n" << LLt << std::endl;
108  std::cout << "GT M =\n" << M << std::endl;
109  }
110  if (areEqual && verbose) {
111  std::cout << "Test " << title << " succeeded" << std::endl;
112  }
113 
114  return areEqual;
115 }
116 
117 
118 TEST_CASE("3 x 3 input", "[cholesky]")
119 {
120  vpMatrix M(3, 3);
121  M[0][0] = 4; M[0][1] = 12; M[0][2] = -16;
122  M[1][0] = 12; M[1][1] = 37; M[1][2] = -43;
123  M[2][0] = -16; M[2][1] = -43; M[2][2] = 98;
124 
125  vpMatrix gtL(3, 3, 0.);
126  gtL[0][0] = 2; // M[0][1] = 0; M[0][2] = 0;
127  gtL[1][0] = 6; gtL[1][1] = 1; // M[1][2] = 0;
128  gtL[2][0] = -8; gtL[2][1] = 5; gtL[2][2] = 3;
129 
130 #if defined(VISP_HAVE_LAPACK)
131  SECTION("LAPACK")
132  {
133  vpMatrix L = M.choleskyByLapack();
134  CHECK(testCholeskyDecomposition(M, gtL, L, "Test 1 Cholesky's decomposition using Lapack", true));
135  }
136 #endif
137 
138 #if defined(VISP_HAVE_EIGEN3)
139  SECTION("EIGEN3")
140  {
141  vpMatrix L = M.choleskyByEigen3();
142  CHECK(testCholeskyDecomposition(M, gtL, L, "Test Cholesky's decomposition using Eigen3", true));
143  }
144 #endif
145 
146 #if defined(VISP_HAVE_OPENCV)
147  SECTION("OPENCV")
148  {
149  vpMatrix L = M.choleskyByOpenCV();
150  CHECK(testCholeskyDecomposition(M, gtL, L, "Test Cholesky's decomposition using OpenCV", true));
151  }
152 #endif
153 }
154 
155 TEST_CASE("4 x 4 input", "[cholesky]")
156 {
157  vpMatrix M(4, 4);
158  M[0][0] = 4.16; M[0][1] = -3.12; M[0][2] = 0.56; M[0][3] = -0.10;
159  M[1][0] = -3.12; M[1][1] = 5.03; M[1][2] = -0.83; M[1][3] = 1.18;
160  M[2][0] = 0.56; M[2][1] = -0.83; M[2][2] = 0.76; M[2][3] = 0.34;
161  M[3][0] = -0.10; M[3][1] = 1.18; M[3][2] = 0.34; M[3][3] = 1.18;
162 
163  vpMatrix gtL(4, 4, 0.);
164  gtL[0][0] = 2.0396;
165  gtL[1][0] = -1.5297; gtL[1][1] = 1.6401;
166  gtL[2][0] = 0.2746; gtL[2][1] = -0.2500; gtL[2][2] = 0.7887;
167  gtL[3][0] = -0.0490; gtL[3][1] = 0.6737; gtL[3][2] = 0.6617; gtL[3][3] = 0.5347;
168 
169 #if defined(VISP_HAVE_LAPACK)
170  SECTION("LAPACK")
171  {
172  vpMatrix L = M.choleskyByLapack();
173  CHECK(testCholeskyDecomposition(M, gtL, L, "Test 1 Cholesky's decomposition using Lapack", true));
174  }
175 #endif
176 
177 #if defined(VISP_HAVE_EIGEN3)
178  SECTION("EIGEN3")
179  {
180  vpMatrix L = M.choleskyByEigen3();
181  CHECK(testCholeskyDecomposition(M, gtL, L, "Test Cholesky's decomposition using Eigen3", true));
182  }
183 #endif
184 
185 #if defined(VISP_HAVE_OPENCV)
186  SECTION("OPENCV")
187  {
188  vpMatrix L = M.choleskyByOpenCV();
189  CHECK(testCholeskyDecomposition(M, gtL, L, "Test Cholesky's decomposition using OpenCV", true));
190  }
191 #endif
192 }
193 
194 int main(int argc, char *argv[])
195 {
196  Catch::Session session; // There must be exactly one instance
197 
198  // Let Catch (using Clara) parse the command line
199  session.applyCommandLine(argc, argv);
200 
201  int numFailed = session.run();
202 
203  // numFailed is clamped to 255 as some unices only use the lower 8 bits.
204  // This clamping has already been applied, so just return it here
205  // You can also do any post run clean-up here
206  return numFailed;
207 }
208 #else
209 #include <iostream>
210 
211 int main()
212 {
213  std::cout << "Test ignored: Catch2 is not available" << std::endl;
214  return EXIT_SUCCESS;
215 }
216 #endif
unsigned int getCols() const
Definition: vpArray2D.h:337
unsigned int getRows() const
Definition: vpArray2D.h:347
static bool equal(double x, double y, double threshold=0.001)
Definition: vpMath.h:459
Implementation of a matrix and operations on matrices.
Definition: vpMatrix.h:169
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