visibility_based_preconditioner.cc 23 KB

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