blob: a8f0852b76d230f8bdd0b9d9ee9ed6d5b5c1825a [file] [log] [blame]
Keir Mierle8ebb0732012-04-30 23:09:08 -07001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2010, 2011, 2012 Google Inc. All rights reserved.
3// http://code.google.com/p/ceres-solver/
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_SCHUR_COMPLEMENT_SOLVER_H_
32#define CERES_INTERNAL_SCHUR_COMPLEMENT_SOLVER_H_
33
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -070034#include <set>
35#include <utility>
Keir Mierle8ebb0732012-04-30 23:09:08 -070036#include "ceres/block_random_access_matrix.h"
37#include "ceres/block_sparse_matrix.h"
38#include "ceres/block_structure.h"
39#include "ceres/linear_solver.h"
40#include "ceres/schur_eliminator.h"
41#include "ceres/suitesparse.h"
42#include "ceres/internal/scoped_ptr.h"
43#include "ceres/types.h"
44
45namespace ceres {
46namespace internal {
47
48class BlockSparseMatrixBase;
49
50// Base class for Schur complement based linear least squares
51// solvers. It assumes that the input linear system Ax = b can be
52// partitioned into
53//
54// E y + F z = b
55//
56// Where x = [y;z] is a partition of the variables. The paritioning
57// of the variables is such that, E'E is a block diagonal
58// matrix. Further, the rows of A are ordered so that for every
59// variable block in y, all the rows containing that variable block
60// occur as a vertically contiguous block. i.e the matrix A looks like
61//
62// E F
63// A = [ y1 0 0 0 | z1 0 0 0 z5]
64// [ y1 0 0 0 | z1 z2 0 0 0]
65// [ 0 y2 0 0 | 0 0 z3 0 0]
66// [ 0 0 y3 0 | z1 z2 z3 z4 z5]
67// [ 0 0 y3 0 | z1 0 0 0 z5]
68// [ 0 0 0 y4 | 0 0 0 0 z5]
69// [ 0 0 0 y4 | 0 z2 0 0 0]
70// [ 0 0 0 y4 | 0 0 0 0 0]
71// [ 0 0 0 0 | z1 0 0 0 0]
72// [ 0 0 0 0 | 0 0 z3 z4 z5]
73//
74// This structure should be reflected in the corresponding
75// CompressedRowBlockStructure object associated with A. The linear
76// system Ax = b should either be well posed or the array D below
77// should be non-null and the diagonal matrix corresponding to it
78// should be non-singular.
79//
80// SchurComplementSolver has two sub-classes.
81//
82// DenseSchurComplementSolver: For problems where the Schur complement
83// matrix is small and dense, or if CHOLMOD/SuiteSparse is not
84// installed. For structure from motion problems, this is solver can
85// be used for problems with upto a few hundred cameras.
86//
87// SparseSchurComplementSolver: For problems where the Schur
88// complement matrix is large and sparse. It requires that
89// CHOLMOD/SuiteSparse be installed, as it uses CHOLMOD to find a
90// sparse Cholesky factorization of the Schur complement. This solver
91// can be used for solving structure from motion problems with tens of
92// thousands of cameras, though depending on the exact sparsity
93// structure, it maybe better to use an iterative solver.
94//
95// The two solvers can be instantiated by calling
96// LinearSolver::CreateLinearSolver with LinearSolver::Options::type
97// set to DENSE_SCHUR and SPARSE_SCHUR
98// respectively. LinearSolver::Options::num_eliminate_blocks should be
99// at least 1.
100class SchurComplementSolver : public BlockSparseMatrixBaseSolver {
101 public:
102 explicit SchurComplementSolver(const LinearSolver::Options& options)
103 : options_(options) {
104 CHECK_GT(options.num_eliminate_blocks, 0);
105 }
106
107 // LinearSolver methods
108 virtual ~SchurComplementSolver() {}
109 virtual LinearSolver::Summary SolveImpl(
110 BlockSparseMatrixBase* A,
111 const double* b,
112 const LinearSolver::PerSolveOptions& per_solve_options,
113 double* x);
114
115 protected:
116 const LinearSolver::Options& options() const { return options_; }
117
118 const BlockRandomAccessMatrix* lhs() const { return lhs_.get(); }
119 void set_lhs(BlockRandomAccessMatrix* lhs) { lhs_.reset(lhs); }
120 const double* rhs() const { return rhs_.get(); }
121 void set_rhs(double* rhs) { rhs_.reset(rhs); }
122
123 private:
124 virtual void InitStorage(const CompressedRowBlockStructure* bs) = 0;
125 virtual bool SolveReducedLinearSystem(double* solution) = 0;
126
127 LinearSolver::Options options_;
128
129 scoped_ptr<SchurEliminatorBase> eliminator_;
130 scoped_ptr<BlockRandomAccessMatrix> lhs_;
131 scoped_array<double> rhs_;
132
Sameer Agarwal237d6592012-05-30 20:34:49 -0700133 CERES_DISALLOW_COPY_AND_ASSIGN(SchurComplementSolver);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700134};
135
136// Dense Cholesky factorization based solver.
137class DenseSchurComplementSolver : public SchurComplementSolver {
138 public:
139 explicit DenseSchurComplementSolver(const LinearSolver::Options& options)
140 : SchurComplementSolver(options) {}
141 virtual ~DenseSchurComplementSolver() {}
142
143 private:
144 virtual void InitStorage(const CompressedRowBlockStructure* bs);
145 virtual bool SolveReducedLinearSystem(double* solution);
146
Sameer Agarwal237d6592012-05-30 20:34:49 -0700147 CERES_DISALLOW_COPY_AND_ASSIGN(DenseSchurComplementSolver);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700148};
149
Sameer Agarwalb0518732012-05-29 00:27:57 -0700150
Keir Mierle8ebb0732012-04-30 23:09:08 -0700151// Sparse Cholesky factorization based solver.
152class SparseSchurComplementSolver : public SchurComplementSolver {
153 public:
154 explicit SparseSchurComplementSolver(const LinearSolver::Options& options);
155 virtual ~SparseSchurComplementSolver();
156
157 private:
158 virtual void InitStorage(const CompressedRowBlockStructure* bs);
159 virtual bool SolveReducedLinearSystem(double* solution);
Sameer Agarwalb0518732012-05-29 00:27:57 -0700160 bool SolveReducedLinearSystemUsingSuiteSparse(double* solution);
161 bool SolveReducedLinearSystemUsingCXSparse(double* solution);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700162
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700163 // Size of the blocks in the Schur complement.
164 vector<int> blocks_;
165
Sameer Agarwalb0518732012-05-29 00:27:57 -0700166#ifndef CERES_NO_SUITESPARSE
Keir Mierle8ebb0732012-04-30 23:09:08 -0700167 SuiteSparse ss_;
168 // Symbolic factorization of the reduced linear system. Precomputed
Sameer Agarwala9d8ef82012-05-14 02:28:05 -0700169 // once and reused in subsequent calls.
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700170 cholmod_factor* factor_;
Sameer Agarwalb0518732012-05-29 00:27:57 -0700171#endif // CERES_NO_SUITESPARSE
Sameer Agarwal237d6592012-05-30 20:34:49 -0700172 CERES_DISALLOW_COPY_AND_ASSIGN(SparseSchurComplementSolver);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700173};
Keir Mierle8ebb0732012-04-30 23:09:08 -0700174
175} // namespace internal
176} // namespace ceres
177
178#endif // CERES_INTERNAL_SCHUR_COMPLEMENT_SOLVER_H_