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 (at) google.com (Keir Mierle) 30 // sameeragarwal (at) google.com (Sameer Agarwal) 31 32 #include "ceres/solver.h" 33 34 #include <vector> 35 #include "ceres/problem.h" 36 #include "ceres/problem_impl.h" 37 #include "ceres/program.h" 38 #include "ceres/solver_impl.h" 39 #include "ceres/stringprintf.h" 40 #include "ceres/wall_time.h" 41 42 namespace ceres { 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 58 59 Solver::Options::~Options() { 60 delete linear_solver_ordering; 61 delete inner_iteration_ordering; 62 } 63 64 Solver::~Solver() {} 65 66 void Solver::Solve(const Solver::Options& options, 67 Problem* problem, 68 Solver::Summary* summary) { 69 double start_time_seconds = internal::WallTimeInSeconds(); 70 internal::ProblemImpl* problem_impl = 71 CHECK_NOTNULL(problem)->problem_impl_.get(); 72 internal::SolverImpl::Solve(options, problem_impl, summary); 73 summary->total_time_in_seconds = 74 internal::WallTimeInSeconds() - start_time_seconds; 75 } 76 77 void Solve(const Solver::Options& options, 78 Problem* problem, 79 Solver::Summary* summary) { 80 Solver solver; 81 solver.Solve(options, problem, summary); 82 } 83 84 Solver::Summary::Summary() 85 // Invalid values for most fields, to ensure that we are not 86 // accidentally reporting default values. 87 : minimizer_type(TRUST_REGION), 88 termination_type(DID_NOT_RUN), 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 num_inner_iteration_steps(-1), 95 preprocessor_time_in_seconds(-1.0), 96 minimizer_time_in_seconds(-1.0), 97 postprocessor_time_in_seconds(-1.0), 98 total_time_in_seconds(-1.0), 99 linear_solver_time_in_seconds(-1.0), 100 residual_evaluation_time_in_seconds(-1.0), 101 jacobian_evaluation_time_in_seconds(-1.0), 102 inner_iteration_time_in_seconds(-1.0), 103 num_parameter_blocks(-1), 104 num_parameters(-1), 105 num_effective_parameters(-1), 106 num_residual_blocks(-1), 107 num_residuals(-1), 108 num_parameter_blocks_reduced(-1), 109 num_parameters_reduced(-1), 110 num_effective_parameters_reduced(-1), 111 num_residual_blocks_reduced(-1), 112 num_residuals_reduced(-1), 113 num_threads_given(-1), 114 num_threads_used(-1), 115 num_linear_solver_threads_given(-1), 116 num_linear_solver_threads_used(-1), 117 linear_solver_type_given(SPARSE_NORMAL_CHOLESKY), 118 linear_solver_type_used(SPARSE_NORMAL_CHOLESKY), 119 inner_iterations_given(false), 120 inner_iterations_used(false), 121 preconditioner_type(IDENTITY), 122 trust_region_strategy_type(LEVENBERG_MARQUARDT), 123 dense_linear_algebra_library_type(EIGEN), 124 sparse_linear_algebra_library_type(SUITE_SPARSE), 125 line_search_direction_type(LBFGS), 126 line_search_type(ARMIJO) { 127 } 128 129 string Solver::Summary::BriefReport() const { 130 string report = "Ceres Solver Report: "; 131 if (termination_type == DID_NOT_RUN) { 132 CHECK(!error.empty()) 133 << "Solver terminated with DID_NOT_RUN but the solver did not " 134 << "return a reason. This is a Ceres error. Please report this " 135 << "to the Ceres team"; 136 return report + "Termination: DID_NOT_RUN, because " + error; 137 } 138 139 internal::StringAppendF(&report, "Iterations: %d", 140 num_successful_steps + num_unsuccessful_steps); 141 internal::StringAppendF(&report, ", Initial cost: %e", initial_cost); 142 143 // If the solver failed or was aborted, then the final_cost has no 144 // meaning. 145 if (termination_type != NUMERICAL_FAILURE && 146 termination_type != USER_ABORT) { 147 internal::StringAppendF(&report, ", Final cost: %e", final_cost); 148 } 149 150 internal::StringAppendF(&report, ", Termination: %s.", 151 SolverTerminationTypeToString(termination_type)); 152 return report; 153 }; 154 155 using internal::StringAppendF; 156 using internal::StringPrintf; 157 158 string Solver::Summary::FullReport() const { 159 string report = 160 "\n" 161 "Ceres Solver Report\n" 162 "-------------------\n"; 163 164 if (termination_type == DID_NOT_RUN) { 165 StringAppendF(&report, " Original\n"); 166 StringAppendF(&report, "Parameter blocks % 10d\n", num_parameter_blocks); 167 StringAppendF(&report, "Parameters % 10d\n", num_parameters); 168 if (num_effective_parameters != num_parameters) { 169 StringAppendF(&report, "Effective parameters% 10d\n", num_parameters); 170 } 171 172 StringAppendF(&report, "Residual blocks % 10d\n", 173 num_residual_blocks); 174 StringAppendF(&report, "Residuals % 10d\n\n", 175 num_residuals); 176 } else { 177 StringAppendF(&report, "%45s %21s\n", "Original", "Reduced"); 178 StringAppendF(&report, "Parameter blocks % 25d% 25d\n", 179 num_parameter_blocks, num_parameter_blocks_reduced); 180 StringAppendF(&report, "Parameters % 25d% 25d\n", 181 num_parameters, num_parameters_reduced); 182 if (num_effective_parameters_reduced != num_parameters_reduced) { 183 StringAppendF(&report, "Effective parameters% 25d% 25d\n", 184 num_effective_parameters, num_effective_parameters_reduced); 185 } 186 StringAppendF(&report, "Residual blocks % 25d% 25d\n", 187 num_residual_blocks, num_residual_blocks_reduced); 188 StringAppendF(&report, "Residual % 25d% 25d\n", 189 num_residuals, num_residuals_reduced); 190 } 191 192 if (minimizer_type == TRUST_REGION) { 193 // TRUST_SEARCH HEADER 194 StringAppendF(&report, "\nMinimizer %19s\n", 195 "TRUST_REGION"); 196 197 if (linear_solver_type_used == DENSE_NORMAL_CHOLESKY || 198 linear_solver_type_used == DENSE_SCHUR || 199 linear_solver_type_used == DENSE_QR) { 200 StringAppendF(&report, "\nDense linear algebra library %15s\n", 201 DenseLinearAlgebraLibraryTypeToString( 202 dense_linear_algebra_library_type)); 203 } 204 205 if (linear_solver_type_used == SPARSE_NORMAL_CHOLESKY || 206 linear_solver_type_used == SPARSE_SCHUR || 207 (linear_solver_type_used == ITERATIVE_SCHUR && 208 (preconditioner_type == CLUSTER_JACOBI || 209 preconditioner_type == CLUSTER_TRIDIAGONAL))) { 210 StringAppendF(&report, "\nSparse linear algebra library %15s\n", 211 SparseLinearAlgebraLibraryTypeToString( 212 sparse_linear_algebra_library_type)); 213 } 214 215 StringAppendF(&report, "Trust region strategy %19s", 216 TrustRegionStrategyTypeToString( 217 trust_region_strategy_type)); 218 if (trust_region_strategy_type == DOGLEG) { 219 if (dogleg_type == TRADITIONAL_DOGLEG) { 220 StringAppendF(&report, " (TRADITIONAL)"); 221 } else { 222 StringAppendF(&report, " (SUBSPACE)"); 223 } 224 } 225 StringAppendF(&report, "\n"); 226 StringAppendF(&report, "\n"); 227 228 StringAppendF(&report, "%45s %21s\n", "Given", "Used"); 229 StringAppendF(&report, "Linear solver %25s%25s\n", 230 LinearSolverTypeToString(linear_solver_type_given), 231 LinearSolverTypeToString(linear_solver_type_used)); 232 233 if (linear_solver_type_given == CGNR || 234 linear_solver_type_given == ITERATIVE_SCHUR) { 235 StringAppendF(&report, "Preconditioner %25s%25s\n", 236 PreconditionerTypeToString(preconditioner_type), 237 PreconditionerTypeToString(preconditioner_type)); 238 } 239 240 StringAppendF(&report, "Threads % 25d% 25d\n", 241 num_threads_given, num_threads_used); 242 StringAppendF(&report, "Linear solver threads % 23d% 25d\n", 243 num_linear_solver_threads_given, 244 num_linear_solver_threads_used); 245 246 if (IsSchurType(linear_solver_type_used)) { 247 string given; 248 StringifyOrdering(linear_solver_ordering_given, &given); 249 string used; 250 StringifyOrdering(linear_solver_ordering_used, &used); 251 StringAppendF(&report, 252 "Linear solver ordering %22s %24s\n", 253 given.c_str(), 254 used.c_str()); 255 } 256 257 if (inner_iterations_given) { 258 StringAppendF(&report, 259 "Use inner iterations %20s %20s\n", 260 inner_iterations_given ? "True" : "False", 261 inner_iterations_used ? "True" : "False"); 262 } 263 264 if (inner_iterations_used) { 265 string given; 266 StringifyOrdering(inner_iteration_ordering_given, &given); 267 string used; 268 StringifyOrdering(inner_iteration_ordering_used, &used); 269 StringAppendF(&report, 270 "Inner iteration ordering %20s %24s\n", 271 given.c_str(), 272 used.c_str()); 273 } 274 } else { 275 // LINE_SEARCH HEADER 276 StringAppendF(&report, "\nMinimizer %19s\n", "LINE_SEARCH"); 277 278 279 string line_search_direction_string; 280 if (line_search_direction_type == LBFGS) { 281 line_search_direction_string = StringPrintf("LBFGS (%d)", max_lbfgs_rank); 282 } else if (line_search_direction_type == NONLINEAR_CONJUGATE_GRADIENT) { 283 line_search_direction_string = 284 NonlinearConjugateGradientTypeToString( 285 nonlinear_conjugate_gradient_type); 286 } else { 287 line_search_direction_string = 288 LineSearchDirectionTypeToString(line_search_direction_type); 289 } 290 291 StringAppendF(&report, "Line search direction %19s\n", 292 line_search_direction_string.c_str()); 293 294 const string line_search_type_string = 295 StringPrintf("%s %s", 296 LineSearchInterpolationTypeToString( 297 line_search_interpolation_type), 298 LineSearchTypeToString(line_search_type)); 299 StringAppendF(&report, "Line search type %19s\n", 300 line_search_type_string.c_str()); 301 StringAppendF(&report, "\n"); 302 303 StringAppendF(&report, "%45s %21s\n", "Given", "Used"); 304 StringAppendF(&report, "Threads % 25d% 25d\n", 305 num_threads_given, num_threads_used); 306 } 307 308 if (termination_type == DID_NOT_RUN) { 309 CHECK(!error.empty()) 310 << "Solver terminated with DID_NOT_RUN but the solver did not " 311 << "return a reason. This is a Ceres error. Please report this " 312 << "to the Ceres team"; 313 StringAppendF(&report, "Termination: %20s\n", 314 "DID_NOT_RUN"); 315 StringAppendF(&report, "Reason: %s\n", error.c_str()); 316 return report; 317 } 318 319 StringAppendF(&report, "\nCost:\n"); 320 StringAppendF(&report, "Initial % 30e\n", initial_cost); 321 if (termination_type != NUMERICAL_FAILURE && 322 termination_type != USER_ABORT) { 323 StringAppendF(&report, "Final % 30e\n", final_cost); 324 StringAppendF(&report, "Change % 30e\n", 325 initial_cost - final_cost); 326 } 327 328 StringAppendF(&report, "\nMinimizer iterations % 16d\n", 329 num_successful_steps + num_unsuccessful_steps); 330 331 // Successful/Unsuccessful steps only matter in the case of the 332 // trust region solver. Line search terminates when it encounters 333 // the first unsuccessful step. 334 if (minimizer_type == TRUST_REGION) { 335 StringAppendF(&report, "Successful steps % 14d\n", 336 num_successful_steps); 337 StringAppendF(&report, "Unsuccessful steps % 14d\n", 338 num_unsuccessful_steps); 339 } 340 if (inner_iterations_used) { 341 StringAppendF(&report, "Steps with inner iterations % 14d\n", 342 num_inner_iteration_steps); 343 } 344 345 StringAppendF(&report, "\nTime (in seconds):\n"); 346 StringAppendF(&report, "Preprocessor %25.3f\n", 347 preprocessor_time_in_seconds); 348 349 StringAppendF(&report, "\n Residual evaluation %23.3f\n", 350 residual_evaluation_time_in_seconds); 351 StringAppendF(&report, " Jacobian evaluation %23.3f\n", 352 jacobian_evaluation_time_in_seconds); 353 354 if (minimizer_type == TRUST_REGION) { 355 StringAppendF(&report, " Linear solver %23.3f\n", 356 linear_solver_time_in_seconds); 357 } 358 359 if (inner_iterations_used) { 360 StringAppendF(&report, " Inner iterations %23.3f\n", 361 inner_iteration_time_in_seconds); 362 } 363 364 StringAppendF(&report, "Minimizer %25.3f\n\n", 365 minimizer_time_in_seconds); 366 367 StringAppendF(&report, "Postprocessor %24.3f\n", 368 postprocessor_time_in_seconds); 369 370 StringAppendF(&report, "Total %25.3f\n\n", 371 total_time_in_seconds); 372 373 StringAppendF(&report, "Termination: %25s\n", 374 SolverTerminationTypeToString(termination_type)); 375 return report; 376 }; 377 378 } // namespace ceres 379