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 | |
| 59 | Solver::~Solver() {} |
| 60 | |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 61 | void Solver::Solve(const Solver::Options& options, |
| 62 | Problem* problem, |
| 63 | Solver::Summary* summary) { |
Petter Strandmark | 76533b3 | 2012-09-04 22:08:50 -0700 | [diff] [blame] | 64 | double start_time_seconds = internal::WallTimeInSeconds(); |
Keir Mierle | 6196cba | 2012-06-18 11:03:40 -0700 | [diff] [blame] | 65 | internal::ProblemImpl* problem_impl = |
| 66 | CHECK_NOTNULL(problem)->problem_impl_.get(); |
| 67 | internal::SolverImpl::Solve(options, problem_impl, summary); |
Petter Strandmark | 76533b3 | 2012-09-04 22:08:50 -0700 | [diff] [blame] | 68 | summary->total_time_in_seconds = |
| 69 | internal::WallTimeInSeconds() - start_time_seconds; |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 70 | } |
| 71 | |
| 72 | void Solve(const Solver::Options& options, |
| 73 | Problem* problem, |
| 74 | Solver::Summary* summary) { |
Keir Mierle | 6196cba | 2012-06-18 11:03:40 -0700 | [diff] [blame] | 75 | Solver solver; |
| 76 | solver.Solve(options, problem, summary); |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 77 | } |
| 78 | |
| 79 | Solver::Summary::Summary() |
| 80 | // Invalid values for most fields, to ensure that we are not |
| 81 | // accidentally reporting default values. |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 82 | : minimizer_type(TRUST_REGION), |
Sameer Agarwal | dcee120 | 2013-12-07 21:48:56 -0800 | [diff] [blame] | 83 | termination_type(FAILURE), |
| 84 | message("ceres::Solve was not called."), |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 85 | initial_cost(-1.0), |
| 86 | final_cost(-1.0), |
| 87 | fixed_cost(-1.0), |
| 88 | num_successful_steps(-1), |
| 89 | num_unsuccessful_steps(-1), |
Sameer Agarwal | 9f9488b | 2013-05-23 09:40:40 -0700 | [diff] [blame] | 90 | num_inner_iteration_steps(-1), |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 91 | preprocessor_time_in_seconds(-1.0), |
| 92 | minimizer_time_in_seconds(-1.0), |
Sameer Agarwal | fa01519 | 2012-06-11 14:21:42 -0700 | [diff] [blame] | 93 | postprocessor_time_in_seconds(-1.0), |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 94 | total_time_in_seconds(-1.0), |
Sameer Agarwal | 42a84b8 | 2013-02-01 12:22:53 -0800 | [diff] [blame] | 95 | linear_solver_time_in_seconds(-1.0), |
| 96 | residual_evaluation_time_in_seconds(-1.0), |
| 97 | jacobian_evaluation_time_in_seconds(-1.0), |
Sameer Agarwal | 9f9488b | 2013-05-23 09:40:40 -0700 | [diff] [blame] | 98 | inner_iteration_time_in_seconds(-1.0), |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 99 | num_parameter_blocks(-1), |
| 100 | num_parameters(-1), |
Keir Mierle | ba94421 | 2013-02-25 12:46:44 -0800 | [diff] [blame] | 101 | num_effective_parameters(-1), |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 102 | num_residual_blocks(-1), |
| 103 | num_residuals(-1), |
| 104 | num_parameter_blocks_reduced(-1), |
| 105 | num_parameters_reduced(-1), |
Keir Mierle | ba94421 | 2013-02-25 12:46:44 -0800 | [diff] [blame] | 106 | num_effective_parameters_reduced(-1), |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 107 | num_residual_blocks_reduced(-1), |
| 108 | num_residuals_reduced(-1), |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 109 | num_threads_given(-1), |
| 110 | num_threads_used(-1), |
| 111 | num_linear_solver_threads_given(-1), |
| 112 | num_linear_solver_threads_used(-1), |
| 113 | linear_solver_type_given(SPARSE_NORMAL_CHOLESKY), |
| 114 | linear_solver_type_used(SPARSE_NORMAL_CHOLESKY), |
Sameer Agarwal | 9f9488b | 2013-05-23 09:40:40 -0700 | [diff] [blame] | 115 | inner_iterations_given(false), |
| 116 | inner_iterations_used(false), |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 117 | preconditioner_type(IDENTITY), |
Sameer Agarwal | f06b9fa | 2013-10-27 21:38:13 -0700 | [diff] [blame] | 118 | visibility_clustering_type(CANONICAL_VIEWS), |
Sameer Agarwal | 97fb6d9 | 2012-06-17 10:08:19 -0700 | [diff] [blame] | 119 | trust_region_strategy_type(LEVENBERG_MARQUARDT), |
Sameer Agarwal | 367b65e | 2013-08-09 10:35:37 -0700 | [diff] [blame] | 120 | dense_linear_algebra_library_type(EIGEN), |
| 121 | sparse_linear_algebra_library_type(SUITE_SPARSE), |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 122 | line_search_direction_type(LBFGS), |
Alex Stewart | 7124c34 | 2013-11-07 16:10:02 +0000 | [diff] [blame] | 123 | line_search_type(ARMIJO), |
| 124 | line_search_interpolation_type(BISECTION), |
| 125 | nonlinear_conjugate_gradient_type(FLETCHER_REEVES), |
| 126 | max_lbfgs_rank(-1) { |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 127 | } |
| 128 | |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 129 | using internal::StringAppendF; |
Sameer Agarwal | 4f010b2 | 2013-07-01 08:01:01 -0700 | [diff] [blame] | 130 | using internal::StringPrintf; |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 131 | |
Sameer Agarwal | dcee120 | 2013-12-07 21:48:56 -0800 | [diff] [blame] | 132 | string Solver::Summary::BriefReport() const { |
| 133 | return StringPrintf("Ceres Solver Report: " |
| 134 | "Iterations: %d, " |
| 135 | "Initial cost: %e, " |
| 136 | "Final cost: %e, " |
| 137 | "Termination: %s", |
| 138 | num_successful_steps + num_unsuccessful_steps, |
| 139 | initial_cost, |
| 140 | final_cost, |
| 141 | TerminationTypeToString(termination_type)); |
| 142 | }; |
| 143 | |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 144 | string Solver::Summary::FullReport() const { |
| 145 | string report = |
| 146 | "\n" |
| 147 | "Ceres Solver Report\n" |
| 148 | "-------------------\n"; |
| 149 | |
Sameer Agarwal | dcee120 | 2013-12-07 21:48:56 -0800 | [diff] [blame] | 150 | StringAppendF(&report, "%45s %21s\n", "Original", "Reduced"); |
| 151 | StringAppendF(&report, "Parameter blocks % 25d% 25d\n", |
| 152 | num_parameter_blocks, num_parameter_blocks_reduced); |
| 153 | StringAppendF(&report, "Parameters % 25d% 25d\n", |
| 154 | num_parameters, num_parameters_reduced); |
| 155 | if (num_effective_parameters_reduced != num_parameters_reduced) { |
| 156 | StringAppendF(&report, "Effective parameters% 25d% 25d\n", |
| 157 | num_effective_parameters, num_effective_parameters_reduced); |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 158 | } |
Sameer Agarwal | dcee120 | 2013-12-07 21:48:56 -0800 | [diff] [blame] | 159 | StringAppendF(&report, "Residual blocks % 25d% 25d\n", |
| 160 | num_residual_blocks, num_residual_blocks_reduced); |
| 161 | StringAppendF(&report, "Residual % 25d% 25d\n", |
| 162 | num_residuals, num_residuals_reduced); |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 163 | |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 164 | if (minimizer_type == TRUST_REGION) { |
Sameer Agarwal | 9f9488b | 2013-05-23 09:40:40 -0700 | [diff] [blame] | 165 | // TRUST_SEARCH HEADER |
Sameer Agarwal | 509f68c | 2013-02-20 01:39:03 -0800 | [diff] [blame] | 166 | StringAppendF(&report, "\nMinimizer %19s\n", |
| 167 | "TRUST_REGION"); |
Sameer Agarwal | 367b65e | 2013-08-09 10:35:37 -0700 | [diff] [blame] | 168 | |
| 169 | if (linear_solver_type_used == DENSE_NORMAL_CHOLESKY || |
| 170 | linear_solver_type_used == DENSE_SCHUR || |
| 171 | linear_solver_type_used == DENSE_QR) { |
| 172 | StringAppendF(&report, "\nDense linear algebra library %15s\n", |
| 173 | DenseLinearAlgebraLibraryTypeToString( |
| 174 | dense_linear_algebra_library_type)); |
| 175 | } |
| 176 | |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 177 | if (linear_solver_type_used == SPARSE_NORMAL_CHOLESKY || |
| 178 | linear_solver_type_used == SPARSE_SCHUR || |
| 179 | (linear_solver_type_used == ITERATIVE_SCHUR && |
Sameer Agarwal | 290b975 | 2013-02-17 16:50:37 -0800 | [diff] [blame] | 180 | (preconditioner_type == CLUSTER_JACOBI || |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 181 | preconditioner_type == CLUSTER_TRIDIAGONAL))) { |
Sameer Agarwal | 9f9488b | 2013-05-23 09:40:40 -0700 | [diff] [blame] | 182 | StringAppendF(&report, "\nSparse linear algebra library %15s\n", |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 183 | SparseLinearAlgebraLibraryTypeToString( |
Sameer Agarwal | 367b65e | 2013-08-09 10:35:37 -0700 | [diff] [blame] | 184 | sparse_linear_algebra_library_type)); |
Sameer Agarwal | 1e28920 | 2012-08-21 18:00:54 -0700 | [diff] [blame] | 185 | } |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 186 | |
Sameer Agarwal | 9f9488b | 2013-05-23 09:40:40 -0700 | [diff] [blame] | 187 | StringAppendF(&report, "Trust region strategy %19s", |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 188 | TrustRegionStrategyTypeToString( |
| 189 | trust_region_strategy_type)); |
| 190 | if (trust_region_strategy_type == DOGLEG) { |
| 191 | if (dogleg_type == TRADITIONAL_DOGLEG) { |
| 192 | StringAppendF(&report, " (TRADITIONAL)"); |
| 193 | } else { |
| 194 | StringAppendF(&report, " (SUBSPACE)"); |
| 195 | } |
| 196 | } |
| 197 | StringAppendF(&report, "\n"); |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 198 | StringAppendF(&report, "\n"); |
Sameer Agarwal | 1e28920 | 2012-08-21 18:00:54 -0700 | [diff] [blame] | 199 | |
Sameer Agarwal | 931c309 | 2013-02-25 09:46:21 -0800 | [diff] [blame] | 200 | StringAppendF(&report, "%45s %21s\n", "Given", "Used"); |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 201 | StringAppendF(&report, "Linear solver %25s%25s\n", |
| 202 | LinearSolverTypeToString(linear_solver_type_given), |
| 203 | LinearSolverTypeToString(linear_solver_type_used)); |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 204 | |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 205 | if (linear_solver_type_given == CGNR || |
| 206 | linear_solver_type_given == ITERATIVE_SCHUR) { |
| 207 | StringAppendF(&report, "Preconditioner %25s%25s\n", |
| 208 | PreconditionerTypeToString(preconditioner_type), |
| 209 | PreconditionerTypeToString(preconditioner_type)); |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 210 | } |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 211 | |
Sameer Agarwal | f06b9fa | 2013-10-27 21:38:13 -0700 | [diff] [blame] | 212 | if (preconditioner_type == CLUSTER_JACOBI || |
| 213 | preconditioner_type == CLUSTER_TRIDIAGONAL) { |
| 214 | StringAppendF(&report, "Visibility clustering%24s%25s\n", |
Sameer Agarwal | 9ba0b35 | 2013-11-05 13:04:56 -0800 | [diff] [blame] | 215 | VisibilityClusteringTypeToString( |
| 216 | visibility_clustering_type), |
| 217 | VisibilityClusteringTypeToString( |
| 218 | visibility_clustering_type)); |
Sameer Agarwal | f06b9fa | 2013-10-27 21:38:13 -0700 | [diff] [blame] | 219 | } |
Sameer Agarwal | 9f9488b | 2013-05-23 09:40:40 -0700 | [diff] [blame] | 220 | StringAppendF(&report, "Threads % 25d% 25d\n", |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 221 | num_threads_given, num_threads_used); |
| 222 | StringAppendF(&report, "Linear solver threads % 23d% 25d\n", |
| 223 | num_linear_solver_threads_given, |
| 224 | num_linear_solver_threads_used); |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 225 | |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 226 | if (IsSchurType(linear_solver_type_used)) { |
| 227 | string given; |
| 228 | StringifyOrdering(linear_solver_ordering_given, &given); |
| 229 | string used; |
| 230 | StringifyOrdering(linear_solver_ordering_used, &used); |
| 231 | StringAppendF(&report, |
| 232 | "Linear solver ordering %22s %24s\n", |
| 233 | given.c_str(), |
| 234 | used.c_str()); |
| 235 | } |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 236 | |
Sameer Agarwal | 9f9488b | 2013-05-23 09:40:40 -0700 | [diff] [blame] | 237 | if (inner_iterations_given) { |
| 238 | StringAppendF(&report, |
| 239 | "Use inner iterations %20s %20s\n", |
| 240 | inner_iterations_given ? "True" : "False", |
| 241 | inner_iterations_used ? "True" : "False"); |
| 242 | } |
| 243 | |
| 244 | if (inner_iterations_used) { |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 245 | string given; |
| 246 | StringifyOrdering(inner_iteration_ordering_given, &given); |
| 247 | string used; |
| 248 | StringifyOrdering(inner_iteration_ordering_used, &used); |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 249 | StringAppendF(&report, |
| 250 | "Inner iteration ordering %20s %24s\n", |
| 251 | given.c_str(), |
| 252 | used.c_str()); |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 253 | } |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 254 | } else { |
Sameer Agarwal | 9f9488b | 2013-05-23 09:40:40 -0700 | [diff] [blame] | 255 | // LINE_SEARCH HEADER |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 256 | StringAppendF(&report, "\nMinimizer %19s\n", "LINE_SEARCH"); |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 257 | |
Sameer Agarwal | 4f010b2 | 2013-07-01 08:01:01 -0700 | [diff] [blame] | 258 | |
| 259 | string line_search_direction_string; |
| 260 | if (line_search_direction_type == LBFGS) { |
| 261 | line_search_direction_string = StringPrintf("LBFGS (%d)", max_lbfgs_rank); |
| 262 | } else if (line_search_direction_type == NONLINEAR_CONJUGATE_GRADIENT) { |
| 263 | line_search_direction_string = |
| 264 | NonlinearConjugateGradientTypeToString( |
| 265 | nonlinear_conjugate_gradient_type); |
| 266 | } else { |
| 267 | line_search_direction_string = |
| 268 | LineSearchDirectionTypeToString(line_search_direction_type); |
| 269 | } |
| 270 | |
| 271 | StringAppendF(&report, "Line search direction %19s\n", |
| 272 | line_search_direction_string.c_str()); |
| 273 | |
| 274 | const string line_search_type_string = |
| 275 | StringPrintf("%s %s", |
Sameer Agarwal | 7a8f797 | 2013-07-03 09:03:55 -0700 | [diff] [blame] | 276 | LineSearchInterpolationTypeToString( |
| 277 | line_search_interpolation_type), |
Sameer Agarwal | 4f010b2 | 2013-07-01 08:01:01 -0700 | [diff] [blame] | 278 | LineSearchTypeToString(line_search_type)); |
| 279 | StringAppendF(&report, "Line search type %19s\n", |
| 280 | line_search_type_string.c_str()); |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 281 | StringAppendF(&report, "\n"); |
| 282 | |
Sameer Agarwal | beb4505 | 2013-02-22 13:37:01 -0800 | [diff] [blame] | 283 | StringAppendF(&report, "%45s %21s\n", "Given", "Used"); |
Sameer Agarwal | 9f9488b | 2013-05-23 09:40:40 -0700 | [diff] [blame] | 284 | StringAppendF(&report, "Threads % 25d% 25d\n", |
Sameer Agarwal | d010de5 | 2013-02-15 14:26:56 -0800 | [diff] [blame] | 285 | num_threads_given, num_threads_used); |
Sameer Agarwal | 977be7c | 2013-01-26 16:01:54 -0800 | [diff] [blame] | 286 | } |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 287 | |
Sameer Agarwal | 9f9488b | 2013-05-23 09:40:40 -0700 | [diff] [blame] | 288 | StringAppendF(&report, "\nCost:\n"); |
| 289 | StringAppendF(&report, "Initial % 30e\n", initial_cost); |
Sameer Agarwal | dcee120 | 2013-12-07 21:48:56 -0800 | [diff] [blame] | 290 | if (termination_type != FAILURE && |
| 291 | termination_type != USER_FAILURE) { |
Sameer Agarwal | 9f9488b | 2013-05-23 09:40:40 -0700 | [diff] [blame] | 292 | StringAppendF(&report, "Final % 30e\n", final_cost); |
| 293 | StringAppendF(&report, "Change % 30e\n", |
| 294 | initial_cost - final_cost); |
| 295 | } |
| 296 | |
| 297 | StringAppendF(&report, "\nMinimizer iterations % 16d\n", |
| 298 | num_successful_steps + num_unsuccessful_steps); |
Sameer Agarwal | 4f010b2 | 2013-07-01 08:01:01 -0700 | [diff] [blame] | 299 | |
| 300 | // Successful/Unsuccessful steps only matter in the case of the |
| 301 | // trust region solver. Line search terminates when it encounters |
| 302 | // the first unsuccessful step. |
| 303 | if (minimizer_type == TRUST_REGION) { |
| 304 | StringAppendF(&report, "Successful steps % 14d\n", |
| 305 | num_successful_steps); |
| 306 | StringAppendF(&report, "Unsuccessful steps % 14d\n", |
| 307 | num_unsuccessful_steps); |
| 308 | } |
Sameer Agarwal | 9f9488b | 2013-05-23 09:40:40 -0700 | [diff] [blame] | 309 | if (inner_iterations_used) { |
| 310 | StringAppendF(&report, "Steps with inner iterations % 14d\n", |
| 311 | num_inner_iteration_steps); |
| 312 | } |
| 313 | |
| 314 | StringAppendF(&report, "\nTime (in seconds):\n"); |
| 315 | StringAppendF(&report, "Preprocessor %25.3f\n", |
| 316 | preprocessor_time_in_seconds); |
| 317 | |
| 318 | StringAppendF(&report, "\n Residual evaluation %23.3f\n", |
| 319 | residual_evaluation_time_in_seconds); |
| 320 | StringAppendF(&report, " Jacobian evaluation %23.3f\n", |
| 321 | jacobian_evaluation_time_in_seconds); |
| 322 | |
| 323 | if (minimizer_type == TRUST_REGION) { |
| 324 | StringAppendF(&report, " Linear solver %23.3f\n", |
| 325 | linear_solver_time_in_seconds); |
| 326 | } |
| 327 | |
| 328 | if (inner_iterations_used) { |
| 329 | StringAppendF(&report, " Inner iterations %23.3f\n", |
| 330 | inner_iteration_time_in_seconds); |
| 331 | } |
| 332 | |
| 333 | StringAppendF(&report, "Minimizer %25.3f\n\n", |
| 334 | minimizer_time_in_seconds); |
| 335 | |
| 336 | StringAppendF(&report, "Postprocessor %24.3f\n", |
| 337 | postprocessor_time_in_seconds); |
| 338 | |
| 339 | StringAppendF(&report, "Total %25.3f\n\n", |
| 340 | total_time_in_seconds); |
| 341 | |
Sameer Agarwal | 1da9292 | 2014-02-23 16:10:33 -0800 | [diff] [blame] | 342 | StringAppendF(&report, "Termination: %25s (%s)\n", |
| 343 | TerminationTypeToString(termination_type), message.c_str()); |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 344 | return report; |
| 345 | }; |
| 346 | |
Sameer Agarwal | dcee120 | 2013-12-07 21:48:56 -0800 | [diff] [blame] | 347 | bool Solver::Summary::IsSolutionUsable() const { |
| 348 | return (termination_type == CONVERGENCE || |
| 349 | termination_type == NO_CONVERGENCE || |
| 350 | termination_type == USER_SUCCESS); |
| 351 | } |
| 352 | |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 353 | } // namespace ceres |