blob: e83892af017a37ac9a1d9fb95de9782d2b308cc4 [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2015 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// An iterative solver for solving the Schur complement/reduced camera
32// linear system that arise in SfM problems.
33
34#ifndef CERES_INTERNAL_IMPLICIT_SCHUR_COMPLEMENT_H_
35#define CERES_INTERNAL_IMPLICIT_SCHUR_COMPLEMENT_H_
36
37#include <memory>
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080038
39#include "ceres/internal/eigen.h"
40#include "ceres/internal/port.h"
Austin Schuh70cc9552019-01-21 19:46:48 -080041#include "ceres/linear_operator.h"
42#include "ceres/linear_solver.h"
43#include "ceres/partitioned_matrix_view.h"
Austin Schuh70cc9552019-01-21 19:46:48 -080044#include "ceres/types.h"
45
46namespace ceres {
47namespace internal {
48
49class BlockSparseMatrix;
50
51// This class implements various linear algebraic operations related
52// to the Schur complement without explicitly forming it.
53//
54//
55// Given a reactangular linear system Ax = b, where
56//
57// A = [E F]
58//
59// The normal equations are given by
60//
61// A'Ax = A'b
62//
63// |E'E E'F||y| = |E'b|
64// |F'E F'F||z| |F'b|
65//
66// and the Schur complement system is given by
67//
68// [F'F - F'E (E'E)^-1 E'F] z = F'b - F'E (E'E)^-1 E'b
69//
70// Now if we wish to solve Ax = b in the least squares sense, one way
71// is to form this Schur complement system and solve it using
72// Preconditioned Conjugate Gradients.
73//
74// The key operation in a conjugate gradient solver is the evaluation of the
75// matrix vector product with the Schur complement
76//
77// S = F'F - F'E (E'E)^-1 E'F
78//
79// It is straightforward to see that matrix vector products with S can
80// be evaluated without storing S in memory. Instead, given (E'E)^-1
81// (which for our purposes is an easily inverted block diagonal
82// matrix), it can be done in terms of matrix vector products with E,
83// F and (E'E)^-1. This class implements this functionality and other
84// auxilliary bits needed to implement a CG solver on the Schur
85// complement using the PartitionedMatrixView object.
86//
87// THREAD SAFETY: This class is nqot thread safe. In particular, the
88// RightMultiply (and the LeftMultiply) methods are not thread safe as
89// they depend on mutable arrays used for the temporaries needed to
90// compute the product y += Sx;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080091class CERES_EXPORT_INTERNAL ImplicitSchurComplement : public LinearOperator {
Austin Schuh70cc9552019-01-21 19:46:48 -080092 public:
93 // num_eliminate_blocks is the number of E blocks in the matrix
94 // A.
95 //
96 // preconditioner indicates whether the inverse of the matrix F'F
97 // should be computed or not as a preconditioner for the Schur
98 // Complement.
99 //
100 // TODO(sameeragarwal): Get rid of the two bools below and replace
101 // them with enums.
102 explicit ImplicitSchurComplement(const LinearSolver::Options& options);
103 virtual ~ImplicitSchurComplement();
104
105 // Initialize the Schur complement for a linear least squares
106 // problem of the form
107 //
108 // |A | x = |b|
109 // |diag(D)| |0|
110 //
111 // If D is null, then it is treated as a zero dimensional matrix. It
112 // is important that the matrix A have a BlockStructure object
113 // associated with it and has a block structure that is compatible
114 // with the SchurComplement solver.
115 void Init(const BlockSparseMatrix& A, const double* D, const double* b);
116
117 // y += Sx, where S is the Schur complement.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800118 void RightMultiply(const double* x, double* y) const final;
Austin Schuh70cc9552019-01-21 19:46:48 -0800119
120 // The Schur complement is a symmetric positive definite matrix,
121 // thus the left and right multiply operators are the same.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800122 void LeftMultiply(const double* x, double* y) const final {
Austin Schuh70cc9552019-01-21 19:46:48 -0800123 RightMultiply(x, y);
124 }
125
126 // y = (E'E)^-1 (E'b - E'F x). Given an estimate of the solution to
127 // the Schur complement system, this method computes the value of
128 // the e_block variables that were eliminated to form the Schur
129 // complement.
130 void BackSubstitute(const double* x, double* y);
131
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800132 int num_rows() const final { return A_->num_cols_f(); }
133 int num_cols() const final { return A_->num_cols_f(); }
134 const Vector& rhs() const { return rhs_; }
Austin Schuh70cc9552019-01-21 19:46:48 -0800135
136 const BlockSparseMatrix* block_diagonal_EtE_inverse() const {
137 return block_diagonal_EtE_inverse_.get();
138 }
139
140 const BlockSparseMatrix* block_diagonal_FtF_inverse() const {
141 return block_diagonal_FtF_inverse_.get();
142 }
143
144 private:
145 void AddDiagonalAndInvert(const double* D, BlockSparseMatrix* matrix);
146 void UpdateRhs();
147
148 const LinearSolver::Options& options_;
149
150 std::unique_ptr<PartitionedMatrixViewBase> A_;
151 const double* D_;
152 const double* b_;
153
154 std::unique_ptr<BlockSparseMatrix> block_diagonal_EtE_inverse_;
155 std::unique_ptr<BlockSparseMatrix> block_diagonal_FtF_inverse_;
156
157 Vector rhs_;
158
159 // Temporary storage vectors used to implement RightMultiply.
160 mutable Vector tmp_rows_;
161 mutable Vector tmp_e_cols_;
162 mutable Vector tmp_e_cols_2_;
163 mutable Vector tmp_f_cols_;
164};
165
166} // namespace internal
167} // namespace ceres
168
169#endif // CERES_INTERNAL_IMPLICIT_SCHUR_COMPLEMENT_H_