reorder_program.cc 25 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644
  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/reorder_program.h"
  31. #include <algorithm>
  32. #include <map>
  33. #include <memory>
  34. #include <numeric>
  35. #include <set>
  36. #include <string>
  37. #include <vector>
  38. #include "Eigen/SparseCore"
  39. #include "ceres/internal/config.h"
  40. #include "ceres/internal/export.h"
  41. #include "ceres/ordered_groups.h"
  42. #include "ceres/parameter_block.h"
  43. #include "ceres/parameter_block_ordering.h"
  44. #include "ceres/problem_impl.h"
  45. #include "ceres/program.h"
  46. #include "ceres/residual_block.h"
  47. #include "ceres/solver.h"
  48. #include "ceres/suitesparse.h"
  49. #include "ceres/triplet_sparse_matrix.h"
  50. #include "ceres/types.h"
  51. #ifdef CERES_USE_EIGEN_SPARSE
  52. #ifndef CERES_NO_EIGEN_METIS
  53. #include <iostream> // Need this because MetisSupport refers to std::cerr.
  54. #include "Eigen/MetisSupport"
  55. #endif
  56. #include "Eigen/OrderingMethods"
  57. #endif
  58. #include "glog/logging.h"
  59. namespace ceres::internal {
  60. namespace {
  61. // Find the minimum index of any parameter block to the given
  62. // residual. Parameter blocks that have indices greater than
  63. // size_of_first_elimination_group are considered to have an index
  64. // equal to size_of_first_elimination_group.
  65. static int MinParameterBlock(const ResidualBlock* residual_block,
  66. int size_of_first_elimination_group) {
  67. int min_parameter_block_position = size_of_first_elimination_group;
  68. for (int i = 0; i < residual_block->NumParameterBlocks(); ++i) {
  69. ParameterBlock* parameter_block = residual_block->parameter_blocks()[i];
  70. if (!parameter_block->IsConstant()) {
  71. CHECK_NE(parameter_block->index(), -1)
  72. << "Did you forget to call Program::SetParameterOffsetsAndIndex()? "
  73. << "This is a Ceres bug; please contact the developers!";
  74. min_parameter_block_position =
  75. std::min(parameter_block->index(), min_parameter_block_position);
  76. }
  77. }
  78. return min_parameter_block_position;
  79. }
  80. Eigen::SparseMatrix<int> CreateBlockJacobian(
  81. const TripletSparseMatrix& block_jacobian_transpose) {
  82. using SparseMatrix = Eigen::SparseMatrix<int>;
  83. using Triplet = Eigen::Triplet<int>;
  84. const int* rows = block_jacobian_transpose.rows();
  85. const int* cols = block_jacobian_transpose.cols();
  86. int num_nonzeros = block_jacobian_transpose.num_nonzeros();
  87. std::vector<Triplet> triplets;
  88. triplets.reserve(num_nonzeros);
  89. for (int i = 0; i < num_nonzeros; ++i) {
  90. triplets.emplace_back(cols[i], rows[i], 1);
  91. }
  92. SparseMatrix block_jacobian(block_jacobian_transpose.num_cols(),
  93. block_jacobian_transpose.num_rows());
  94. block_jacobian.setFromTriplets(triplets.begin(), triplets.end());
  95. return block_jacobian;
  96. }
  97. void OrderingForSparseNormalCholeskyUsingSuiteSparse(
  98. const LinearSolverOrderingType linear_solver_ordering_type,
  99. const TripletSparseMatrix& tsm_block_jacobian_transpose,
  100. const std::vector<ParameterBlock*>& parameter_blocks,
  101. const ParameterBlockOrdering& parameter_block_ordering,
  102. int* ordering) {
  103. #ifdef CERES_NO_SUITESPARSE
  104. // "Void"ing values to avoid compiler warnings about unused parameters
  105. (void)linear_solver_ordering_type;
  106. (void)tsm_block_jacobian_transpose;
  107. (void)parameter_blocks;
  108. (void)parameter_block_ordering;
  109. (void)ordering;
  110. LOG(FATAL) << "Congratulations, you found a Ceres bug! "
  111. << "Please report this error to the developers.";
  112. #else
  113. SuiteSparse ss;
  114. cholmod_sparse* block_jacobian_transpose = ss.CreateSparseMatrix(
  115. const_cast<TripletSparseMatrix*>(&tsm_block_jacobian_transpose));
  116. if (linear_solver_ordering_type == ceres::AMD) {
  117. if (parameter_block_ordering.NumGroups() <= 1) {
  118. // The user did not supply a useful ordering so just go ahead
  119. // and use AMD.
  120. ss.Ordering(block_jacobian_transpose, OrderingType::AMD, ordering);
  121. } else {
  122. // The user supplied an ordering, so use CAMD.
  123. std::vector<int> constraints;
  124. constraints.reserve(parameter_blocks.size());
  125. for (auto* parameter_block : parameter_blocks) {
  126. constraints.push_back(parameter_block_ordering.GroupId(
  127. parameter_block->mutable_user_state()));
  128. }
  129. // Renumber the entries of constraints to be contiguous integers
  130. // as CAMD requires that the group ids be in the range [0,
  131. // parameter_blocks.size() - 1].
  132. MapValuesToContiguousRange(constraints.size(), constraints.data());
  133. ss.ConstrainedApproximateMinimumDegreeOrdering(
  134. block_jacobian_transpose, constraints.data(), ordering);
  135. }
  136. } else if (linear_solver_ordering_type == ceres::NESDIS) {
  137. // If nested dissection is chosen as an ordering algorithm, then
  138. // ignore any user provided linear_solver_ordering.
  139. CHECK(SuiteSparse::IsNestedDissectionAvailable())
  140. << "Congratulations, you found a Ceres bug! "
  141. << "Please report this error to the developers.";
  142. ss.Ordering(block_jacobian_transpose, OrderingType::NESDIS, ordering);
  143. } else {
  144. LOG(FATAL) << "Congratulations, you found a Ceres bug! "
  145. << "Please report this error to the developers.";
  146. }
  147. ss.Free(block_jacobian_transpose);
  148. #endif // CERES_NO_SUITESPARSE
  149. }
  150. void OrderingForSparseNormalCholeskyUsingEigenSparse(
  151. const LinearSolverOrderingType linear_solver_ordering_type,
  152. const TripletSparseMatrix& tsm_block_jacobian_transpose,
  153. int* ordering) {
  154. #ifndef CERES_USE_EIGEN_SPARSE
  155. LOG(FATAL) << "SPARSE_NORMAL_CHOLESKY cannot be used with EIGEN_SPARSE "
  156. "because Ceres was not built with support for "
  157. "Eigen's SimplicialLDLT decomposition. "
  158. "This requires enabling building with -DEIGENSPARSE=ON.";
  159. #else
  160. // TODO(sameeragarwal): This conversion from a TripletSparseMatrix
  161. // to a Eigen::Triplet matrix is unfortunate, but unavoidable for
  162. // now. It is not a significant performance penalty in the grand
  163. // scheme of things. The right thing to do here would be to get a
  164. // compressed row sparse matrix representation of the jacobian and
  165. // go from there. But that is a project for another day.
  166. using SparseMatrix = Eigen::SparseMatrix<int>;
  167. const SparseMatrix block_jacobian =
  168. CreateBlockJacobian(tsm_block_jacobian_transpose);
  169. const SparseMatrix block_hessian =
  170. block_jacobian.transpose() * block_jacobian;
  171. Eigen::PermutationMatrix<Eigen::Dynamic, Eigen::Dynamic, int> perm;
  172. if (linear_solver_ordering_type == ceres::AMD) {
  173. Eigen::AMDOrdering<int> amd_ordering;
  174. amd_ordering(block_hessian, perm);
  175. } else {
  176. #ifndef CERES_NO_EIGEN_METIS
  177. Eigen::MetisOrdering<int> metis_ordering;
  178. metis_ordering(block_hessian, perm);
  179. #else
  180. perm.setIdentity(block_hessian.rows());
  181. #endif
  182. }
  183. for (int i = 0; i < block_hessian.rows(); ++i) {
  184. ordering[i] = perm.indices()[i];
  185. }
  186. #endif // CERES_USE_EIGEN_SPARSE
  187. }
  188. } // namespace
  189. bool ApplyOrdering(const ProblemImpl::ParameterMap& parameter_map,
  190. const ParameterBlockOrdering& ordering,
  191. Program* program,
  192. std::string* error) {
  193. const int num_parameter_blocks = program->NumParameterBlocks();
  194. if (ordering.NumElements() != num_parameter_blocks) {
  195. *error = StringPrintf(
  196. "User specified ordering does not have the same "
  197. "number of parameters as the problem. The problem"
  198. "has %d blocks while the ordering has %d blocks.",
  199. num_parameter_blocks,
  200. ordering.NumElements());
  201. return false;
  202. }
  203. std::vector<ParameterBlock*>* parameter_blocks =
  204. program->mutable_parameter_blocks();
  205. parameter_blocks->clear();
  206. // TODO(sameeragarwal): Investigate whether this should be a set or an
  207. // unordered_set.
  208. const std::map<int, std::set<double*>>& groups = ordering.group_to_elements();
  209. for (const auto& p : groups) {
  210. const std::set<double*>& group = p.second;
  211. for (double* parameter_block_ptr : group) {
  212. auto it = parameter_map.find(parameter_block_ptr);
  213. if (it == parameter_map.end()) {
  214. *error = StringPrintf(
  215. "User specified ordering contains a pointer "
  216. "to a double that is not a parameter block in "
  217. "the problem. The invalid double is in group: %d",
  218. p.first);
  219. return false;
  220. }
  221. parameter_blocks->push_back(it->second);
  222. }
  223. }
  224. return true;
  225. }
  226. bool LexicographicallyOrderResidualBlocks(
  227. const int size_of_first_elimination_group,
  228. Program* program,
  229. std::string* /*error*/) {
  230. CHECK_GE(size_of_first_elimination_group, 1)
  231. << "Congratulations, you found a Ceres bug! Please report this error "
  232. << "to the developers.";
  233. // Create a histogram of the number of residuals for each E block. There is an
  234. // extra bucket at the end to catch all non-eliminated F blocks.
  235. std::vector<int> residual_blocks_per_e_block(size_of_first_elimination_group +
  236. 1);
  237. std::vector<ResidualBlock*>* residual_blocks =
  238. program->mutable_residual_blocks();
  239. std::vector<int> min_position_per_residual(residual_blocks->size());
  240. for (int i = 0; i < residual_blocks->size(); ++i) {
  241. ResidualBlock* residual_block = (*residual_blocks)[i];
  242. int position =
  243. MinParameterBlock(residual_block, size_of_first_elimination_group);
  244. min_position_per_residual[i] = position;
  245. DCHECK_LE(position, size_of_first_elimination_group);
  246. residual_blocks_per_e_block[position]++;
  247. }
  248. // Run a cumulative sum on the histogram, to obtain offsets to the start of
  249. // each histogram bucket (where each bucket is for the residuals for that
  250. // E-block).
  251. std::vector<int> offsets(size_of_first_elimination_group + 1);
  252. std::partial_sum(residual_blocks_per_e_block.begin(),
  253. residual_blocks_per_e_block.end(),
  254. offsets.begin());
  255. CHECK_EQ(offsets.back(), residual_blocks->size())
  256. << "Congratulations, you found a Ceres bug! Please report this error "
  257. << "to the developers.";
  258. CHECK(find(residual_blocks_per_e_block.begin(),
  259. residual_blocks_per_e_block.end() - 1,
  260. 0) == residual_blocks_per_e_block.end() - 1)
  261. << "Congratulations, you found a Ceres bug! Please report this error "
  262. << "to the developers.";
  263. // Fill in each bucket with the residual blocks for its corresponding E block.
  264. // Each bucket is individually filled from the back of the bucket to the front
  265. // of the bucket. The filling order among the buckets is dictated by the
  266. // residual blocks. This loop uses the offsets as counters; subtracting one
  267. // from each offset as a residual block is placed in the bucket. When the
  268. // filling is finished, the offset pointers should have shifted down one
  269. // entry (this is verified below).
  270. std::vector<ResidualBlock*> reordered_residual_blocks(
  271. (*residual_blocks).size(), static_cast<ResidualBlock*>(nullptr));
  272. for (int i = 0; i < residual_blocks->size(); ++i) {
  273. int bucket = min_position_per_residual[i];
  274. // Decrement the cursor, which should now point at the next empty position.
  275. offsets[bucket]--;
  276. // Sanity.
  277. CHECK(reordered_residual_blocks[offsets[bucket]] == nullptr)
  278. << "Congratulations, you found a Ceres bug! Please report this error "
  279. << "to the developers.";
  280. reordered_residual_blocks[offsets[bucket]] = (*residual_blocks)[i];
  281. }
  282. // Sanity check #1: The difference in bucket offsets should match the
  283. // histogram sizes.
  284. for (int i = 0; i < size_of_first_elimination_group; ++i) {
  285. CHECK_EQ(residual_blocks_per_e_block[i], offsets[i + 1] - offsets[i])
  286. << "Congratulations, you found a Ceres bug! Please report this error "
  287. << "to the developers.";
  288. }
  289. // Sanity check #2: No nullptr's left behind.
  290. for (auto* residual_block : reordered_residual_blocks) {
  291. CHECK(residual_block != nullptr)
  292. << "Congratulations, you found a Ceres bug! Please report this error "
  293. << "to the developers.";
  294. }
  295. // Now that the residuals are collected by E block, swap them in place.
  296. swap(*program->mutable_residual_blocks(), reordered_residual_blocks);
  297. return true;
  298. }
  299. // Pre-order the columns corresponding to the Schur complement if
  300. // possible.
  301. static void ReorderSchurComplementColumnsUsingSuiteSparse(
  302. const ParameterBlockOrdering& parameter_block_ordering, Program* program) {
  303. #ifdef CERES_NO_SUITESPARSE
  304. // "Void"ing values to avoid compiler warnings about unused parameters
  305. (void)parameter_block_ordering;
  306. (void)program;
  307. #else
  308. SuiteSparse ss;
  309. std::vector<int> constraints;
  310. std::vector<ParameterBlock*>& parameter_blocks =
  311. *(program->mutable_parameter_blocks());
  312. for (auto* parameter_block : parameter_blocks) {
  313. constraints.push_back(parameter_block_ordering.GroupId(
  314. parameter_block->mutable_user_state()));
  315. }
  316. // Renumber the entries of constraints to be contiguous integers as
  317. // CAMD requires that the group ids be in the range [0,
  318. // parameter_blocks.size() - 1].
  319. MapValuesToContiguousRange(constraints.size(), constraints.data());
  320. // Compute a block sparse presentation of J'.
  321. std::unique_ptr<TripletSparseMatrix> tsm_block_jacobian_transpose(
  322. program->CreateJacobianBlockSparsityTranspose());
  323. cholmod_sparse* block_jacobian_transpose =
  324. ss.CreateSparseMatrix(tsm_block_jacobian_transpose.get());
  325. std::vector<int> ordering(parameter_blocks.size(), 0);
  326. ss.ConstrainedApproximateMinimumDegreeOrdering(
  327. block_jacobian_transpose, constraints.data(), ordering.data());
  328. ss.Free(block_jacobian_transpose);
  329. const std::vector<ParameterBlock*> parameter_blocks_copy(parameter_blocks);
  330. for (int i = 0; i < program->NumParameterBlocks(); ++i) {
  331. parameter_blocks[i] = parameter_blocks_copy[ordering[i]];
  332. }
  333. program->SetParameterOffsetsAndIndex();
  334. #endif
  335. }
  336. static void ReorderSchurComplementColumnsUsingEigen(
  337. LinearSolverOrderingType ordering_type,
  338. const int size_of_first_elimination_group,
  339. const ProblemImpl::ParameterMap& /*parameter_map*/,
  340. Program* program) {
  341. #if defined(CERES_USE_EIGEN_SPARSE)
  342. std::unique_ptr<TripletSparseMatrix> tsm_block_jacobian_transpose(
  343. program->CreateJacobianBlockSparsityTranspose());
  344. using SparseMatrix = Eigen::SparseMatrix<int>;
  345. const SparseMatrix block_jacobian =
  346. CreateBlockJacobian(*tsm_block_jacobian_transpose);
  347. const int num_rows = block_jacobian.rows();
  348. const int num_cols = block_jacobian.cols();
  349. // Vertically partition the jacobian in parameter blocks of type E
  350. // and F.
  351. const SparseMatrix E =
  352. block_jacobian.block(0, 0, num_rows, size_of_first_elimination_group);
  353. const SparseMatrix F =
  354. block_jacobian.block(0,
  355. size_of_first_elimination_group,
  356. num_rows,
  357. num_cols - size_of_first_elimination_group);
  358. // Block sparsity pattern of the schur complement.
  359. const SparseMatrix block_schur_complement =
  360. F.transpose() * F - F.transpose() * E * E.transpose() * F;
  361. Eigen::PermutationMatrix<Eigen::Dynamic, Eigen::Dynamic, int> perm;
  362. if (ordering_type == ceres::AMD) {
  363. Eigen::AMDOrdering<int> amd_ordering;
  364. amd_ordering(block_schur_complement, perm);
  365. } else {
  366. #ifndef CERES_NO_EIGEN_METIS
  367. Eigen::MetisOrdering<int> metis_ordering;
  368. metis_ordering(block_schur_complement, perm);
  369. #else
  370. perm.setIdentity(block_schur_complement.rows());
  371. #endif
  372. }
  373. const std::vector<ParameterBlock*>& parameter_blocks =
  374. program->parameter_blocks();
  375. std::vector<ParameterBlock*> ordering(num_cols);
  376. // The ordering of the first size_of_first_elimination_group does
  377. // not matter, so we preserve the existing ordering.
  378. for (int i = 0; i < size_of_first_elimination_group; ++i) {
  379. ordering[i] = parameter_blocks[i];
  380. }
  381. // For the rest of the blocks, use the ordering computed using AMD.
  382. for (int i = 0; i < block_schur_complement.cols(); ++i) {
  383. ordering[size_of_first_elimination_group + i] =
  384. parameter_blocks[size_of_first_elimination_group + perm.indices()[i]];
  385. }
  386. swap(*program->mutable_parameter_blocks(), ordering);
  387. program->SetParameterOffsetsAndIndex();
  388. #endif
  389. }
  390. bool ReorderProgramForSchurTypeLinearSolver(
  391. const LinearSolverType linear_solver_type,
  392. const SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type,
  393. const LinearSolverOrderingType linear_solver_ordering_type,
  394. const ProblemImpl::ParameterMap& parameter_map,
  395. ParameterBlockOrdering* parameter_block_ordering,
  396. Program* program,
  397. std::string* error) {
  398. if (parameter_block_ordering->NumElements() !=
  399. program->NumParameterBlocks()) {
  400. *error = StringPrintf(
  401. "The program has %d parameter blocks, but the parameter block "
  402. "ordering has %d parameter blocks.",
  403. program->NumParameterBlocks(),
  404. parameter_block_ordering->NumElements());
  405. return false;
  406. }
  407. if (parameter_block_ordering->NumGroups() == 1) {
  408. // If the user supplied an parameter_block_ordering with just one
  409. // group, it is equivalent to the user supplying nullptr as an
  410. // parameter_block_ordering. Ceres is completely free to choose the
  411. // parameter block ordering as it sees fit. For Schur type solvers,
  412. // this means that the user wishes for Ceres to identify the
  413. // e_blocks, which we do by computing a maximal independent set.
  414. std::vector<ParameterBlock*> schur_ordering;
  415. const int size_of_first_elimination_group =
  416. ComputeStableSchurOrdering(*program, &schur_ordering);
  417. CHECK_EQ(schur_ordering.size(), program->NumParameterBlocks())
  418. << "Congratulations, you found a Ceres bug! Please report this error "
  419. << "to the developers.";
  420. // Update the parameter_block_ordering object.
  421. for (int i = 0; i < schur_ordering.size(); ++i) {
  422. double* parameter_block = schur_ordering[i]->mutable_user_state();
  423. const int group_id = (i < size_of_first_elimination_group) ? 0 : 1;
  424. parameter_block_ordering->AddElementToGroup(parameter_block, group_id);
  425. }
  426. // We could call ApplyOrdering but this is cheaper and
  427. // simpler.
  428. swap(*program->mutable_parameter_blocks(), schur_ordering);
  429. } else {
  430. // The user provided an ordering with more than one elimination
  431. // group.
  432. // Verify that the first elimination group is an independent set.
  433. // TODO(sameeragarwal): Investigate if this should be a set or an
  434. // unordered_set.
  435. const std::set<double*>& first_elimination_group =
  436. parameter_block_ordering->group_to_elements().begin()->second;
  437. if (!program->IsParameterBlockSetIndependent(first_elimination_group)) {
  438. *error = StringPrintf(
  439. "The first elimination group in the parameter block "
  440. "ordering of size %zd is not an independent set",
  441. first_elimination_group.size());
  442. return false;
  443. }
  444. if (!ApplyOrdering(
  445. parameter_map, *parameter_block_ordering, program, error)) {
  446. return false;
  447. }
  448. }
  449. program->SetParameterOffsetsAndIndex();
  450. const int size_of_first_elimination_group =
  451. parameter_block_ordering->group_to_elements().begin()->second.size();
  452. if (linear_solver_type == SPARSE_SCHUR) {
  453. if (sparse_linear_algebra_library_type == SUITE_SPARSE &&
  454. linear_solver_ordering_type == ceres::AMD) {
  455. // Preordering support for schur complement only works with AMD
  456. // for now, since we are using CAMD.
  457. //
  458. // TODO(sameeragarwal): It maybe worth adding pre-ordering support for
  459. // nested dissection too.
  460. ReorderSchurComplementColumnsUsingSuiteSparse(*parameter_block_ordering,
  461. program);
  462. } else if (sparse_linear_algebra_library_type == EIGEN_SPARSE) {
  463. ReorderSchurComplementColumnsUsingEigen(linear_solver_ordering_type,
  464. size_of_first_elimination_group,
  465. parameter_map,
  466. program);
  467. }
  468. }
  469. // Schur type solvers also require that their residual blocks be
  470. // lexicographically ordered.
  471. return LexicographicallyOrderResidualBlocks(
  472. size_of_first_elimination_group, program, error);
  473. }
  474. bool ReorderProgramForSparseCholesky(
  475. const SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type,
  476. const LinearSolverOrderingType linear_solver_ordering_type,
  477. const ParameterBlockOrdering& parameter_block_ordering,
  478. int start_row_block,
  479. Program* program,
  480. std::string* error) {
  481. if (parameter_block_ordering.NumElements() != program->NumParameterBlocks()) {
  482. *error = StringPrintf(
  483. "The program has %d parameter blocks, but the parameter block "
  484. "ordering has %d parameter blocks.",
  485. program->NumParameterBlocks(),
  486. parameter_block_ordering.NumElements());
  487. return false;
  488. }
  489. // Compute a block sparse presentation of J'.
  490. std::unique_ptr<TripletSparseMatrix> tsm_block_jacobian_transpose(
  491. program->CreateJacobianBlockSparsityTranspose(start_row_block));
  492. std::vector<int> ordering(program->NumParameterBlocks(), 0);
  493. std::vector<ParameterBlock*>& parameter_blocks =
  494. *(program->mutable_parameter_blocks());
  495. if (sparse_linear_algebra_library_type == SUITE_SPARSE) {
  496. OrderingForSparseNormalCholeskyUsingSuiteSparse(
  497. linear_solver_ordering_type,
  498. *tsm_block_jacobian_transpose,
  499. parameter_blocks,
  500. parameter_block_ordering,
  501. ordering.data());
  502. } else if (sparse_linear_algebra_library_type == ACCELERATE_SPARSE) {
  503. // Accelerate does not provide a function to perform reordering without
  504. // performing a full symbolic factorisation. As such, we have nothing
  505. // to gain from trying to reorder the problem here, as it will happen
  506. // in AppleAccelerateCholesky::Factorize() (once) and reordering here
  507. // would involve performing two symbolic factorisations instead of one
  508. // which would have a negative overall impact on performance.
  509. return true;
  510. } else if (sparse_linear_algebra_library_type == EIGEN_SPARSE) {
  511. OrderingForSparseNormalCholeskyUsingEigenSparse(
  512. linear_solver_ordering_type,
  513. *tsm_block_jacobian_transpose,
  514. ordering.data());
  515. }
  516. // Apply ordering.
  517. const std::vector<ParameterBlock*> parameter_blocks_copy(parameter_blocks);
  518. for (int i = 0; i < program->NumParameterBlocks(); ++i) {
  519. parameter_blocks[i] = parameter_blocks_copy[ordering[i]];
  520. }
  521. program->SetParameterOffsetsAndIndex();
  522. return true;
  523. }
  524. int ReorderResidualBlocksByPartition(
  525. const std::unordered_set<ResidualBlockId>& bottom_residual_blocks,
  526. Program* program) {
  527. auto residual_blocks = program->mutable_residual_blocks();
  528. auto it = std::partition(residual_blocks->begin(),
  529. residual_blocks->end(),
  530. [&bottom_residual_blocks](ResidualBlock* r) {
  531. return bottom_residual_blocks.count(r) == 0;
  532. });
  533. return it - residual_blocks->begin();
  534. }
  535. bool AreJacobianColumnsOrdered(
  536. const LinearSolverType linear_solver_type,
  537. const PreconditionerType preconditioner_type,
  538. const SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type,
  539. const LinearSolverOrderingType linear_solver_ordering_type) {
  540. if (sparse_linear_algebra_library_type == SUITE_SPARSE) {
  541. if (linear_solver_type == SPARSE_NORMAL_CHOLESKY ||
  542. (linear_solver_type == CGNR && preconditioner_type == SUBSET)) {
  543. return true;
  544. }
  545. if (linear_solver_type == SPARSE_SCHUR &&
  546. linear_solver_ordering_type == ceres::AMD) {
  547. return true;
  548. }
  549. return false;
  550. }
  551. if (sparse_linear_algebra_library_type == ceres::EIGEN_SPARSE) {
  552. if (linear_solver_type == SPARSE_NORMAL_CHOLESKY ||
  553. linear_solver_type == SPARSE_SCHUR ||
  554. (linear_solver_type == CGNR && preconditioner_type == SUBSET)) {
  555. return true;
  556. }
  557. return false;
  558. }
  559. if (sparse_linear_algebra_library_type == ceres::ACCELERATE_SPARSE) {
  560. // Apple's accelerate framework does not allow direct access to
  561. // ordering algorithms, so jacobian columns are never pre-ordered.
  562. return false;
  563. }
  564. return false;
  565. }
  566. } // namespace ceres::internal