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/
5
// Redistribution and use in source and binary forms, with or without
6
// modification, are permitted provided that the following conditions are met:
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.
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.
29
// Author: sameeragarwal@google.com (Sameer Agarwal)
31
#include "ceres/visibility_based_preconditioner.h"
40
#include <glog/logging.h>
41
#include "Eigen/Dense"
42
#include "ceres/block_random_access_sparse_matrix.h"
43
#include "ceres/block_sparse_matrix.h"
44
#include "ceres/canonical_views_clustering.h"
45
#include "ceres/collections_port.h"
46
#include "ceres/detect_structure.h"
47
#include "ceres/graph.h"
48
#include "ceres/graph_algorithms.h"
49
#include "ceres/linear_solver.h"
50
#include "ceres/schur_eliminator.h"
51
#include "ceres/visibility.h"
52
#include "ceres/internal/scoped_ptr.h"
57
// TODO(sameeragarwal): Currently these are magic weights for the
58
// preconditioner construction. Move these higher up into the Options
59
// struct and provide some guidelines for choosing them.
61
// This will require some more work on the clustering algorithm and
62
// possibly some more refactoring of the code.
63
static const double kSizePenaltyWeight = 3.0;
64
static const double kSimilarityPenaltyWeight = 0.0;
66
#ifndef CERES_NO_SUITESPARSE
67
VisibilityBasedPreconditioner::VisibilityBasedPreconditioner(
68
const CompressedRowBlockStructure& bs,
69
const LinearSolver::Options& options)
74
CHECK_GT(options_.num_eliminate_blocks, 0);
75
CHECK(options_.preconditioner_type == SCHUR_JACOBI ||
76
options_.preconditioner_type == CLUSTER_JACOBI ||
77
options_.preconditioner_type == CLUSTER_TRIDIAGONAL)
78
<< "Unknown preconditioner type: " << options_.preconditioner_type;
79
num_blocks_ = bs.cols.size() - options_.num_eliminate_blocks;
80
CHECK_GT(num_blocks_, 0)
81
<< "Jacobian should have atleast 1 f_block for "
82
<< "visibility based preconditioning.";
84
// Vector of camera block sizes
85
block_size_.resize(num_blocks_);
86
for (int i = 0; i < num_blocks_; ++i) {
87
block_size_[i] = bs.cols[i + options_.num_eliminate_blocks].size;
90
const time_t start_time = time(NULL);
91
switch (options_.preconditioner_type) {
93
ComputeSchurJacobiSparsity(bs);
96
ComputeClusterJacobiSparsity(bs);
98
case CLUSTER_TRIDIAGONAL:
99
ComputeClusterTridiagonalSparsity(bs);
102
LOG(FATAL) << "Unknown preconditioner type";
104
const time_t structure_time = time(NULL);
106
const time_t storage_time = time(NULL);
108
const time_t eliminator_time = time(NULL);
110
// Allocate temporary storage for a vector used during
112
tmp_rhs_ = CHECK_NOTNULL(ss_.CreateDenseVector(NULL,
115
const time_t init_time = time(NULL);
116
VLOG(2) << "init time: "
117
<< init_time - start_time
118
<< " structure time: " << structure_time - start_time
119
<< " storage time:" << storage_time - structure_time
120
<< " eliminator time: " << eliminator_time - storage_time;
123
VisibilityBasedPreconditioner::~VisibilityBasedPreconditioner() {
124
if (factor_ != NULL) {
128
if (tmp_rhs_ != NULL) {
134
// Determine the sparsity structure of the SCHUR_JACOBI
135
// preconditioner. SCHUR_JACOBI is an extreme case of a visibility
136
// based preconditioner where each camera block corresponds to a
137
// cluster and there is no interaction between clusters.
138
void VisibilityBasedPreconditioner::ComputeSchurJacobiSparsity(
139
const CompressedRowBlockStructure& bs) {
140
num_clusters_ = num_blocks_;
141
cluster_membership_.resize(num_blocks_);
142
cluster_pairs_.clear();
144
// Each camea block is a member of its own cluster and the only
145
// cluster pairs are the self edges (i,i).
146
for (int i = 0; i < num_clusters_; ++i) {
147
cluster_membership_[i] = i;
148
cluster_pairs_.insert(make_pair(i, i));
152
// Determine the sparsity structure of the CLUSTER_JACOBI
153
// preconditioner. It clusters cameras using their scene
154
// visibility. The clusters form the diagonal blocks of the
155
// preconditioner matrix.
156
void VisibilityBasedPreconditioner::ComputeClusterJacobiSparsity(
157
const CompressedRowBlockStructure& bs) {
158
vector<set<int> > visibility;
159
ComputeVisibility(bs, options_.num_eliminate_blocks, &visibility);
160
CHECK_EQ(num_blocks_, visibility.size());
161
ClusterCameras(visibility);
162
cluster_pairs_.clear();
163
for (int i = 0; i < num_clusters_; ++i) {
164
cluster_pairs_.insert(make_pair(i, i));
168
// Determine the sparsity structure of the CLUSTER_TRIDIAGONAL
169
// preconditioner. It clusters cameras using using the scene
170
// visibility and then finds the strongly interacting pairs of
171
// clusters by constructing another graph with the clusters as
172
// vertices and approximating it with a degree-2 maximum spanning
173
// forest. The set of edges in this forest are the cluster pairs.
174
void VisibilityBasedPreconditioner::ComputeClusterTridiagonalSparsity(
175
const CompressedRowBlockStructure& bs) {
176
vector<set<int> > visibility;
177
ComputeVisibility(bs, options_.num_eliminate_blocks, &visibility);
178
CHECK_EQ(num_blocks_, visibility.size());
179
ClusterCameras(visibility);
181
// Construct a weighted graph on the set of clusters, where the
182
// edges are the number of 3D points/e_blocks visible in both the
183
// clusters at the ends of the edge. Return an approximate degree-2
184
// maximum spanning forest of this graph.
185
vector<set<int> > cluster_visibility;
186
ComputeClusterVisibility(visibility, &cluster_visibility);
187
scoped_ptr<Graph<int> > cluster_graph(
188
CHECK_NOTNULL(CreateClusterGraph(cluster_visibility)));
189
scoped_ptr<Graph<int> > forest(
190
CHECK_NOTNULL(Degree2MaximumSpanningForest(*cluster_graph)));
191
ForestToClusterPairs(*forest, &cluster_pairs_);
194
// Allocate storage for the preconditioner matrix.
195
void VisibilityBasedPreconditioner::InitStorage(
196
const CompressedRowBlockStructure& bs) {
197
ComputeBlockPairsInPreconditioner(bs);
198
m_.reset(new BlockRandomAccessSparseMatrix(block_size_, block_pairs_));
201
// Call the canonical views algorithm and cluster the cameras based on
202
// their visibility sets. The visibility set of a camera is the set of
203
// e_blocks/3D points in the scene that are seen by it.
205
// The cluster_membership_ vector is updated to indicate cluster
206
// memberships for each camera block.
207
void VisibilityBasedPreconditioner::ClusterCameras(
208
const vector<set<int> >& visibility) {
209
scoped_ptr<Graph<int> > schur_complement_graph(
210
CHECK_NOTNULL(CreateSchurComplementGraph(visibility)));
212
CanonicalViewsClusteringOptions options;
213
options.size_penalty_weight = kSizePenaltyWeight;
214
options.similarity_penalty_weight = kSimilarityPenaltyWeight;
217
HashMap<int, int> membership;
218
ComputeCanonicalViewsClustering(*schur_complement_graph,
222
num_clusters_ = centers.size();
223
CHECK_GT(num_clusters_, 0);
224
VLOG(2) << "num_clusters: " << num_clusters_;
225
FlattenMembershipMap(membership, &cluster_membership_);
228
// Compute the block sparsity structure of the Schur complement
229
// matrix. For each pair of cameras contributing a non-zero cell to
230
// the schur complement, determine if that cell is present in the
231
// preconditioner or not.
233
// A pair of cameras contribute a cell to the preconditioner if they
234
// are part of the same cluster or if the the two clusters that they
235
// belong have an edge connecting them in the degree-2 maximum
238
// For example, a camera pair (i,j) where i belonges to cluster1 and
239
// j belongs to cluster2 (assume that cluster1 < cluster2).
241
// The cell corresponding to (i,j) is present in the preconditioner
242
// if cluster1 == cluster2 or the pair (cluster1, cluster2) were
243
// connected by an edge in the degree-2 maximum spanning forest.
245
// Since we have already expanded the forest into a set of camera
246
// pairs/edges, including self edges, the check can be reduced to
247
// checking membership of (cluster1, cluster2) in cluster_pairs_.
248
void VisibilityBasedPreconditioner::ComputeBlockPairsInPreconditioner(
249
const CompressedRowBlockStructure& bs) {
250
block_pairs_.clear();
251
for (int i = 0; i < num_blocks_; ++i) {
252
block_pairs_.insert(make_pair(i, i));
256
set<pair<int, int> > skipped_pairs;
257
const int num_row_blocks = bs.rows.size();
258
const int num_eliminate_blocks = options_.num_eliminate_blocks;
260
// Iterate over each row of the matrix. The block structure of the
261
// matrix is assumed to be sorted in order of the e_blocks/point
262
// blocks. Thus all row blocks containing an e_block/point occur
263
// contiguously. Further, if present, an e_block is always the first
264
// parameter block in each row block. These structural assumptions
265
// are common to all Schur complement based solvers in Ceres.
267
// For each e_block/point block we identify the set of cameras
268
// seeing it. The cross product of this set with itself is the set
269
// of non-zero cells contibuted by this e_block.
271
// The time complexity of this is O(nm^2) where, n is the number of
272
// 3d points and m is the maximum number of cameras seeing any
273
// point, which for most scenes is a fairly small number.
274
while (r < num_row_blocks) {
275
int e_block_id = bs.rows[r].cells.front().block_id;
276
if (e_block_id >= num_eliminate_blocks) {
277
// Skip the rows whose first block is an f_block.
282
for (; r < num_row_blocks; ++r) {
283
const CompressedRow& row = bs.rows[r];
284
if (row.cells.front().block_id != e_block_id) {
288
// Iterate over the blocks in the row, ignoring the first block
289
// since it is the one to be eliminated and adding the rest to
290
// the list of f_blocks associated with this e_block.
291
for (int c = 1; c < row.cells.size(); ++c) {
292
const Cell& cell = row.cells[c];
293
const int f_block_id = cell.block_id - num_eliminate_blocks;
294
CHECK_GE(f_block_id, 0);
295
f_blocks.insert(f_block_id);
299
for (set<int>::const_iterator block1 = f_blocks.begin();
300
block1 != f_blocks.end();
302
set<int>::const_iterator block2 = block1;
304
for (; block2 != f_blocks.end(); ++block2) {
305
if (IsBlockPairInPreconditioner(*block1, *block2)) {
306
block_pairs_.insert(make_pair(*block1, *block2));
308
skipped_pairs.insert(make_pair(*block1, *block2));
314
// The remaining rows which do not contain any e_blocks.
315
for (; r < num_row_blocks; ++r) {
316
const CompressedRow& row = bs.rows[r];
317
CHECK_GE(row.cells.front().block_id, num_eliminate_blocks);
318
for (int i = 0; i < row.cells.size(); ++i) {
319
const int block1 = row.cells[i].block_id - num_eliminate_blocks;
320
for (int j = 0; j < row.cells.size(); ++j) {
321
const int block2 = row.cells[j].block_id - num_eliminate_blocks;
322
if (block1 <= block2) {
323
if (IsBlockPairInPreconditioner(block1, block2)) {
324
block_pairs_.insert(make_pair(block1, block2));
326
skipped_pairs.insert(make_pair(block1, block2));
333
VLOG(1) << "Block pair stats: "
334
<< block_pairs_.size() << " included "
335
<< skipped_pairs.size() << " excluded";
338
// Initialize the SchurEliminator.
339
void VisibilityBasedPreconditioner::InitEliminator(
340
const CompressedRowBlockStructure& bs) {
341
LinearSolver::Options eliminator_options;
342
eliminator_options.num_eliminate_blocks = options_.num_eliminate_blocks;
343
eliminator_options.num_threads = options_.num_threads;
345
DetectStructure(bs, options_.num_eliminate_blocks,
346
&eliminator_options.row_block_size,
347
&eliminator_options.e_block_size,
348
&eliminator_options.f_block_size);
350
eliminator_.reset(SchurEliminatorBase::Create(eliminator_options));
351
eliminator_->Init(options_.num_eliminate_blocks, &bs);
354
// Update the values of the preconditioner matrix and factorize it.
355
bool VisibilityBasedPreconditioner::Update(const BlockSparseMatrixBase& A,
357
const time_t start_time = time(NULL);
358
const int num_rows = m_->num_rows();
359
CHECK_GT(num_rows, 0);
361
// We need a dummy rhs vector and a dummy b vector since the Schur
362
// eliminator combines the computation of the reduced camera matrix
363
// with the computation of the right hand side of that linear
366
// TODO(sameeragarwal): Perhaps its worth refactoring the
367
// SchurEliminator::Eliminate function to allow NULL for the rhs. As
368
// of now it does not seem to be worth the effort.
369
Vector rhs = Vector::Zero(m_->num_rows());
370
Vector b = Vector::Zero(A.num_rows());
372
// Compute a subset of the entries of the Schur complement.
373
eliminator_->Eliminate(&A, b.data(), D, m_.get(), rhs.data());
375
// Try factorizing the matrix. For SCHUR_JACOBI and CLUSTER_JACOBI,
376
// this should always succeed modulo some numerical/conditioning
377
// problems. For CLUSTER_TRIDIAGONAL, in general the preconditioner
378
// matrix as constructed is not positive definite. However, we will
379
// go ahead and try factorizing it. If it works, great, otherwise we
380
// scale all the cells in the preconditioner corresponding to the
381
// edges in the degree-2 forest and that guarantees positive
382
// definiteness. The proof of this fact can be found in Lemma 1 in
383
// "Visibility Based Preconditioning for Bundle Adjustment".
385
// Doing the factorization like this saves us matrix mass when
386
// scaling is not needed, which is quite often in our experience.
387
bool status = Factorize();
389
// The scaling only affects the tri-diagonal case, since
390
// ScaleOffDiagonalBlocks only pays attenion to the cells that
391
// belong to the edges of the degree-2 forest. In the SCHUR_JACOBI
392
// and the CLUSTER_JACOBI cases, the preconditioner is guaranteed to
393
// be positive semidefinite.
394
if (!status && options_.preconditioner_type == CLUSTER_TRIDIAGONAL) {
395
VLOG(1) << "Unscaled factorization failed. Retrying with off-diagonal "
397
ScaleOffDiagonalCells();
398
status = Factorize();
401
VLOG(2) << "Compute time: " << time(NULL) - start_time;
405
// Consider the preconditioner matrix as meta-block matrix, whose
406
// blocks correspond to the clusters. Then cluster pairs corresponding
407
// to edges in the degree-2 forest are off diagonal entries of this
408
// matrix. Scaling these off-diagonal entries by 1/2 forces this
409
// matrix to be positive definite.
410
void VisibilityBasedPreconditioner::ScaleOffDiagonalCells() {
411
for (set< pair<int, int> >::const_iterator it = block_pairs_.begin();
412
it != block_pairs_.end();
414
const int block1 = it->first;
415
const int block2 = it->second;
416
if (!IsBlockPairOffDiagonal(block1, block2)) {
420
int r, c, row_stride, col_stride;
421
CellInfo* cell_info = m_->GetCell(block1, block2,
423
&row_stride, &col_stride);
424
CHECK(cell_info != NULL)
425
<< "Cell missing for block pair (" << block1 << "," << block2 << ")"
426
<< " cluster pair (" << cluster_membership_[block1]
427
<< " " << cluster_membership_[block2] << ")";
429
// Ah the magic of tri-diagonal matrices and diagonal
430
// dominance. See Lemma 1 in "Visibility Based Preconditioning
431
// For Bundle Adjustment".
432
MatrixRef m(cell_info->values, row_stride, col_stride);
433
m.block(r, c, block_size_[block1], block_size_[block2]) *= 0.5;
437
// Compute the sparse Cholesky factorization of the preconditioner
439
bool VisibilityBasedPreconditioner::Factorize() {
440
// Extract the TripletSparseMatrix that is used for actually storing
441
// S and convert it into a cholmod_sparse object.
442
cholmod_sparse* lhs = ss_.CreateSparseMatrix(
443
down_cast<BlockRandomAccessSparseMatrix*>(
444
m_.get())->mutable_matrix());
446
// The matrix is symmetric, and the upper triangular part of the
447
// matrix contains the values.
450
// Symbolic factorization is computed if we don't already have one
452
if (factor_ == NULL) {
453
factor_ = ss_.AnalyzeCholesky(lhs);
456
bool status = ss_.Cholesky(lhs, factor_);
461
void VisibilityBasedPreconditioner::RightMultiply(const double* x,
465
SuiteSparse* ss = const_cast<SuiteSparse*>(&ss_);
467
const int num_rows = m_->num_rows();
468
memcpy(CHECK_NOTNULL(tmp_rhs_)->x, x, m_->num_rows() * sizeof(*x));
469
cholmod_dense* solution = CHECK_NOTNULL(ss->Solve(factor_, tmp_rhs_));
470
memcpy(y, solution->x, sizeof(*y) * num_rows);
474
int VisibilityBasedPreconditioner::num_rows() const {
475
return m_->num_rows();
478
// Classify camera/f_block pairs as in and out of the preconditioner,
479
// based on whether the cluster pair that they belong to is in the
480
// preconditioner or not.
481
bool VisibilityBasedPreconditioner::IsBlockPairInPreconditioner(
483
const int block2) const {
484
int cluster1 = cluster_membership_[block1];
485
int cluster2 = cluster_membership_[block2];
486
if (cluster1 > cluster2) {
487
std::swap(cluster1, cluster2);
489
return (cluster_pairs_.count(make_pair(cluster1, cluster2)) > 0);
492
bool VisibilityBasedPreconditioner::IsBlockPairOffDiagonal(
494
const int block2) const {
495
return (cluster_membership_[block1] != cluster_membership_[block2]);
498
// Convert a graph into a list of edges that includes self edges for
500
void VisibilityBasedPreconditioner::ForestToClusterPairs(
501
const Graph<int>& forest,
502
HashSet<pair<int, int> >* cluster_pairs) const {
503
CHECK_NOTNULL(cluster_pairs)->clear();
504
const HashSet<int>& vertices = forest.vertices();
505
CHECK_EQ(vertices.size(), num_clusters_);
507
// Add all the cluster pairs corresponding to the edges in the
509
for (HashSet<int>::const_iterator it1 = vertices.begin();
510
it1 != vertices.end();
512
const int cluster1 = *it1;
513
cluster_pairs->insert(make_pair(cluster1, cluster1));
514
const HashSet<int>& neighbors = forest.Neighbors(cluster1);
515
for (HashSet<int>::const_iterator it2 = neighbors.begin();
516
it2 != neighbors.end();
518
const int cluster2 = *it2;
519
if (cluster1 < cluster2) {
520
cluster_pairs->insert(make_pair(cluster1, cluster2));
526
// The visibilty set of a cluster is the union of the visibilty sets
527
// of all its cameras. In other words, the set of points visible to
528
// any camera in the cluster.
529
void VisibilityBasedPreconditioner::ComputeClusterVisibility(
530
const vector<set<int> >& visibility,
531
vector<set<int> >* cluster_visibility) const {
532
CHECK_NOTNULL(cluster_visibility)->resize(0);
533
cluster_visibility->resize(num_clusters_);
534
for (int i = 0; i < num_blocks_; ++i) {
535
const int cluster_id = cluster_membership_[i];
536
(*cluster_visibility)[cluster_id].insert(visibility[i].begin(),
537
visibility[i].end());
541
// Construct a graph whose vertices are the clusters, and the edge
542
// weights are the number of 3D points visible to cameras in both the
544
Graph<int>* VisibilityBasedPreconditioner::CreateClusterGraph(
545
const vector<set<int> >& cluster_visibility) const {
546
Graph<int>* cluster_graph = new Graph<int>;
548
for (int i = 0; i < num_clusters_; ++i) {
549
cluster_graph->AddVertex(i);
552
for (int i = 0; i < num_clusters_; ++i) {
553
const set<int>& cluster_i = cluster_visibility[i];
554
for (int j = i+1; j < num_clusters_; ++j) {
555
vector<int> intersection;
556
const set<int>& cluster_j = cluster_visibility[j];
557
set_intersection(cluster_i.begin(), cluster_i.end(),
558
cluster_j.begin(), cluster_j.end(),
559
back_inserter(intersection));
561
if (intersection.size() > 0) {
562
// Clusters interact strongly when they share a large number
563
// of 3D points. The degree-2 maximum spanning forest
564
// alorithm, iterates on the edges in decreasing order of
565
// their weight, which is the number of points shared by the
566
// two cameras that it connects.
567
cluster_graph->AddEdge(i, j, intersection.size());
571
return cluster_graph;
574
// Canonical views clustering returns a HashMap from vertices to
575
// cluster ids. Convert this into a flat array for quick lookup. It is
576
// possible that some of the vertices may not be associated with any
577
// cluster. In that case, randomly assign them to one of the clusters.
578
void VisibilityBasedPreconditioner::FlattenMembershipMap(
579
const HashMap<int, int>& membership_map,
580
vector<int>* membership_vector) const {
581
CHECK_NOTNULL(membership_vector)->resize(0);
582
membership_vector->resize(num_blocks_, -1);
583
// Iterate over the cluster membership map and update the
584
// cluster_membership_ vector assigning arbitrary cluster ids to
585
// the few cameras that have not been clustered.
586
for (HashMap<int, int>::const_iterator it = membership_map.begin();
587
it != membership_map.end();
589
const int camera_id = it->first;
590
int cluster_id = it->second;
592
// If the view was not clustered, randomly assign it to one of the
593
// clusters. This preserves the mathematical correctness of the
594
// preconditioner. If there are too many views which are not
595
// clustered, it may lead to some quality degradation though.
597
// TODO(sameeragarwal): Check if a large number of views have not
598
// been clustered and deal with it?
599
if (cluster_id == -1) {
600
cluster_id = camera_id % num_clusters_;
603
membership_vector->at(camera_id) = cluster_id;
607
#endif // CERES_NO_SUITESPARSE
609
} // namespace internal