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@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/internal/port.h" 39 40#include "ceres/block_random_access_matrix.h" 41#include "ceres/block_sparse_matrix.h" 42#include "ceres/block_structure.h" 43#include "ceres/cxsparse.h" 44#include "ceres/linear_solver.h" 45#include "ceres/schur_eliminator.h" 46#include "ceres/suitesparse.h" 47#include "ceres/internal/scoped_ptr.h" 48#include "ceres/types.h" 49 50#ifdef CERES_USE_EIGEN_SPARSE 51#include "Eigen/SparseCholesky" 52#endif 53 54namespace ceres { 55namespace internal { 56 57class BlockSparseMatrix; 58 59// Base class for Schur complement based linear least squares 60// solvers. It assumes that the input linear system Ax = b can be 61// partitioned into 62// 63// E y + F z = b 64// 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 70// 71// E F 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] 82// 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. 88// 89// SchurComplementSolver has two sub-classes. 90// 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. 95// 96// SparseSchurComplementSolver: For problems where the Schur 97// complement matrix is large and sparse. It requires that 98// CHOLMOD/SuiteSparse be installed, as it uses CHOLMOD to find a 99// sparse Cholesky factorization of the Schur complement. This solver 100// can be used for solving structure from motion problems with tens of 101// thousands of cameras, though depending on the exact sparsity 102// structure, it maybe better to use an iterative solver. 103// 104// The two solvers can be instantiated by calling 105// LinearSolver::CreateLinearSolver with LinearSolver::Options::type 106// set to DENSE_SCHUR and SPARSE_SCHUR 107// respectively. LinearSolver::Options::elimination_groups[0] should be 108// at least 1. 109class SchurComplementSolver : public BlockSparseMatrixSolver { 110 public: 111 explicit SchurComplementSolver(const LinearSolver::Options& options) 112 : options_(options) { 113 CHECK_GT(options.elimination_groups.size(), 1); 114 CHECK_GT(options.elimination_groups[0], 0); 115 } 116 117 // LinearSolver methods 118 virtual ~SchurComplementSolver() {} 119 virtual LinearSolver::Summary SolveImpl( 120 BlockSparseMatrix* A, 121 const double* b, 122 const LinearSolver::PerSolveOptions& per_solve_options, 123 double* x); 124 125 protected: 126 const LinearSolver::Options& options() const { return options_; } 127 128 const BlockRandomAccessMatrix* lhs() const { return lhs_.get(); } 129 void set_lhs(BlockRandomAccessMatrix* lhs) { lhs_.reset(lhs); } 130 const double* rhs() const { return rhs_.get(); } 131 void set_rhs(double* rhs) { rhs_.reset(rhs); } 132 133 private: 134 virtual void InitStorage(const CompressedRowBlockStructure* bs) = 0; 135 virtual LinearSolver::Summary SolveReducedLinearSystem( 136 double* solution) = 0; 137 138 LinearSolver::Options options_; 139 140 scoped_ptr<SchurEliminatorBase> eliminator_; 141 scoped_ptr<BlockRandomAccessMatrix> lhs_; 142 scoped_array<double> rhs_; 143 144 CERES_DISALLOW_COPY_AND_ASSIGN(SchurComplementSolver); 145}; 146 147// Dense Cholesky factorization based solver. 148class DenseSchurComplementSolver : public SchurComplementSolver { 149 public: 150 explicit DenseSchurComplementSolver(const LinearSolver::Options& options) 151 : SchurComplementSolver(options) {} 152 virtual ~DenseSchurComplementSolver() {} 153 154 private: 155 virtual void InitStorage(const CompressedRowBlockStructure* bs); 156 virtual LinearSolver::Summary SolveReducedLinearSystem( 157 double* solution); 158 159 CERES_DISALLOW_COPY_AND_ASSIGN(DenseSchurComplementSolver); 160}; 161 162// Sparse Cholesky factorization based solver. 163class SparseSchurComplementSolver : public SchurComplementSolver { 164 public: 165 explicit SparseSchurComplementSolver(const LinearSolver::Options& options); 166 virtual ~SparseSchurComplementSolver(); 167 168 private: 169 virtual void InitStorage(const CompressedRowBlockStructure* bs); 170 virtual LinearSolver::Summary SolveReducedLinearSystem( 171 double* solution); 172 LinearSolver::Summary SolveReducedLinearSystemUsingSuiteSparse( 173 double* solution); 174 LinearSolver::Summary SolveReducedLinearSystemUsingCXSparse( 175 double* solution); 176 LinearSolver::Summary SolveReducedLinearSystemUsingEigen( 177 double* solution); 178 179 // Size of the blocks in the Schur complement. 180 vector<int> blocks_; 181 182 SuiteSparse ss_; 183 // Symbolic factorization of the reduced linear system. Precomputed 184 // once and reused in subsequent calls. 185 cholmod_factor* factor_; 186 187 CXSparse cxsparse_; 188 // Cached factorization 189 cs_dis* cxsparse_factor_; 190 191#ifdef CERES_USE_EIGEN_SPARSE 192 typedef Eigen::SimplicialLDLT<Eigen::SparseMatrix<double> > SimplicialLDLT; 193 scoped_ptr<SimplicialLDLT> simplicial_ldlt_; 194#endif 195 196 CERES_DISALLOW_COPY_AND_ASSIGN(SparseSchurComplementSolver); 197}; 198 199} // namespace internal 200} // namespace ceres 201 202#endif // CERES_INTERNAL_SCHUR_COMPLEMENT_SOLVER_H_ 203