Imported Upstream version ceres 1.13.0
[platform/upstream/ceres-solver.git] / internal / ceres / invert_psd_matrix.h
1 // 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 #ifndef CERES_INTERNAL_INVERT_PSD_MATRIX_H_
32 #define CERES_INTERNAL_INVERT_PSD_MATRIX_H_
33
34 #include "ceres/internal/eigen.h"
35 #include "glog/logging.h"
36 #include "Eigen/Dense"
37
38 namespace ceres {
39 namespace internal {
40
41 // Helper routine to compute the inverse or pseudo-inverse of a
42 // symmetric positive semi-definite matrix.
43 //
44 // assume_full_rank controls whether a Cholesky factorization or an
45 // Singular Value Decomposition is used to compute the inverse and the
46 // pseudo-inverse respectively.
47 //
48 // The template parameter kSize can either be Eigen::Dynamic or a
49 // positive integer equal to the number of rows of m.
50 template <int kSize>
51 typename EigenTypes<kSize, kSize>::Matrix InvertPSDMatrix(
52     const bool assume_full_rank,
53     const typename EigenTypes<kSize, kSize>::Matrix& m) {
54   const int size = m.rows();
55
56   // If the matrix can be assumed to be full rank, then just use the
57   // Cholesky factorization to invert it.
58   if (assume_full_rank) {
59     return m.template selfadjointView<Eigen::Upper>().llt().solve(
60         Matrix::Identity(size, size));
61   }
62
63   Eigen::JacobiSVD<Matrix> svd(m, Eigen::ComputeThinU | Eigen::ComputeThinV);
64   const double tolerance =
65       std::numeric_limits<double>::epsilon() * size * svd.singularValues()(0);
66
67   return svd.matrixV() *
68          (svd.singularValues().array() > tolerance)
69              .select(svd.singularValues().array().inverse(), 0)
70              .matrix()
71              .asDiagonal() *
72          svd.matrixU().adjoint();
73 }
74
75 }  // namespace internal
76 }  // namespace ceres
77
78 #endif // CERES_INTERNAL_INVERT_PSD_MATRIX_H_