Home | History | Annotate | Download | only in ceres
      1 // 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 (at) google.com (Sameer Agarwal)
     30 
     31 #ifndef CERES_INTERNAL_SCHUR_COMPLEMENT_SOLVER_H_
     32 #define CERES_INTERNAL_SCHUR_COMPLEMENT_SOLVER_H_
     33 
     34 #include <set>
     35 #include <utility>
     36 #include <vector>
     37 
     38 #include "ceres/block_random_access_matrix.h"
     39 #include "ceres/block_sparse_matrix.h"
     40 #include "ceres/block_structure.h"
     41 #include "ceres/cxsparse.h"
     42 #include "ceres/linear_solver.h"
     43 #include "ceres/schur_eliminator.h"
     44 #include "ceres/suitesparse.h"
     45 #include "ceres/internal/scoped_ptr.h"
     46 #include "ceres/types.h"
     47 
     48 namespace ceres {
     49 namespace internal {
     50 
     51 class BlockSparseMatrix;
     52 
     53 // Base class for Schur complement based linear least squares
     54 // solvers. It assumes that the input linear system Ax = b can be
     55 // partitioned into
     56 //
     57 //  E y + F z = b
     58 //
     59 // Where x = [y;z] is a partition of the variables.  The paritioning
     60 // of the variables is such that, E'E is a block diagonal
     61 // matrix. Further, the rows of A are ordered so that for every
     62 // variable block in y, all the rows containing that variable block
     63 // occur as a vertically contiguous block. i.e the matrix A looks like
     64 //
     65 //              E                 F
     66 //  A = [ y1   0   0   0 |  z1    0    0   0    z5]
     67 //      [ y1   0   0   0 |  z1   z2    0   0     0]
     68 //      [  0  y2   0   0 |   0    0   z3   0     0]
     69 //      [  0   0  y3   0 |  z1   z2   z3  z4    z5]
     70 //      [  0   0  y3   0 |  z1    0    0   0    z5]
     71 //      [  0   0   0  y4 |   0    0    0   0    z5]
     72 //      [  0   0   0  y4 |   0   z2    0   0     0]
     73 //      [  0   0   0  y4 |   0    0    0   0     0]
     74 //      [  0   0   0   0 |  z1    0    0   0     0]
     75 //      [  0   0   0   0 |   0    0   z3  z4    z5]
     76 //
     77 // This structure should be reflected in the corresponding
     78 // CompressedRowBlockStructure object associated with A. The linear
     79 // system Ax = b should either be well posed or the array D below
     80 // should be non-null and the diagonal matrix corresponding to it
     81 // should be non-singular.
     82 //
     83 // SchurComplementSolver has two sub-classes.
     84 //
     85 // DenseSchurComplementSolver: For problems where the Schur complement
     86 // matrix is small and dense, or if CHOLMOD/SuiteSparse is not
     87 // installed. For structure from motion problems, this is solver can
     88 // be used for problems with upto a few hundred cameras.
     89 //
     90 // SparseSchurComplementSolver: For problems where the Schur
     91 // complement matrix is large and sparse. It requires that
     92 // CHOLMOD/SuiteSparse be installed, as it uses CHOLMOD to find a
     93 // sparse Cholesky factorization of the Schur complement. This solver
     94 // can be used for solving structure from motion problems with tens of
     95 // thousands of cameras, though depending on the exact sparsity
     96 // structure, it maybe better to use an iterative solver.
     97 //
     98 // The two solvers can be instantiated by calling
     99 // LinearSolver::CreateLinearSolver with LinearSolver::Options::type
    100 // set to DENSE_SCHUR and SPARSE_SCHUR
    101 // respectively. LinearSolver::Options::elimination_groups[0] should be
    102 // at least 1.
    103 class SchurComplementSolver : public BlockSparseMatrixSolver {
    104  public:
    105   explicit SchurComplementSolver(const LinearSolver::Options& options)
    106       : options_(options) {
    107     CHECK_GT(options.elimination_groups.size(), 1);
    108     CHECK_GT(options.elimination_groups[0], 0);
    109   }
    110 
    111   // LinearSolver methods
    112   virtual ~SchurComplementSolver() {}
    113   virtual LinearSolver::Summary SolveImpl(
    114       BlockSparseMatrix* A,
    115       const double* b,
    116       const LinearSolver::PerSolveOptions& per_solve_options,
    117       double* x);
    118 
    119  protected:
    120   const LinearSolver::Options& options() const { return options_; }
    121 
    122   const BlockRandomAccessMatrix* lhs() const { return lhs_.get(); }
    123   void set_lhs(BlockRandomAccessMatrix* lhs) { lhs_.reset(lhs); }
    124   const double* rhs() const { return rhs_.get(); }
    125   void set_rhs(double* rhs) { rhs_.reset(rhs); }
    126 
    127  private:
    128   virtual void InitStorage(const CompressedRowBlockStructure* bs) = 0;
    129   virtual bool SolveReducedLinearSystem(double* solution) = 0;
    130 
    131   LinearSolver::Options options_;
    132 
    133   scoped_ptr<SchurEliminatorBase> eliminator_;
    134   scoped_ptr<BlockRandomAccessMatrix> lhs_;
    135   scoped_array<double> rhs_;
    136 
    137   CERES_DISALLOW_COPY_AND_ASSIGN(SchurComplementSolver);
    138 };
    139 
    140 // Dense Cholesky factorization based solver.
    141 class DenseSchurComplementSolver : public SchurComplementSolver {
    142  public:
    143   explicit DenseSchurComplementSolver(const LinearSolver::Options& options)
    144       : SchurComplementSolver(options) {}
    145   virtual ~DenseSchurComplementSolver() {}
    146 
    147  private:
    148   virtual void InitStorage(const CompressedRowBlockStructure* bs);
    149   virtual bool SolveReducedLinearSystem(double* solution);
    150 
    151   CERES_DISALLOW_COPY_AND_ASSIGN(DenseSchurComplementSolver);
    152 };
    153 
    154 #if !defined(CERES_NO_SUITESPARSE) || !defined(CERES_NO_CXSPARE)
    155 // Sparse Cholesky factorization based solver.
    156 class SparseSchurComplementSolver : public SchurComplementSolver {
    157  public:
    158   explicit SparseSchurComplementSolver(const LinearSolver::Options& options);
    159   virtual ~SparseSchurComplementSolver();
    160 
    161  private:
    162   virtual void InitStorage(const CompressedRowBlockStructure* bs);
    163   virtual bool SolveReducedLinearSystem(double* solution);
    164   bool SolveReducedLinearSystemUsingSuiteSparse(double* solution);
    165   bool SolveReducedLinearSystemUsingCXSparse(double* solution);
    166 
    167   // Size of the blocks in the Schur complement.
    168   vector<int> blocks_;
    169 
    170   SuiteSparse ss_;
    171   // Symbolic factorization of the reduced linear system. Precomputed
    172   // once and reused in subsequent calls.
    173   cholmod_factor* factor_;
    174 
    175   CXSparse cxsparse_;
    176   // Cached factorization
    177   cs_dis* cxsparse_factor_;
    178   CERES_DISALLOW_COPY_AND_ASSIGN(SparseSchurComplementSolver);
    179 };
    180 
    181 #endif  // !defined(CERES_NO_SUITESPARSE) || !defined(CERES_NO_CXSPARE)
    182 }  // namespace internal
    183 }  // namespace ceres
    184 
    185 #endif  // CERES_INTERNAL_SCHUR_COMPLEMENT_SOLVER_H_
    186