• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
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 // Preconditioners for linear systems that arise in Structure from
32 // Motion problems. VisibilityBasedPreconditioner implements three
33 // preconditioners:
34 //
35 //  SCHUR_JACOBI
36 //  CLUSTER_JACOBI
37 //  CLUSTER_TRIDIAGONAL
38 //
39 // Detailed descriptions of these preconditions beyond what is
40 // documented here can be found in
41 //
42 // Bundle Adjustment in the Large
43 // S. Agarwal, N. Snavely, S. Seitz & R. Szeliski, ECCV 2010
44 // http://www.cs.washington.edu/homes/sagarwal/bal.pdf
45 //
46 // Visibility Based Preconditioning for Bundle Adjustment
47 // A. Kushal & S. Agarwal, submitted to CVPR 2012
48 // http://www.cs.washington.edu/homes/sagarwal/vbp.pdf
49 //
50 // The three preconditioners share enough code that its most efficient
51 // to implement them as part of the same code base.
52 
53 #ifndef CERES_INTERNAL_VISIBILITY_BASED_PRECONDITIONER_H_
54 #define CERES_INTERNAL_VISIBILITY_BASED_PRECONDITIONER_H_
55 
56 #include <set>
57 #include <vector>
58 #include <utility>
59 #include "ceres/collections_port.h"
60 #include "ceres/graph.h"
61 #include "ceres/linear_solver.h"
62 #include "ceres/linear_operator.h"
63 #include "ceres/suitesparse.h"
64 #include "ceres/internal/macros.h"
65 #include "ceres/internal/scoped_ptr.h"
66 
67 namespace ceres {
68 namespace internal {
69 
70 class BlockRandomAccessSparseMatrix;
71 class BlockSparseMatrixBase;
72 struct CompressedRowBlockStructure;
73 class SchurEliminatorBase;
74 
75 // This class implements three preconditioners for Structure from
76 // Motion/Bundle Adjustment problems. The name
77 // VisibilityBasedPreconditioner comes from the fact that the sparsity
78 // structure of the preconditioner matrix is determined by analyzing
79 // the visibility structure of the scene, i.e. which cameras see which
80 // points.
81 //
82 // Strictly speaking, SCHUR_JACOBI is not a visibility based
83 // preconditioner but it is an extreme case of CLUSTER_JACOBI, where
84 // every cluster contains exactly one camera block. Treating it as a
85 // special case of CLUSTER_JACOBI makes it easy to implement as part
86 // of the same code base with no significant loss of performance.
87 //
88 // In the following, we will only discuss CLUSTER_JACOBI and
89 // CLUSTER_TRIDIAGONAL.
90 //
91 // The key idea of visibility based preconditioning is to identify
92 // cameras that we expect have strong interactions, and then using the
93 // entries in the Schur complement matrix corresponding to these
94 // camera pairs as an approximation to the full Schur complement.
95 //
96 // CLUSTER_JACOBI identifies these camera pairs by clustering cameras,
97 // and considering all non-zero camera pairs within each cluster. The
98 // clustering in the current implementation is done using the
99 // Canonical Views algorithm of Simon et al. (see
100 // canonical_views_clustering.h). For the purposes of clustering, the
101 // similarity or the degree of interaction between a pair of cameras
102 // is measured by counting the number of points visible in both the
103 // cameras. Thus the name VisibilityBasedPreconditioner. Further, if we
104 // were to permute the parameter blocks such that all the cameras in
105 // the same cluster occur contiguously, the preconditioner matrix will
106 // be a block diagonal matrix with blocks corresponding to the
107 // clusters. Thus in analogy with the Jacobi preconditioner we refer
108 // to this as the CLUSTER_JACOBI preconditioner.
109 //
110 // CLUSTER_TRIDIAGONAL adds more mass to the CLUSTER_JACOBI
111 // preconditioner by considering the interaction between clusters and
112 // identifying strong interactions between cluster pairs. This is done
113 // by constructing a weighted graph on the clusters, with the weight
114 // on the edges connecting two clusters proportional to the number of
115 // 3D points visible to cameras in both the clusters. A degree-2
116 // maximum spanning forest is identified in this graph and the camera
117 // pairs contained in the edges of this forest are added to the
118 // preconditioner. The detailed reasoning for this construction is
119 // explained in the paper mentioned above.
120 //
121 // Degree-2 spanning trees and forests have the property that they
122 // correspond to tri-diagonal matrices. Thus there exist a permutation
123 // of the camera blocks under which the CLUSTER_TRIDIAGONAL
124 // preconditioner matrix is a block tridiagonal matrix, and thus the
125 // name for the preconditioner.
126 //
127 // Thread Safety: This class is NOT thread safe.
128 //
129 // Example usage:
130 //
131 //   LinearSolver::Options options;
132 //   options.preconditioner_type = CLUSTER_JACOBI;
133 //   options.num_eliminate_blocks = num_points;
134 //   VisibilityBasedPreconditioner preconditioner(
135 //      *A.block_structure(), options);
136 //   preconditioner.Update(A, NULL);
137 //   preconditioner.RightMultiply(x, y);
138 //
139 
140 #ifndef CERES_NO_SUITESPARSE
141 class VisibilityBasedPreconditioner : public LinearOperator {
142  public:
143   // Initialize the symbolic structure of the preconditioner. bs is
144   // the block structure of the linear system to be solved. It is used
145   // to determine the sparsity structure of the preconditioner matrix.
146   //
147   // It has the same structural requirement as other Schur complement
148   // based solvers. Please see schur_eliminator.h for more details.
149   //
150   // LinearSolver::Options::num_eliminate_blocks should be set to the
151   // number of e_blocks in the block structure.
152   //
153   // TODO(sameeragarwal): The use of LinearSolver::Options should
154   // ultimately be replaced with Preconditioner::Options and some sort
155   // of preconditioner factory along the lines of
156   // LinearSolver::CreateLinearSolver. I will wait to do this till I
157   // create a general purpose block Jacobi preconditioner for general
158   // sparse problems along with a CGLS solver.
159   VisibilityBasedPreconditioner(const CompressedRowBlockStructure& bs,
160                                 const LinearSolver::Options& options);
161   virtual ~VisibilityBasedPreconditioner();
162 
163   // Update the numerical value of the preconditioner for the linear
164   // system:
165   //
166   //  |   A   | x = |b|
167   //  |diag(D)|     |0|
168   //
169   // for some vector b. It is important that the matrix A have the
170   // same block structure as the one used to construct this object.
171   //
172   // D can be NULL, in which case its interpreted as a diagonal matrix
173   // of size zero.
174   bool Update(const BlockSparseMatrixBase& A, const double* D);
175 
176 
177   // LinearOperator interface. Since the operator is symmetric,
178   // LeftMultiply and num_cols are just calls to RightMultiply and
179   // num_rows respectively. Update() must be called before
180   // RightMultiply can be called.
181   virtual void RightMultiply(const double* x, double* y) const;
LeftMultiply(const double * x,double * y)182   virtual void LeftMultiply(const double* x, double* y) const {
183     RightMultiply(x, y);
184   }
185   virtual int num_rows() const;
num_cols()186   virtual int num_cols() const { return num_rows(); }
187 
188   friend class VisibilityBasedPreconditionerTest;
189  private:
190   void ComputeSchurJacobiSparsity(const CompressedRowBlockStructure& bs);
191   void ComputeClusterJacobiSparsity(const CompressedRowBlockStructure& bs);
192   void ComputeClusterTridiagonalSparsity(const CompressedRowBlockStructure& bs);
193   void InitStorage(const CompressedRowBlockStructure& bs);
194   void InitEliminator(const CompressedRowBlockStructure& bs);
195   bool Factorize();
196   void ScaleOffDiagonalCells();
197 
198   void ClusterCameras(const vector< set<int> >& visibility);
199   void FlattenMembershipMap(const HashMap<int, int>& membership_map,
200                             vector<int>* membership_vector) const;
201   void ComputeClusterVisibility(const vector<set<int> >& visibility,
202                                 vector<set<int> >* cluster_visibility) const;
203   Graph<int>* CreateClusterGraph(const vector<set<int> >& visibility) const;
204   void ForestToClusterPairs(const Graph<int>& forest,
205                             HashSet<pair<int, int> >* cluster_pairs) const;
206   void ComputeBlockPairsInPreconditioner(const CompressedRowBlockStructure& bs);
207   bool IsBlockPairInPreconditioner(int block1, int block2) const;
208   bool IsBlockPairOffDiagonal(int block1, int block2) const;
209 
210   LinearSolver::Options options_;
211 
212   // Number of parameter blocks in the schur complement.
213   int num_blocks_;
214   int num_clusters_;
215 
216   // Sizes of the blocks in the schur complement.
217   vector<int> block_size_;
218 
219   // Mapping from cameras to clusters.
220   vector<int> cluster_membership_;
221 
222   // Non-zero camera pairs from the schur complement matrix that are
223   // present in the preconditioner, sorted by row (first element of
224   // each pair), then column (second).
225   set<pair<int, int> > block_pairs_;
226 
227   // Set of cluster pairs (including self pairs (i,i)) in the
228   // preconditioner.
229   HashSet<pair<int, int> > cluster_pairs_;
230   scoped_ptr<SchurEliminatorBase> eliminator_;
231 
232   // Preconditioner matrix.
233   scoped_ptr<BlockRandomAccessSparseMatrix> m_;
234 
235   // RightMultiply is a const method for LinearOperators. It is
236   // implemented using CHOLMOD's sparse triangular matrix solve
237   // function. This however requires non-const access to the
238   // SuiteSparse context object, even though it does not result in any
239   // of the state of the preconditioner being modified.
240   SuiteSparse ss_;
241 
242   // Symbolic and numeric factorization of the preconditioner.
243   cholmod_factor* factor_;
244 
245   // Temporary vector used by RightMultiply.
246   cholmod_dense* tmp_rhs_;
247   CERES_DISALLOW_COPY_AND_ASSIGN(VisibilityBasedPreconditioner);
248 };
249 #else  // SuiteSparse
250 // If SuiteSparse is not compiled in, the preconditioner is not
251 // available.
252 class VisibilityBasedPreconditioner : public LinearOperator {
253  public:
VisibilityBasedPreconditioner(const CompressedRowBlockStructure & bs,const LinearSolver::Options & options)254   VisibilityBasedPreconditioner(const CompressedRowBlockStructure& bs,
255                                 const LinearSolver::Options& options) {
256     LOG(FATAL) << "Visibility based preconditioning is not available. Please "
257         "build Ceres with SuiteSparse.";
258   }
~VisibilityBasedPreconditioner()259   virtual ~VisibilityBasedPreconditioner() {}
RightMultiply(const double * x,double * y)260   virtual void RightMultiply(const double* x, double* y) const {}
LeftMultiply(const double * x,double * y)261   virtual void LeftMultiply(const double* x, double* y) const {}
num_rows()262   virtual int num_rows() const { return -1; }
num_cols()263   virtual int num_cols() const { return -1; }
Update(const BlockSparseMatrixBase & A,const double * D)264   bool Update(const BlockSparseMatrixBase& A, const double* D) {
265     return false;
266   }
267 };
268 #endif  // CERES_NO_SUITESPARSE
269 
270 }  // namespace internal
271 }  // namespace ceres
272 
273 #endif  // CERES_INTERNAL_VISIBILITY_BASED_PRECONDITIONER_H_
274