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