blob: dfde1221b61921d7991f5086ca5dcc5f61895208 [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2015 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: keir@google.com (Keir Mierle)
30// sameeragarwal@google.com (Sameer Agarwal)
31
32#include "ceres/solver.h"
33
34#include <algorithm>
35#include <memory>
36#include <sstream> // NOLINT
37#include <vector>
38
39#include "ceres/casts.h"
40#include "ceres/context.h"
41#include "ceres/context_impl.h"
42#include "ceres/detect_structure.h"
43#include "ceres/gradient_checking_cost_function.h"
44#include "ceres/internal/port.h"
45#include "ceres/parameter_block_ordering.h"
46#include "ceres/preprocessor.h"
47#include "ceres/problem.h"
48#include "ceres/problem_impl.h"
49#include "ceres/program.h"
50#include "ceres/schur_templates.h"
51#include "ceres/solver_utils.h"
52#include "ceres/stringprintf.h"
53#include "ceres/types.h"
54#include "ceres/wall_time.h"
55
56namespace ceres {
57namespace {
58
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080059using internal::StringAppendF;
60using internal::StringPrintf;
Austin Schuh70cc9552019-01-21 19:46:48 -080061using std::map;
62using std::string;
63using std::vector;
Austin Schuh70cc9552019-01-21 19:46:48 -080064
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080065#define OPTION_OP(x, y, OP) \
66 if (!(options.x OP y)) { \
67 std::stringstream ss; \
68 ss << "Invalid configuration. "; \
69 ss << string("Solver::Options::" #x " = ") << options.x << ". "; \
70 ss << "Violated constraint: "; \
71 ss << string("Solver::Options::" #x " " #OP " " #y); \
72 *error = ss.str(); \
73 return false; \
Austin Schuh70cc9552019-01-21 19:46:48 -080074 }
75
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080076#define OPTION_OP_OPTION(x, y, OP) \
77 if (!(options.x OP options.y)) { \
78 std::stringstream ss; \
79 ss << "Invalid configuration. "; \
80 ss << string("Solver::Options::" #x " = ") << options.x << ". "; \
81 ss << string("Solver::Options::" #y " = ") << options.y << ". "; \
82 ss << "Violated constraint: "; \
83 ss << string("Solver::Options::" #x); \
84 ss << string(#OP " Solver::Options::" #y "."); \
85 *error = ss.str(); \
86 return false; \
Austin Schuh70cc9552019-01-21 19:46:48 -080087 }
88
89#define OPTION_GE(x, y) OPTION_OP(x, y, >=);
90#define OPTION_GT(x, y) OPTION_OP(x, y, >);
91#define OPTION_LE(x, y) OPTION_OP(x, y, <=);
92#define OPTION_LT(x, y) OPTION_OP(x, y, <);
93#define OPTION_LE_OPTION(x, y) OPTION_OP_OPTION(x, y, <=)
94#define OPTION_LT_OPTION(x, y) OPTION_OP_OPTION(x, y, <)
95
96bool CommonOptionsAreValid(const Solver::Options& options, string* error) {
97 OPTION_GE(max_num_iterations, 0);
98 OPTION_GE(max_solver_time_in_seconds, 0.0);
99 OPTION_GE(function_tolerance, 0.0);
100 OPTION_GE(gradient_tolerance, 0.0);
101 OPTION_GE(parameter_tolerance, 0.0);
102 OPTION_GT(num_threads, 0);
103 if (options.check_gradients) {
104 OPTION_GT(gradient_check_relative_precision, 0.0);
105 OPTION_GT(gradient_check_numeric_derivative_relative_step_size, 0.0);
106 }
107 return true;
108}
109
110bool TrustRegionOptionsAreValid(const Solver::Options& options, string* error) {
111 OPTION_GT(initial_trust_region_radius, 0.0);
112 OPTION_GT(min_trust_region_radius, 0.0);
113 OPTION_GT(max_trust_region_radius, 0.0);
114 OPTION_LE_OPTION(min_trust_region_radius, max_trust_region_radius);
115 OPTION_LE_OPTION(min_trust_region_radius, initial_trust_region_radius);
116 OPTION_LE_OPTION(initial_trust_region_radius, max_trust_region_radius);
117 OPTION_GE(min_relative_decrease, 0.0);
118 OPTION_GE(min_lm_diagonal, 0.0);
119 OPTION_GE(max_lm_diagonal, 0.0);
120 OPTION_LE_OPTION(min_lm_diagonal, max_lm_diagonal);
121 OPTION_GE(max_num_consecutive_invalid_steps, 0);
122 OPTION_GT(eta, 0.0);
123 OPTION_GE(min_linear_solver_iterations, 0);
124 OPTION_GE(max_linear_solver_iterations, 1);
125 OPTION_LE_OPTION(min_linear_solver_iterations, max_linear_solver_iterations);
126
127 if (options.use_inner_iterations) {
128 OPTION_GE(inner_iteration_tolerance, 0.0);
129 }
130
Austin Schuh70cc9552019-01-21 19:46:48 -0800131 if (options.use_nonmonotonic_steps) {
132 OPTION_GT(max_consecutive_nonmonotonic_steps, 0);
133 }
134
135 if (options.linear_solver_type == ITERATIVE_SCHUR &&
136 options.use_explicit_schur_complement &&
137 options.preconditioner_type != SCHUR_JACOBI) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800138 *error =
139 "use_explicit_schur_complement only supports "
Austin Schuh70cc9552019-01-21 19:46:48 -0800140 "SCHUR_JACOBI as the preconditioner.";
141 return false;
142 }
143
144 if (options.dense_linear_algebra_library_type == LAPACK &&
145 !IsDenseLinearAlgebraLibraryTypeAvailable(LAPACK) &&
146 (options.linear_solver_type == DENSE_NORMAL_CHOLESKY ||
147 options.linear_solver_type == DENSE_QR ||
148 options.linear_solver_type == DENSE_SCHUR)) {
149 *error = StringPrintf(
150 "Can't use %s with "
151 "Solver::Options::dense_linear_algebra_library_type = LAPACK "
152 "because LAPACK was not enabled when Ceres was built.",
153 LinearSolverTypeToString(options.linear_solver_type));
154 return false;
155 }
156
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800157 {
158 const char* sparse_linear_algebra_library_name =
159 SparseLinearAlgebraLibraryTypeToString(
160 options.sparse_linear_algebra_library_type);
Austin Schuh70cc9552019-01-21 19:46:48 -0800161 const char* name = nullptr;
162 if (options.linear_solver_type == SPARSE_NORMAL_CHOLESKY ||
163 options.linear_solver_type == SPARSE_SCHUR) {
164 name = LinearSolverTypeToString(options.linear_solver_type);
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800165 } else if ((options.linear_solver_type == ITERATIVE_SCHUR &&
166 (options.preconditioner_type == CLUSTER_JACOBI ||
167 options.preconditioner_type == CLUSTER_TRIDIAGONAL)) ||
168 (options.linear_solver_type == CGNR &&
169 options.preconditioner_type == SUBSET)) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800170 name = PreconditionerTypeToString(options.preconditioner_type);
171 }
172
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800173 if (name) {
174 if (options.sparse_linear_algebra_library_type == NO_SPARSE) {
175 *error = StringPrintf(
176 "Can't use %s with "
177 "Solver::Options::sparse_linear_algebra_library_type = %s.",
178 name,
179 sparse_linear_algebra_library_name);
180 return false;
181 } else if (!IsSparseLinearAlgebraLibraryTypeAvailable(
182 options.sparse_linear_algebra_library_type)) {
183 *error = StringPrintf(
184 "Can't use %s with "
185 "Solver::Options::sparse_linear_algebra_library_type = %s, "
186 "because support was not enabled when Ceres Solver was built.",
187 name,
188 sparse_linear_algebra_library_name);
189 return false;
190 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800191 }
192 }
193
194 if (options.trust_region_strategy_type == DOGLEG) {
195 if (options.linear_solver_type == ITERATIVE_SCHUR ||
196 options.linear_solver_type == CGNR) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800197 *error =
198 "DOGLEG only supports exact factorization based linear "
Austin Schuh70cc9552019-01-21 19:46:48 -0800199 "solvers. If you want to use an iterative solver please "
200 "use LEVENBERG_MARQUARDT as the trust_region_strategy_type";
201 return false;
202 }
203 }
204
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800205 if (!options.trust_region_minimizer_iterations_to_dump.empty() &&
Austin Schuh70cc9552019-01-21 19:46:48 -0800206 options.trust_region_problem_dump_format_type != CONSOLE &&
207 options.trust_region_problem_dump_directory.empty()) {
208 *error = "Solver::Options::trust_region_problem_dump_directory is empty.";
209 return false;
210 }
211
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800212 if (options.dynamic_sparsity) {
213 if (options.linear_solver_type != SPARSE_NORMAL_CHOLESKY) {
214 *error =
215 "Dynamic sparsity is only supported with SPARSE_NORMAL_CHOLESKY.";
216 return false;
217 }
218 if (options.sparse_linear_algebra_library_type == ACCELERATE_SPARSE) {
219 *error =
220 "ACCELERATE_SPARSE is not currently supported with dynamic sparsity.";
221 return false;
222 }
223 }
224
225 if (options.linear_solver_type == CGNR &&
226 options.preconditioner_type == SUBSET &&
227 options.residual_blocks_for_subset_preconditioner.empty()) {
228 *error =
229 "When using SUBSET preconditioner, "
230 "Solver::Options::residual_blocks_for_subset_preconditioner cannot be "
231 "empty";
Austin Schuh70cc9552019-01-21 19:46:48 -0800232 return false;
233 }
234
235 return true;
236}
237
238bool LineSearchOptionsAreValid(const Solver::Options& options, string* error) {
239 OPTION_GT(max_lbfgs_rank, 0);
240 OPTION_GT(min_line_search_step_size, 0.0);
241 OPTION_GT(max_line_search_step_contraction, 0.0);
242 OPTION_LT(max_line_search_step_contraction, 1.0);
243 OPTION_LT_OPTION(max_line_search_step_contraction,
244 min_line_search_step_contraction);
245 OPTION_LE(min_line_search_step_contraction, 1.0);
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800246 OPTION_GE(max_num_line_search_step_size_iterations,
247 (options.minimizer_type == ceres::TRUST_REGION ? 0 : 1));
Austin Schuh70cc9552019-01-21 19:46:48 -0800248 OPTION_GT(line_search_sufficient_function_decrease, 0.0);
249 OPTION_LT_OPTION(line_search_sufficient_function_decrease,
250 line_search_sufficient_curvature_decrease);
251 OPTION_LT(line_search_sufficient_curvature_decrease, 1.0);
252 OPTION_GT(max_line_search_step_expansion, 1.0);
253
254 if ((options.line_search_direction_type == ceres::BFGS ||
255 options.line_search_direction_type == ceres::LBFGS) &&
256 options.line_search_type != ceres::WOLFE) {
257 *error =
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800258 string("Invalid configuration: Solver::Options::line_search_type = ") +
259 string(LineSearchTypeToString(options.line_search_type)) +
260 string(
261 ". When using (L)BFGS, "
262 "Solver::Options::line_search_type must be set to WOLFE.");
Austin Schuh70cc9552019-01-21 19:46:48 -0800263 return false;
264 }
265
266 // Warn user if they have requested BISECTION interpolation, but constraints
267 // on max/min step size change during line search prevent bisection scaling
268 // from occurring. Warn only, as this is likely a user mistake, but one which
269 // does not prevent us from continuing.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800270 if (options.line_search_interpolation_type == ceres::BISECTION &&
271 (options.max_line_search_step_contraction > 0.5 ||
272 options.min_line_search_step_contraction < 0.5)) {
273 LOG(WARNING)
274 << "Line search interpolation type is BISECTION, but specified "
275 << "max_line_search_step_contraction: "
276 << options.max_line_search_step_contraction << ", and "
277 << "min_line_search_step_contraction: "
278 << options.min_line_search_step_contraction
279 << ", prevent bisection (0.5) scaling, continuing with solve "
280 "regardless.";
281 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800282 return true;
283}
284
285#undef OPTION_OP
286#undef OPTION_OP_OPTION
287#undef OPTION_GT
288#undef OPTION_GE
289#undef OPTION_LE
290#undef OPTION_LT
291#undef OPTION_LE_OPTION
292#undef OPTION_LT_OPTION
293
294void StringifyOrdering(const vector<int>& ordering, string* report) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800295 if (ordering.empty()) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800296 internal::StringAppendF(report, "AUTOMATIC");
297 return;
298 }
299
300 for (int i = 0; i < ordering.size() - 1; ++i) {
301 internal::StringAppendF(report, "%d,", ordering[i]);
302 }
303 internal::StringAppendF(report, "%d", ordering.back());
304}
305
306void SummarizeGivenProgram(const internal::Program& program,
307 Solver::Summary* summary) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800308 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800309 summary->num_parameter_blocks = program.NumParameterBlocks();
310 summary->num_parameters = program.NumParameters();
311 summary->num_effective_parameters = program.NumEffectiveParameters();
312 summary->num_residual_blocks = program.NumResidualBlocks();
313 summary->num_residuals = program.NumResiduals();
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800314 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800315}
316
317void SummarizeReducedProgram(const internal::Program& program,
318 Solver::Summary* summary) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800319 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800320 summary->num_parameter_blocks_reduced = program.NumParameterBlocks();
321 summary->num_parameters_reduced = program.NumParameters();
322 summary->num_effective_parameters_reduced = program.NumEffectiveParameters();
323 summary->num_residual_blocks_reduced = program.NumResidualBlocks();
324 summary->num_residuals_reduced = program.NumResiduals();
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800325 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800326}
327
328void PreSolveSummarize(const Solver::Options& options,
329 const internal::ProblemImpl* problem,
330 Solver::Summary* summary) {
331 SummarizeGivenProgram(problem->program(), summary);
332 internal::OrderingToGroupSizes(options.linear_solver_ordering.get(),
333 &(summary->linear_solver_ordering_given));
334 internal::OrderingToGroupSizes(options.inner_iteration_ordering.get(),
335 &(summary->inner_iteration_ordering_given));
336
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800337 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800338 summary->dense_linear_algebra_library_type = options.dense_linear_algebra_library_type; // NOLINT
339 summary->dogleg_type = options.dogleg_type;
340 summary->inner_iteration_time_in_seconds = 0.0;
341 summary->num_line_search_steps = 0;
342 summary->line_search_cost_evaluation_time_in_seconds = 0.0;
343 summary->line_search_gradient_evaluation_time_in_seconds = 0.0;
344 summary->line_search_polynomial_minimization_time_in_seconds = 0.0;
345 summary->line_search_total_time_in_seconds = 0.0;
346 summary->inner_iterations_given = options.use_inner_iterations;
347 summary->line_search_direction_type = options.line_search_direction_type; // NOLINT
348 summary->line_search_interpolation_type = options.line_search_interpolation_type; // NOLINT
349 summary->line_search_type = options.line_search_type;
350 summary->linear_solver_type_given = options.linear_solver_type;
351 summary->max_lbfgs_rank = options.max_lbfgs_rank;
352 summary->minimizer_type = options.minimizer_type;
353 summary->nonlinear_conjugate_gradient_type = options.nonlinear_conjugate_gradient_type; // NOLINT
354 summary->num_threads_given = options.num_threads;
355 summary->preconditioner_type_given = options.preconditioner_type;
356 summary->sparse_linear_algebra_library_type = options.sparse_linear_algebra_library_type; // NOLINT
357 summary->trust_region_strategy_type = options.trust_region_strategy_type; // NOLINT
358 summary->visibility_clustering_type = options.visibility_clustering_type; // NOLINT
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800359 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800360}
361
362void PostSolveSummarize(const internal::PreprocessedProblem& pp,
363 Solver::Summary* summary) {
364 internal::OrderingToGroupSizes(pp.options.linear_solver_ordering.get(),
365 &(summary->linear_solver_ordering_used));
366 internal::OrderingToGroupSizes(pp.options.inner_iteration_ordering.get(),
367 &(summary->inner_iteration_ordering_used));
368
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800369 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800370 summary->inner_iterations_used = pp.inner_iteration_minimizer.get() != NULL; // NOLINT
371 summary->linear_solver_type_used = pp.linear_solver_options.type;
372 summary->num_threads_used = pp.options.num_threads;
373 summary->preconditioner_type_used = pp.options.preconditioner_type;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800374 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800375
376 internal::SetSummaryFinalCost(summary);
377
378 if (pp.reduced_program.get() != NULL) {
379 SummarizeReducedProgram(*pp.reduced_program, summary);
380 }
381
382 using internal::CallStatistics;
383
384 // It is possible that no evaluator was created. This would be the
385 // case if the preprocessor failed, or if the reduced problem did
386 // not contain any parameter blocks. Thus, only extract the
387 // evaluator statistics if one exists.
388 if (pp.evaluator.get() != NULL) {
389 const map<string, CallStatistics>& evaluator_statistics =
390 pp.evaluator->Statistics();
391 {
392 const CallStatistics& call_stats = FindWithDefault(
393 evaluator_statistics, "Evaluator::Residual", CallStatistics());
394
395 summary->residual_evaluation_time_in_seconds = call_stats.time;
396 summary->num_residual_evaluations = call_stats.calls;
397 }
398 {
399 const CallStatistics& call_stats = FindWithDefault(
400 evaluator_statistics, "Evaluator::Jacobian", CallStatistics());
401
402 summary->jacobian_evaluation_time_in_seconds = call_stats.time;
403 summary->num_jacobian_evaluations = call_stats.calls;
404 }
405 }
406
407 // Again, like the evaluator, there may or may not be a linear
408 // solver from which we can extract run time statistics. In
409 // particular the line search solver does not use a linear solver.
410 if (pp.linear_solver.get() != NULL) {
411 const map<string, CallStatistics>& linear_solver_statistics =
412 pp.linear_solver->Statistics();
413 const CallStatistics& call_stats = FindWithDefault(
414 linear_solver_statistics, "LinearSolver::Solve", CallStatistics());
415 summary->num_linear_solves = call_stats.calls;
416 summary->linear_solver_time_in_seconds = call_stats.time;
417 }
418}
419
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800420void Minimize(internal::PreprocessedProblem* pp, Solver::Summary* summary) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800421 using internal::Minimizer;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800422 using internal::Program;
Austin Schuh70cc9552019-01-21 19:46:48 -0800423
424 Program* program = pp->reduced_program.get();
425 if (pp->reduced_program->NumParameterBlocks() == 0) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800426 summary->message =
427 "Function tolerance reached. "
Austin Schuh70cc9552019-01-21 19:46:48 -0800428 "No non-constant parameter blocks found.";
429 summary->termination_type = CONVERGENCE;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800430 if (pp->options.logging_type != SILENT) {
431 VLOG(1) << summary->message;
432 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800433 summary->initial_cost = summary->fixed_cost;
434 summary->final_cost = summary->fixed_cost;
435 return;
436 }
437
438 const Vector original_reduced_parameters = pp->reduced_parameters;
439 std::unique_ptr<Minimizer> minimizer(
440 Minimizer::Create(pp->options.minimizer_type));
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800441 minimizer->Minimize(
442 pp->minimizer_options, pp->reduced_parameters.data(), summary);
Austin Schuh70cc9552019-01-21 19:46:48 -0800443
444 program->StateVectorToParameterBlocks(
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800445 summary->IsSolutionUsable() ? pp->reduced_parameters.data()
446 : original_reduced_parameters.data());
Austin Schuh70cc9552019-01-21 19:46:48 -0800447 program->CopyParameterBlockStateToUserState();
448}
449
450std::string SchurStructureToString(const int row_block_size,
451 const int e_block_size,
452 const int f_block_size) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800453 const std::string row = (row_block_size == Eigen::Dynamic)
454 ? "d"
455 : internal::StringPrintf("%d", row_block_size);
Austin Schuh70cc9552019-01-21 19:46:48 -0800456
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800457 const std::string e = (e_block_size == Eigen::Dynamic)
458 ? "d"
459 : internal::StringPrintf("%d", e_block_size);
Austin Schuh70cc9552019-01-21 19:46:48 -0800460
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800461 const std::string f = (f_block_size == Eigen::Dynamic)
462 ? "d"
463 : internal::StringPrintf("%d", f_block_size);
Austin Schuh70cc9552019-01-21 19:46:48 -0800464
465 return internal::StringPrintf("%s,%s,%s", row.c_str(), e.c_str(), f.c_str());
466}
467
468} // namespace
469
470bool Solver::Options::IsValid(string* error) const {
471 if (!CommonOptionsAreValid(*this, error)) {
472 return false;
473 }
474
475 if (minimizer_type == TRUST_REGION &&
476 !TrustRegionOptionsAreValid(*this, error)) {
477 return false;
478 }
479
480 // We do not know if the problem is bounds constrained or not, if it
481 // is then the trust region solver will also use the line search
482 // solver to do a projection onto the box constraints, so make sure
483 // that the line search options are checked independent of what
484 // minimizer algorithm is being used.
485 return LineSearchOptionsAreValid(*this, error);
486}
487
488Solver::~Solver() {}
489
490void Solver::Solve(const Solver::Options& options,
491 Problem* problem,
492 Solver::Summary* summary) {
493 using internal::PreprocessedProblem;
494 using internal::Preprocessor;
495 using internal::ProblemImpl;
496 using internal::Program;
497 using internal::WallTimeInSeconds;
498
499 CHECK(problem != nullptr);
500 CHECK(summary != nullptr);
501
502 double start_time = WallTimeInSeconds();
503 *summary = Summary();
504 if (!options.IsValid(&summary->message)) {
505 LOG(ERROR) << "Terminating: " << summary->message;
506 return;
507 }
508
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800509 ProblemImpl* problem_impl = problem->impl_.get();
Austin Schuh70cc9552019-01-21 19:46:48 -0800510 Program* program = problem_impl->mutable_program();
511 PreSolveSummarize(options, problem_impl, summary);
512
513 // If gradient_checking is enabled, wrap all cost functions in a
514 // gradient checker and install a callback that terminates if any gradient
515 // error is detected.
516 std::unique_ptr<internal::ProblemImpl> gradient_checking_problem;
517 internal::GradientCheckingIterationCallback gradient_checking_callback;
518 Solver::Options modified_options = options;
519 if (options.check_gradients) {
520 modified_options.callbacks.push_back(&gradient_checking_callback);
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800521 gradient_checking_problem.reset(CreateGradientCheckingProblemImpl(
522 problem_impl,
523 options.gradient_check_numeric_derivative_relative_step_size,
524 options.gradient_check_relative_precision,
525 &gradient_checking_callback));
Austin Schuh70cc9552019-01-21 19:46:48 -0800526 problem_impl = gradient_checking_problem.get();
527 program = problem_impl->mutable_program();
528 }
529
530 // Make sure that all the parameter blocks states are set to the
531 // values provided by the user.
532 program->SetParameterBlockStatePtrsToUserStatePtrs();
533
534 // The main thread also does work so we only need to launch num_threads - 1.
535 problem_impl->context()->EnsureMinimumThreads(options.num_threads - 1);
536
537 std::unique_ptr<Preprocessor> preprocessor(
538 Preprocessor::Create(modified_options.minimizer_type));
539 PreprocessedProblem pp;
540
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800541 const bool status =
542 preprocessor->Preprocess(modified_options, problem_impl, &pp);
Austin Schuh70cc9552019-01-21 19:46:48 -0800543
544 // We check the linear_solver_options.type rather than
545 // modified_options.linear_solver_type because, depending on the
546 // lack of a Schur structure, the preprocessor may change the linear
547 // solver type.
548 if (IsSchurType(pp.linear_solver_options.type)) {
549 // TODO(sameeragarwal): We can likely eliminate the duplicate call
550 // to DetectStructure here and inside the linear solver, by
551 // calling this in the preprocessor.
552 int row_block_size;
553 int e_block_size;
554 int f_block_size;
555 DetectStructure(*static_cast<internal::BlockSparseMatrix*>(
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800556 pp.minimizer_options.jacobian.get())
557 ->block_structure(),
Austin Schuh70cc9552019-01-21 19:46:48 -0800558 pp.linear_solver_options.elimination_groups[0],
559 &row_block_size,
560 &e_block_size,
561 &f_block_size);
562 summary->schur_structure_given =
563 SchurStructureToString(row_block_size, e_block_size, f_block_size);
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800564 internal::GetBestSchurTemplateSpecialization(
565 &row_block_size, &e_block_size, &f_block_size);
Austin Schuh70cc9552019-01-21 19:46:48 -0800566 summary->schur_structure_used =
567 SchurStructureToString(row_block_size, e_block_size, f_block_size);
568 }
569
570 summary->fixed_cost = pp.fixed_cost;
571 summary->preprocessor_time_in_seconds = WallTimeInSeconds() - start_time;
572
573 if (status) {
574 const double minimizer_start_time = WallTimeInSeconds();
575 Minimize(&pp, summary);
576 summary->minimizer_time_in_seconds =
577 WallTimeInSeconds() - minimizer_start_time;
578 } else {
579 summary->message = pp.error;
580 }
581
582 const double postprocessor_start_time = WallTimeInSeconds();
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800583 problem_impl = problem->impl_.get();
Austin Schuh70cc9552019-01-21 19:46:48 -0800584 program = problem_impl->mutable_program();
585 // On exit, ensure that the parameter blocks again point at the user
586 // provided values and the parameter blocks are numbered according
587 // to their position in the original user provided program.
588 program->SetParameterBlockStatePtrsToUserStatePtrs();
589 program->SetParameterOffsetsAndIndex();
590 PostSolveSummarize(pp, summary);
591 summary->postprocessor_time_in_seconds =
592 WallTimeInSeconds() - postprocessor_start_time;
593
594 // If the gradient checker reported an error, we want to report FAILURE
595 // instead of USER_FAILURE and provide the error log.
596 if (gradient_checking_callback.gradient_error_detected()) {
597 summary->termination_type = FAILURE;
598 summary->message = gradient_checking_callback.error_log();
599 }
600
601 summary->total_time_in_seconds = WallTimeInSeconds() - start_time;
602}
603
604void Solve(const Solver::Options& options,
605 Problem* problem,
606 Solver::Summary* summary) {
607 Solver solver;
608 solver.Solve(options, problem, summary);
609}
610
611string Solver::Summary::BriefReport() const {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800612 return StringPrintf(
613 "Ceres Solver Report: "
614 "Iterations: %d, "
615 "Initial cost: %e, "
616 "Final cost: %e, "
617 "Termination: %s",
618 num_successful_steps + num_unsuccessful_steps,
619 initial_cost,
620 final_cost,
621 TerminationTypeToString(termination_type));
Austin Schuh70cc9552019-01-21 19:46:48 -0800622}
623
624string Solver::Summary::FullReport() const {
625 using internal::VersionString;
626
627 string report = string("\nSolver Summary (v " + VersionString() + ")\n\n");
628
629 StringAppendF(&report, "%45s %21s\n", "Original", "Reduced");
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800630 StringAppendF(&report,
631 "Parameter blocks % 25d% 25d\n",
632 num_parameter_blocks,
633 num_parameter_blocks_reduced);
634 StringAppendF(&report,
635 "Parameters % 25d% 25d\n",
636 num_parameters,
637 num_parameters_reduced);
Austin Schuh70cc9552019-01-21 19:46:48 -0800638 if (num_effective_parameters_reduced != num_parameters_reduced) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800639 StringAppendF(&report,
640 "Effective parameters% 25d% 25d\n",
641 num_effective_parameters,
642 num_effective_parameters_reduced);
Austin Schuh70cc9552019-01-21 19:46:48 -0800643 }
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800644 StringAppendF(&report,
645 "Residual blocks % 25d% 25d\n",
646 num_residual_blocks,
647 num_residual_blocks_reduced);
648 StringAppendF(&report,
649 "Residuals % 25d% 25d\n",
650 num_residuals,
651 num_residuals_reduced);
Austin Schuh70cc9552019-01-21 19:46:48 -0800652
653 if (minimizer_type == TRUST_REGION) {
654 // TRUST_SEARCH HEADER
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800655 StringAppendF(
656 &report, "\nMinimizer %19s\n", "TRUST_REGION");
Austin Schuh70cc9552019-01-21 19:46:48 -0800657
658 if (linear_solver_type_used == DENSE_NORMAL_CHOLESKY ||
659 linear_solver_type_used == DENSE_SCHUR ||
660 linear_solver_type_used == DENSE_QR) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800661 StringAppendF(&report,
662 "\nDense linear algebra library %15s\n",
Austin Schuh70cc9552019-01-21 19:46:48 -0800663 DenseLinearAlgebraLibraryTypeToString(
664 dense_linear_algebra_library_type));
665 }
666
667 if (linear_solver_type_used == SPARSE_NORMAL_CHOLESKY ||
668 linear_solver_type_used == SPARSE_SCHUR ||
669 (linear_solver_type_used == ITERATIVE_SCHUR &&
670 (preconditioner_type_used == CLUSTER_JACOBI ||
671 preconditioner_type_used == CLUSTER_TRIDIAGONAL))) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800672 StringAppendF(&report,
673 "\nSparse linear algebra library %15s\n",
Austin Schuh70cc9552019-01-21 19:46:48 -0800674 SparseLinearAlgebraLibraryTypeToString(
675 sparse_linear_algebra_library_type));
676 }
677
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800678 StringAppendF(&report,
679 "Trust region strategy %19s",
680 TrustRegionStrategyTypeToString(trust_region_strategy_type));
Austin Schuh70cc9552019-01-21 19:46:48 -0800681 if (trust_region_strategy_type == DOGLEG) {
682 if (dogleg_type == TRADITIONAL_DOGLEG) {
683 StringAppendF(&report, " (TRADITIONAL)");
684 } else {
685 StringAppendF(&report, " (SUBSPACE)");
686 }
687 }
688 StringAppendF(&report, "\n");
689 StringAppendF(&report, "\n");
690
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800691 StringAppendF(&report, "%45s %21s\n", "Given", "Used");
692 StringAppendF(&report,
693 "Linear solver %25s%25s\n",
Austin Schuh70cc9552019-01-21 19:46:48 -0800694 LinearSolverTypeToString(linear_solver_type_given),
695 LinearSolverTypeToString(linear_solver_type_used));
696
697 if (linear_solver_type_given == CGNR ||
698 linear_solver_type_given == ITERATIVE_SCHUR) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800699 StringAppendF(&report,
700 "Preconditioner %25s%25s\n",
Austin Schuh70cc9552019-01-21 19:46:48 -0800701 PreconditionerTypeToString(preconditioner_type_given),
702 PreconditionerTypeToString(preconditioner_type_used));
703 }
704
705 if (preconditioner_type_used == CLUSTER_JACOBI ||
706 preconditioner_type_used == CLUSTER_TRIDIAGONAL) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800707 StringAppendF(
708 &report,
709 "Visibility clustering%24s%25s\n",
710 VisibilityClusteringTypeToString(visibility_clustering_type),
711 VisibilityClusteringTypeToString(visibility_clustering_type));
Austin Schuh70cc9552019-01-21 19:46:48 -0800712 }
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800713 StringAppendF(&report,
714 "Threads % 25d% 25d\n",
715 num_threads_given,
716 num_threads_used);
Austin Schuh70cc9552019-01-21 19:46:48 -0800717
718 string given;
719 StringifyOrdering(linear_solver_ordering_given, &given);
720 string used;
721 StringifyOrdering(linear_solver_ordering_used, &used);
722 StringAppendF(&report,
723 "Linear solver ordering %22s %24s\n",
724 given.c_str(),
725 used.c_str());
726 if (IsSchurType(linear_solver_type_used)) {
727 StringAppendF(&report,
728 "Schur structure %22s %24s\n",
729 schur_structure_given.c_str(),
730 schur_structure_used.c_str());
731 }
732
733 if (inner_iterations_given) {
734 StringAppendF(&report,
735 "Use inner iterations %20s %20s\n",
736 inner_iterations_given ? "True" : "False",
737 inner_iterations_used ? "True" : "False");
738 }
739
740 if (inner_iterations_used) {
741 string given;
742 StringifyOrdering(inner_iteration_ordering_given, &given);
743 string used;
744 StringifyOrdering(inner_iteration_ordering_used, &used);
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800745 StringAppendF(&report,
746 "Inner iteration ordering %20s %24s\n",
747 given.c_str(),
748 used.c_str());
Austin Schuh70cc9552019-01-21 19:46:48 -0800749 }
750 } else {
751 // LINE_SEARCH HEADER
752 StringAppendF(&report, "\nMinimizer %19s\n", "LINE_SEARCH");
753
Austin Schuh70cc9552019-01-21 19:46:48 -0800754 string line_search_direction_string;
755 if (line_search_direction_type == LBFGS) {
756 line_search_direction_string = StringPrintf("LBFGS (%d)", max_lbfgs_rank);
757 } else if (line_search_direction_type == NONLINEAR_CONJUGATE_GRADIENT) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800758 line_search_direction_string = NonlinearConjugateGradientTypeToString(
759 nonlinear_conjugate_gradient_type);
Austin Schuh70cc9552019-01-21 19:46:48 -0800760 } else {
761 line_search_direction_string =
762 LineSearchDirectionTypeToString(line_search_direction_type);
763 }
764
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800765 StringAppendF(&report,
766 "Line search direction %19s\n",
Austin Schuh70cc9552019-01-21 19:46:48 -0800767 line_search_direction_string.c_str());
768
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800769 const string line_search_type_string = StringPrintf(
770 "%s %s",
771 LineSearchInterpolationTypeToString(line_search_interpolation_type),
772 LineSearchTypeToString(line_search_type));
773 StringAppendF(&report,
774 "Line search type %19s\n",
Austin Schuh70cc9552019-01-21 19:46:48 -0800775 line_search_type_string.c_str());
776 StringAppendF(&report, "\n");
777
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800778 StringAppendF(&report, "%45s %21s\n", "Given", "Used");
779 StringAppendF(&report,
780 "Threads % 25d% 25d\n",
781 num_threads_given,
782 num_threads_used);
Austin Schuh70cc9552019-01-21 19:46:48 -0800783 }
784
785 StringAppendF(&report, "\nCost:\n");
786 StringAppendF(&report, "Initial % 30e\n", initial_cost);
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800787 if (termination_type != FAILURE && termination_type != USER_FAILURE) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800788 StringAppendF(&report, "Final % 30e\n", final_cost);
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800789 StringAppendF(&report, "Change % 30e\n", initial_cost - final_cost);
Austin Schuh70cc9552019-01-21 19:46:48 -0800790 }
791
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800792 StringAppendF(&report,
793 "\nMinimizer iterations % 16d\n",
Austin Schuh70cc9552019-01-21 19:46:48 -0800794 num_successful_steps + num_unsuccessful_steps);
795
796 // Successful/Unsuccessful steps only matter in the case of the
797 // trust region solver. Line search terminates when it encounters
798 // the first unsuccessful step.
799 if (minimizer_type == TRUST_REGION) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800800 StringAppendF(&report,
801 "Successful steps % 14d\n",
Austin Schuh70cc9552019-01-21 19:46:48 -0800802 num_successful_steps);
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800803 StringAppendF(&report,
804 "Unsuccessful steps % 14d\n",
Austin Schuh70cc9552019-01-21 19:46:48 -0800805 num_unsuccessful_steps);
806 }
807 if (inner_iterations_used) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800808 StringAppendF(&report,
809 "Steps with inner iterations % 14d\n",
Austin Schuh70cc9552019-01-21 19:46:48 -0800810 num_inner_iteration_steps);
811 }
812
813 const bool line_search_used =
814 (minimizer_type == LINE_SEARCH ||
815 (minimizer_type == TRUST_REGION && is_constrained));
816
817 if (line_search_used) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800818 StringAppendF(&report,
819 "Line search steps % 14d\n",
Austin Schuh70cc9552019-01-21 19:46:48 -0800820 num_line_search_steps);
821 }
822
823 StringAppendF(&report, "\nTime (in seconds):\n");
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800824 StringAppendF(
825 &report, "Preprocessor %25.6f\n", preprocessor_time_in_seconds);
Austin Schuh70cc9552019-01-21 19:46:48 -0800826
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800827 StringAppendF(&report,
828 "\n Residual only evaluation %18.6f (%d)\n",
829 residual_evaluation_time_in_seconds,
830 num_residual_evaluations);
Austin Schuh70cc9552019-01-21 19:46:48 -0800831 if (line_search_used) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800832 StringAppendF(&report,
833 " Line search cost evaluation %10.6f\n",
Austin Schuh70cc9552019-01-21 19:46:48 -0800834 line_search_cost_evaluation_time_in_seconds);
835 }
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800836 StringAppendF(&report,
837 " Jacobian & residual evaluation %12.6f (%d)\n",
838 jacobian_evaluation_time_in_seconds,
839 num_jacobian_evaluations);
Austin Schuh70cc9552019-01-21 19:46:48 -0800840 if (line_search_used) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800841 StringAppendF(&report,
842 " Line search gradient evaluation %6.6f\n",
Austin Schuh70cc9552019-01-21 19:46:48 -0800843 line_search_gradient_evaluation_time_in_seconds);
844 }
845
846 if (minimizer_type == TRUST_REGION) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800847 StringAppendF(&report,
848 " Linear solver %23.6f (%d)\n",
849 linear_solver_time_in_seconds,
850 num_linear_solves);
Austin Schuh70cc9552019-01-21 19:46:48 -0800851 }
852
853 if (inner_iterations_used) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800854 StringAppendF(&report,
855 " Inner iterations %23.6f\n",
Austin Schuh70cc9552019-01-21 19:46:48 -0800856 inner_iteration_time_in_seconds);
857 }
858
859 if (line_search_used) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800860 StringAppendF(&report,
861 " Line search polynomial minimization %.6f\n",
Austin Schuh70cc9552019-01-21 19:46:48 -0800862 line_search_polynomial_minimization_time_in_seconds);
863 }
864
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800865 StringAppendF(
866 &report, "Minimizer %25.6f\n\n", minimizer_time_in_seconds);
Austin Schuh70cc9552019-01-21 19:46:48 -0800867
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800868 StringAppendF(
869 &report, "Postprocessor %24.6f\n", postprocessor_time_in_seconds);
Austin Schuh70cc9552019-01-21 19:46:48 -0800870
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800871 StringAppendF(
872 &report, "Total %25.6f\n\n", total_time_in_seconds);
Austin Schuh70cc9552019-01-21 19:46:48 -0800873
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800874 StringAppendF(&report,
875 "Termination: %25s (%s)\n",
876 TerminationTypeToString(termination_type),
877 message.c_str());
Austin Schuh70cc9552019-01-21 19:46:48 -0800878 return report;
879}
880
881bool Solver::Summary::IsSolutionUsable() const {
882 return internal::IsSolutionUsable(*this);
883}
884
885} // namespace ceres