blob: c382f9538da452180f14ad95aa444c5d3604cb69 [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>
Petter Strandmark1e3cbd92012-08-29 09:39:56 -070036
Keir Mierle8ebb0732012-04-30 23:09:08 -070037#include "ceres/block_random_access_matrix.h"
38#include "ceres/block_sparse_matrix.h"
39#include "ceres/block_structure.h"
Petter Strandmark1e3cbd92012-08-29 09:39:56 -070040#include "ceres/cxsparse.h"
Keir Mierle8ebb0732012-04-30 23:09:08 -070041#include "ceres/linear_solver.h"
42#include "ceres/schur_eliminator.h"
43#include "ceres/suitesparse.h"
44#include "ceres/internal/scoped_ptr.h"
45#include "ceres/types.h"
46
47namespace ceres {
48namespace internal {
49
50class BlockSparseMatrixBase;
51
52// Base class for Schur complement based linear least squares
53// solvers. It assumes that the input linear system Ax = b can be
54// partitioned into
55//
56// E y + F z = b
57//
58// Where x = [y;z] is a partition of the variables. The paritioning
59// of the variables is such that, E'E is a block diagonal
60// matrix. Further, the rows of A are ordered so that for every
61// variable block in y, all the rows containing that variable block
62// occur as a vertically contiguous block. i.e the matrix A looks like
63//
64// E F
65// A = [ y1 0 0 0 | z1 0 0 0 z5]
66// [ y1 0 0 0 | z1 z2 0 0 0]
67// [ 0 y2 0 0 | 0 0 z3 0 0]
68// [ 0 0 y3 0 | z1 z2 z3 z4 z5]
69// [ 0 0 y3 0 | z1 0 0 0 z5]
70// [ 0 0 0 y4 | 0 0 0 0 z5]
71// [ 0 0 0 y4 | 0 z2 0 0 0]
72// [ 0 0 0 y4 | 0 0 0 0 0]
73// [ 0 0 0 0 | z1 0 0 0 0]
74// [ 0 0 0 0 | 0 0 z3 z4 z5]
75//
76// This structure should be reflected in the corresponding
77// CompressedRowBlockStructure object associated with A. The linear
78// system Ax = b should either be well posed or the array D below
79// should be non-null and the diagonal matrix corresponding to it
80// should be non-singular.
81//
82// SchurComplementSolver has two sub-classes.
83//
84// DenseSchurComplementSolver: For problems where the Schur complement
85// matrix is small and dense, or if CHOLMOD/SuiteSparse is not
86// installed. For structure from motion problems, this is solver can
87// be used for problems with upto a few hundred cameras.
88//
89// SparseSchurComplementSolver: For problems where the Schur
90// complement matrix is large and sparse. It requires that
91// CHOLMOD/SuiteSparse be installed, as it uses CHOLMOD to find a
92// sparse Cholesky factorization of the Schur complement. This solver
93// can be used for solving structure from motion problems with tens of
94// thousands of cameras, though depending on the exact sparsity
95// structure, it maybe better to use an iterative solver.
96//
97// The two solvers can be instantiated by calling
98// LinearSolver::CreateLinearSolver with LinearSolver::Options::type
99// set to DENSE_SCHUR and SPARSE_SCHUR
Sameer Agarwal0c52f1e2012-09-17 11:30:14 -0700100// respectively. LinearSolver::Options::elimination_groups[0] should be
Keir Mierle8ebb0732012-04-30 23:09:08 -0700101// at least 1.
102class SchurComplementSolver : public BlockSparseMatrixBaseSolver {
103 public:
104 explicit SchurComplementSolver(const LinearSolver::Options& options)
105 : options_(options) {
Sameer Agarwal0c52f1e2012-09-17 11:30:14 -0700106 CHECK_GT(options.elimination_groups.size(), 1);
107 CHECK_GT(options.elimination_groups[0], 0);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700108 }
109
110 // LinearSolver methods
111 virtual ~SchurComplementSolver() {}
112 virtual LinearSolver::Summary SolveImpl(
113 BlockSparseMatrixBase* A,
114 const double* b,
115 const LinearSolver::PerSolveOptions& per_solve_options,
116 double* x);
117
118 protected:
119 const LinearSolver::Options& options() const { return options_; }
120
121 const BlockRandomAccessMatrix* lhs() const { return lhs_.get(); }
122 void set_lhs(BlockRandomAccessMatrix* lhs) { lhs_.reset(lhs); }
123 const double* rhs() const { return rhs_.get(); }
124 void set_rhs(double* rhs) { rhs_.reset(rhs); }
125
126 private:
127 virtual void InitStorage(const CompressedRowBlockStructure* bs) = 0;
128 virtual bool SolveReducedLinearSystem(double* solution) = 0;
129
130 LinearSolver::Options options_;
131
132 scoped_ptr<SchurEliminatorBase> eliminator_;
133 scoped_ptr<BlockRandomAccessMatrix> lhs_;
134 scoped_array<double> rhs_;
135
Sameer Agarwal237d6592012-05-30 20:34:49 -0700136 CERES_DISALLOW_COPY_AND_ASSIGN(SchurComplementSolver);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700137};
138
139// Dense Cholesky factorization based solver.
140class DenseSchurComplementSolver : public SchurComplementSolver {
141 public:
142 explicit DenseSchurComplementSolver(const LinearSolver::Options& options)
143 : SchurComplementSolver(options) {}
144 virtual ~DenseSchurComplementSolver() {}
145
146 private:
147 virtual void InitStorage(const CompressedRowBlockStructure* bs);
148 virtual bool SolveReducedLinearSystem(double* solution);
149
Sameer Agarwal237d6592012-05-30 20:34:49 -0700150 CERES_DISALLOW_COPY_AND_ASSIGN(DenseSchurComplementSolver);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700151};
152
Sameer Agarwalb0518732012-05-29 00:27:57 -0700153
Keir Mierle8ebb0732012-04-30 23:09:08 -0700154// Sparse Cholesky factorization based solver.
155class SparseSchurComplementSolver : public SchurComplementSolver {
156 public:
157 explicit SparseSchurComplementSolver(const LinearSolver::Options& options);
158 virtual ~SparseSchurComplementSolver();
159
160 private:
161 virtual void InitStorage(const CompressedRowBlockStructure* bs);
162 virtual bool SolveReducedLinearSystem(double* solution);
Sameer Agarwalb0518732012-05-29 00:27:57 -0700163 bool SolveReducedLinearSystemUsingSuiteSparse(double* solution);
164 bool SolveReducedLinearSystemUsingCXSparse(double* solution);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700165
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700166 // Size of the blocks in the Schur complement.
167 vector<int> blocks_;
168
Sameer Agarwalb0518732012-05-29 00:27:57 -0700169#ifndef CERES_NO_SUITESPARSE
Keir Mierle8ebb0732012-04-30 23:09:08 -0700170 SuiteSparse ss_;
171 // Symbolic factorization of the reduced linear system. Precomputed
Sameer Agarwala9d8ef82012-05-14 02:28:05 -0700172 // once and reused in subsequent calls.
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700173 cholmod_factor* factor_;
Sameer Agarwalb0518732012-05-29 00:27:57 -0700174#endif // CERES_NO_SUITESPARSE
Petter Strandmark1e3cbd92012-08-29 09:39:56 -0700175
176#ifndef CERES_NO_CXSPARSE
177 CXSparse cxsparse_;
178 // Cached factorization
179 cs_dis* cxsparse_factor_;
180#endif // CERES_NO_CXSPARSE
Sameer Agarwal237d6592012-05-30 20:34:49 -0700181 CERES_DISALLOW_COPY_AND_ASSIGN(SparseSchurComplementSolver);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700182};
Keir Mierle8ebb0732012-04-30 23:09:08 -0700183
184} // namespace internal
185} // namespace ceres
186
187#endif // CERES_INTERNAL_SCHUR_COMPLEMENT_SOLVER_H_