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"
39
#include "Eigen/Dense"
40
#include "ceres/block_random_access_sparse_matrix.h"
41
#include "ceres/block_sparse_matrix.h"
42
#include "ceres/canonical_views_clustering.h"
43
#include "ceres/collections_port.h"
44
#include "ceres/detect_structure.h"
45
#include "ceres/graph.h"
46
#include "ceres/graph_algorithms.h"
47
#include "ceres/internal/scoped_ptr.h"
48
#include "ceres/linear_solver.h"
49
#include "ceres/schur_eliminator.h"
50
#include "ceres/visibility.h"
51
#include "glog/logging.h"
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.
60
// This will require some more work on the clustering algorithm and
61
// possibly some more refactoring of the code.
62
static const double kSizePenaltyWeight = 3.0;
63
static const double kSimilarityPenaltyWeight = 0.0;
65
#ifndef CERES_NO_SUITESPARSE
66
VisibilityBasedPreconditioner::VisibilityBasedPreconditioner(
67
const CompressedRowBlockStructure& bs,
68
const LinearSolver::Options& options)
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.";
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;
89
const time_t start_time = time(NULL);
90
switch (options_.preconditioner_type) {
92
ComputeSchurJacobiSparsity(bs);
95
ComputeClusterJacobiSparsity(bs);
97
case CLUSTER_TRIDIAGONAL:
98
ComputeClusterTridiagonalSparsity(bs);
101
LOG(FATAL) << "Unknown preconditioner type";
103
const time_t structure_time = time(NULL);
105
const time_t storage_time = time(NULL);
107
const time_t eliminator_time = time(NULL);
109
// Allocate temporary storage for a vector used during
111
tmp_rhs_ = CHECK_NOTNULL(ss_.CreateDenseVector(NULL,
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;
122
VisibilityBasedPreconditioner::~VisibilityBasedPreconditioner() {
123
if (factor_ != NULL) {
127
if (tmp_rhs_ != NULL) {
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.
137
void VisibilityBasedPreconditioner::ComputeSchurJacobiSparsity(
138
const CompressedRowBlockStructure& bs) {
139
num_clusters_ = num_blocks_;
140
cluster_membership_.resize(num_blocks_);
141
cluster_pairs_.clear();
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));
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.
155
void 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));
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.
173
void 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);
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_);
193
// Allocate storage for the preconditioner matrix.
194
void VisibilityBasedPreconditioner::InitStorage(
195
const CompressedRowBlockStructure& bs) {
196
ComputeBlockPairsInPreconditioner(bs);
197
m_.reset(new BlockRandomAccessSparseMatrix(block_size_, block_pairs_));
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.
204
// The cluster_membership_ vector is updated to indicate cluster
205
// memberships for each camera block.
206
void VisibilityBasedPreconditioner::ClusterCameras(
207
const vector<set<int> >& visibility) {
208
scoped_ptr<Graph<int> > schur_complement_graph(
209
CHECK_NOTNULL(CreateSchurComplementGraph(visibility)));
211
CanonicalViewsClusteringOptions options;
212
options.size_penalty_weight = kSizePenaltyWeight;
213
options.similarity_penalty_weight = kSimilarityPenaltyWeight;
216
HashMap<int, int> membership;
217
ComputeCanonicalViewsClustering(*schur_complement_graph,
221
num_clusters_ = centers.size();
222
CHECK_GT(num_clusters_, 0);
223
VLOG(2) << "num_clusters: " << num_clusters_;
224
FlattenMembershipMap(membership, &cluster_membership_);
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.
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
237
// For example, a camera pair (i,j) where i belonges to cluster1 and
238
// j belongs to cluster2 (assume that cluster1 < cluster2).
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.
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_.
247
void 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));
255
const int num_row_blocks = bs.rows.size();
256
const int num_eliminate_blocks = options_.num_eliminate_blocks;
258
// Iterate over each row of the matrix. The block structure of the
259
// matrix is assumed to be sorted in order of the e_blocks/point
260
// blocks. Thus all row blocks containing an e_block/point occur
261
// contiguously. Further, if present, an e_block is always the first
262
// parameter block in each row block. These structural assumptions
263
// are common to all Schur complement based solvers in Ceres.
265
// For each e_block/point block we identify the set of cameras
266
// seeing it. The cross product of this set with itself is the set
267
// of non-zero cells contibuted by this e_block.
269
// The time complexity of this is O(nm^2) where, n is the number of
270
// 3d points and m is the maximum number of cameras seeing any
271
// point, which for most scenes is a fairly small number.
272
while (r < num_row_blocks) {
273
int e_block_id = bs.rows[r].cells.front().block_id;
274
if (e_block_id >= num_eliminate_blocks) {
275
// Skip the rows whose first block is an f_block.
280
for (; r < num_row_blocks; ++r) {
281
const CompressedRow& row = bs.rows[r];
282
if (row.cells.front().block_id != e_block_id) {
286
// Iterate over the blocks in the row, ignoring the first block
287
// since it is the one to be eliminated and adding the rest to
288
// the list of f_blocks associated with this e_block.
289
for (int c = 1; c < row.cells.size(); ++c) {
290
const Cell& cell = row.cells[c];
291
const int f_block_id = cell.block_id - num_eliminate_blocks;
292
CHECK_GE(f_block_id, 0);
293
f_blocks.insert(f_block_id);
297
for (set<int>::const_iterator block1 = f_blocks.begin();
298
block1 != f_blocks.end();
300
set<int>::const_iterator block2 = block1;
302
for (; block2 != f_blocks.end(); ++block2) {
303
if (IsBlockPairInPreconditioner(*block1, *block2)) {
304
block_pairs_.insert(make_pair(*block1, *block2));
310
// The remaining rows which do not contain any e_blocks.
311
for (; r < num_row_blocks; ++r) {
312
const CompressedRow& row = bs.rows[r];
313
CHECK_GE(row.cells.front().block_id, num_eliminate_blocks);
314
for (int i = 0; i < row.cells.size(); ++i) {
315
const int block1 = row.cells[i].block_id - num_eliminate_blocks;
316
for (int j = 0; j < row.cells.size(); ++j) {
317
const int block2 = row.cells[j].block_id - num_eliminate_blocks;
318
if (block1 <= block2) {
319
if (IsBlockPairInPreconditioner(block1, block2)) {
320
block_pairs_.insert(make_pair(block1, block2));
327
VLOG(1) << "Block pair stats: " << block_pairs_.size();
330
// Initialize the SchurEliminator.
331
void VisibilityBasedPreconditioner::InitEliminator(
332
const CompressedRowBlockStructure& bs) {
333
LinearSolver::Options eliminator_options;
334
eliminator_options.num_eliminate_blocks = options_.num_eliminate_blocks;
335
eliminator_options.num_threads = options_.num_threads;
337
DetectStructure(bs, options_.num_eliminate_blocks,
338
&eliminator_options.row_block_size,
339
&eliminator_options.e_block_size,
340
&eliminator_options.f_block_size);
342
eliminator_.reset(SchurEliminatorBase::Create(eliminator_options));
343
eliminator_->Init(options_.num_eliminate_blocks, &bs);
346
// Update the values of the preconditioner matrix and factorize it.
347
bool VisibilityBasedPreconditioner::Update(const BlockSparseMatrixBase& A,
349
const time_t start_time = time(NULL);
350
const int num_rows = m_->num_rows();
351
CHECK_GT(num_rows, 0);
353
// We need a dummy rhs vector and a dummy b vector since the Schur
354
// eliminator combines the computation of the reduced camera matrix
355
// with the computation of the right hand side of that linear
358
// TODO(sameeragarwal): Perhaps its worth refactoring the
359
// SchurEliminator::Eliminate function to allow NULL for the rhs. As
360
// of now it does not seem to be worth the effort.
361
Vector rhs = Vector::Zero(m_->num_rows());
362
Vector b = Vector::Zero(A.num_rows());
364
// Compute a subset of the entries of the Schur complement.
365
eliminator_->Eliminate(&A, b.data(), D, m_.get(), rhs.data());
367
// Try factorizing the matrix. For SCHUR_JACOBI and CLUSTER_JACOBI,
368
// this should always succeed modulo some numerical/conditioning
369
// problems. For CLUSTER_TRIDIAGONAL, in general the preconditioner
370
// matrix as constructed is not positive definite. However, we will
371
// go ahead and try factorizing it. If it works, great, otherwise we
372
// scale all the cells in the preconditioner corresponding to the
373
// edges in the degree-2 forest and that guarantees positive
374
// definiteness. The proof of this fact can be found in Lemma 1 in
375
// "Visibility Based Preconditioning for Bundle Adjustment".
377
// Doing the factorization like this saves us matrix mass when
378
// scaling is not needed, which is quite often in our experience.
379
bool status = Factorize();
381
// The scaling only affects the tri-diagonal case, since
382
// ScaleOffDiagonalBlocks only pays attenion to the cells that
383
// belong to the edges of the degree-2 forest. In the SCHUR_JACOBI
384
// and the CLUSTER_JACOBI cases, the preconditioner is guaranteed to
385
// be positive semidefinite.
386
if (!status && options_.preconditioner_type == CLUSTER_TRIDIAGONAL) {
387
VLOG(1) << "Unscaled factorization failed. Retrying with off-diagonal "
389
ScaleOffDiagonalCells();
390
status = Factorize();
393
VLOG(2) << "Compute time: " << time(NULL) - start_time;
397
// Consider the preconditioner matrix as meta-block matrix, whose
398
// blocks correspond to the clusters. Then cluster pairs corresponding
399
// to edges in the degree-2 forest are off diagonal entries of this
400
// matrix. Scaling these off-diagonal entries by 1/2 forces this
401
// matrix to be positive definite.
402
void VisibilityBasedPreconditioner::ScaleOffDiagonalCells() {
403
for (set< pair<int, int> >::const_iterator it = block_pairs_.begin();
404
it != block_pairs_.end();
406
const int block1 = it->first;
407
const int block2 = it->second;
408
if (!IsBlockPairOffDiagonal(block1, block2)) {
412
int r, c, row_stride, col_stride;
413
CellInfo* cell_info = m_->GetCell(block1, block2,
415
&row_stride, &col_stride);
416
CHECK(cell_info != NULL)
417
<< "Cell missing for block pair (" << block1 << "," << block2 << ")"
418
<< " cluster pair (" << cluster_membership_[block1]
419
<< " " << cluster_membership_[block2] << ")";
421
// Ah the magic of tri-diagonal matrices and diagonal
422
// dominance. See Lemma 1 in "Visibility Based Preconditioning
423
// For Bundle Adjustment".
424
MatrixRef m(cell_info->values, row_stride, col_stride);
425
m.block(r, c, block_size_[block1], block_size_[block2]) *= 0.5;
429
// Compute the sparse Cholesky factorization of the preconditioner
431
bool VisibilityBasedPreconditioner::Factorize() {
432
// Extract the TripletSparseMatrix that is used for actually storing
433
// S and convert it into a cholmod_sparse object.
434
cholmod_sparse* lhs = ss_.CreateSparseMatrix(
435
down_cast<BlockRandomAccessSparseMatrix*>(
436
m_.get())->mutable_matrix());
438
// The matrix is symmetric, and the upper triangular part of the
439
// matrix contains the values.
442
// Symbolic factorization is computed if we don't already have one handy.
443
if (factor_ == NULL) {
444
if (options_.use_block_amd) {
445
factor_ = ss_.BlockAnalyzeCholesky(lhs, block_size_, block_size_);
447
factor_ = ss_.AnalyzeCholesky(lhs);
451
cholmod_print_common("Symbolic Analysis", ss_.mutable_cc());
455
CHECK_NOTNULL(factor_);
457
bool status = ss_.Cholesky(lhs, factor_);
462
void VisibilityBasedPreconditioner::RightMultiply(const double* x,
466
SuiteSparse* ss = const_cast<SuiteSparse*>(&ss_);
468
const int num_rows = m_->num_rows();
469
memcpy(CHECK_NOTNULL(tmp_rhs_)->x, x, m_->num_rows() * sizeof(*x));
470
cholmod_dense* solution = CHECK_NOTNULL(ss->Solve(factor_, tmp_rhs_));
471
memcpy(y, solution->x, sizeof(*y) * num_rows);
475
int VisibilityBasedPreconditioner::num_rows() const {
476
return m_->num_rows();
479
// Classify camera/f_block pairs as in and out of the preconditioner,
480
// based on whether the cluster pair that they belong to is in the
481
// preconditioner or not.
482
bool VisibilityBasedPreconditioner::IsBlockPairInPreconditioner(
484
const int block2) const {
485
int cluster1 = cluster_membership_[block1];
486
int cluster2 = cluster_membership_[block2];
487
if (cluster1 > cluster2) {
488
std::swap(cluster1, cluster2);
490
return (cluster_pairs_.count(make_pair(cluster1, cluster2)) > 0);
493
bool VisibilityBasedPreconditioner::IsBlockPairOffDiagonal(
495
const int block2) const {
496
return (cluster_membership_[block1] != cluster_membership_[block2]);
499
// Convert a graph into a list of edges that includes self edges for
501
void VisibilityBasedPreconditioner::ForestToClusterPairs(
502
const Graph<int>& forest,
503
HashSet<pair<int, int> >* cluster_pairs) const {
504
CHECK_NOTNULL(cluster_pairs)->clear();
505
const HashSet<int>& vertices = forest.vertices();
506
CHECK_EQ(vertices.size(), num_clusters_);
508
// Add all the cluster pairs corresponding to the edges in the
510
for (HashSet<int>::const_iterator it1 = vertices.begin();
511
it1 != vertices.end();
513
const int cluster1 = *it1;
514
cluster_pairs->insert(make_pair(cluster1, cluster1));
515
const HashSet<int>& neighbors = forest.Neighbors(cluster1);
516
for (HashSet<int>::const_iterator it2 = neighbors.begin();
517
it2 != neighbors.end();
519
const int cluster2 = *it2;
520
if (cluster1 < cluster2) {
521
cluster_pairs->insert(make_pair(cluster1, cluster2));
527
// The visibilty set of a cluster is the union of the visibilty sets
528
// of all its cameras. In other words, the set of points visible to
529
// any camera in the cluster.
530
void VisibilityBasedPreconditioner::ComputeClusterVisibility(
531
const vector<set<int> >& visibility,
532
vector<set<int> >* cluster_visibility) const {
533
CHECK_NOTNULL(cluster_visibility)->resize(0);
534
cluster_visibility->resize(num_clusters_);
535
for (int i = 0; i < num_blocks_; ++i) {
536
const int cluster_id = cluster_membership_[i];
537
(*cluster_visibility)[cluster_id].insert(visibility[i].begin(),
538
visibility[i].end());
542
// Construct a graph whose vertices are the clusters, and the edge
543
// weights are the number of 3D points visible to cameras in both the
545
Graph<int>* VisibilityBasedPreconditioner::CreateClusterGraph(
546
const vector<set<int> >& cluster_visibility) const {
547
Graph<int>* cluster_graph = new Graph<int>;
549
for (int i = 0; i < num_clusters_; ++i) {
550
cluster_graph->AddVertex(i);
553
for (int i = 0; i < num_clusters_; ++i) {
554
const set<int>& cluster_i = cluster_visibility[i];
555
for (int j = i+1; j < num_clusters_; ++j) {
556
vector<int> intersection;
557
const set<int>& cluster_j = cluster_visibility[j];
558
set_intersection(cluster_i.begin(), cluster_i.end(),
559
cluster_j.begin(), cluster_j.end(),
560
back_inserter(intersection));
562
if (intersection.size() > 0) {
563
// Clusters interact strongly when they share a large number
564
// of 3D points. The degree-2 maximum spanning forest
565
// alorithm, iterates on the edges in decreasing order of
566
// their weight, which is the number of points shared by the
567
// two cameras that it connects.
568
cluster_graph->AddEdge(i, j, intersection.size());
572
return cluster_graph;
575
// Canonical views clustering returns a HashMap from vertices to
576
// cluster ids. Convert this into a flat array for quick lookup. It is
577
// possible that some of the vertices may not be associated with any
578
// cluster. In that case, randomly assign them to one of the clusters.
579
void VisibilityBasedPreconditioner::FlattenMembershipMap(
580
const HashMap<int, int>& membership_map,
581
vector<int>* membership_vector) const {
582
CHECK_NOTNULL(membership_vector)->resize(0);
583
membership_vector->resize(num_blocks_, -1);
584
// Iterate over the cluster membership map and update the
585
// cluster_membership_ vector assigning arbitrary cluster ids to
586
// the few cameras that have not been clustered.
587
for (HashMap<int, int>::const_iterator it = membership_map.begin();
588
it != membership_map.end();
590
const int camera_id = it->first;
591
int cluster_id = it->second;
593
// If the view was not clustered, randomly assign it to one of the
594
// clusters. This preserves the mathematical correctness of the
595
// preconditioner. If there are too many views which are not
596
// clustered, it may lead to some quality degradation though.
598
// TODO(sameeragarwal): Check if a large number of views have not
599
// been clustered and deal with it?
600
if (cluster_id == -1) {
601
cluster_id = camera_id % num_clusters_;
604
membership_vector->at(camera_id) = cluster_id;
608
#endif // CERES_NO_SUITESPARSE
610
} // namespace internal