blob: 5a01006218d26b5e8e62d833e1bc5cf412e1d1b9 [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#include "ceres/visibility_based_preconditioner.h"
32
33#include <algorithm>
34#include <functional>
35#include <iterator>
Keir Mierle8ebb0732012-04-30 23:09:08 -070036#include <set>
37#include <utility>
38#include <vector>
39#include <glog/logging.h>
40#include "Eigen/Dense"
41#include "ceres/block_random_access_sparse_matrix.h"
42#include "ceres/block_sparse_matrix.h"
43#include "ceres/canonical_views_clustering.h"
44#include "ceres/collections_port.h"
45#include "ceres/detect_structure.h"
46#include "ceres/graph.h"
47#include "ceres/graph_algorithms.h"
48#include "ceres/linear_solver.h"
49#include "ceres/schur_eliminator.h"
50#include "ceres/visibility.h"
51#include "ceres/internal/scoped_ptr.h"
52
53namespace ceres {
54namespace internal {
55
56// TODO(sameeragarwal): Currently these are magic weights for the
57// preconditioner construction. Move these higher up into the Options
58// struct and provide some guidelines for choosing them.
59//
60// This will require some more work on the clustering algorithm and
61// possibly some more refactoring of the code.
62static const double kSizePenaltyWeight = 3.0;
63static const double kSimilarityPenaltyWeight = 0.0;
64
65#ifndef CERES_NO_SUITESPARSE
66VisibilityBasedPreconditioner::VisibilityBasedPreconditioner(
67 const CompressedRowBlockStructure& bs,
68 const LinearSolver::Options& options)
69 : options_(options),
70 num_blocks_(0),
71 num_clusters_(0),
72 factor_(NULL) {
73 CHECK_GT(options_.num_eliminate_blocks, 0);
74 CHECK(options_.preconditioner_type == SCHUR_JACOBI ||
75 options_.preconditioner_type == CLUSTER_JACOBI ||
76 options_.preconditioner_type == CLUSTER_TRIDIAGONAL)
77 << "Unknown preconditioner type: " << options_.preconditioner_type;
78 num_blocks_ = bs.cols.size() - options_.num_eliminate_blocks;
79 CHECK_GT(num_blocks_, 0)
80 << "Jacobian should have atleast 1 f_block for "
81 << "visibility based preconditioning.";
82
83 // Vector of camera block sizes
84 block_size_.resize(num_blocks_);
85 for (int i = 0; i < num_blocks_; ++i) {
86 block_size_[i] = bs.cols[i + options_.num_eliminate_blocks].size;
87 }
88
89 const time_t start_time = time(NULL);
90 switch (options_.preconditioner_type) {
91 case SCHUR_JACOBI:
92 ComputeSchurJacobiSparsity(bs);
93 break;
94 case CLUSTER_JACOBI:
95 ComputeClusterJacobiSparsity(bs);
96 break;
97 case CLUSTER_TRIDIAGONAL:
98 ComputeClusterTridiagonalSparsity(bs);
99 break;
100 default:
101 LOG(FATAL) << "Unknown preconditioner type";
102 }
103 const time_t structure_time = time(NULL);
104 InitStorage(bs);
105 const time_t storage_time = time(NULL);
106 InitEliminator(bs);
107 const time_t eliminator_time = time(NULL);
108
109 // Allocate temporary storage for a vector used during
110 // RightMultiply.
111 tmp_rhs_ = CHECK_NOTNULL(ss_.CreateDenseVector(NULL,
112 m_->num_rows(),
113 m_->num_rows()));
114 const time_t init_time = time(NULL);
115 VLOG(2) << "init time: "
116 << init_time - start_time
117 << " structure time: " << structure_time - start_time
118 << " storage time:" << storage_time - structure_time
119 << " eliminator time: " << eliminator_time - storage_time;
120}
121
122VisibilityBasedPreconditioner::~VisibilityBasedPreconditioner() {
123 if (factor_ != NULL) {
124 ss_.Free(factor_);
125 factor_ = NULL;
126 }
127 if (tmp_rhs_ != NULL) {
128 ss_.Free(tmp_rhs_);
129 tmp_rhs_ = NULL;
130 }
131}
132
133// Determine the sparsity structure of the SCHUR_JACOBI
134// preconditioner. SCHUR_JACOBI is an extreme case of a visibility
135// based preconditioner where each camera block corresponds to a
136// cluster and there is no interaction between clusters.
137void VisibilityBasedPreconditioner::ComputeSchurJacobiSparsity(
138 const CompressedRowBlockStructure& bs) {
139 num_clusters_ = num_blocks_;
140 cluster_membership_.resize(num_blocks_);
141 cluster_pairs_.clear();
142
143 // Each camea block is a member of its own cluster and the only
144 // cluster pairs are the self edges (i,i).
145 for (int i = 0; i < num_clusters_; ++i) {
146 cluster_membership_[i] = i;
147 cluster_pairs_.insert(make_pair(i, i));
148 }
149}
150
151// Determine the sparsity structure of the CLUSTER_JACOBI
152// preconditioner. It clusters cameras using their scene
153// visibility. The clusters form the diagonal blocks of the
154// preconditioner matrix.
155void VisibilityBasedPreconditioner::ComputeClusterJacobiSparsity(
156 const CompressedRowBlockStructure& bs) {
157 vector<set<int> > visibility;
158 ComputeVisibility(bs, options_.num_eliminate_blocks, &visibility);
159 CHECK_EQ(num_blocks_, visibility.size());
160 ClusterCameras(visibility);
161 cluster_pairs_.clear();
162 for (int i = 0; i < num_clusters_; ++i) {
163 cluster_pairs_.insert(make_pair(i, i));
164 }
165}
166
167// Determine the sparsity structure of the CLUSTER_TRIDIAGONAL
168// preconditioner. It clusters cameras using using the scene
169// visibility and then finds the strongly interacting pairs of
170// clusters by constructing another graph with the clusters as
171// vertices and approximating it with a degree-2 maximum spanning
172// forest. The set of edges in this forest are the cluster pairs.
173void VisibilityBasedPreconditioner::ComputeClusterTridiagonalSparsity(
174 const CompressedRowBlockStructure& bs) {
175 vector<set<int> > visibility;
176 ComputeVisibility(bs, options_.num_eliminate_blocks, &visibility);
177 CHECK_EQ(num_blocks_, visibility.size());
178 ClusterCameras(visibility);
179
180 // Construct a weighted graph on the set of clusters, where the
181 // edges are the number of 3D points/e_blocks visible in both the
182 // clusters at the ends of the edge. Return an approximate degree-2
183 // maximum spanning forest of this graph.
184 vector<set<int> > cluster_visibility;
185 ComputeClusterVisibility(visibility, &cluster_visibility);
186 scoped_ptr<Graph<int> > cluster_graph(
187 CHECK_NOTNULL(CreateClusterGraph(cluster_visibility)));
188 scoped_ptr<Graph<int> > forest(
189 CHECK_NOTNULL(Degree2MaximumSpanningForest(*cluster_graph)));
190 ForestToClusterPairs(*forest, &cluster_pairs_);
191}
192
193// Allocate storage for the preconditioner matrix.
194void VisibilityBasedPreconditioner::InitStorage(
195 const CompressedRowBlockStructure& bs) {
196 ComputeBlockPairsInPreconditioner(bs);
197 m_.reset(new BlockRandomAccessSparseMatrix(block_size_, block_pairs_));
198}
199
200// Call the canonical views algorithm and cluster the cameras based on
201// their visibility sets. The visibility set of a camera is the set of
202// e_blocks/3D points in the scene that are seen by it.
203//
204// The cluster_membership_ vector is updated to indicate cluster
205// memberships for each camera block.
206void VisibilityBasedPreconditioner::ClusterCameras(
207 const vector<set<int> >& visibility) {
208 scoped_ptr<Graph<int> > schur_complement_graph(
209 CHECK_NOTNULL(CreateSchurComplementGraph(visibility)));
210
211 CanonicalViewsClusteringOptions options;
212 options.size_penalty_weight = kSizePenaltyWeight;
213 options.similarity_penalty_weight = kSimilarityPenaltyWeight;
214
215 vector<int> centers;
216 HashMap<int, int> membership;
217 ComputeCanonicalViewsClustering(*schur_complement_graph,
218 options,
219 &centers,
220 &membership);
221 num_clusters_ = centers.size();
222 CHECK_GT(num_clusters_, 0);
223 VLOG(2) << "num_clusters: " << num_clusters_;
224 FlattenMembershipMap(membership, &cluster_membership_);
225}
226
227// Compute the block sparsity structure of the Schur complement
228// matrix. For each pair of cameras contributing a non-zero cell to
229// the schur complement, determine if that cell is present in the
230// preconditioner or not.
231//
232// A pair of cameras contribute a cell to the preconditioner if they
233// are part of the same cluster or if the the two clusters that they
234// belong have an edge connecting them in the degree-2 maximum
235// spanning forest.
236//
237// For example, a camera pair (i,j) where i belonges to cluster1 and
238// j belongs to cluster2 (assume that cluster1 < cluster2).
239//
240// The cell corresponding to (i,j) is present in the preconditioner
241// if cluster1 == cluster2 or the pair (cluster1, cluster2) were
242// connected by an edge in the degree-2 maximum spanning forest.
243//
244// Since we have already expanded the forest into a set of camera
245// pairs/edges, including self edges, the check can be reduced to
246// checking membership of (cluster1, cluster2) in cluster_pairs_.
247void VisibilityBasedPreconditioner::ComputeBlockPairsInPreconditioner(
248 const CompressedRowBlockStructure& bs) {
249 block_pairs_.clear();
250 for (int i = 0; i < num_blocks_; ++i) {
251 block_pairs_.insert(make_pair(i, i));
252 }
253
254 int r = 0;
255 set<pair<int, int> > skipped_pairs;
256 const int num_row_blocks = bs.rows.size();
257 const int num_eliminate_blocks = options_.num_eliminate_blocks;
258
259 // Iterate over each row of the matrix. The block structure of the
260 // matrix is assumed to be sorted in order of the e_blocks/point
261 // blocks. Thus all row blocks containing an e_block/point occur
262 // contiguously. Further, if present, an e_block is always the first
263 // parameter block in each row block. These structural assumptions
264 // are common to all Schur complement based solvers in Ceres.
265 //
266 // For each e_block/point block we identify the set of cameras
267 // seeing it. The cross product of this set with itself is the set
268 // of non-zero cells contibuted by this e_block.
269 //
270 // The time complexity of this is O(nm^2) where, n is the number of
271 // 3d points and m is the maximum number of cameras seeing any
272 // point, which for most scenes is a fairly small number.
273 while (r < num_row_blocks) {
274 int e_block_id = bs.rows[r].cells.front().block_id;
275 if (e_block_id >= num_eliminate_blocks) {
276 // Skip the rows whose first block is an f_block.
277 break;
278 }
279
280 set<int> f_blocks;
281 for (; r < num_row_blocks; ++r) {
282 const CompressedRow& row = bs.rows[r];
283 if (row.cells.front().block_id != e_block_id) {
284 break;
285 }
286
287 // Iterate over the blocks in the row, ignoring the first block
288 // since it is the one to be eliminated and adding the rest to
289 // the list of f_blocks associated with this e_block.
290 for (int c = 1; c < row.cells.size(); ++c) {
291 const Cell& cell = row.cells[c];
292 const int f_block_id = cell.block_id - num_eliminate_blocks;
293 CHECK_GE(f_block_id, 0);
294 f_blocks.insert(f_block_id);
295 }
296 }
297
298 for (set<int>::const_iterator block1 = f_blocks.begin();
299 block1 != f_blocks.end();
300 ++block1) {
301 set<int>::const_iterator block2 = block1;
302 ++block2;
303 for (; block2 != f_blocks.end(); ++block2) {
304 if (IsBlockPairInPreconditioner(*block1, *block2)) {
305 block_pairs_.insert(make_pair(*block1, *block2));
306 } else {
307 skipped_pairs.insert(make_pair(*block1, *block2));
308 }
309 }
310 }
311 }
312
313 // The remaining rows which do not contain any e_blocks.
314 for (; r < num_row_blocks; ++r) {
315 const CompressedRow& row = bs.rows[r];
316 CHECK_GE(row.cells.front().block_id, num_eliminate_blocks);
317 for (int i = 0; i < row.cells.size(); ++i) {
318 const int block1 = row.cells[i].block_id - num_eliminate_blocks;
319 for (int j = 0; j < row.cells.size(); ++j) {
320 const int block2 = row.cells[j].block_id - num_eliminate_blocks;
321 if (block1 <= block2) {
322 if (IsBlockPairInPreconditioner(block1, block2)) {
323 block_pairs_.insert(make_pair(block1, block2));
324 } else {
325 skipped_pairs.insert(make_pair(block1, block2));
326 }
327 }
328 }
329 }
330 }
331
332 VLOG(1) << "Block pair stats: "
333 << block_pairs_.size() << " included "
334 << skipped_pairs.size() << " excluded";
335}
336
337// Initialize the SchurEliminator.
338void VisibilityBasedPreconditioner::InitEliminator(
339 const CompressedRowBlockStructure& bs) {
340 LinearSolver::Options eliminator_options;
341 eliminator_options.num_eliminate_blocks = options_.num_eliminate_blocks;
342 eliminator_options.num_threads = options_.num_threads;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700343
344 DetectStructure(bs, options_.num_eliminate_blocks,
345 &eliminator_options.row_block_size,
346 &eliminator_options.e_block_size,
347 &eliminator_options.f_block_size);
348
349 eliminator_.reset(SchurEliminatorBase::Create(eliminator_options));
350 eliminator_->Init(options_.num_eliminate_blocks, &bs);
351}
352
Sameer Agarwala9d8ef82012-05-14 02:28:05 -0700353// Update the values of the preconditioner matrix and factorize it.
354bool VisibilityBasedPreconditioner::Update(const BlockSparseMatrixBase& A,
355 const double* D) {
Keir Mierle8ebb0732012-04-30 23:09:08 -0700356 const time_t start_time = time(NULL);
357 const int num_rows = m_->num_rows();
358 CHECK_GT(num_rows, 0);
359
360 // We need a dummy rhs vector and a dummy b vector since the Schur
361 // eliminator combines the computation of the reduced camera matrix
362 // with the computation of the right hand side of that linear
363 // system.
364 //
365 // TODO(sameeragarwal): Perhaps its worth refactoring the
366 // SchurEliminator::Eliminate function to allow NULL for the rhs. As
367 // of now it does not seem to be worth the effort.
368 Vector rhs = Vector::Zero(m_->num_rows());
369 Vector b = Vector::Zero(A.num_rows());
370
371 // Compute a subset of the entries of the Schur complement.
372 eliminator_->Eliminate(&A, b.data(), D, m_.get(), rhs.data());
373
374 // Try factorizing the matrix. For SCHUR_JACOBI and CLUSTER_JACOBI,
375 // this should always succeed modulo some numerical/conditioning
376 // problems. For CLUSTER_TRIDIAGONAL, in general the preconditioner
377 // matrix as constructed is not positive definite. However, we will
378 // go ahead and try factorizing it. If it works, great, otherwise we
379 // scale all the cells in the preconditioner corresponding to the
380 // edges in the degree-2 forest and that guarantees positive
381 // definiteness. The proof of this fact can be found in Lemma 1 in
382 // "Visibility Based Preconditioning for Bundle Adjustment".
383 //
384 // Doing the factorization like this saves us matrix mass when
385 // scaling is not needed, which is quite often in our experience.
386 bool status = Factorize();
387
388 // The scaling only affects the tri-diagonal case, since
389 // ScaleOffDiagonalBlocks only pays attenion to the cells that
390 // belong to the edges of the degree-2 forest. In the SCHUR_JACOBI
391 // and the CLUSTER_JACOBI cases, the preconditioner is guaranteed to
392 // be positive semidefinite.
393 if (!status && options_.preconditioner_type == CLUSTER_TRIDIAGONAL) {
394 VLOG(1) << "Unscaled factorization failed. Retrying with off-diagonal "
395 << "scaling";
396 ScaleOffDiagonalCells();
397 status = Factorize();
398 }
399
400 VLOG(2) << "Compute time: " << time(NULL) - start_time;
401 return status;
402}
403
404// Consider the preconditioner matrix as meta-block matrix, whose
405// blocks correspond to the clusters. Then cluster pairs corresponding
406// to edges in the degree-2 forest are off diagonal entries of this
407// matrix. Scaling these off-diagonal entries by 1/2 forces this
408// matrix to be positive definite.
409void VisibilityBasedPreconditioner::ScaleOffDiagonalCells() {
410 for (set< pair<int, int> >::const_iterator it = block_pairs_.begin();
411 it != block_pairs_.end();
412 ++it) {
413 const int block1 = it->first;
414 const int block2 = it->second;
415 if (!IsBlockPairOffDiagonal(block1, block2)) {
416 continue;
417 }
418
419 int r, c, row_stride, col_stride;
420 CellInfo* cell_info = m_->GetCell(block1, block2,
421 &r, &c,
422 &row_stride, &col_stride);
423 CHECK(cell_info != NULL)
424 << "Cell missing for block pair (" << block1 << "," << block2 << ")"
425 << " cluster pair (" << cluster_membership_[block1]
426 << " " << cluster_membership_[block2] << ")";
427
428 // Ah the magic of tri-diagonal matrices and diagonal
429 // dominance. See Lemma 1 in "Visibility Based Preconditioning
430 // For Bundle Adjustment".
431 MatrixRef m(cell_info->values, row_stride, col_stride);
432 m.block(r, c, block_size_[block1], block_size_[block2]) *= 0.5;
433 }
434}
435
436// Compute the sparse Cholesky factorization of the preconditioner
437// matrix.
438bool VisibilityBasedPreconditioner::Factorize() {
439 // Extract the TripletSparseMatrix that is used for actually storing
440 // S and convert it into a cholmod_sparse object.
441 cholmod_sparse* lhs = ss_.CreateSparseMatrix(
442 down_cast<BlockRandomAccessSparseMatrix*>(
443 m_.get())->mutable_matrix());
444
445 // The matrix is symmetric, and the upper triangular part of the
446 // matrix contains the values.
447 lhs->stype = 1;
448
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700449 // Symbolic factorization is computed if we don't already have one handy.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700450 if (factor_ == NULL) {
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700451 if (options_.use_block_amd) {
452 factor_ = ss_.BlockAnalyzeCholesky(lhs, block_size_, block_size_);
453 } else {
454 factor_ = ss_.AnalyzeCholesky(lhs);
455 }
Keir Mierle8ebb0732012-04-30 23:09:08 -0700456 }
457
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700458 if (VLOG_IS_ON(2)) {
459 cholmod_print_common("Symbolic Analysis", ss_.mutable_cc());
460 }
461
462 CHECK_NOTNULL(factor_);
463
Keir Mierle8ebb0732012-04-30 23:09:08 -0700464 bool status = ss_.Cholesky(lhs, factor_);
465 ss_.Free(lhs);
466 return status;
467}
468
469void VisibilityBasedPreconditioner::RightMultiply(const double* x,
470 double* y) const {
471 CHECK_NOTNULL(x);
472 CHECK_NOTNULL(y);
473 SuiteSparse* ss = const_cast<SuiteSparse*>(&ss_);
474
475 const int num_rows = m_->num_rows();
476 memcpy(CHECK_NOTNULL(tmp_rhs_)->x, x, m_->num_rows() * sizeof(*x));
477 cholmod_dense* solution = CHECK_NOTNULL(ss->Solve(factor_, tmp_rhs_));
478 memcpy(y, solution->x, sizeof(*y) * num_rows);
479 ss->Free(solution);
480}
481
482int VisibilityBasedPreconditioner::num_rows() const {
483 return m_->num_rows();
484}
485
486// Classify camera/f_block pairs as in and out of the preconditioner,
487// based on whether the cluster pair that they belong to is in the
488// preconditioner or not.
489bool VisibilityBasedPreconditioner::IsBlockPairInPreconditioner(
490 const int block1,
491 const int block2) const {
492 int cluster1 = cluster_membership_[block1];
493 int cluster2 = cluster_membership_[block2];
494 if (cluster1 > cluster2) {
495 std::swap(cluster1, cluster2);
496 }
497 return (cluster_pairs_.count(make_pair(cluster1, cluster2)) > 0);
498}
499
500bool VisibilityBasedPreconditioner::IsBlockPairOffDiagonal(
501 const int block1,
502 const int block2) const {
503 return (cluster_membership_[block1] != cluster_membership_[block2]);
504}
505
506// Convert a graph into a list of edges that includes self edges for
507// each vertex.
508void VisibilityBasedPreconditioner::ForestToClusterPairs(
509 const Graph<int>& forest,
510 HashSet<pair<int, int> >* cluster_pairs) const {
511 CHECK_NOTNULL(cluster_pairs)->clear();
512 const HashSet<int>& vertices = forest.vertices();
513 CHECK_EQ(vertices.size(), num_clusters_);
514
515 // Add all the cluster pairs corresponding to the edges in the
516 // forest.
517 for (HashSet<int>::const_iterator it1 = vertices.begin();
518 it1 != vertices.end();
519 ++it1) {
520 const int cluster1 = *it1;
521 cluster_pairs->insert(make_pair(cluster1, cluster1));
522 const HashSet<int>& neighbors = forest.Neighbors(cluster1);
523 for (HashSet<int>::const_iterator it2 = neighbors.begin();
524 it2 != neighbors.end();
525 ++it2) {
526 const int cluster2 = *it2;
527 if (cluster1 < cluster2) {
528 cluster_pairs->insert(make_pair(cluster1, cluster2));
529 }
530 }
531 }
532}
533
534// The visibilty set of a cluster is the union of the visibilty sets
535// of all its cameras. In other words, the set of points visible to
536// any camera in the cluster.
537void VisibilityBasedPreconditioner::ComputeClusterVisibility(
538 const vector<set<int> >& visibility,
539 vector<set<int> >* cluster_visibility) const {
540 CHECK_NOTNULL(cluster_visibility)->resize(0);
541 cluster_visibility->resize(num_clusters_);
542 for (int i = 0; i < num_blocks_; ++i) {
543 const int cluster_id = cluster_membership_[i];
544 (*cluster_visibility)[cluster_id].insert(visibility[i].begin(),
545 visibility[i].end());
546 }
547}
548
549// Construct a graph whose vertices are the clusters, and the edge
550// weights are the number of 3D points visible to cameras in both the
551// vertices.
552Graph<int>* VisibilityBasedPreconditioner::CreateClusterGraph(
553 const vector<set<int> >& cluster_visibility) const {
554 Graph<int>* cluster_graph = new Graph<int>;
555
556 for (int i = 0; i < num_clusters_; ++i) {
557 cluster_graph->AddVertex(i);
558 }
559
560 for (int i = 0; i < num_clusters_; ++i) {
561 const set<int>& cluster_i = cluster_visibility[i];
562 for (int j = i+1; j < num_clusters_; ++j) {
563 vector<int> intersection;
564 const set<int>& cluster_j = cluster_visibility[j];
565 set_intersection(cluster_i.begin(), cluster_i.end(),
566 cluster_j.begin(), cluster_j.end(),
567 back_inserter(intersection));
568
569 if (intersection.size() > 0) {
570 // Clusters interact strongly when they share a large number
571 // of 3D points. The degree-2 maximum spanning forest
572 // alorithm, iterates on the edges in decreasing order of
573 // their weight, which is the number of points shared by the
574 // two cameras that it connects.
575 cluster_graph->AddEdge(i, j, intersection.size());
576 }
577 }
578 }
579 return cluster_graph;
580}
581
582// Canonical views clustering returns a HashMap from vertices to
583// cluster ids. Convert this into a flat array for quick lookup. It is
584// possible that some of the vertices may not be associated with any
585// cluster. In that case, randomly assign them to one of the clusters.
586void VisibilityBasedPreconditioner::FlattenMembershipMap(
587 const HashMap<int, int>& membership_map,
588 vector<int>* membership_vector) const {
589 CHECK_NOTNULL(membership_vector)->resize(0);
590 membership_vector->resize(num_blocks_, -1);
591 // Iterate over the cluster membership map and update the
592 // cluster_membership_ vector assigning arbitrary cluster ids to
593 // the few cameras that have not been clustered.
594 for (HashMap<int, int>::const_iterator it = membership_map.begin();
595 it != membership_map.end();
596 ++it) {
597 const int camera_id = it->first;
598 int cluster_id = it->second;
599
600 // If the view was not clustered, randomly assign it to one of the
601 // clusters. This preserves the mathematical correctness of the
602 // preconditioner. If there are too many views which are not
603 // clustered, it may lead to some quality degradation though.
604 //
605 // TODO(sameeragarwal): Check if a large number of views have not
606 // been clustered and deal with it?
607 if (cluster_id == -1) {
608 cluster_id = camera_id % num_clusters_;
609 }
610
611 membership_vector->at(camera_id) = cluster_id;
612 }
613}
614
615#endif // CERES_NO_SUITESPARSE
616
617} // namespace internal
618} // namespace ceres