blob: 279eeab27129141d03956e16a868e6c398ab37e4 [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2017 Google Inc. All rights reserved.
3// http://ceres-solver.org/
4//
5// Redistribution and use in source and binary forms, with or without
6// modification, are permitted provided that the following conditions are met:
7//
8// * Redistributions of source code must retain the above copyright notice,
9// this list of conditions and the following disclaimer.
10// * Redistributions in binary form must reproduce the above copyright notice,
11// this list of conditions and the following disclaimer in the documentation
12// and/or other materials provided with the distribution.
13// * Neither the name of Google Inc. nor the names of its contributors may be
14// used to endorse or promote products derived from this software without
15// specific prior written permission.
16//
17// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
18// AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
19// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
20// ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
21// LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
22// CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
23// SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
24// INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
25// CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
26// ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
27// POSSIBILITY OF SUCH DAMAGE.
28//
29// Author: sameeragarwal@google.com (Sameer Agarwal)
30
31#include "ceres/invert_psd_matrix.h"
32
33#include "ceres/internal/eigen.h"
34#include "gtest/gtest.h"
35
36namespace ceres {
37namespace internal {
38
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080039static constexpr bool kFullRank = true;
40static constexpr bool kRankDeficient = false;
Austin Schuh70cc9552019-01-21 19:46:48 -080041
42template <int kSize>
43typename EigenTypes<kSize, kSize>::Matrix RandomPSDMatrixWithEigenValues(
44 const typename EigenTypes<kSize>::Vector& eigenvalues) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080045 typename EigenTypes<kSize, kSize>::Matrix m(eigenvalues.rows(),
46 eigenvalues.rows());
Austin Schuh70cc9552019-01-21 19:46:48 -080047 m.setRandom();
48 Eigen::SelfAdjointEigenSolver<typename EigenTypes<kSize, kSize>::Matrix> es(
49 m);
50 return es.eigenvectors() * eigenvalues.asDiagonal() *
51 es.eigenvectors().transpose();
52}
53
54TEST(InvertPSDMatrix, Identity3x3) {
55 const Matrix m = Matrix::Identity(3, 3);
56 const Matrix inverse_m = InvertPSDMatrix<3>(kFullRank, m);
57 EXPECT_NEAR((inverse_m - m).norm() / m.norm(),
58 0.0,
59 std::numeric_limits<double>::epsilon());
60}
61
62TEST(InvertPSDMatrix, FullRank5x5) {
63 EigenTypes<5>::Vector eigenvalues;
64 eigenvalues.setRandom();
65 eigenvalues = eigenvalues.array().abs().matrix();
66 const Matrix m = RandomPSDMatrixWithEigenValues<5>(eigenvalues);
67 const Matrix inverse_m = InvertPSDMatrix<5>(kFullRank, m);
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080068 EXPECT_NEAR((m * inverse_m - Matrix::Identity(5, 5)).norm() / 5.0,
69 0.0,
70 10 * std::numeric_limits<double>::epsilon());
Austin Schuh70cc9552019-01-21 19:46:48 -080071}
72
73TEST(InvertPSDMatrix, RankDeficient5x5) {
74 EigenTypes<5>::Vector eigenvalues;
75 eigenvalues.setRandom();
76 eigenvalues = eigenvalues.array().abs().matrix();
77 eigenvalues(3) = 0.0;
78 const Matrix m = RandomPSDMatrixWithEigenValues<5>(eigenvalues);
79 const Matrix inverse_m = InvertPSDMatrix<5>(kRankDeficient, m);
80 Matrix pseudo_identity = Matrix::Identity(5, 5);
81 pseudo_identity(3, 3) = 0.0;
82 EXPECT_NEAR((m * inverse_m * m - m).norm() / m.norm(),
83 0.0,
84 10 * std::numeric_limits<double>::epsilon());
85}
86
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080087TEST(InvertPSDMatrix, DynamicFullRank5x5) {
88 EigenTypes<Eigen::Dynamic>::Vector eigenvalues(5);
89 eigenvalues.setRandom();
90 eigenvalues = eigenvalues.array().abs().matrix();
91 const Matrix m = RandomPSDMatrixWithEigenValues<Eigen::Dynamic>(eigenvalues);
92 const Matrix inverse_m = InvertPSDMatrix<Eigen::Dynamic>(kFullRank, m);
93 EXPECT_NEAR((m * inverse_m - Matrix::Identity(5, 5)).norm() / 5.0,
94 0.0,
95 10 * std::numeric_limits<double>::epsilon());
96}
97
98TEST(InvertPSDMatrix, DynamicRankDeficient5x5) {
99 EigenTypes<Eigen::Dynamic>::Vector eigenvalues(5);
100 eigenvalues.setRandom();
101 eigenvalues = eigenvalues.array().abs().matrix();
102 eigenvalues(3) = 0.0;
103 const Matrix m = RandomPSDMatrixWithEigenValues<Eigen::Dynamic>(eigenvalues);
104 const Matrix inverse_m = InvertPSDMatrix<Eigen::Dynamic>(kRankDeficient, m);
105 Matrix pseudo_identity = Matrix::Identity(5, 5);
106 pseudo_identity(3, 3) = 0.0;
107 EXPECT_NEAR((m * inverse_m * m - m).norm() / m.norm(),
108 0.0,
109 10 * std::numeric_limits<double>::epsilon());
110}
111
Austin Schuh70cc9552019-01-21 19:46:48 -0800112} // namespace internal
113} // namespace ceres