blob: 4fb6e4ddca3f04b7bb88f4d38f2ab035a77b8f90 [file] [log] [blame]
Keir Mierle8ebb0732012-04-30 23:09:08 -07001// 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: sameeragarwal@google.com (Sameer Agarwal)
30
31#ifndef CERES_PUBLIC_SOLVER_H_
32#define CERES_PUBLIC_SOLVER_H_
33
34#include <cmath>
35#include <string>
36#include <vector>
Sameer Agarwal4997cbc2012-07-02 12:44:34 -070037#include "ceres/crs_matrix.h"
Keir Mierle8ebb0732012-04-30 23:09:08 -070038#include "ceres/internal/macros.h"
39#include "ceres/internal/port.h"
Sameer Agarwal4997cbc2012-07-02 12:44:34 -070040#include "ceres/iteration_callback.h"
Sameer Agarwal2c94eed2012-10-01 16:34:37 -070041#include "ceres/ordered_groups.h"
Keir Mierle8ebb0732012-04-30 23:09:08 -070042#include "ceres/types.h"
43
44namespace ceres {
45
46class Problem;
47
48// Interface for non-linear least squares solvers.
Björn Piltz5d7eed82014-04-23 22:13:37 +020049class CERES_EXPORT Solver {
Keir Mierle8ebb0732012-04-30 23:09:08 -070050 public:
51 virtual ~Solver();
52
53 // The options structure contains, not surprisingly, options that control how
54 // the solver operates. The defaults should be suitable for a wide range of
55 // problems; however, better performance is often obtainable with tweaking.
56 //
57 // The constants are defined inside types.h
Björn Piltz5d7eed82014-04-23 22:13:37 +020058 struct CERES_EXPORT Options {
Keir Mierle8ebb0732012-04-30 23:09:08 -070059 // Default constructor that sets up a generic sparse problem.
60 Options() {
Sameer Agarwalf4d01642012-11-26 12:55:58 -080061 minimizer_type = TRUST_REGION;
Sameer Agarwald010de52013-02-15 14:26:56 -080062 line_search_direction_type = LBFGS;
Alex Stewart9aa0e3c2013-07-05 20:22:37 +010063 line_search_type = WOLFE;
Sameer Agarwalf4d01642012-11-26 12:55:58 -080064 nonlinear_conjugate_gradient_type = FLETCHER_REEVES;
Sameer Agarwalaed99612012-11-29 10:33:19 -080065 max_lbfgs_rank = 20;
Alex Stewart9aa0e3c2013-07-05 20:22:37 +010066 use_approximate_eigenvalue_bfgs_scaling = false;
Sameer Agarwal09244012013-06-30 14:33:23 -070067 line_search_interpolation_type = CUBIC;
68 min_line_search_step_size = 1e-9;
Alex Stewart9aa0e3c2013-07-05 20:22:37 +010069 line_search_sufficient_function_decrease = 1e-4;
70 max_line_search_step_contraction = 1e-3;
71 min_line_search_step_contraction = 0.6;
72 max_num_line_search_step_size_iterations = 20;
73 max_num_line_search_direction_restarts = 5;
74 line_search_sufficient_curvature_decrease = 0.9;
75 max_line_search_step_expansion = 10.0;
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -070076 trust_region_strategy_type = LEVENBERG_MARQUARDT;
Markus Moll51cf7cb2012-08-20 20:10:20 +020077 dogleg_type = TRADITIONAL_DOGLEG;
Sameer Agarwala8f87d72012-08-08 10:38:31 -070078 use_nonmonotonic_steps = false;
79 max_consecutive_nonmonotonic_steps = 5;
Keir Mierle8ebb0732012-04-30 23:09:08 -070080 max_num_iterations = 50;
Sameer Agarwalfa015192012-06-11 14:21:42 -070081 max_solver_time_in_seconds = 1e9;
Keir Mierle8ebb0732012-04-30 23:09:08 -070082 num_threads = 1;
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -070083 initial_trust_region_radius = 1e4;
84 max_trust_region_radius = 1e16;
85 min_trust_region_radius = 1e-32;
Keir Mierle8ebb0732012-04-30 23:09:08 -070086 min_relative_decrease = 1e-3;
Sameer Agarwaleeedd2e2013-07-07 23:04:31 -070087 min_lm_diagonal = 1e-6;
88 max_lm_diagonal = 1e32;
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -070089 max_num_consecutive_invalid_steps = 5;
Keir Mierle8ebb0732012-04-30 23:09:08 -070090 function_tolerance = 1e-6;
91 gradient_tolerance = 1e-10;
92 parameter_tolerance = 1e-8;
Sameer Agarwalb0518732012-05-29 00:27:57 -070093
94#if defined(CERES_NO_SUITESPARSE) && defined(CERES_NO_CXSPARSE)
Keir Mierle8ebb0732012-04-30 23:09:08 -070095 linear_solver_type = DENSE_QR;
Sameer Agarwalb0518732012-05-29 00:27:57 -070096#else
97 linear_solver_type = SPARSE_NORMAL_CHOLESKY;
98#endif
99
Sameer Agarwal97fb6d92012-06-17 10:08:19 -0700100 preconditioner_type = JACOBI;
Sameer Agarwalf06b9fa2013-10-27 21:38:13 -0700101 visibility_clustering_type = CANONICAL_VIEWS;
Sameer Agarwald61b68a2013-08-16 17:02:56 -0700102 dense_linear_algebra_library_type = EIGEN;
Sameer Agarwal367b65e2013-08-09 10:35:37 -0700103 sparse_linear_algebra_library_type = SUITE_SPARSE;
Sameer Agarwalb0518732012-05-29 00:27:57 -0700104#if defined(CERES_NO_SUITESPARSE) && !defined(CERES_NO_CXSPARSE)
Sameer Agarwal367b65e2013-08-09 10:35:37 -0700105 sparse_linear_algebra_library_type = CX_SPARSE;
Sameer Agarwalb0518732012-05-29 00:27:57 -0700106#endif
107
Sameer Agarwald61b68a2013-08-16 17:02:56 -0700108
Sameer Agarwal97fb6d92012-06-17 10:08:19 -0700109 num_linear_solver_threads = 1;
Sameer Agarwal9189f4e2013-04-19 17:09:49 -0700110 use_postordering = false;
Richard Stebbing32530782014-04-26 07:42:23 +0100111 dynamic_sparsity = false;
Sameer Agarwaleeedd2e2013-07-07 23:04:31 -0700112 min_linear_solver_iterations = 1;
113 max_linear_solver_iterations = 500;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700114 eta = 1e-1;
115 jacobi_scaling = true;
Sameer Agarwal09396322013-05-28 22:29:36 -0700116 use_inner_iterations = false;
117 inner_iteration_tolerance = 1e-3;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700118 logging_type = PER_MINIMIZER_ITERATION;
119 minimizer_progress_to_stdout = false;
Sameer Agarwalc4a32912013-06-13 22:00:48 -0700120 trust_region_problem_dump_directory = "/tmp";
121 trust_region_problem_dump_format_type = TEXTFILE;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700122 check_gradients = false;
123 gradient_check_relative_precision = 1e-8;
124 numeric_derivative_relative_step_size = 1e-6;
125 update_state_every_iteration = false;
126 }
127
128 // Minimizer options ----------------------------------------
129
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800130 // Ceres supports the two major families of optimization strategies -
131 // Trust Region and Line Search.
132 //
133 // 1. The line search approach first finds a descent direction
134 // along which the objective function will be reduced and then
135 // computes a step size that decides how far should move along
136 // that direction. The descent direction can be computed by
137 // various methods, such as gradient descent, Newton's method and
138 // Quasi-Newton method. The step size can be determined either
139 // exactly or inexactly.
140 //
141 // 2. The trust region approach approximates the objective
142 // function using using a model function (often a quadratic) over
143 // a subset of the search space known as the trust region. If the
144 // model function succeeds in minimizing the true objective
145 // function the trust region is expanded; conversely, otherwise it
146 // is contracted and the model optimization problem is solved
147 // again.
148 //
149 // Trust region methods are in some sense dual to line search methods:
150 // trust region methods first choose a step size (the size of the
151 // trust region) and then a step direction while line search methods
152 // first choose a step direction and then a step size.
153 MinimizerType minimizer_type;
154
155 LineSearchDirectionType line_search_direction_type;
156 LineSearchType line_search_type;
157 NonlinearConjugateGradientType nonlinear_conjugate_gradient_type;
158
Sameer Agarwalaed99612012-11-29 10:33:19 -0800159 // The LBFGS hessian approximation is a low rank approximation to
160 // the inverse of the Hessian matrix. The rank of the
161 // approximation determines (linearly) the space and time
162 // complexity of using the approximation. Higher the rank, the
Sameer Agarwal2293cb52012-11-29 16:00:18 -0800163 // better is the quality of the approximation. The increase in
164 // quality is however is bounded for a number of reasons.
165 //
166 // 1. The method only uses secant information and not actual
167 // derivatives.
168 //
169 // 2. The Hessian approximation is constrained to be positive
170 // definite.
171 //
172 // So increasing this rank to a large number will cost time and
173 // space complexity without the corresponding increase in solution
174 // quality. There are no hard and fast rules for choosing the
175 // maximum rank. The best choice usually requires some problem
176 // specific experimentation.
177 //
178 // For more theoretical and implementation details of the LBFGS
179 // method, please see:
Sameer Agarwalaed99612012-11-29 10:33:19 -0800180 //
181 // Nocedal, J. (1980). "Updating Quasi-Newton Matrices with
182 // Limited Storage". Mathematics of Computation 35 (151): 773–782.
183 int max_lbfgs_rank;
184
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100185 // As part of the (L)BFGS update step (BFGS) / right-multiply step (L-BFGS),
186 // the initial inverse Hessian approximation is taken to be the Identity.
187 // However, Oren showed that using instead I * \gamma, where \gamma is
188 // chosen to approximate an eigenvalue of the true inverse Hessian can
189 // result in improved convergence in a wide variety of cases. Setting
190 // use_approximate_eigenvalue_bfgs_scaling to true enables this scaling.
191 //
192 // It is important to note that approximate eigenvalue scaling does not
193 // always improve convergence, and that it can in fact significantly degrade
194 // performance for certain classes of problem, which is why it is disabled
195 // by default. In particular it can degrade performance when the
196 // sensitivity of the problem to different parameters varies significantly,
197 // as in this case a single scalar factor fails to capture this variation
198 // and detrimentally downscales parts of the jacobian approximation which
199 // correspond to low-sensitivity parameters. It can also reduce the
200 // robustness of the solution to errors in the jacobians.
201 //
202 // Oren S.S., Self-scaling variable metric (SSVM) algorithms
203 // Part II: Implementation and experiments, Management Science,
204 // 20(5), 863-874, 1974.
205 bool use_approximate_eigenvalue_bfgs_scaling;
206
Sameer Agarwal09244012013-06-30 14:33:23 -0700207 // Degree of the polynomial used to approximate the objective
208 // function. Valid values are BISECTION, QUADRATIC and CUBIC.
209 //
210 // BISECTION corresponds to pure backtracking search with no
211 // interpolation.
212 LineSearchInterpolationType line_search_interpolation_type;
213
214 // If during the line search, the step_size falls below this
215 // value, it is truncated to zero.
216 double min_line_search_step_size;
217
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100218 // Line search parameters.
Sameer Agarwal09244012013-06-30 14:33:23 -0700219
220 // Solving the line search problem exactly is computationally
221 // prohibitive. Fortunately, line search based optimization
222 // algorithms can still guarantee convergence if instead of an
223 // exact solution, the line search algorithm returns a solution
224 // which decreases the value of the objective function
225 // sufficiently. More precisely, we are looking for a step_size
226 // s.t.
227 //
228 // f(step_size) <= f(0) + sufficient_decrease * f'(0) * step_size
229 //
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100230 double line_search_sufficient_function_decrease;
Sameer Agarwal09244012013-06-30 14:33:23 -0700231
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100232 // In each iteration of the line search,
Sameer Agarwal09244012013-06-30 14:33:23 -0700233 //
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100234 // new_step_size >= max_line_search_step_contraction * step_size
Sameer Agarwal09244012013-06-30 14:33:23 -0700235 //
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100236 // Note that by definition, for contraction:
237 //
238 // 0 < max_step_contraction < min_step_contraction < 1
239 //
240 double max_line_search_step_contraction;
Sameer Agarwal09244012013-06-30 14:33:23 -0700241
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100242 // In each iteration of the line search,
Sameer Agarwal09244012013-06-30 14:33:23 -0700243 //
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100244 // new_step_size <= min_line_search_step_contraction * step_size
Sameer Agarwal09244012013-06-30 14:33:23 -0700245 //
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100246 // Note that by definition, for contraction:
247 //
248 // 0 < max_step_contraction < min_step_contraction < 1
249 //
250 double min_line_search_step_contraction;
251
252 // Maximum number of trial step size iterations during each line search,
253 // if a step size satisfying the search conditions cannot be found within
254 // this number of trials, the line search will terminate.
255 int max_num_line_search_step_size_iterations;
256
257 // Maximum number of restarts of the line search direction algorithm before
258 // terminating the optimization. Restarts of the line search direction
259 // algorithm occur when the current algorithm fails to produce a new descent
260 // direction. This typically indicates a numerical failure, or a breakdown
261 // in the validity of the approximations used.
262 int max_num_line_search_direction_restarts;
263
264 // The strong Wolfe conditions consist of the Armijo sufficient
265 // decrease condition, and an additional requirement that the
266 // step-size be chosen s.t. the _magnitude_ ('strong' Wolfe
267 // conditions) of the gradient along the search direction
268 // decreases sufficiently. Precisely, this second condition
269 // is that we seek a step_size s.t.
270 //
271 // |f'(step_size)| <= sufficient_curvature_decrease * |f'(0)|
272 //
273 // Where f() is the line search objective and f'() is the derivative
274 // of f w.r.t step_size (d f / d step_size).
275 double line_search_sufficient_curvature_decrease;
276
277 // During the bracketing phase of the Wolfe search, the step size is
278 // increased until either a point satisfying the Wolfe conditions is
279 // found, or an upper bound for a bracket containing a point satisfying
280 // the conditions is found. Precisely, at each iteration of the
281 // expansion:
282 //
283 // new_step_size <= max_step_expansion * step_size.
284 //
285 // By definition for expansion, max_step_expansion > 1.0.
286 double max_line_search_step_expansion;
Sameer Agarwal09244012013-06-30 14:33:23 -0700287
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700288 TrustRegionStrategyType trust_region_strategy_type;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700289
Markus Moll51cf7cb2012-08-20 20:10:20 +0200290 // Type of dogleg strategy to use.
291 DoglegType dogleg_type;
292
Sameer Agarwala8f87d72012-08-08 10:38:31 -0700293 // The classical trust region methods are descent methods, in that
294 // they only accept a point if it strictly reduces the value of
295 // the objective function.
296 //
297 // Relaxing this requirement allows the algorithm to be more
298 // efficient in the long term at the cost of some local increase
299 // in the value of the objective function.
300 //
301 // This is because allowing for non-decreasing objective function
302 // values in a princpled manner allows the algorithm to "jump over
303 // boulders" as the method is not restricted to move into narrow
304 // valleys while preserving its convergence properties.
305 //
306 // Setting use_nonmonotonic_steps to true enables the
307 // non-monotonic trust region algorithm as described by Conn,
308 // Gould & Toint in "Trust Region Methods", Section 10.1.
309 //
310 // The parameter max_consecutive_nonmonotonic_steps controls the
311 // window size used by the step selection algorithm to accept
312 // non-monotonic steps.
313 //
314 // Even though the value of the objective function may be larger
315 // than the minimum value encountered over the course of the
316 // optimization, the final parameters returned to the user are the
317 // ones corresponding to the minimum cost over all iterations.
318 bool use_nonmonotonic_steps;
319 int max_consecutive_nonmonotonic_steps;
320
Keir Mierle8ebb0732012-04-30 23:09:08 -0700321 // Maximum number of iterations for the minimizer to run for.
322 int max_num_iterations;
323
324 // Maximum time for which the minimizer should run for.
Sameer Agarwalfa015192012-06-11 14:21:42 -0700325 double max_solver_time_in_seconds;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700326
327 // Number of threads used by Ceres for evaluating the cost and
328 // jacobians.
329 int num_threads;
330
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700331 // Trust region minimizer settings.
332 double initial_trust_region_radius;
333 double max_trust_region_radius;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700334
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700335 // Minimizer terminates when the trust region radius becomes
336 // smaller than this value.
337 double min_trust_region_radius;
Sameer Agarwal835911e2012-05-14 12:41:10 -0700338
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700339 // Lower bound for the relative decrease before a step is
340 // accepted.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700341 double min_relative_decrease;
342
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700343 // For the Levenberg-Marquadt algorithm, the scaled diagonal of
344 // the normal equations J'J is used to control the size of the
345 // trust region. Extremely small and large values along the
346 // diagonal can make this regularization scheme
Sameer Agarwaleeedd2e2013-07-07 23:04:31 -0700347 // fail. max_lm_diagonal and min_lm_diagonal, clamp the values of
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700348 // diag(J'J) from above and below. In the normal course of
349 // operation, the user should not have to modify these parameters.
Sameer Agarwaleeedd2e2013-07-07 23:04:31 -0700350 double min_lm_diagonal;
351 double max_lm_diagonal;
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700352
353 // Sometimes due to numerical conditioning problems or linear
354 // solver flakiness, the trust region strategy may return a
355 // numerically invalid step that can be fixed by reducing the
356 // trust region size. So the TrustRegionMinimizer allows for a few
357 // successive invalid steps before it declares NUMERICAL_FAILURE.
358 int max_num_consecutive_invalid_steps;
359
Keir Mierle8ebb0732012-04-30 23:09:08 -0700360 // Minimizer terminates when
361 //
362 // (new_cost - old_cost) < function_tolerance * old_cost;
363 //
364 double function_tolerance;
365
366 // Minimizer terminates when
367 //
368 // max_i |gradient_i| < gradient_tolerance * max_i|initial_gradient_i|
369 //
370 // This value should typically be 1e-4 * function_tolerance.
371 double gradient_tolerance;
372
373 // Minimizer terminates when
374 //
375 // |step|_2 <= parameter_tolerance * ( |x|_2 + parameter_tolerance)
376 //
377 double parameter_tolerance;
378
379 // Linear least squares solver options -------------------------------------
380
381 LinearSolverType linear_solver_type;
382
383 // Type of preconditioner to use with the iterative linear solvers.
384 PreconditionerType preconditioner_type;
385
Sameer Agarwal9ba0b352013-11-05 13:04:56 -0800386 // Type of clustering algorithm to use for visibility based
Sameer Agarwalf06b9fa2013-10-27 21:38:13 -0700387 // preconditioning. This option is used only when the
388 // preconditioner_type is CLUSTER_JACOBI or CLUSTER_TRIDIAGONAL.
389 VisibilityClusteringType visibility_clustering_type;
390
Sameer Agarwal367b65e2013-08-09 10:35:37 -0700391 // Ceres supports using multiple dense linear algebra libraries
392 // for dense matrix factorizations. Currently EIGEN and LAPACK are
393 // the valid choices. EIGEN is always available, LAPACK refers to
394 // the system BLAS + LAPACK library which may or may not be
395 // available.
396 //
397 // This setting affects the DENSE_QR, DENSE_NORMAL_CHOLESKY and
398 // DENSE_SCHUR solvers. For small to moderate sized probem EIGEN
399 // is a fine choice but for large problems, an optimized LAPACK +
400 // BLAS implementation can make a substantial difference in
401 // performance.
402 DenseLinearAlgebraLibraryType dense_linear_algebra_library_type;
Sameer Agarwalb0518732012-05-29 00:27:57 -0700403
Sameer Agarwald61b68a2013-08-16 17:02:56 -0700404 // Ceres supports using multiple sparse linear algebra libraries
405 // for sparse matrix ordering and factorizations. Currently,
406 // SUITE_SPARSE and CX_SPARSE are the valid choices, depending on
407 // whether they are linked into Ceres at build time.
408 SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type;
409
Keir Mierle8ebb0732012-04-30 23:09:08 -0700410 // Number of threads used by Ceres to solve the Newton
411 // step. Currently only the SPARSE_SCHUR solver is capable of
412 // using this setting.
413 int num_linear_solver_threads;
414
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700415 // The order in which variables are eliminated in a linear solver
416 // can have a significant of impact on the efficiency and accuracy
417 // of the method. e.g., when doing sparse Cholesky factorization,
418 // there are matrices for which a good ordering will give a
419 // Cholesky factor with O(n) storage, where as a bad ordering will
420 // result in an completely dense factor.
421 //
422 // Ceres allows the user to provide varying amounts of hints to
423 // the solver about the variable elimination ordering to use. This
424 // can range from no hints, where the solver is free to decide the
425 // best possible ordering based on the user's choices like the
426 // linear solver being used, to an exact order in which the
427 // variables should be eliminated, and a variety of possibilities
428 // in between.
429 //
Sameer Agarwal2c94eed2012-10-01 16:34:37 -0700430 // Instances of the ParameterBlockOrdering class are used to
431 // communicate this information to Ceres.
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700432 //
Sameer Agarwal2c94eed2012-10-01 16:34:37 -0700433 // Formally an ordering is an ordered partitioning of the
434 // parameter blocks, i.e, each parameter block belongs to exactly
435 // one group, and each group has a unique non-negative integer
436 // associated with it, that determines its order in the set of
437 // groups.
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700438 //
439 // Given such an ordering, Ceres ensures that the parameter blocks in
440 // the lowest numbered group are eliminated first, and then the
441 // parmeter blocks in the next lowest numbered group and so on. Within
442 // each group, Ceres is free to order the parameter blocks as it
443 // chooses.
444 //
Sameer Agarwal65625f72012-09-17 12:06:57 -0700445 // If NULL, then all parameter blocks are assumed to be in the
446 // same group and the solver is free to decide the best
Sameer Agarwal2c94eed2012-10-01 16:34:37 -0700447 // ordering.
448 //
449 // e.g. Consider the linear system
450 //
451 // x + y = 3
452 // 2x + 3y = 7
453 //
454 // There are two ways in which it can be solved. First eliminating x
455 // from the two equations, solving for y and then back substituting
456 // for x, or first eliminating y, solving for x and back substituting
457 // for y. The user can construct three orderings here.
458 //
459 // {0: x}, {1: y} - eliminate x first.
460 // {0: y}, {1: x} - eliminate y first.
461 // {0: x, y} - Solver gets to decide the elimination order.
462 //
463 // Thus, to have Ceres determine the ordering automatically using
464 // heuristics, put all the variables in group 0 and to control the
465 // ordering for every variable, create groups 0..N-1, one per
466 // variable, in the desired order.
467 //
468 // Bundle Adjustment
469 // -----------------
470 //
471 // A particular case of interest is bundle adjustment, where the user
472 // has two options. The default is to not specify an ordering at all,
473 // the solver will see that the user wants to use a Schur type solver
474 // and figure out the right elimination ordering.
475 //
476 // But if the user already knows what parameter blocks are points and
477 // what are cameras, they can save preprocessing time by partitioning
478 // the parameter blocks into two groups, one for the points and one
479 // for the cameras, where the group containing the points has an id
480 // smaller than the group containing cameras.
Sameer Agarwalbb05be32014-04-13 14:22:19 -0700481 shared_ptr<ParameterBlockOrdering> linear_solver_ordering;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700482
Sameer Agarwal9189f4e2013-04-19 17:09:49 -0700483 // Sparse Cholesky factorization algorithms use a fill-reducing
484 // ordering to permute the columns of the Jacobian matrix. There
485 // are two ways of doing this.
486
487 // 1. Compute the Jacobian matrix in some order and then have the
488 // factorization algorithm permute the columns of the Jacobian.
489
490 // 2. Compute the Jacobian with its columns already permuted.
491
492 // The first option incurs a significant memory penalty. The
493 // factorization algorithm has to make a copy of the permuted
494 // Jacobian matrix, thus Ceres pre-permutes the columns of the
495 // Jacobian matrix and generally speaking, there is no performance
496 // penalty for doing so.
497
498 // In some rare cases, it is worth using a more complicated
499 // reordering algorithm which has slightly better runtime
500 // performance at the expense of an extra copy of the Jacobian
Sameer Agarwalcbdeb792013-04-22 10:18:18 -0700501 // matrix. Setting use_postordering to true enables this tradeoff.
Sameer Agarwal9189f4e2013-04-19 17:09:49 -0700502 bool use_postordering;
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700503
Richard Stebbing32530782014-04-26 07:42:23 +0100504 // Some non-linear least squares problems are symbolically dense but
505 // numerically sparse. i.e. at any given state only a small number
506 // of jacobian entries are non-zero, but the position and number of
507 // non-zeros is different depending on the state. For these problems
508 // it can be useful to factorize the sparse jacobian at each solver
509 // iteration instead of including all of the zero entries in a single
510 // general factorization.
511 //
512 // If your problem does not have this property (or you do not know),
513 // then it is probably best to keep this false, otherwise it will
514 // likely lead to worse performance.
515
516 // This settings affects the SPARSE_NORMAL_CHOLESKY solver.
517 bool dynamic_sparsity;
518
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700519 // Some non-linear least squares problems have additional
520 // structure in the way the parameter blocks interact that it is
521 // beneficial to modify the way the trust region step is computed.
522 //
523 // e.g., consider the following regression problem
524 //
525 // y = a_1 exp(b_1 x) + a_2 exp(b_3 x^2 + c_1)
526 //
527 // Given a set of pairs{(x_i, y_i)}, the user wishes to estimate
528 // a_1, a_2, b_1, b_2, and c_1.
529 //
530 // Notice here that the expression on the left is linear in a_1
531 // and a_2, and given any value for b_1, b_2 and c_1, it is
532 // possible to use linear regression to estimate the optimal
533 // values of a_1 and a_2. Indeed, its possible to analytically
534 // eliminate the variables a_1 and a_2 from the problem all
535 // together. Problems like these are known as separable least
536 // squares problem and the most famous algorithm for solving them
537 // is the Variable Projection algorithm invented by Golub &
538 // Pereyra.
539 //
540 // Similar structure can be found in the matrix factorization with
541 // missing data problem. There the corresponding algorithm is
542 // known as Wiberg's algorithm.
543 //
544 // Ruhe & Wedin (Algorithms for Separable Nonlinear Least Squares
545 // Problems, SIAM Reviews, 22(3), 1980) present an analyis of
546 // various algorithms for solving separable non-linear least
547 // squares problems and refer to "Variable Projection" as
548 // Algorithm I in their paper.
549 //
550 // Implementing Variable Projection is tedious and expensive, and
551 // they present a simpler algorithm, which they refer to as
552 // Algorithm II, where once the Newton/Trust Region step has been
553 // computed for the whole problem (a_1, a_2, b_1, b_2, c_1) and
554 // additional optimization step is performed to estimate a_1 and
555 // a_2 exactly.
556 //
557 // This idea can be generalized to cases where the residual is not
558 // linear in a_1 and a_2, i.e., Solve for the trust region step
559 // for the full problem, and then use it as the starting point to
560 // further optimize just a_1 and a_2. For the linear case, this
561 // amounts to doing a single linear least squares solve. For
562 // non-linear problems, any method for solving the a_1 and a_2
563 // optimization problems will do. The only constraint on a_1 and
564 // a_2 is that they do not co-occur in any residual block.
565 //
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700566 // This idea can be further generalized, by not just optimizing
567 // (a_1, a_2), but decomposing the graph corresponding to the
568 // Hessian matrix's sparsity structure in a collection of
569 // non-overlapping independent sets and optimizing each of them.
570 //
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700571 // Setting "use_inner_iterations" to true enables the use of this
572 // non-linear generalization of Ruhe & Wedin's Algorithm II. This
573 // version of Ceres has a higher iteration complexity, but also
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700574 // displays better convergence behaviour per iteration. Setting
575 // Solver::Options::num_threads to the maximum number possible is
576 // highly recommended.
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700577 bool use_inner_iterations;
578
579 // If inner_iterations is true, then the user has two choices.
580 //
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700581 // 1. Let the solver heuristically decide which parameter blocks
582 // to optimize in each inner iteration. To do this leave
583 // Solver::Options::inner_iteration_ordering untouched.
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700584 //
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700585 // 2. Specify a collection of of ordered independent sets. Where
586 // the lower numbered groups are optimized before the higher
Sameer Agarwal09396322013-05-28 22:29:36 -0700587 // number groups. Each group must be an independent set. Not
588 // all parameter blocks need to be present in the ordering.
Sameer Agarwalbb05be32014-04-13 14:22:19 -0700589 shared_ptr<ParameterBlockOrdering> inner_iteration_ordering;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700590
Sameer Agarwald4cb94b2013-05-22 09:13:27 -0700591 // Generally speaking, inner iterations make significant progress
592 // in the early stages of the solve and then their contribution
593 // drops down sharply, at which point the time spent doing inner
594 // iterations is not worth it.
595 //
Sameer Agarwal09396322013-05-28 22:29:36 -0700596 // Once the relative decrease in the objective function due to
597 // inner iterations drops below inner_iteration_tolerance, the use
598 // of inner iterations in subsequent trust region minimizer
599 // iterations is disabled.
Sameer Agarwald4cb94b2013-05-22 09:13:27 -0700600 double inner_iteration_tolerance;
601
Keir Mierle8ebb0732012-04-30 23:09:08 -0700602 // Minimum number of iterations for which the linear solver should
603 // run, even if the convergence criterion is satisfied.
Sameer Agarwaleeedd2e2013-07-07 23:04:31 -0700604 int min_linear_solver_iterations;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700605
606 // Maximum number of iterations for which the linear solver should
607 // run. If the solver does not converge in less than
Sameer Agarwaleeedd2e2013-07-07 23:04:31 -0700608 // max_linear_solver_iterations, then it returns MAX_ITERATIONS,
609 // as its termination type.
610 int max_linear_solver_iterations;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700611
612 // Forcing sequence parameter. The truncated Newton solver uses
613 // this number to control the relative accuracy with which the
614 // Newton step is computed.
615 //
616 // This constant is passed to ConjugateGradientsSolver which uses
617 // it to terminate the iterations when
618 //
619 // (Q_i - Q_{i-1})/Q_i < eta/i
620 double eta;
621
622 // Normalize the jacobian using Jacobi scaling before calling
623 // the linear least squares solver.
624 bool jacobi_scaling;
625
626 // Logging options ---------------------------------------------------------
627
628 LoggingType logging_type;
629
630 // By default the Minimizer progress is logged to VLOG(1), which
631 // is sent to STDERR depending on the vlog level. If this flag is
632 // set to true, and logging_type is not SILENT, the logging output
633 // is sent to STDOUT.
634 bool minimizer_progress_to_stdout;
635
Sameer Agarwalc4a32912013-06-13 22:00:48 -0700636 // List of iterations at which the minimizer should dump the trust
637 // region problem. Useful for testing and benchmarking. If empty
638 // (default), no problems are dumped.
639 vector<int> trust_region_minimizer_iterations_to_dump;
640
641 // Directory to which the problems should be written to. Should be
642 // non-empty if trust_region_minimizer_iterations_to_dump is
643 // non-empty and trust_region_problem_dump_format_type is not
644 // CONSOLE.
645 string trust_region_problem_dump_directory;
646 DumpFormatType trust_region_problem_dump_format_type;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700647
Keir Mierle8ebb0732012-04-30 23:09:08 -0700648 // Finite differences options ----------------------------------------------
649
650 // Check all jacobians computed by each residual block with finite
651 // differences. This is expensive since it involves computing the
652 // derivative by normal means (e.g. user specified, autodiff,
653 // etc), then also computing it using finite differences. The
654 // results are compared, and if they differ substantially, details
655 // are printed to the log.
656 bool check_gradients;
657
658 // Relative precision to check for in the gradient checker. If the
659 // relative difference between an element in a jacobian exceeds
660 // this number, then the jacobian for that cost term is dumped.
661 double gradient_check_relative_precision;
662
663 // Relative shift used for taking numeric derivatives. For finite
664 // differencing, each dimension is evaluated at slightly shifted
665 // values; for the case of central difference, this is what gets
666 // evaluated:
667 //
668 // delta = numeric_derivative_relative_step_size;
669 // f_initial = f(x)
670 // f_forward = f((1 + delta) * x)
671 // f_backward = f((1 - delta) * x)
672 //
673 // The finite differencing is done along each dimension. The
674 // reason to use a relative (rather than absolute) step size is
675 // that this way, numeric differentation works for functions where
676 // the arguments are typically large (e.g. 1e9) and when the
677 // values are small (e.g. 1e-5). It is possible to construct
678 // "torture cases" which break this finite difference heuristic,
679 // but they do not come up often in practice.
680 //
681 // TODO(keir): Pick a smarter number than the default above! In
682 // theory a good choice is sqrt(eps) * x, which for doubles means
683 // about 1e-8 * x. However, I have found this number too
684 // optimistic. This number should be exposed for users to change.
685 double numeric_derivative_relative_step_size;
686
687 // If true, the user's parameter blocks are updated at the end of
688 // every Minimizer iteration, otherwise they are updated when the
689 // Minimizer terminates. This is useful if, for example, the user
690 // wishes to visualize the state of the optimization every
691 // iteration.
692 bool update_state_every_iteration;
693
694 // Callbacks that are executed at the end of each iteration of the
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700695 // Minimizer. An iteration may terminate midway, either due to
696 // numerical failures or because one of the convergence tests has
697 // been satisfied. In this case none of the callbacks are
698 // executed.
699
700 // Callbacks are executed in the order that they are specified in
701 // this vector. By default, parameter blocks are updated only at
702 // the end of the optimization, i.e when the Minimizer
703 // terminates. This behaviour is controlled by
Keir Mierle8ebb0732012-04-30 23:09:08 -0700704 // update_state_every_variable. If the user wishes to have access
705 // to the update parameter blocks when his/her callbacks are
706 // executed, then set update_state_every_iteration to true.
707 //
708 // The solver does NOT take ownership of these pointers.
709 vector<IterationCallback*> callbacks;
Sameer Agarwal1b7f3b52012-08-09 21:46:19 -0700710
711 // If non-empty, a summary of the execution of the solver is
712 // recorded to this file.
713 string solver_log;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700714 };
715
Björn Piltz5d7eed82014-04-23 22:13:37 +0200716 struct CERES_EXPORT Summary {
Keir Mierle8ebb0732012-04-30 23:09:08 -0700717 Summary();
718
719 // A brief one line description of the state of the solver after
720 // termination.
721 string BriefReport() const;
722
723 // A full multiline description of the state of the solver after
724 // termination.
725 string FullReport() const;
726
Sameer Agarwaldcee1202013-12-07 21:48:56 -0800727 bool IsSolutionUsable() const;
728
Keir Mierle8ebb0732012-04-30 23:09:08 -0700729 // Minimizer summary -------------------------------------------------
Sameer Agarwald010de52013-02-15 14:26:56 -0800730 MinimizerType minimizer_type;
731
Sameer Agarwaldcee1202013-12-07 21:48:56 -0800732 TerminationType termination_type;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700733
Sameer Agarwaldcee1202013-12-07 21:48:56 -0800734 // Reason why the solver terminated.
735 string message;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700736
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700737 // Cost of the problem (value of the objective function) before
738 // the optimization.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700739 double initial_cost;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700740
741 // Cost of the problem (value of the objective function) after the
742 // optimization.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700743 double final_cost;
744
745 // The part of the total cost that comes from residual blocks that
746 // were held fixed by the preprocessor because all the parameter
747 // blocks that they depend on were fixed.
748 double fixed_cost;
749
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700750 // IterationSummary for each minimizer iteration in order.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700751 vector<IterationSummary> iterations;
752
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700753 // Number of minimizer iterations in which the step was
754 // accepted. Unless use_non_monotonic_steps is true this is also
755 // the number of steps in which the objective function value/cost
756 // went down.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700757 int num_successful_steps;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700758
759 // Number of minimizer iterations in which the step was rejected
760 // either because it did not reduce the cost enough or the step
761 // was not numerically valid.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700762 int num_unsuccessful_steps;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700763
764 // Number of times inner iterations were performed.
Sameer Agarwal9f9488b2013-05-23 09:40:40 -0700765 int num_inner_iteration_steps;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700766
Sameer Agarwalf0b071b2013-05-31 13:22:51 -0700767 // All times reported below are wall times.
768
Sameer Agarwalfa015192012-06-11 14:21:42 -0700769 // When the user calls Solve, before the actual optimization
770 // occurs, Ceres performs a number of preprocessing steps. These
771 // include error checks, memory allocations, and reorderings. This
772 // time is accounted for as preprocessing time.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700773 double preprocessor_time_in_seconds;
Sameer Agarwalfa015192012-06-11 14:21:42 -0700774
775 // Time spent in the TrustRegionMinimizer.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700776 double minimizer_time_in_seconds;
Sameer Agarwalfa015192012-06-11 14:21:42 -0700777
778 // After the Minimizer is finished, some time is spent in
779 // re-evaluating residuals etc. This time is accounted for in the
780 // postprocessor time.
781 double postprocessor_time_in_seconds;
782
783 // Some total of all time spent inside Ceres when Solve is called.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700784 double total_time_in_seconds;
785
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700786 // Time (in seconds) spent in the linear solver computing the
787 // trust region step.
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800788 double linear_solver_time_in_seconds;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700789
790 // Time (in seconds) spent evaluating the residual vector.
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800791 double residual_evaluation_time_in_seconds;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700792
793 // Time (in seconds) spent evaluating the jacobian matrix.
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800794 double jacobian_evaluation_time_in_seconds;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700795
796 // Time (in seconds) spent doing inner iterations.
Sameer Agarwal9f9488b2013-05-23 09:40:40 -0700797 double inner_iteration_time_in_seconds;
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800798
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700799 // Number of parameter blocks in the problem.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700800 int num_parameter_blocks;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700801
802 // Number of parameters in the probem.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700803 int num_parameters;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700804
805 // Dimension of the tangent space of the problem (or the number of
806 // columns in the Jacobian for the problem). This is different
807 // from num_parameters if a parameter block is associated with a
808 // LocalParameterization
Keir Mierleba944212013-02-25 12:46:44 -0800809 int num_effective_parameters;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700810
811 // Number of residual blocks in the problem.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700812 int num_residual_blocks;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700813
814 // Number of residuals in the problem.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700815 int num_residuals;
816
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700817 // Number of parameter blocks in the problem after the inactive
818 // and constant parameter blocks have been removed. A parameter
819 // block is inactive if no residual block refers to it.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700820 int num_parameter_blocks_reduced;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700821
822 // Number of parameters in the reduced problem.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700823 int num_parameters_reduced;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700824
825 // Dimension of the tangent space of the reduced problem (or the
826 // number of columns in the Jacobian for the reduced
827 // problem). This is different from num_parameters_reduced if a
828 // parameter block in the reduced problem is associated with a
829 // LocalParameterization.
Keir Mierleba944212013-02-25 12:46:44 -0800830 int num_effective_parameters_reduced;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700831
832 // Number of residual blocks in the reduced problem.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700833 int num_residual_blocks_reduced;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700834
835 // Number of residuals in the reduced problem.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700836 int num_residuals_reduced;
837
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700838 // Number of threads specified by the user for Jacobian and
839 // residual evaluation.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700840 int num_threads_given;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700841
842 // Number of threads actually used by the solver for Jacobian and
843 // residual evaluation. This number is not equal to
844 // num_threads_given if OpenMP is not available.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700845 int num_threads_used;
846
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700847 // Number of threads specified by the user for solving the trust
848 // region problem.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700849 int num_linear_solver_threads_given;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700850
851 // Number of threads actually used by the solver for solving the
852 // trust region problem. This number is not equal to
853 // num_threads_given if OpenMP is not available.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700854 int num_linear_solver_threads_used;
855
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700856 // Type of the linear solver requested by the user.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700857 LinearSolverType linear_solver_type_given;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700858
859 // Type of the linear solver actually used. This may be different
860 // from linear_solver_type_given if Ceres determines that the
861 // problem structure is not compatible with the linear solver
862 // requested or if the linear solver requested by the user is not
863 // available, e.g. The user requested SPARSE_NORMAL_CHOLESKY but
864 // no sparse linear algebra library was available.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700865 LinearSolverType linear_solver_type_used;
866
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700867 // Size of the elimination groups given by the user as hints to
868 // the linear solver.
Sameer Agarwal977be7c2013-01-26 16:01:54 -0800869 vector<int> linear_solver_ordering_given;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700870
871 // Size of the parameter groups used by the solver when ordering
872 // the columns of the Jacobian. This maybe different from
873 // linear_solver_ordering_given if the user left
874 // linear_solver_ordering_given blank and asked for an automatic
875 // ordering, or if the problem contains some constant or inactive
876 // parameter blocks.
Sameer Agarwal977be7c2013-01-26 16:01:54 -0800877 vector<int> linear_solver_ordering_used;
878
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700879 // True if the user asked for inner iterations to be used as part
880 // of the optimization.
Sameer Agarwal9f9488b2013-05-23 09:40:40 -0700881 bool inner_iterations_given;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700882
883 // True if the user asked for inner iterations to be used as part
884 // of the optimization and the problem structure was such that
885 // they were actually performed. e.g., in a problem with just one
886 // parameter block, inner iterations are not performed.
Sameer Agarwal9f9488b2013-05-23 09:40:40 -0700887 bool inner_iterations_used;
888
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700889 // Size of the parameter groups given by the user for performing
890 // inner iterations.
Sameer Agarwal9f9488b2013-05-23 09:40:40 -0700891 vector<int> inner_iteration_ordering_given;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700892
893 // Size of the parameter groups given used by the solver for
894 // performing inner iterations. This maybe different from
895 // inner_iteration_ordering_given if the user left
896 // inner_iteration_ordering_given blank and asked for an automatic
897 // ordering, or if the problem contains some constant or inactive
898 // parameter blocks.
Sameer Agarwal9f9488b2013-05-23 09:40:40 -0700899 vector<int> inner_iteration_ordering_used;
900
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700901 // Type of preconditioner used for solving the trust region
902 // step. Only meaningful when an iterative linear solver is used.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700903 PreconditionerType preconditioner_type;
Sameer Agarwal97fb6d92012-06-17 10:08:19 -0700904
Sameer Agarwalf06b9fa2013-10-27 21:38:13 -0700905 // Type of clustering algorithm used for visibility based
906 // preconditioning. Only meaningful when the preconditioner_type
907 // is CLUSTER_JACOBI or CLUSTER_TRIDIAGONAL.
908 VisibilityClusteringType visibility_clustering_type;
909
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700910 // Type of trust region strategy.
Sameer Agarwal97fb6d92012-06-17 10:08:19 -0700911 TrustRegionStrategyType trust_region_strategy_type;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700912
913 // Type of dogleg strategy used for solving the trust region
914 // problem.
Sameer Agarwal1e289202012-08-21 18:00:54 -0700915 DoglegType dogleg_type;
Sameer Agarwald010de52013-02-15 14:26:56 -0800916
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700917 // Type of the dense linear algebra library used.
Sameer Agarwal367b65e2013-08-09 10:35:37 -0700918 DenseLinearAlgebraLibraryType dense_linear_algebra_library_type;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700919
920 // Type of the sparse linear algebra library used.
Sameer Agarwal367b65e2013-08-09 10:35:37 -0700921 SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type;
Sameer Agarwal977be7c2013-01-26 16:01:54 -0800922
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700923 // Type of line search direction used.
Sameer Agarwald010de52013-02-15 14:26:56 -0800924 LineSearchDirectionType line_search_direction_type;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700925
926 // Type of the line search algorithm used.
Sameer Agarwald010de52013-02-15 14:26:56 -0800927 LineSearchType line_search_type;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700928
929 // When performing line search, the degree of the polynomial used
930 // to approximate the objective function.
Sameer Agarwal4f010b22013-07-01 08:01:01 -0700931 LineSearchInterpolationType line_search_interpolation_type;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700932
933 // If the line search direction is NONLINEAR_CONJUGATE_GRADIENT,
934 // then this indicates the particular variant of non-linear
935 // conjugate gradient used.
Sameer Agarwal67ccb732013-07-03 06:28:34 -0700936 NonlinearConjugateGradientType nonlinear_conjugate_gradient_type;
Sameer Agarwal09244012013-06-30 14:33:23 -0700937
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700938 // If the type of the line search direction is LBFGS, then this
939 // indicates the rank of the Hessian approximation.
Sameer Agarwald010de52013-02-15 14:26:56 -0800940 int max_lbfgs_rank;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700941 };
942
943 // Once a least squares problem has been built, this function takes
944 // the problem and optimizes it based on the values of the options
945 // parameters. Upon return, a detailed summary of the work performed
946 // by the preprocessor, the non-linear minmizer and the linear
947 // solver are reported in the summary object.
948 virtual void Solve(const Options& options,
949 Problem* problem,
950 Solver::Summary* summary);
951};
952
953// Helper function which avoids going through the interface.
Björn Piltz5d7eed82014-04-23 22:13:37 +0200954CERES_EXPORT void Solve(const Solver::Options& options,
Keir Mierle8ebb0732012-04-30 23:09:08 -0700955 Problem* problem,
956 Solver::Summary* summary);
957
958} // namespace ceres
959
960#endif // CERES_PUBLIC_SOLVER_H_