1 // Ceres Solver - A fast non-linear least squares minimizer
2 // Copyright 2015 Google Inc. All rights reserved.
3 // http://ceres-solver.org/
5 // Redistribution and use in source and binary forms, with or without
6 // modification, are permitted provided that the following conditions are met:
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.
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.
29 // Author: sameeragarwal@google.com (Sameer Agarwal)
31 #ifndef CERES_INTERNAL_SCHUR_COMPLEMENT_SOLVER_H_
32 #define CERES_INTERNAL_SCHUR_COMPLEMENT_SOLVER_H_
38 #include "ceres/block_random_access_diagonal_matrix.h"
39 #include "ceres/block_random_access_matrix.h"
40 #include "ceres/block_sparse_matrix.h"
41 #include "ceres/block_structure.h"
42 #include "ceres/internal/port.h"
43 #include "ceres/internal/scoped_ptr.h"
44 #include "ceres/linear_solver.h"
45 #include "ceres/schur_eliminator.h"
46 #include "ceres/types.h"
48 #ifdef CERES_USE_EIGEN_SPARSE
49 #include "Eigen/SparseCholesky"
50 #include "Eigen/OrderingMethods"
56 class BlockSparseMatrix;
59 // Base class for Schur complement based linear least squares
60 // solvers. It assumes that the input linear system Ax = b can be
65 // Where x = [y;z] is a partition of the variables. The paritioning
66 // of the variables is such that, E'E is a block diagonal
67 // matrix. Further, the rows of A are ordered so that for every
68 // variable block in y, all the rows containing that variable block
69 // occur as a vertically contiguous block. i.e the matrix A looks like
72 // A = [ y1 0 0 0 | z1 0 0 0 z5]
73 // [ y1 0 0 0 | z1 z2 0 0 0]
74 // [ 0 y2 0 0 | 0 0 z3 0 0]
75 // [ 0 0 y3 0 | z1 z2 z3 z4 z5]
76 // [ 0 0 y3 0 | z1 0 0 0 z5]
77 // [ 0 0 0 y4 | 0 0 0 0 z5]
78 // [ 0 0 0 y4 | 0 z2 0 0 0]
79 // [ 0 0 0 y4 | 0 0 0 0 0]
80 // [ 0 0 0 0 | z1 0 0 0 0]
81 // [ 0 0 0 0 | 0 0 z3 z4 z5]
83 // This structure should be reflected in the corresponding
84 // CompressedRowBlockStructure object associated with A. The linear
85 // system Ax = b should either be well posed or the array D below
86 // should be non-null and the diagonal matrix corresponding to it
87 // should be non-singular.
89 // SchurComplementSolver has two sub-classes.
91 // DenseSchurComplementSolver: For problems where the Schur complement
92 // matrix is small and dense, or if CHOLMOD/SuiteSparse is not
93 // installed. For structure from motion problems, this is solver can
94 // be used for problems with upto a few hundred cameras.
96 // SparseSchurComplementSolver: For problems where the Schur
97 // complement matrix is large and sparse. It requires that Ceres be
98 // build with at least one sparse linear algebra library, as it
99 // computes a sparse Cholesky factorization of the Schur complement.
101 // This solver can be used for solving structure from motion problems
102 // with tens of thousands of cameras, though depending on the exact
103 // sparsity structure, it maybe better to use an iterative solver.
105 // The two solvers can be instantiated by calling
106 // LinearSolver::CreateLinearSolver with LinearSolver::Options::type
107 // set to DENSE_SCHUR and SPARSE_SCHUR
108 // respectively. LinearSolver::Options::elimination_groups[0] should
110 class SchurComplementSolver : public BlockSparseMatrixSolver {
112 explicit SchurComplementSolver(const LinearSolver::Options& options)
113 : options_(options) {
114 CHECK_GT(options.elimination_groups.size(), 1);
115 CHECK_GT(options.elimination_groups[0], 0);
118 // LinearSolver methods
119 virtual ~SchurComplementSolver() {}
120 virtual LinearSolver::Summary SolveImpl(
121 BlockSparseMatrix* A,
123 const LinearSolver::PerSolveOptions& per_solve_options,
127 const LinearSolver::Options& options() const { return options_; }
129 const BlockRandomAccessMatrix* lhs() const { return lhs_.get(); }
130 void set_lhs(BlockRandomAccessMatrix* lhs) { lhs_.reset(lhs); }
131 const double* rhs() const { return rhs_.get(); }
132 void set_rhs(double* rhs) { rhs_.reset(rhs); }
135 virtual void InitStorage(const CompressedRowBlockStructure* bs) = 0;
136 virtual LinearSolver::Summary SolveReducedLinearSystem(
137 const LinearSolver::PerSolveOptions& per_solve_options,
138 double* solution) = 0;
140 LinearSolver::Options options_;
142 scoped_ptr<SchurEliminatorBase> eliminator_;
143 scoped_ptr<BlockRandomAccessMatrix> lhs_;
144 scoped_array<double> rhs_;
146 CERES_DISALLOW_COPY_AND_ASSIGN(SchurComplementSolver);
149 // Dense Cholesky factorization based solver.
150 class DenseSchurComplementSolver : public SchurComplementSolver {
152 explicit DenseSchurComplementSolver(const LinearSolver::Options& options)
153 : SchurComplementSolver(options) {}
154 virtual ~DenseSchurComplementSolver() {}
157 virtual void InitStorage(const CompressedRowBlockStructure* bs);
158 virtual LinearSolver::Summary SolveReducedLinearSystem(
159 const LinearSolver::PerSolveOptions& per_solve_options,
162 CERES_DISALLOW_COPY_AND_ASSIGN(DenseSchurComplementSolver);
165 // Sparse Cholesky factorization based solver.
166 class SparseSchurComplementSolver : public SchurComplementSolver {
168 explicit SparseSchurComplementSolver(const LinearSolver::Options& options);
169 virtual ~SparseSchurComplementSolver();
172 virtual void InitStorage(const CompressedRowBlockStructure* bs);
173 virtual LinearSolver::Summary SolveReducedLinearSystem(
174 const LinearSolver::PerSolveOptions& per_solve_options,
176 LinearSolver::Summary SolveReducedLinearSystemUsingConjugateGradients(
177 const LinearSolver::PerSolveOptions& per_solve_options,
180 // Size of the blocks in the Schur complement.
181 std::vector<int> blocks_;
182 scoped_ptr<SparseCholesky> sparse_cholesky_;
183 scoped_ptr<BlockRandomAccessDiagonalMatrix> preconditioner_;
184 CERES_DISALLOW_COPY_AND_ASSIGN(SparseSchurComplementSolver);
187 } // namespace internal
190 #endif // CERES_INTERNAL_SCHUR_COMPLEMENT_SOLVER_H_