Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 1 | // Ceres Solver - A fast non-linear least squares minimizer |
| 2 | // Copyright 2010, 2011, 2012 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/solver.h" |
| 33 | |
| 34 | #include <vector> |
Keir Mierle | 6196cba | 2012-06-18 11:03:40 -0700 | [diff] [blame] | 35 | #include "ceres/problem.h" |
| 36 | #include "ceres/problem_impl.h" |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 37 | #include "ceres/program.h" |
| 38 | #include "ceres/solver_impl.h" |
| 39 | #include "ceres/stringprintf.h" |
Petter Strandmark | 76533b3 | 2012-09-04 22:08:50 -0700 | [diff] [blame] | 40 | #include "ceres/wall_time.h" |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 41 | |
| 42 | namespace ceres { |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 43 | namespace { |
| 44 | |
| 45 | void StringifyOrdering(const vector<int>& ordering, string* report) { |
| 46 | if (ordering.size() == 0) { |
| 47 | internal::StringAppendF(report, "AUTOMATIC"); |
| 48 | return; |
| 49 | } |
| 50 | |
| 51 | for (int i = 0; i < ordering.size() - 1; ++i) { |
| 52 | internal::StringAppendF(report, "%d, ", ordering[i]); |
| 53 | } |
| 54 | internal::StringAppendF(report, "%d", ordering.back()); |
| 55 | } |
| 56 | |
| 57 | } // namespace |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 58 | |
Sameer Agarwal | ba8d967 | 2012-10-02 00:48:57 -0700 | [diff] [blame] | 59 | Solver::Options::~Options() { |
Sameer Agarwal | 68b32a9 | 2012-10-06 23:10:51 -0700 | [diff] [blame] | 60 | delete linear_solver_ordering; |
Sameer Agarwal | ba8d967 | 2012-10-02 00:48:57 -0700 | [diff] [blame] | 61 | delete inner_iteration_ordering; |
| 62 | } |
| 63 | |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 64 | Solver::~Solver() {} |
| 65 | |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 66 | void Solver::Solve(const Solver::Options& options, |
| 67 | Problem* problem, |
| 68 | Solver::Summary* summary) { |
Petter Strandmark | 76533b3 | 2012-09-04 22:08:50 -0700 | [diff] [blame] | 69 | double start_time_seconds = internal::WallTimeInSeconds(); |
Keir Mierle | 6196cba | 2012-06-18 11:03:40 -0700 | [diff] [blame] | 70 | internal::ProblemImpl* problem_impl = |
| 71 | CHECK_NOTNULL(problem)->problem_impl_.get(); |
| 72 | internal::SolverImpl::Solve(options, problem_impl, summary); |
Petter Strandmark | 76533b3 | 2012-09-04 22:08:50 -0700 | [diff] [blame] | 73 | summary->total_time_in_seconds = |
| 74 | internal::WallTimeInSeconds() - start_time_seconds; |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 75 | } |
| 76 | |
| 77 | void Solve(const Solver::Options& options, |
| 78 | Problem* problem, |
| 79 | Solver::Summary* summary) { |
Keir Mierle | 6196cba | 2012-06-18 11:03:40 -0700 | [diff] [blame] | 80 | Solver solver; |
| 81 | solver.Solve(options, problem, summary); |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 82 | } |
| 83 | |
| 84 | Solver::Summary::Summary() |
| 85 | // Invalid values for most fields, to ensure that we are not |
| 86 | // accidentally reporting default values. |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 87 | : minimizer_type(TRUST_REGION), |
| 88 | termination_type(DID_NOT_RUN), |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 89 | initial_cost(-1.0), |
| 90 | final_cost(-1.0), |
| 91 | fixed_cost(-1.0), |
| 92 | num_successful_steps(-1), |
| 93 | num_unsuccessful_steps(-1), |
| 94 | preprocessor_time_in_seconds(-1.0), |
| 95 | minimizer_time_in_seconds(-1.0), |
Sameer Agarwal | fa01519 | 2012-06-11 14:21:42 -0700 | [diff] [blame] | 96 | postprocessor_time_in_seconds(-1.0), |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 97 | total_time_in_seconds(-1.0), |
Sameer Agarwal | 42a84b8 | 2013-02-01 12:22:53 -0800 | [diff] [blame] | 98 | linear_solver_time_in_seconds(-1.0), |
| 99 | residual_evaluation_time_in_seconds(-1.0), |
| 100 | jacobian_evaluation_time_in_seconds(-1.0), |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 101 | num_parameter_blocks(-1), |
| 102 | num_parameters(-1), |
Keir Mierle | ba94421 | 2013-02-25 12:46:44 -0800 | [diff] [blame] | 103 | num_effective_parameters(-1), |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 104 | num_residual_blocks(-1), |
| 105 | num_residuals(-1), |
| 106 | num_parameter_blocks_reduced(-1), |
| 107 | num_parameters_reduced(-1), |
Keir Mierle | ba94421 | 2013-02-25 12:46:44 -0800 | [diff] [blame] | 108 | num_effective_parameters_reduced(-1), |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 109 | num_residual_blocks_reduced(-1), |
| 110 | num_residuals_reduced(-1), |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 111 | num_threads_given(-1), |
| 112 | num_threads_used(-1), |
| 113 | num_linear_solver_threads_given(-1), |
| 114 | num_linear_solver_threads_used(-1), |
| 115 | linear_solver_type_given(SPARSE_NORMAL_CHOLESKY), |
| 116 | linear_solver_type_used(SPARSE_NORMAL_CHOLESKY), |
| 117 | preconditioner_type(IDENTITY), |
Sameer Agarwal | 97fb6d9 | 2012-06-17 10:08:19 -0700 | [diff] [blame] | 118 | trust_region_strategy_type(LEVENBERG_MARQUARDT), |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 119 | inner_iterations(false), |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 120 | sparse_linear_algebra_library(SUITE_SPARSE), |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 121 | line_search_direction_type(LBFGS), |
| 122 | line_search_type(ARMIJO) { |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 123 | } |
| 124 | |
| 125 | string Solver::Summary::BriefReport() const { |
| 126 | string report = "Ceres Solver Report: "; |
| 127 | if (termination_type == DID_NOT_RUN) { |
| 128 | CHECK(!error.empty()) |
| 129 | << "Solver terminated with DID_NOT_RUN but the solver did not " |
| 130 | << "return a reason. This is a Ceres error. Please report this " |
| 131 | << "to the Ceres team"; |
| 132 | return report + "Termination: DID_NOT_RUN, because " + error; |
| 133 | } |
| 134 | |
| 135 | internal::StringAppendF(&report, "Iterations: %d", |
| 136 | num_successful_steps + num_unsuccessful_steps); |
| 137 | internal::StringAppendF(&report, ", Initial cost: %e", initial_cost); |
| 138 | |
| 139 | // If the solver failed or was aborted, then the final_cost has no |
| 140 | // meaning. |
| 141 | if (termination_type != NUMERICAL_FAILURE && |
| 142 | termination_type != USER_ABORT) { |
| 143 | internal::StringAppendF(&report, ", Final cost: %e", final_cost); |
| 144 | } |
| 145 | |
| 146 | internal::StringAppendF(&report, ", Termination: %s.", |
| 147 | SolverTerminationTypeToString(termination_type)); |
| 148 | return report; |
| 149 | }; |
| 150 | |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 151 | using internal::StringAppendF; |
| 152 | |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 153 | string Solver::Summary::FullReport() const { |
| 154 | string report = |
| 155 | "\n" |
| 156 | "Ceres Solver Report\n" |
| 157 | "-------------------\n"; |
| 158 | |
| 159 | if (termination_type == DID_NOT_RUN) { |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 160 | StringAppendF(&report, " Original\n"); |
Keir Mierle | ba94421 | 2013-02-25 12:46:44 -0800 | [diff] [blame] | 161 | StringAppendF(&report, "Parameter blocks % 10d\n", num_parameter_blocks); |
| 162 | StringAppendF(&report, "Parameters % 10d\n", num_parameters); |
| 163 | if (num_effective_parameters != num_parameters) { |
| 164 | StringAppendF(&report, "Effective parameters% 10d\n", num_parameters); |
| 165 | } |
| 166 | |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 167 | StringAppendF(&report, "Residual blocks % 10d\n", |
| 168 | num_residual_blocks); |
| 169 | StringAppendF(&report, "Residuals % 10d\n\n", |
| 170 | num_residuals); |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 171 | } else { |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 172 | StringAppendF(&report, "%45s %21s\n", "Original", "Reduced"); |
| 173 | StringAppendF(&report, "Parameter blocks % 25d% 25d\n", |
| 174 | num_parameter_blocks, num_parameter_blocks_reduced); |
| 175 | StringAppendF(&report, "Parameters % 25d% 25d\n", |
| 176 | num_parameters, num_parameters_reduced); |
Keir Mierle | ba94421 | 2013-02-25 12:46:44 -0800 | [diff] [blame] | 177 | if (num_effective_parameters_reduced != num_parameters_reduced) { |
| 178 | StringAppendF(&report, "Effective parameters% 25d% 25d\n", |
| 179 | num_effective_parameters, num_effective_parameters_reduced); |
| 180 | } |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 181 | StringAppendF(&report, "Residual blocks % 25d% 25d\n", |
| 182 | num_residual_blocks, num_residual_blocks_reduced); |
| 183 | StringAppendF(&report, "Residual % 25d% 25d\n", |
| 184 | num_residuals, num_residuals_reduced); |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 185 | } |
| 186 | |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 187 | // TODO(sameeragarwal): Refactor this into separate functions. |
Sameer Agarwal | 97fb6d9 | 2012-06-17 10:08:19 -0700 | [diff] [blame] | 188 | |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 189 | if (minimizer_type == TRUST_REGION) { |
Sameer Agarwal | 509f68c | 2013-02-20 01:39:03 -0800 | [diff] [blame] | 190 | StringAppendF(&report, "\nMinimizer %19s\n", |
| 191 | "TRUST_REGION"); |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 192 | if (linear_solver_type_used == SPARSE_NORMAL_CHOLESKY || |
| 193 | linear_solver_type_used == SPARSE_SCHUR || |
| 194 | (linear_solver_type_used == ITERATIVE_SCHUR && |
Sameer Agarwal | 290b975 | 2013-02-17 16:50:37 -0800 | [diff] [blame] | 195 | (preconditioner_type == CLUSTER_JACOBI || |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 196 | preconditioner_type == CLUSTER_TRIDIAGONAL))) { |
| 197 | StringAppendF(&report, "\nSparse Linear Algebra Library %15s\n", |
| 198 | SparseLinearAlgebraLibraryTypeToString( |
| 199 | sparse_linear_algebra_library)); |
Sameer Agarwal | 1e28920 | 2012-08-21 18:00:54 -0700 | [diff] [blame] | 200 | } |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 201 | |
| 202 | StringAppendF(&report, "Trust Region Strategy %19s", |
| 203 | TrustRegionStrategyTypeToString( |
| 204 | trust_region_strategy_type)); |
| 205 | if (trust_region_strategy_type == DOGLEG) { |
| 206 | if (dogleg_type == TRADITIONAL_DOGLEG) { |
| 207 | StringAppendF(&report, " (TRADITIONAL)"); |
| 208 | } else { |
| 209 | StringAppendF(&report, " (SUBSPACE)"); |
| 210 | } |
| 211 | } |
| 212 | StringAppendF(&report, "\n"); |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 213 | StringAppendF(&report, "\n"); |
Sameer Agarwal | 1e28920 | 2012-08-21 18:00:54 -0700 | [diff] [blame] | 214 | |
Sameer Agarwal | 931c309 | 2013-02-25 09:46:21 -0800 | [diff] [blame] | 215 | StringAppendF(&report, "%45s %21s\n", "Given", "Used"); |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 216 | StringAppendF(&report, "Linear solver %25s%25s\n", |
| 217 | LinearSolverTypeToString(linear_solver_type_given), |
| 218 | LinearSolverTypeToString(linear_solver_type_used)); |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 219 | |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 220 | if (linear_solver_type_given == CGNR || |
| 221 | linear_solver_type_given == ITERATIVE_SCHUR) { |
| 222 | StringAppendF(&report, "Preconditioner %25s%25s\n", |
| 223 | PreconditionerTypeToString(preconditioner_type), |
| 224 | PreconditionerTypeToString(preconditioner_type)); |
| 225 | } else { |
| 226 | StringAppendF(&report, "Preconditioner %25s%25s\n", |
| 227 | "N/A", "N/A"); |
| 228 | } |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 229 | |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 230 | StringAppendF(&report, "Threads: % 25d% 25d\n", |
| 231 | num_threads_given, num_threads_used); |
| 232 | StringAppendF(&report, "Linear solver threads % 23d% 25d\n", |
| 233 | num_linear_solver_threads_given, |
| 234 | num_linear_solver_threads_used); |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 235 | |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 236 | if (IsSchurType(linear_solver_type_used)) { |
| 237 | string given; |
| 238 | StringifyOrdering(linear_solver_ordering_given, &given); |
| 239 | string used; |
| 240 | StringifyOrdering(linear_solver_ordering_used, &used); |
| 241 | StringAppendF(&report, |
| 242 | "Linear solver ordering %22s %24s\n", |
| 243 | given.c_str(), |
| 244 | used.c_str()); |
| 245 | } |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 246 | |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 247 | if (inner_iterations) { |
| 248 | string given; |
| 249 | StringifyOrdering(inner_iteration_ordering_given, &given); |
| 250 | string used; |
| 251 | StringifyOrdering(inner_iteration_ordering_used, &used); |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 252 | StringAppendF(&report, |
| 253 | "Inner iteration ordering %20s %24s\n", |
| 254 | given.c_str(), |
| 255 | used.c_str()); |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 256 | } |
| 257 | |
| 258 | if (termination_type == DID_NOT_RUN) { |
| 259 | CHECK(!error.empty()) |
| 260 | << "Solver terminated with DID_NOT_RUN but the solver did not " |
| 261 | << "return a reason. This is a Ceres error. Please report this " |
| 262 | << "to the Ceres team"; |
| 263 | StringAppendF(&report, "Termination: %20s\n", |
| 264 | "DID_NOT_RUN"); |
| 265 | StringAppendF(&report, "Reason: %s\n", error.c_str()); |
| 266 | return report; |
| 267 | } |
| 268 | |
| 269 | StringAppendF(&report, "\nCost:\n"); |
| 270 | StringAppendF(&report, "Initial % 30e\n", initial_cost); |
Sameer Agarwal | 509f68c | 2013-02-20 01:39:03 -0800 | [diff] [blame] | 271 | if (termination_type != NUMERICAL_FAILURE && |
| 272 | termination_type != USER_ABORT) { |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 273 | StringAppendF(&report, "Final % 30e\n", final_cost); |
| 274 | StringAppendF(&report, "Change % 30e\n", |
| 275 | initial_cost - final_cost); |
| 276 | } |
| 277 | |
| 278 | StringAppendF(&report, "\nNumber of iterations:\n"); |
| 279 | StringAppendF(&report, "Successful % 20d\n", |
| 280 | num_successful_steps); |
| 281 | StringAppendF(&report, "Unsuccessful % 20d\n", |
| 282 | num_unsuccessful_steps); |
| 283 | StringAppendF(&report, "Total % 20d\n", |
| 284 | num_successful_steps + num_unsuccessful_steps); |
| 285 | |
| 286 | StringAppendF(&report, "\nTime (in seconds):\n"); |
| 287 | StringAppendF(&report, "Preprocessor %25.3f\n", |
| 288 | preprocessor_time_in_seconds); |
| 289 | StringAppendF(&report, "\n Residual Evaluations %22.3f\n", |
| 290 | residual_evaluation_time_in_seconds); |
| 291 | StringAppendF(&report, " Jacobian Evaluations %22.3f\n", |
| 292 | jacobian_evaluation_time_in_seconds); |
| 293 | StringAppendF(&report, " Linear Solver %23.3f\n", |
| 294 | linear_solver_time_in_seconds); |
| 295 | StringAppendF(&report, "Minimizer %25.3f\n\n", |
| 296 | minimizer_time_in_seconds); |
| 297 | |
| 298 | StringAppendF(&report, "Postprocessor %24.3f\n", |
| 299 | postprocessor_time_in_seconds); |
| 300 | |
| 301 | StringAppendF(&report, "Total %25.3f\n\n", |
| 302 | total_time_in_seconds); |
| 303 | |
| 304 | StringAppendF(&report, "Termination: %25s\n", |
| 305 | SolverTerminationTypeToString(termination_type)); |
| 306 | } else { |
| 307 | // LINE_SEARCH |
| 308 | StringAppendF(&report, "\nMinimizer %19s\n", "LINE_SEARCH"); |
| 309 | if (line_search_direction_type == LBFGS) { |
Sameer Agarwal | 931c309 | 2013-02-25 09:46:21 -0800 | [diff] [blame] | 310 | StringAppendF(&report, "Line search direction %19s(%d)\n", |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 311 | LineSearchDirectionTypeToString(line_search_direction_type), |
| 312 | max_lbfgs_rank); |
| 313 | } else { |
Sameer Agarwal | 931c309 | 2013-02-25 09:46:21 -0800 | [diff] [blame] | 314 | StringAppendF(&report, "Line search direction %19s\n", |
Sameer Agarwal | 509f68c | 2013-02-20 01:39:03 -0800 | [diff] [blame] | 315 | LineSearchDirectionTypeToString( |
| 316 | line_search_direction_type)); |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 317 | } |
Sameer Agarwal | 931c309 | 2013-02-25 09:46:21 -0800 | [diff] [blame] | 318 | StringAppendF(&report, "Line search type %19s\n", |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 319 | LineSearchTypeToString(line_search_type)); |
| 320 | |
| 321 | StringAppendF(&report, "\n"); |
| 322 | |
Sameer Agarwal | beb4505 | 2013-02-22 13:37:01 -0800 | [diff] [blame] | 323 | StringAppendF(&report, "%45s %21s\n", "Given", "Used"); |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 324 | StringAppendF(&report, "Threads: % 25d% 25d\n", |
| 325 | num_threads_given, num_threads_used); |
| 326 | |
| 327 | if (termination_type == DID_NOT_RUN) { |
| 328 | CHECK(!error.empty()) |
| 329 | << "Solver terminated with DID_NOT_RUN but the solver did not " |
| 330 | << "return a reason. This is a Ceres error. Please report this " |
| 331 | << "to the Ceres team"; |
| 332 | StringAppendF(&report, "Termination: %20s\n", |
| 333 | "DID_NOT_RUN"); |
| 334 | StringAppendF(&report, "Reason: %s\n", error.c_str()); |
| 335 | return report; |
| 336 | } |
| 337 | |
| 338 | StringAppendF(&report, "\nCost:\n"); |
| 339 | StringAppendF(&report, "Initial % 30e\n", initial_cost); |
Sameer Agarwal | 509f68c | 2013-02-20 01:39:03 -0800 | [diff] [blame] | 340 | if (termination_type != NUMERICAL_FAILURE && |
| 341 | termination_type != USER_ABORT) { |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 342 | StringAppendF(&report, "Final % 30e\n", final_cost); |
| 343 | StringAppendF(&report, "Change % 30e\n", |
| 344 | initial_cost - final_cost); |
| 345 | } |
| 346 | |
Sameer Agarwal | 36dc14d | 2013-02-25 10:33:10 -0800 | [diff] [blame] | 347 | StringAppendF(&report, "\nNumber of iterations: % 20d\n", |
| 348 | static_cast<int>(iterations.size() - 1)); |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 349 | |
| 350 | StringAppendF(&report, "\nTime (in seconds):\n"); |
| 351 | StringAppendF(&report, "Preprocessor %25.3f\n", |
| 352 | preprocessor_time_in_seconds); |
| 353 | StringAppendF(&report, "\n Residual Evaluations %22.3f\n", |
| 354 | residual_evaluation_time_in_seconds); |
| 355 | StringAppendF(&report, " Jacobian Evaluations %22.3f\n", |
| 356 | jacobian_evaluation_time_in_seconds); |
| 357 | StringAppendF(&report, "Minimizer %25.3f\n\n", |
| 358 | minimizer_time_in_seconds); |
| 359 | |
| 360 | StringAppendF(&report, "Postprocessor %24.3f\n", |
| 361 | postprocessor_time_in_seconds); |
| 362 | |
| 363 | StringAppendF(&report, "Total %25.3f\n\n", |
| 364 | total_time_in_seconds); |
| 365 | |
| 366 | StringAppendF(&report, "Termination: %25s\n", |
| 367 | SolverTerminationTypeToString(termination_type)); |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 368 | } |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 369 | |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 370 | return report; |
| 371 | }; |
| 372 | |
| 373 | } // namespace ceres |