• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1 // Ceres Solver - A fast non-linear least squares minimizer
2 // Copyright 2014 Google Inc. All rights reserved.
3 // http://code.google.com/p/ceres-solver/
4 //
5 // Redistribution and use in source and binary forms, with or without
6 // modification, are permitted provided that the following conditions are met:
7 //
8 // * Redistributions of source code must retain the above copyright notice,
9 //   this list of conditions and the following disclaimer.
10 // * Redistributions in binary form must reproduce the above copyright notice,
11 //   this list of conditions and the following disclaimer in the documentation
12 //   and/or other materials provided with the distribution.
13 // * Neither the name of Google Inc. nor the names of its contributors may be
14 //   used to endorse or promote products derived from this software without
15 //   specific prior written permission.
16 //
17 // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
18 // AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
19 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
20 // ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
21 // LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
22 // CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
23 // SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
24 // INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
25 // CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
26 // ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
27 // POSSIBILITY OF SUCH DAMAGE.
28 //
29 // Author: keir@google.com (Keir Mierle)
30 //         sameeragarwal@google.com (Sameer Agarwal)
31 
32 #include "ceres/internal/port.h"
33 #include "ceres/solver.h"
34 
35 #include <sstream>   // NOLINT
36 #include <vector>
37 
38 #include "ceres/problem.h"
39 #include "ceres/problem_impl.h"
40 #include "ceres/program.h"
41 #include "ceres/solver_impl.h"
42 #include "ceres/stringprintf.h"
43 #include "ceres/types.h"
44 #include "ceres/version.h"
45 #include "ceres/wall_time.h"
46 
47 namespace ceres {
48 namespace {
49 
50 #define OPTION_OP(x, y, OP)                                             \
51   if (!(options.x OP y)) {                                              \
52     std::stringstream ss;                                               \
53     ss << "Invalid configuration. ";                                    \
54     ss << string("Solver::Options::" #x " = ") << options.x << ". ";    \
55     ss << "Violated constraint: ";                                      \
56     ss << string("Solver::Options::" #x " " #OP " "#y);                 \
57     *error = ss.str();                                                  \
58     return false;                                                       \
59   }
60 
61 #define OPTION_OP_OPTION(x, y, OP)                                      \
62   if (!(options.x OP options.y)) {                                      \
63     std::stringstream ss;                                               \
64     ss << "Invalid configuration. ";                                    \
65     ss << string("Solver::Options::" #x " = ") << options.x << ". ";    \
66     ss << string("Solver::Options::" #y " = ") << options.y << ". ";    \
67     ss << "Violated constraint: ";                                      \
68     ss << string("Solver::Options::" #x );                              \
69     ss << string(#OP " Solver::Options::" #y ".");                      \
70     *error = ss.str();                                                  \
71     return false;                                                       \
72   }
73 
74 #define OPTION_GE(x, y) OPTION_OP(x, y, >=);
75 #define OPTION_GT(x, y) OPTION_OP(x, y, >);
76 #define OPTION_LE(x, y) OPTION_OP(x, y, <=);
77 #define OPTION_LT(x, y) OPTION_OP(x, y, <);
78 #define OPTION_LE_OPTION(x, y) OPTION_OP_OPTION(x ,y, <=)
79 #define OPTION_LT_OPTION(x, y) OPTION_OP_OPTION(x ,y, <)
80 
CommonOptionsAreValid(const Solver::Options & options,string * error)81 bool CommonOptionsAreValid(const Solver::Options& options, string* error) {
82   OPTION_GE(max_num_iterations, 0);
83   OPTION_GE(max_solver_time_in_seconds, 0.0);
84   OPTION_GE(function_tolerance, 0.0);
85   OPTION_GE(gradient_tolerance, 0.0);
86   OPTION_GE(parameter_tolerance, 0.0);
87   OPTION_GT(num_threads, 0);
88   OPTION_GT(num_linear_solver_threads, 0);
89   if (options.check_gradients) {
90     OPTION_GT(gradient_check_relative_precision, 0.0);
91     OPTION_GT(numeric_derivative_relative_step_size, 0.0);
92   }
93   return true;
94 }
95 
TrustRegionOptionsAreValid(const Solver::Options & options,string * error)96 bool TrustRegionOptionsAreValid(const Solver::Options& options, string* error) {
97   OPTION_GT(initial_trust_region_radius, 0.0);
98   OPTION_GT(min_trust_region_radius, 0.0);
99   OPTION_GT(max_trust_region_radius, 0.0);
100   OPTION_LE_OPTION(min_trust_region_radius, max_trust_region_radius);
101   OPTION_LE_OPTION(min_trust_region_radius, initial_trust_region_radius);
102   OPTION_LE_OPTION(initial_trust_region_radius, max_trust_region_radius);
103   OPTION_GE(min_relative_decrease, 0.0);
104   OPTION_GE(min_lm_diagonal, 0.0);
105   OPTION_GE(max_lm_diagonal, 0.0);
106   OPTION_LE_OPTION(min_lm_diagonal, max_lm_diagonal);
107   OPTION_GE(max_num_consecutive_invalid_steps, 0);
108   OPTION_GT(eta, 0.0);
109   OPTION_GE(min_linear_solver_iterations, 1);
110   OPTION_GE(max_linear_solver_iterations, 1);
111   OPTION_LE_OPTION(min_linear_solver_iterations, max_linear_solver_iterations);
112 
113   if (options.use_inner_iterations) {
114     OPTION_GE(inner_iteration_tolerance, 0.0);
115   }
116 
117   if (options.use_nonmonotonic_steps) {
118     OPTION_GT(max_consecutive_nonmonotonic_steps, 0);
119   }
120 
121   if (options.preconditioner_type == CLUSTER_JACOBI &&
122       options.sparse_linear_algebra_library_type != SUITE_SPARSE) {
123     *error =  "CLUSTER_JACOBI requires "
124         "Solver::Options::sparse_linear_algebra_library_type to be "
125         "SUITE_SPARSE";
126     return false;
127   }
128 
129   if (options.preconditioner_type == CLUSTER_TRIDIAGONAL &&
130       options.sparse_linear_algebra_library_type != SUITE_SPARSE) {
131     *error =  "CLUSTER_TRIDIAGONAL requires "
132         "Solver::Options::sparse_linear_algebra_library_type to be "
133         "SUITE_SPARSE";
134     return false;
135   }
136 
137 #ifdef CERES_NO_LAPACK
138   if (options.dense_linear_algebra_library_type == LAPACK) {
139     if (options.linear_solver_type == DENSE_NORMAL_CHOLESKY) {
140       *error = "Can't use DENSE_NORMAL_CHOLESKY with LAPACK because "
141           "LAPACK was not enabled when Ceres was built.";
142       return false;
143     }
144 
145     if (options.linear_solver_type == DENSE_QR) {
146       *error = "Can't use DENSE_QR with LAPACK because "
147           "LAPACK was not enabled when Ceres was built.";
148       return false;
149     }
150 
151     if (options.linear_solver_type == DENSE_SCHUR) {
152       *error = "Can't use DENSE_SCHUR with LAPACK because "
153           "LAPACK was not enabled when Ceres was built.";
154       return false;
155     }
156   }
157 #endif
158 
159 #ifdef CERES_NO_SUITESPARSE
160   if (options.sparse_linear_algebra_library_type == SUITE_SPARSE) {
161     if (options.linear_solver_type == SPARSE_NORMAL_CHOLESKY) {
162       *error = "Can't use SPARSE_NORMAL_CHOLESKY with SUITESPARSE because "
163              "SuiteSparse was not enabled when Ceres was built.";
164       return false;
165     }
166 
167     if (options.linear_solver_type == SPARSE_SCHUR) {
168       *error = "Can't use SPARSE_SCHUR with SUITESPARSE because "
169           "SuiteSparse was not enabled when Ceres was built.";
170       return false;
171     }
172 
173     if (options.preconditioner_type == CLUSTER_JACOBI) {
174       *error =  "CLUSTER_JACOBI preconditioner not supported. "
175           "SuiteSparse was not enabled when Ceres was built.";
176       return false;
177     }
178 
179     if (options.preconditioner_type == CLUSTER_TRIDIAGONAL) {
180       *error =  "CLUSTER_TRIDIAGONAL preconditioner not supported. "
181           "SuiteSparse was not enabled when Ceres was built.";
182     return false;
183     }
184   }
185 #endif
186 
187 #ifdef CERES_NO_CXSPARSE
188   if (options.sparse_linear_algebra_library_type == CX_SPARSE) {
189     if (options.linear_solver_type == SPARSE_NORMAL_CHOLESKY) {
190       *error = "Can't use SPARSE_NORMAL_CHOLESKY with CX_SPARSE because "
191              "CXSparse was not enabled when Ceres was built.";
192       return false;
193     }
194 
195     if (options.linear_solver_type == SPARSE_SCHUR) {
196       *error = "Can't use SPARSE_SCHUR with CX_SPARSE because "
197           "CXSparse was not enabled when Ceres was built.";
198       return false;
199     }
200   }
201 #endif
202 
203   if (options.trust_region_strategy_type == DOGLEG) {
204     if (options.linear_solver_type == ITERATIVE_SCHUR ||
205         options.linear_solver_type == CGNR) {
206       *error = "DOGLEG only supports exact factorization based linear "
207           "solvers. If you want to use an iterative solver please "
208           "use LEVENBERG_MARQUARDT as the trust_region_strategy_type";
209       return false;
210     }
211   }
212 
213   if (options.trust_region_minimizer_iterations_to_dump.size() > 0 &&
214       options.trust_region_problem_dump_format_type != CONSOLE &&
215       options.trust_region_problem_dump_directory.empty()) {
216     *error = "Solver::Options::trust_region_problem_dump_directory is empty.";
217     return false;
218   }
219 
220   if (options.dynamic_sparsity &&
221       options.linear_solver_type != SPARSE_NORMAL_CHOLESKY) {
222     *error = "Dynamic sparsity is only supported with SPARSE_NORMAL_CHOLESKY.";
223     return false;
224   }
225 
226   return true;
227 }
228 
LineSearchOptionsAreValid(const Solver::Options & options,string * error)229 bool LineSearchOptionsAreValid(const Solver::Options& options, string* error) {
230   OPTION_GT(max_lbfgs_rank, 0);
231   OPTION_GT(min_line_search_step_size, 0.0);
232   OPTION_GT(max_line_search_step_contraction, 0.0);
233   OPTION_LT(max_line_search_step_contraction, 1.0);
234   OPTION_LT_OPTION(max_line_search_step_contraction,
235                    min_line_search_step_contraction);
236   OPTION_LE(min_line_search_step_contraction, 1.0);
237   OPTION_GT(max_num_line_search_step_size_iterations, 0);
238   OPTION_GT(line_search_sufficient_function_decrease, 0.0);
239   OPTION_LT_OPTION(line_search_sufficient_function_decrease,
240                    line_search_sufficient_curvature_decrease);
241   OPTION_LT(line_search_sufficient_curvature_decrease, 1.0);
242   OPTION_GT(max_line_search_step_expansion, 1.0);
243 
244   if ((options.line_search_direction_type == ceres::BFGS ||
245        options.line_search_direction_type == ceres::LBFGS) &&
246       options.line_search_type != ceres::WOLFE) {
247 
248     *error =
249         string("Invalid configuration: Solver::Options::line_search_type = ")
250         + string(LineSearchTypeToString(options.line_search_type))
251         + string(". When using (L)BFGS, "
252                  "Solver::Options::line_search_type must be set to WOLFE.");
253     return false;
254   }
255 
256   // Warn user if they have requested BISECTION interpolation, but constraints
257   // on max/min step size change during line search prevent bisection scaling
258   // from occurring. Warn only, as this is likely a user mistake, but one which
259   // does not prevent us from continuing.
260   LOG_IF(WARNING,
261          (options.line_search_interpolation_type == ceres::BISECTION &&
262           (options.max_line_search_step_contraction > 0.5 ||
263            options.min_line_search_step_contraction < 0.5)))
264       << "Line search interpolation type is BISECTION, but specified "
265       << "max_line_search_step_contraction: "
266       << options.max_line_search_step_contraction << ", and "
267       << "min_line_search_step_contraction: "
268       << options.min_line_search_step_contraction
269       << ", prevent bisection (0.5) scaling, continuing with solve regardless.";
270 
271   return true;
272 }
273 
274 #undef OPTION_OP
275 #undef OPTION_OP_OPTION
276 #undef OPTION_GT
277 #undef OPTION_GE
278 #undef OPTION_LE
279 #undef OPTION_LT
280 #undef OPTION_LE_OPTION
281 #undef OPTION_LT_OPTION
282 
StringifyOrdering(const vector<int> & ordering,string * report)283 void StringifyOrdering(const vector<int>& ordering, string* report) {
284   if (ordering.size() == 0) {
285     internal::StringAppendF(report, "AUTOMATIC");
286     return;
287   }
288 
289   for (int i = 0; i < ordering.size() - 1; ++i) {
290     internal::StringAppendF(report, "%d, ", ordering[i]);
291   }
292   internal::StringAppendF(report, "%d", ordering.back());
293 }
294 
295 } // namespace
296 
IsValid(string * error) const297 bool Solver::Options::IsValid(string* error) const {
298   if (!CommonOptionsAreValid(*this, error)) {
299     return false;
300   }
301 
302   if (minimizer_type == TRUST_REGION) {
303     return TrustRegionOptionsAreValid(*this, error);
304   }
305 
306   CHECK_EQ(minimizer_type, LINE_SEARCH);
307   return LineSearchOptionsAreValid(*this, error);
308 }
309 
~Solver()310 Solver::~Solver() {}
311 
Solve(const Solver::Options & options,Problem * problem,Solver::Summary * summary)312 void Solver::Solve(const Solver::Options& options,
313                    Problem* problem,
314                    Solver::Summary* summary) {
315   double start_time_seconds = internal::WallTimeInSeconds();
316   CHECK_NOTNULL(problem);
317   CHECK_NOTNULL(summary);
318 
319   *summary = Summary();
320   if (!options.IsValid(&summary->message)) {
321     LOG(ERROR) << "Terminating: " << summary->message;
322     return;
323   }
324 
325   internal::ProblemImpl* problem_impl = problem->problem_impl_.get();
326   internal::SolverImpl::Solve(options, problem_impl, summary);
327   summary->total_time_in_seconds =
328       internal::WallTimeInSeconds() - start_time_seconds;
329 }
330 
Solve(const Solver::Options & options,Problem * problem,Solver::Summary * summary)331 void Solve(const Solver::Options& options,
332            Problem* problem,
333            Solver::Summary* summary) {
334   Solver solver;
335   solver.Solve(options, problem, summary);
336 }
337 
Summary()338 Solver::Summary::Summary()
339     // Invalid values for most fields, to ensure that we are not
340     // accidentally reporting default values.
341     : minimizer_type(TRUST_REGION),
342       termination_type(FAILURE),
343       message("ceres::Solve was not called."),
344       initial_cost(-1.0),
345       final_cost(-1.0),
346       fixed_cost(-1.0),
347       num_successful_steps(-1),
348       num_unsuccessful_steps(-1),
349       num_inner_iteration_steps(-1),
350       preprocessor_time_in_seconds(-1.0),
351       minimizer_time_in_seconds(-1.0),
352       postprocessor_time_in_seconds(-1.0),
353       total_time_in_seconds(-1.0),
354       linear_solver_time_in_seconds(-1.0),
355       residual_evaluation_time_in_seconds(-1.0),
356       jacobian_evaluation_time_in_seconds(-1.0),
357       inner_iteration_time_in_seconds(-1.0),
358       num_parameter_blocks(-1),
359       num_parameters(-1),
360       num_effective_parameters(-1),
361       num_residual_blocks(-1),
362       num_residuals(-1),
363       num_parameter_blocks_reduced(-1),
364       num_parameters_reduced(-1),
365       num_effective_parameters_reduced(-1),
366       num_residual_blocks_reduced(-1),
367       num_residuals_reduced(-1),
368       num_threads_given(-1),
369       num_threads_used(-1),
370       num_linear_solver_threads_given(-1),
371       num_linear_solver_threads_used(-1),
372       linear_solver_type_given(SPARSE_NORMAL_CHOLESKY),
373       linear_solver_type_used(SPARSE_NORMAL_CHOLESKY),
374       inner_iterations_given(false),
375       inner_iterations_used(false),
376       preconditioner_type(IDENTITY),
377       visibility_clustering_type(CANONICAL_VIEWS),
378       trust_region_strategy_type(LEVENBERG_MARQUARDT),
379       dense_linear_algebra_library_type(EIGEN),
380       sparse_linear_algebra_library_type(SUITE_SPARSE),
381       line_search_direction_type(LBFGS),
382       line_search_type(ARMIJO),
383       line_search_interpolation_type(BISECTION),
384       nonlinear_conjugate_gradient_type(FLETCHER_REEVES),
385       max_lbfgs_rank(-1) {
386 }
387 
388 using internal::StringAppendF;
389 using internal::StringPrintf;
390 
BriefReport() const391 string Solver::Summary::BriefReport() const {
392   return StringPrintf("Ceres Solver Report: "
393                       "Iterations: %d, "
394                       "Initial cost: %e, "
395                       "Final cost: %e, "
396                       "Termination: %s",
397                       num_successful_steps + num_unsuccessful_steps,
398                       initial_cost,
399                       final_cost,
400                       TerminationTypeToString(termination_type));
401 };
402 
FullReport() const403 string Solver::Summary::FullReport() const {
404   string report =
405       "\n"
406       "Ceres Solver v" CERES_VERSION_STRING " Solve Report\n"
407       "----------------------------------\n";
408 
409   StringAppendF(&report, "%45s    %21s\n", "Original", "Reduced");
410   StringAppendF(&report, "Parameter blocks    % 25d% 25d\n",
411                 num_parameter_blocks, num_parameter_blocks_reduced);
412   StringAppendF(&report, "Parameters          % 25d% 25d\n",
413                 num_parameters, num_parameters_reduced);
414   if (num_effective_parameters_reduced != num_parameters_reduced) {
415     StringAppendF(&report, "Effective parameters% 25d% 25d\n",
416                   num_effective_parameters, num_effective_parameters_reduced);
417   }
418   StringAppendF(&report, "Residual blocks     % 25d% 25d\n",
419                 num_residual_blocks, num_residual_blocks_reduced);
420   StringAppendF(&report, "Residual            % 25d% 25d\n",
421                 num_residuals, num_residuals_reduced);
422 
423   if (minimizer_type == TRUST_REGION) {
424     // TRUST_SEARCH HEADER
425     StringAppendF(&report, "\nMinimizer                 %19s\n",
426                   "TRUST_REGION");
427 
428     if (linear_solver_type_used == DENSE_NORMAL_CHOLESKY ||
429         linear_solver_type_used == DENSE_SCHUR ||
430         linear_solver_type_used == DENSE_QR) {
431       StringAppendF(&report, "\nDense linear algebra library  %15s\n",
432                     DenseLinearAlgebraLibraryTypeToString(
433                         dense_linear_algebra_library_type));
434     }
435 
436     if (linear_solver_type_used == SPARSE_NORMAL_CHOLESKY ||
437         linear_solver_type_used == SPARSE_SCHUR ||
438         (linear_solver_type_used == ITERATIVE_SCHUR &&
439          (preconditioner_type == CLUSTER_JACOBI ||
440           preconditioner_type == CLUSTER_TRIDIAGONAL))) {
441       StringAppendF(&report, "\nSparse linear algebra library %15s\n",
442                     SparseLinearAlgebraLibraryTypeToString(
443                         sparse_linear_algebra_library_type));
444     }
445 
446     StringAppendF(&report, "Trust region strategy     %19s",
447                   TrustRegionStrategyTypeToString(
448                       trust_region_strategy_type));
449     if (trust_region_strategy_type == DOGLEG) {
450       if (dogleg_type == TRADITIONAL_DOGLEG) {
451         StringAppendF(&report, " (TRADITIONAL)");
452       } else {
453         StringAppendF(&report, " (SUBSPACE)");
454       }
455     }
456     StringAppendF(&report, "\n");
457     StringAppendF(&report, "\n");
458 
459     StringAppendF(&report, "%45s    %21s\n", "Given",  "Used");
460     StringAppendF(&report, "Linear solver       %25s%25s\n",
461                   LinearSolverTypeToString(linear_solver_type_given),
462                   LinearSolverTypeToString(linear_solver_type_used));
463 
464     if (linear_solver_type_given == CGNR ||
465         linear_solver_type_given == ITERATIVE_SCHUR) {
466       StringAppendF(&report, "Preconditioner      %25s%25s\n",
467                     PreconditionerTypeToString(preconditioner_type),
468                     PreconditionerTypeToString(preconditioner_type));
469     }
470 
471     if (preconditioner_type == CLUSTER_JACOBI ||
472         preconditioner_type == CLUSTER_TRIDIAGONAL) {
473       StringAppendF(&report, "Visibility clustering%24s%25s\n",
474                     VisibilityClusteringTypeToString(
475                         visibility_clustering_type),
476                     VisibilityClusteringTypeToString(
477                         visibility_clustering_type));
478     }
479     StringAppendF(&report, "Threads             % 25d% 25d\n",
480                   num_threads_given, num_threads_used);
481     StringAppendF(&report, "Linear solver threads % 23d% 25d\n",
482                   num_linear_solver_threads_given,
483                   num_linear_solver_threads_used);
484 
485     if (IsSchurType(linear_solver_type_used)) {
486       string given;
487       StringifyOrdering(linear_solver_ordering_given, &given);
488       string used;
489       StringifyOrdering(linear_solver_ordering_used, &used);
490       StringAppendF(&report,
491                     "Linear solver ordering %22s %24s\n",
492                     given.c_str(),
493                     used.c_str());
494     }
495 
496     if (inner_iterations_given) {
497       StringAppendF(&report,
498                     "Use inner iterations     %20s     %20s\n",
499                     inner_iterations_given ? "True" : "False",
500                     inner_iterations_used ? "True" : "False");
501     }
502 
503     if (inner_iterations_used) {
504       string given;
505       StringifyOrdering(inner_iteration_ordering_given, &given);
506       string used;
507       StringifyOrdering(inner_iteration_ordering_used, &used);
508     StringAppendF(&report,
509                   "Inner iteration ordering %20s %24s\n",
510                   given.c_str(),
511                   used.c_str());
512     }
513   } else {
514     // LINE_SEARCH HEADER
515     StringAppendF(&report, "\nMinimizer                 %19s\n", "LINE_SEARCH");
516 
517 
518     string line_search_direction_string;
519     if (line_search_direction_type == LBFGS) {
520       line_search_direction_string = StringPrintf("LBFGS (%d)", max_lbfgs_rank);
521     } else if (line_search_direction_type == NONLINEAR_CONJUGATE_GRADIENT) {
522       line_search_direction_string =
523           NonlinearConjugateGradientTypeToString(
524               nonlinear_conjugate_gradient_type);
525     } else {
526       line_search_direction_string =
527           LineSearchDirectionTypeToString(line_search_direction_type);
528     }
529 
530     StringAppendF(&report, "Line search direction     %19s\n",
531                   line_search_direction_string.c_str());
532 
533     const string line_search_type_string =
534         StringPrintf("%s %s",
535                      LineSearchInterpolationTypeToString(
536                          line_search_interpolation_type),
537                      LineSearchTypeToString(line_search_type));
538     StringAppendF(&report, "Line search type          %19s\n",
539                   line_search_type_string.c_str());
540     StringAppendF(&report, "\n");
541 
542     StringAppendF(&report, "%45s    %21s\n", "Given",  "Used");
543     StringAppendF(&report, "Threads             % 25d% 25d\n",
544                   num_threads_given, num_threads_used);
545   }
546 
547   StringAppendF(&report, "\nCost:\n");
548   StringAppendF(&report, "Initial        % 30e\n", initial_cost);
549   if (termination_type != FAILURE &&
550       termination_type != USER_FAILURE) {
551     StringAppendF(&report, "Final          % 30e\n", final_cost);
552     StringAppendF(&report, "Change         % 30e\n",
553                   initial_cost - final_cost);
554   }
555 
556   StringAppendF(&report, "\nMinimizer iterations         % 16d\n",
557                 num_successful_steps + num_unsuccessful_steps);
558 
559   // Successful/Unsuccessful steps only matter in the case of the
560   // trust region solver. Line search terminates when it encounters
561   // the first unsuccessful step.
562   if (minimizer_type == TRUST_REGION) {
563     StringAppendF(&report, "Successful steps               % 14d\n",
564                   num_successful_steps);
565     StringAppendF(&report, "Unsuccessful steps             % 14d\n",
566                   num_unsuccessful_steps);
567   }
568   if (inner_iterations_used) {
569     StringAppendF(&report, "Steps with inner iterations    % 14d\n",
570                   num_inner_iteration_steps);
571   }
572 
573   StringAppendF(&report, "\nTime (in seconds):\n");
574   StringAppendF(&report, "Preprocessor        %25.3f\n",
575                 preprocessor_time_in_seconds);
576 
577   StringAppendF(&report, "\n  Residual evaluation %23.3f\n",
578                 residual_evaluation_time_in_seconds);
579   StringAppendF(&report, "  Jacobian evaluation %23.3f\n",
580                 jacobian_evaluation_time_in_seconds);
581 
582   if (minimizer_type == TRUST_REGION) {
583     StringAppendF(&report, "  Linear solver       %23.3f\n",
584                   linear_solver_time_in_seconds);
585   }
586 
587   if (inner_iterations_used) {
588     StringAppendF(&report, "  Inner iterations    %23.3f\n",
589                   inner_iteration_time_in_seconds);
590   }
591 
592   StringAppendF(&report, "Minimizer           %25.3f\n\n",
593                 minimizer_time_in_seconds);
594 
595   StringAppendF(&report, "Postprocessor        %24.3f\n",
596                 postprocessor_time_in_seconds);
597 
598   StringAppendF(&report, "Total               %25.3f\n\n",
599                 total_time_in_seconds);
600 
601   StringAppendF(&report, "Termination:        %25s (%s)\n",
602                 TerminationTypeToString(termination_type), message.c_str());
603   return report;
604 };
605 
IsSolutionUsable() const606 bool Solver::Summary::IsSolutionUsable() const {
607   return (termination_type == CONVERGENCE ||
608           termination_type == NO_CONVERGENCE ||
609           termination_type == USER_SUCCESS);
610 }
611 
612 }  // namespace ceres
613