blob: 5da9c113b505ff3893752076cf89c833f959c175 [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
39static const bool kFullRank = true;
40static const bool kRankDeficient = false;
41
42template <int kSize>
43typename EigenTypes<kSize, kSize>::Matrix RandomPSDMatrixWithEigenValues(
44 const typename EigenTypes<kSize>::Vector& eigenvalues) {
45 typename EigenTypes<kSize, kSize>::Matrix m;
46 m.setRandom();
47 Eigen::SelfAdjointEigenSolver<typename EigenTypes<kSize, kSize>::Matrix> es(
48 m);
49 return es.eigenvectors() * eigenvalues.asDiagonal() *
50 es.eigenvectors().transpose();
51}
52
53TEST(InvertPSDMatrix, Identity3x3) {
54 const Matrix m = Matrix::Identity(3, 3);
55 const Matrix inverse_m = InvertPSDMatrix<3>(kFullRank, m);
56 EXPECT_NEAR((inverse_m - m).norm() / m.norm(),
57 0.0,
58 std::numeric_limits<double>::epsilon());
59}
60
61TEST(InvertPSDMatrix, FullRank5x5) {
62 EigenTypes<5>::Vector eigenvalues;
63 eigenvalues.setRandom();
64 eigenvalues = eigenvalues.array().abs().matrix();
65 const Matrix m = RandomPSDMatrixWithEigenValues<5>(eigenvalues);
66 const Matrix inverse_m = InvertPSDMatrix<5>(kFullRank, m);
67 EXPECT_NEAR((m * inverse_m - Matrix::Identity(5,5)).norm() / 5.0, 0.0,
68 std::numeric_limits<double>::epsilon());
69}
70
71TEST(InvertPSDMatrix, RankDeficient5x5) {
72 EigenTypes<5>::Vector eigenvalues;
73 eigenvalues.setRandom();
74 eigenvalues = eigenvalues.array().abs().matrix();
75 eigenvalues(3) = 0.0;
76 const Matrix m = RandomPSDMatrixWithEigenValues<5>(eigenvalues);
77 const Matrix inverse_m = InvertPSDMatrix<5>(kRankDeficient, m);
78 Matrix pseudo_identity = Matrix::Identity(5, 5);
79 pseudo_identity(3, 3) = 0.0;
80 EXPECT_NEAR((m * inverse_m * m - m).norm() / m.norm(),
81 0.0,
82 10 * std::numeric_limits<double>::epsilon());
83}
84
85} // namespace internal
86} // namespace ceres