blob: cfeb65e1eecd2242fd2fe477c50dad71c838b97d [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_ITERATIVE_SCHUR_COMPLEMENT_SOLVER_H_
32#define CERES_INTERNAL_ITERATIVE_SCHUR_COMPLEMENT_SOLVER_H_
33
34#include "ceres/linear_solver.h"
35#include "ceres/internal/eigen.h"
36#include "ceres/internal/scoped_ptr.h"
37#include "ceres/types.h"
38
39namespace ceres {
40namespace internal {
41
42class BlockSparseMatrix;
43class BlockSparseMatrixBase;
44class ImplicitSchurComplement;
45class VisibilityBasedPreconditioner;
46
47// This class implements an iterative solver for the linear least
48// squares problems that have a bi-partitte sparsity structure common
49// to Structure from Motion problems.
50//
51// The algorithm used by this solver was developed in a series of
52// papers - "Agarwal et al, Bundle Adjustment in the Large, ECCV 2010"
53// and "Wu et al, Multicore Bundle Adjustment, submitted to CVPR
54// 2011" at the Univeristy of Washington.
55//
56// The key idea is that one can run Conjugate Gradients on the Schur
57// Complement system without explicitly forming the Schur Complement
58// in memory. The heavy lifting for this is done by the
59// ImplicitSchurComplement class. Not forming the Schur complement in
60// memory and factoring it results in substantial savings in time and
61// memory. Further, iterative solvers like this open up the
62// possibility of solving the Newton equations in a non-linear solver
63// only approximately and terminating early, thereby saving even more
64// time.
65//
66// For the curious, running CG on the Schur complement is the same as
67// running CG on the Normal Equations with an SSOR preconditioner. For
68// a proof of this fact and others related to this solver please see
69// the section on Domain Decomposition Methods in Saad's book
70// "Iterative Methods for Sparse Linear Systems".
71class IterativeSchurComplementSolver : public BlockSparseMatrixBaseSolver {
72 public:
73 explicit IterativeSchurComplementSolver(
74 const LinearSolver::Options& options);
75
76 virtual ~IterativeSchurComplementSolver();
77
78 private:
79 virtual LinearSolver::Summary SolveImpl(
80 BlockSparseMatrixBase* A,
81 const double* b,
82 const LinearSolver::PerSolveOptions& options,
83 double* x);
84
85 LinearSolver::Options options_;
86 scoped_ptr<internal::ImplicitSchurComplement> schur_complement_;
87 scoped_ptr<VisibilityBasedPreconditioner> visibility_based_preconditioner_;
88 Vector reduced_linear_system_solution_;
89};
90
91} // namespace internal
92} // namespace ceres
93
94#endif // CERES_INTERNAL_ITERATIVE_SCHUR_COMPLEMENT_SOLVER_H_