trust_region_preprocessor.cc 16 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408
  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/trust_region_preprocessor.h"
  31. #include <numeric>
  32. #include <string>
  33. #include <vector>
  34. #include "ceres/callbacks.h"
  35. #include "ceres/context_impl.h"
  36. #include "ceres/evaluator.h"
  37. #include "ceres/linear_solver.h"
  38. #include "ceres/minimizer.h"
  39. #include "ceres/parameter_block.h"
  40. #include "ceres/preconditioner.h"
  41. #include "ceres/preprocessor.h"
  42. #include "ceres/problem_impl.h"
  43. #include "ceres/program.h"
  44. #include "ceres/reorder_program.h"
  45. #include "ceres/suitesparse.h"
  46. #include "ceres/trust_region_strategy.h"
  47. #include "ceres/wall_time.h"
  48. namespace ceres::internal {
  49. namespace {
  50. std::shared_ptr<ParameterBlockOrdering> CreateDefaultLinearSolverOrdering(
  51. const Program& program) {
  52. std::shared_ptr<ParameterBlockOrdering> ordering =
  53. std::make_shared<ParameterBlockOrdering>();
  54. const std::vector<ParameterBlock*>& parameter_blocks =
  55. program.parameter_blocks();
  56. for (auto* parameter_block : parameter_blocks) {
  57. ordering->AddElementToGroup(
  58. const_cast<double*>(parameter_block->user_state()), 0);
  59. }
  60. return ordering;
  61. }
  62. // Check if all the user supplied values in the parameter blocks are
  63. // sane or not, and if the program is feasible or not.
  64. bool IsProgramValid(const Program& program, std::string* error) {
  65. return (program.ParameterBlocksAreFinite(error) && program.IsFeasible(error));
  66. }
  67. void AlternateLinearSolverAndPreconditionerForSchurTypeLinearSolver(
  68. Solver::Options* options) {
  69. if (!IsSchurType(options->linear_solver_type)) {
  70. return;
  71. }
  72. const LinearSolverType linear_solver_type_given = options->linear_solver_type;
  73. const PreconditionerType preconditioner_type_given =
  74. options->preconditioner_type;
  75. options->linear_solver_type =
  76. LinearSolver::LinearSolverForZeroEBlocks(linear_solver_type_given);
  77. std::string message;
  78. if (linear_solver_type_given == ITERATIVE_SCHUR) {
  79. options->preconditioner_type =
  80. Preconditioner::PreconditionerForZeroEBlocks(preconditioner_type_given);
  81. message =
  82. StringPrintf("No E blocks. Switching from %s(%s) to %s(%s).",
  83. LinearSolverTypeToString(linear_solver_type_given),
  84. PreconditionerTypeToString(preconditioner_type_given),
  85. LinearSolverTypeToString(options->linear_solver_type),
  86. PreconditionerTypeToString(options->preconditioner_type));
  87. } else {
  88. message =
  89. StringPrintf("No E blocks. Switching from %s to %s.",
  90. LinearSolverTypeToString(linear_solver_type_given),
  91. LinearSolverTypeToString(options->linear_solver_type));
  92. }
  93. if (options->logging_type != SILENT) {
  94. VLOG(1) << message;
  95. }
  96. }
  97. // Reorder the program to reduce fill-in and increase cache coherency.
  98. bool ReorderProgram(PreprocessedProblem* pp) {
  99. const Solver::Options& options = pp->options;
  100. if (IsSchurType(options.linear_solver_type)) {
  101. return ReorderProgramForSchurTypeLinearSolver(
  102. options.linear_solver_type,
  103. options.sparse_linear_algebra_library_type,
  104. options.linear_solver_ordering_type,
  105. pp->problem->parameter_map(),
  106. options.linear_solver_ordering.get(),
  107. pp->reduced_program.get(),
  108. &pp->error);
  109. }
  110. if (options.linear_solver_type == SPARSE_NORMAL_CHOLESKY &&
  111. !options.dynamic_sparsity) {
  112. return ReorderProgramForSparseCholesky(
  113. options.sparse_linear_algebra_library_type,
  114. options.linear_solver_ordering_type,
  115. *options.linear_solver_ordering,
  116. 0, /* use all the rows of the jacobian */
  117. pp->reduced_program.get(),
  118. &pp->error);
  119. }
  120. if (options.linear_solver_type == CGNR &&
  121. options.preconditioner_type == SUBSET) {
  122. pp->linear_solver_options.subset_preconditioner_start_row_block =
  123. ReorderResidualBlocksByPartition(
  124. options.residual_blocks_for_subset_preconditioner,
  125. pp->reduced_program.get());
  126. return ReorderProgramForSparseCholesky(
  127. options.sparse_linear_algebra_library_type,
  128. options.linear_solver_ordering_type,
  129. *options.linear_solver_ordering,
  130. pp->linear_solver_options.subset_preconditioner_start_row_block,
  131. pp->reduced_program.get(),
  132. &pp->error);
  133. }
  134. return true;
  135. }
  136. // Configure and create a linear solver object. In doing so, if a
  137. // sparse direct factorization based linear solver is being used, then
  138. // find a fill reducing ordering and reorder the program as needed
  139. // too.
  140. bool SetupLinearSolver(PreprocessedProblem* pp) {
  141. Solver::Options& options = pp->options;
  142. pp->linear_solver_options = LinearSolver::Options();
  143. if (!options.linear_solver_ordering) {
  144. // If the user has not supplied a linear solver ordering, then we
  145. // assume that they are giving all the freedom to us in choosing
  146. // the best possible ordering. This intent can be indicated by
  147. // putting all the parameter blocks in the same elimination group.
  148. options.linear_solver_ordering =
  149. CreateDefaultLinearSolverOrdering(*pp->reduced_program);
  150. } else {
  151. // If the user supplied an ordering, then check if the first
  152. // elimination group is still non-empty after the reduced problem
  153. // has been constructed.
  154. //
  155. // This is important for Schur type linear solvers, where the
  156. // first elimination group is special -- it needs to be an
  157. // independent set.
  158. //
  159. // If the first elimination group is empty, then we cannot use the
  160. // user's requested linear solver (and a preconditioner as the
  161. // case may be) so we must use a different one.
  162. ParameterBlockOrdering* ordering = options.linear_solver_ordering.get();
  163. const int min_group_id = ordering->MinNonZeroGroup();
  164. ordering->Remove(pp->removed_parameter_blocks);
  165. if (IsSchurType(options.linear_solver_type) &&
  166. min_group_id != ordering->MinNonZeroGroup()) {
  167. AlternateLinearSolverAndPreconditionerForSchurTypeLinearSolver(&options);
  168. }
  169. }
  170. // Reorder the program to reduce fill in and improve cache coherency
  171. // of the Jacobian.
  172. if (!ReorderProgram(pp)) {
  173. return false;
  174. }
  175. // Configure the linear solver.
  176. pp->linear_solver_options.min_num_iterations =
  177. options.min_linear_solver_iterations;
  178. pp->linear_solver_options.max_num_iterations =
  179. options.max_linear_solver_iterations;
  180. pp->linear_solver_options.type = options.linear_solver_type;
  181. pp->linear_solver_options.preconditioner_type = options.preconditioner_type;
  182. pp->linear_solver_options.use_spse_initialization =
  183. options.use_spse_initialization;
  184. pp->linear_solver_options.spse_tolerance = options.spse_tolerance;
  185. pp->linear_solver_options.max_num_spse_iterations =
  186. options.max_num_spse_iterations;
  187. pp->linear_solver_options.visibility_clustering_type =
  188. options.visibility_clustering_type;
  189. pp->linear_solver_options.sparse_linear_algebra_library_type =
  190. options.sparse_linear_algebra_library_type;
  191. pp->linear_solver_options.dense_linear_algebra_library_type =
  192. options.dense_linear_algebra_library_type;
  193. pp->linear_solver_options.use_explicit_schur_complement =
  194. options.use_explicit_schur_complement;
  195. pp->linear_solver_options.dynamic_sparsity = options.dynamic_sparsity;
  196. pp->linear_solver_options.use_mixed_precision_solves =
  197. options.use_mixed_precision_solves;
  198. pp->linear_solver_options.max_num_refinement_iterations =
  199. options.max_num_refinement_iterations;
  200. pp->linear_solver_options.num_threads = options.num_threads;
  201. pp->linear_solver_options.context = pp->problem->context();
  202. if (IsSchurType(pp->linear_solver_options.type)) {
  203. OrderingToGroupSizes(options.linear_solver_ordering.get(),
  204. &pp->linear_solver_options.elimination_groups);
  205. // Schur type solvers expect at least two elimination groups. If
  206. // there is only one elimination group, then it is guaranteed that
  207. // this group only contains e_blocks. Thus we add a dummy
  208. // elimination group with zero blocks in it.
  209. if (pp->linear_solver_options.elimination_groups.size() == 1) {
  210. pp->linear_solver_options.elimination_groups.push_back(0);
  211. }
  212. }
  213. if (!options.dynamic_sparsity &&
  214. AreJacobianColumnsOrdered(options.linear_solver_type,
  215. options.preconditioner_type,
  216. options.sparse_linear_algebra_library_type,
  217. options.linear_solver_ordering_type)) {
  218. pp->linear_solver_options.ordering_type = OrderingType::NATURAL;
  219. } else {
  220. if (options.linear_solver_ordering_type == ceres::AMD) {
  221. pp->linear_solver_options.ordering_type = OrderingType::AMD;
  222. } else if (options.linear_solver_ordering_type == ceres::NESDIS) {
  223. pp->linear_solver_options.ordering_type = OrderingType::NESDIS;
  224. } else {
  225. LOG(FATAL) << "Congratulations you have found a bug in Ceres Solver."
  226. << " Please report this to the maintainers. : "
  227. << options.linear_solver_ordering_type;
  228. }
  229. }
  230. pp->linear_solver = LinearSolver::Create(pp->linear_solver_options);
  231. return (pp->linear_solver != nullptr);
  232. }
  233. // Configure and create the evaluator.
  234. bool SetupEvaluator(PreprocessedProblem* pp) {
  235. const Solver::Options& options = pp->options;
  236. pp->evaluator_options = Evaluator::Options();
  237. pp->evaluator_options.linear_solver_type = options.linear_solver_type;
  238. pp->evaluator_options.sparse_linear_algebra_library_type =
  239. options.sparse_linear_algebra_library_type;
  240. pp->evaluator_options.num_eliminate_blocks = 0;
  241. if (IsSchurType(options.linear_solver_type)) {
  242. pp->evaluator_options.num_eliminate_blocks =
  243. options.linear_solver_ordering->group_to_elements()
  244. .begin()
  245. ->second.size();
  246. }
  247. pp->evaluator_options.num_threads = options.num_threads;
  248. pp->evaluator_options.dynamic_sparsity = options.dynamic_sparsity;
  249. pp->evaluator_options.context = pp->problem->context();
  250. pp->evaluator_options.evaluation_callback =
  251. pp->reduced_program->mutable_evaluation_callback();
  252. pp->evaluator = Evaluator::Create(
  253. pp->evaluator_options, pp->reduced_program.get(), &pp->error);
  254. return (pp->evaluator != nullptr);
  255. }
  256. // If the user requested inner iterations, then find an inner
  257. // iteration ordering as needed and configure and create a
  258. // CoordinateDescentMinimizer object to perform the inner iterations.
  259. bool SetupInnerIterationMinimizer(PreprocessedProblem* pp) {
  260. Solver::Options& options = pp->options;
  261. if (!options.use_inner_iterations) {
  262. return true;
  263. }
  264. if (pp->reduced_program->mutable_evaluation_callback()) {
  265. pp->error = "Inner iterations cannot be used with EvaluationCallbacks";
  266. return false;
  267. }
  268. // With just one parameter block, the outer iteration of the trust
  269. // region method and inner iterations are doing exactly the same
  270. // thing, and thus inner iterations are not needed.
  271. if (pp->reduced_program->NumParameterBlocks() == 1) {
  272. LOG(WARNING) << "Reduced problem only contains one parameter block."
  273. << "Disabling inner iterations.";
  274. return true;
  275. }
  276. if (options.inner_iteration_ordering != nullptr) {
  277. // If the user supplied an ordering, then remove the set of
  278. // inactive parameter blocks from it
  279. options.inner_iteration_ordering->Remove(pp->removed_parameter_blocks);
  280. if (options.inner_iteration_ordering->NumElements() == 0) {
  281. LOG(WARNING) << "No remaining elements in the inner iteration ordering.";
  282. return true;
  283. }
  284. // Validate the reduced ordering.
  285. if (!CoordinateDescentMinimizer::IsOrderingValid(
  286. *pp->reduced_program,
  287. *options.inner_iteration_ordering,
  288. &pp->error)) {
  289. return false;
  290. }
  291. } else {
  292. // The user did not supply an ordering, so create one.
  293. options.inner_iteration_ordering =
  294. CoordinateDescentMinimizer::CreateOrdering(*pp->reduced_program);
  295. }
  296. pp->inner_iteration_minimizer =
  297. std::make_unique<CoordinateDescentMinimizer>(pp->problem->context());
  298. return pp->inner_iteration_minimizer->Init(*pp->reduced_program,
  299. pp->problem->parameter_map(),
  300. *options.inner_iteration_ordering,
  301. &pp->error);
  302. }
  303. // Configure and create a TrustRegionMinimizer object.
  304. bool SetupMinimizerOptions(PreprocessedProblem* pp) {
  305. const Solver::Options& options = pp->options;
  306. SetupCommonMinimizerOptions(pp);
  307. pp->minimizer_options.is_constrained =
  308. pp->reduced_program->IsBoundsConstrained();
  309. pp->minimizer_options.jacobian = pp->evaluator->CreateJacobian();
  310. if (pp->minimizer_options.jacobian == nullptr) {
  311. pp->error =
  312. "Unable to create Jacobian matrix. Likely because it is too large.";
  313. return false;
  314. }
  315. pp->minimizer_options.inner_iteration_minimizer =
  316. pp->inner_iteration_minimizer;
  317. TrustRegionStrategy::Options strategy_options;
  318. strategy_options.linear_solver = pp->linear_solver.get();
  319. strategy_options.initial_radius = options.initial_trust_region_radius;
  320. strategy_options.max_radius = options.max_trust_region_radius;
  321. strategy_options.min_lm_diagonal = options.min_lm_diagonal;
  322. strategy_options.max_lm_diagonal = options.max_lm_diagonal;
  323. strategy_options.trust_region_strategy_type =
  324. options.trust_region_strategy_type;
  325. strategy_options.dogleg_type = options.dogleg_type;
  326. strategy_options.context = pp->problem->context();
  327. strategy_options.num_threads = options.num_threads;
  328. pp->minimizer_options.trust_region_strategy =
  329. TrustRegionStrategy::Create(strategy_options);
  330. CHECK(pp->minimizer_options.trust_region_strategy != nullptr);
  331. return true;
  332. }
  333. } // namespace
  334. bool TrustRegionPreprocessor::Preprocess(const Solver::Options& options,
  335. ProblemImpl* problem,
  336. PreprocessedProblem* pp) {
  337. CHECK(pp != nullptr);
  338. pp->options = options;
  339. ChangeNumThreadsIfNeeded(&pp->options);
  340. pp->problem = problem;
  341. Program* program = problem->mutable_program();
  342. if (!IsProgramValid(*program, &pp->error)) {
  343. return false;
  344. }
  345. pp->reduced_program = program->CreateReducedProgram(
  346. &pp->removed_parameter_blocks, &pp->fixed_cost, &pp->error);
  347. if (pp->reduced_program.get() == nullptr) {
  348. return false;
  349. }
  350. if (pp->reduced_program->NumParameterBlocks() == 0) {
  351. // The reduced problem has no parameter or residual blocks. There
  352. // is nothing more to do.
  353. return true;
  354. }
  355. if (!SetupLinearSolver(pp) || !SetupEvaluator(pp) ||
  356. !SetupInnerIterationMinimizer(pp)) {
  357. return false;
  358. }
  359. return SetupMinimizerOptions(pp);
  360. }
  361. } // namespace ceres::internal