blob: 7776c470ebab0fb8ff4d4cac9e716c5c651cc48f [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.
49class Solver {
50 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
58 struct Options {
59 // 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 Agarwal68b32a92012-10-06 23:10:51 -0700110 linear_solver_ordering = NULL;
Sameer Agarwal9189f4e2013-04-19 17:09:49 -0700111 use_postordering = 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;
118 inner_iteration_ordering = NULL;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700119 logging_type = PER_MINIMIZER_ITERATION;
120 minimizer_progress_to_stdout = false;
Sameer Agarwalc4a32912013-06-13 22:00:48 -0700121 trust_region_problem_dump_directory = "/tmp";
122 trust_region_problem_dump_format_type = TEXTFILE;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700123 check_gradients = false;
124 gradient_check_relative_precision = 1e-8;
125 numeric_derivative_relative_step_size = 1e-6;
126 update_state_every_iteration = false;
127 }
128
Sameer Agarwal65625f72012-09-17 12:06:57 -0700129 ~Options();
Keir Mierle8ebb0732012-04-30 23:09:08 -0700130 // Minimizer options ----------------------------------------
131
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800132 // Ceres supports the two major families of optimization strategies -
133 // Trust Region and Line Search.
134 //
135 // 1. The line search approach first finds a descent direction
136 // along which the objective function will be reduced and then
137 // computes a step size that decides how far should move along
138 // that direction. The descent direction can be computed by
139 // various methods, such as gradient descent, Newton's method and
140 // Quasi-Newton method. The step size can be determined either
141 // exactly or inexactly.
142 //
143 // 2. The trust region approach approximates the objective
144 // function using using a model function (often a quadratic) over
145 // a subset of the search space known as the trust region. If the
146 // model function succeeds in minimizing the true objective
147 // function the trust region is expanded; conversely, otherwise it
148 // is contracted and the model optimization problem is solved
149 // again.
150 //
151 // Trust region methods are in some sense dual to line search methods:
152 // trust region methods first choose a step size (the size of the
153 // trust region) and then a step direction while line search methods
154 // first choose a step direction and then a step size.
155 MinimizerType minimizer_type;
156
157 LineSearchDirectionType line_search_direction_type;
158 LineSearchType line_search_type;
159 NonlinearConjugateGradientType nonlinear_conjugate_gradient_type;
160
Sameer Agarwalaed99612012-11-29 10:33:19 -0800161 // The LBFGS hessian approximation is a low rank approximation to
162 // the inverse of the Hessian matrix. The rank of the
163 // approximation determines (linearly) the space and time
164 // complexity of using the approximation. Higher the rank, the
Sameer Agarwal2293cb52012-11-29 16:00:18 -0800165 // better is the quality of the approximation. The increase in
166 // quality is however is bounded for a number of reasons.
167 //
168 // 1. The method only uses secant information and not actual
169 // derivatives.
170 //
171 // 2. The Hessian approximation is constrained to be positive
172 // definite.
173 //
174 // So increasing this rank to a large number will cost time and
175 // space complexity without the corresponding increase in solution
176 // quality. There are no hard and fast rules for choosing the
177 // maximum rank. The best choice usually requires some problem
178 // specific experimentation.
179 //
180 // For more theoretical and implementation details of the LBFGS
181 // method, please see:
Sameer Agarwalaed99612012-11-29 10:33:19 -0800182 //
183 // Nocedal, J. (1980). "Updating Quasi-Newton Matrices with
184 // Limited Storage". Mathematics of Computation 35 (151): 773–782.
185 int max_lbfgs_rank;
186
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100187 // As part of the (L)BFGS update step (BFGS) / right-multiply step (L-BFGS),
188 // the initial inverse Hessian approximation is taken to be the Identity.
189 // However, Oren showed that using instead I * \gamma, where \gamma is
190 // chosen to approximate an eigenvalue of the true inverse Hessian can
191 // result in improved convergence in a wide variety of cases. Setting
192 // use_approximate_eigenvalue_bfgs_scaling to true enables this scaling.
193 //
194 // It is important to note that approximate eigenvalue scaling does not
195 // always improve convergence, and that it can in fact significantly degrade
196 // performance for certain classes of problem, which is why it is disabled
197 // by default. In particular it can degrade performance when the
198 // sensitivity of the problem to different parameters varies significantly,
199 // as in this case a single scalar factor fails to capture this variation
200 // and detrimentally downscales parts of the jacobian approximation which
201 // correspond to low-sensitivity parameters. It can also reduce the
202 // robustness of the solution to errors in the jacobians.
203 //
204 // Oren S.S., Self-scaling variable metric (SSVM) algorithms
205 // Part II: Implementation and experiments, Management Science,
206 // 20(5), 863-874, 1974.
207 bool use_approximate_eigenvalue_bfgs_scaling;
208
Sameer Agarwal09244012013-06-30 14:33:23 -0700209 // Degree of the polynomial used to approximate the objective
210 // function. Valid values are BISECTION, QUADRATIC and CUBIC.
211 //
212 // BISECTION corresponds to pure backtracking search with no
213 // interpolation.
214 LineSearchInterpolationType line_search_interpolation_type;
215
216 // If during the line search, the step_size falls below this
217 // value, it is truncated to zero.
218 double min_line_search_step_size;
219
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100220 // Line search parameters.
Sameer Agarwal09244012013-06-30 14:33:23 -0700221
222 // Solving the line search problem exactly is computationally
223 // prohibitive. Fortunately, line search based optimization
224 // algorithms can still guarantee convergence if instead of an
225 // exact solution, the line search algorithm returns a solution
226 // which decreases the value of the objective function
227 // sufficiently. More precisely, we are looking for a step_size
228 // s.t.
229 //
230 // f(step_size) <= f(0) + sufficient_decrease * f'(0) * step_size
231 //
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100232 double line_search_sufficient_function_decrease;
Sameer Agarwal09244012013-06-30 14:33:23 -0700233
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100234 // In each iteration of the line search,
Sameer Agarwal09244012013-06-30 14:33:23 -0700235 //
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100236 // new_step_size >= max_line_search_step_contraction * step_size
Sameer Agarwal09244012013-06-30 14:33:23 -0700237 //
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100238 // Note that by definition, for contraction:
239 //
240 // 0 < max_step_contraction < min_step_contraction < 1
241 //
242 double max_line_search_step_contraction;
Sameer Agarwal09244012013-06-30 14:33:23 -0700243
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100244 // In each iteration of the line search,
Sameer Agarwal09244012013-06-30 14:33:23 -0700245 //
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100246 // new_step_size <= min_line_search_step_contraction * step_size
Sameer Agarwal09244012013-06-30 14:33:23 -0700247 //
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100248 // Note that by definition, for contraction:
249 //
250 // 0 < max_step_contraction < min_step_contraction < 1
251 //
252 double min_line_search_step_contraction;
253
254 // Maximum number of trial step size iterations during each line search,
255 // if a step size satisfying the search conditions cannot be found within
256 // this number of trials, the line search will terminate.
257 int max_num_line_search_step_size_iterations;
258
259 // Maximum number of restarts of the line search direction algorithm before
260 // terminating the optimization. Restarts of the line search direction
261 // algorithm occur when the current algorithm fails to produce a new descent
262 // direction. This typically indicates a numerical failure, or a breakdown
263 // in the validity of the approximations used.
264 int max_num_line_search_direction_restarts;
265
266 // The strong Wolfe conditions consist of the Armijo sufficient
267 // decrease condition, and an additional requirement that the
268 // step-size be chosen s.t. the _magnitude_ ('strong' Wolfe
269 // conditions) of the gradient along the search direction
270 // decreases sufficiently. Precisely, this second condition
271 // is that we seek a step_size s.t.
272 //
273 // |f'(step_size)| <= sufficient_curvature_decrease * |f'(0)|
274 //
275 // Where f() is the line search objective and f'() is the derivative
276 // of f w.r.t step_size (d f / d step_size).
277 double line_search_sufficient_curvature_decrease;
278
279 // During the bracketing phase of the Wolfe search, the step size is
280 // increased until either a point satisfying the Wolfe conditions is
281 // found, or an upper bound for a bracket containing a point satisfying
282 // the conditions is found. Precisely, at each iteration of the
283 // expansion:
284 //
285 // new_step_size <= max_step_expansion * step_size.
286 //
287 // By definition for expansion, max_step_expansion > 1.0.
288 double max_line_search_step_expansion;
Sameer Agarwal09244012013-06-30 14:33:23 -0700289
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700290 TrustRegionStrategyType trust_region_strategy_type;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700291
Markus Moll51cf7cb2012-08-20 20:10:20 +0200292 // Type of dogleg strategy to use.
293 DoglegType dogleg_type;
294
Sameer Agarwala8f87d72012-08-08 10:38:31 -0700295 // The classical trust region methods are descent methods, in that
296 // they only accept a point if it strictly reduces the value of
297 // the objective function.
298 //
299 // Relaxing this requirement allows the algorithm to be more
300 // efficient in the long term at the cost of some local increase
301 // in the value of the objective function.
302 //
303 // This is because allowing for non-decreasing objective function
304 // values in a princpled manner allows the algorithm to "jump over
305 // boulders" as the method is not restricted to move into narrow
306 // valleys while preserving its convergence properties.
307 //
308 // Setting use_nonmonotonic_steps to true enables the
309 // non-monotonic trust region algorithm as described by Conn,
310 // Gould & Toint in "Trust Region Methods", Section 10.1.
311 //
312 // The parameter max_consecutive_nonmonotonic_steps controls the
313 // window size used by the step selection algorithm to accept
314 // non-monotonic steps.
315 //
316 // Even though the value of the objective function may be larger
317 // than the minimum value encountered over the course of the
318 // optimization, the final parameters returned to the user are the
319 // ones corresponding to the minimum cost over all iterations.
320 bool use_nonmonotonic_steps;
321 int max_consecutive_nonmonotonic_steps;
322
Keir Mierle8ebb0732012-04-30 23:09:08 -0700323 // Maximum number of iterations for the minimizer to run for.
324 int max_num_iterations;
325
326 // Maximum time for which the minimizer should run for.
Sameer Agarwalfa015192012-06-11 14:21:42 -0700327 double max_solver_time_in_seconds;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700328
329 // Number of threads used by Ceres for evaluating the cost and
330 // jacobians.
331 int num_threads;
332
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700333 // Trust region minimizer settings.
334 double initial_trust_region_radius;
335 double max_trust_region_radius;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700336
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700337 // Minimizer terminates when the trust region radius becomes
338 // smaller than this value.
339 double min_trust_region_radius;
Sameer Agarwal835911e2012-05-14 12:41:10 -0700340
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700341 // Lower bound for the relative decrease before a step is
342 // accepted.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700343 double min_relative_decrease;
344
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700345 // For the Levenberg-Marquadt algorithm, the scaled diagonal of
346 // the normal equations J'J is used to control the size of the
347 // trust region. Extremely small and large values along the
348 // diagonal can make this regularization scheme
Sameer Agarwaleeedd2e2013-07-07 23:04:31 -0700349 // fail. max_lm_diagonal and min_lm_diagonal, clamp the values of
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700350 // diag(J'J) from above and below. In the normal course of
351 // operation, the user should not have to modify these parameters.
Sameer Agarwaleeedd2e2013-07-07 23:04:31 -0700352 double min_lm_diagonal;
353 double max_lm_diagonal;
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700354
355 // Sometimes due to numerical conditioning problems or linear
356 // solver flakiness, the trust region strategy may return a
357 // numerically invalid step that can be fixed by reducing the
358 // trust region size. So the TrustRegionMinimizer allows for a few
359 // successive invalid steps before it declares NUMERICAL_FAILURE.
360 int max_num_consecutive_invalid_steps;
361
Keir Mierle8ebb0732012-04-30 23:09:08 -0700362 // Minimizer terminates when
363 //
364 // (new_cost - old_cost) < function_tolerance * old_cost;
365 //
366 double function_tolerance;
367
368 // Minimizer terminates when
369 //
370 // max_i |gradient_i| < gradient_tolerance * max_i|initial_gradient_i|
371 //
372 // This value should typically be 1e-4 * function_tolerance.
373 double gradient_tolerance;
374
375 // Minimizer terminates when
376 //
377 // |step|_2 <= parameter_tolerance * ( |x|_2 + parameter_tolerance)
378 //
379 double parameter_tolerance;
380
381 // Linear least squares solver options -------------------------------------
382
383 LinearSolverType linear_solver_type;
384
385 // Type of preconditioner to use with the iterative linear solvers.
386 PreconditionerType preconditioner_type;
387
Sameer Agarwal9ba0b352013-11-05 13:04:56 -0800388 // Type of clustering algorithm to use for visibility based
Sameer Agarwalf06b9fa2013-10-27 21:38:13 -0700389 // preconditioning. This option is used only when the
390 // preconditioner_type is CLUSTER_JACOBI or CLUSTER_TRIDIAGONAL.
391 VisibilityClusteringType visibility_clustering_type;
392
Sameer Agarwal367b65e2013-08-09 10:35:37 -0700393 // Ceres supports using multiple dense linear algebra libraries
394 // for dense matrix factorizations. Currently EIGEN and LAPACK are
395 // the valid choices. EIGEN is always available, LAPACK refers to
396 // the system BLAS + LAPACK library which may or may not be
397 // available.
398 //
399 // This setting affects the DENSE_QR, DENSE_NORMAL_CHOLESKY and
400 // DENSE_SCHUR solvers. For small to moderate sized probem EIGEN
401 // is a fine choice but for large problems, an optimized LAPACK +
402 // BLAS implementation can make a substantial difference in
403 // performance.
404 DenseLinearAlgebraLibraryType dense_linear_algebra_library_type;
Sameer Agarwalb0518732012-05-29 00:27:57 -0700405
Sameer Agarwald61b68a2013-08-16 17:02:56 -0700406 // Ceres supports using multiple sparse linear algebra libraries
407 // for sparse matrix ordering and factorizations. Currently,
408 // SUITE_SPARSE and CX_SPARSE are the valid choices, depending on
409 // whether they are linked into Ceres at build time.
410 SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type;
411
Keir Mierle8ebb0732012-04-30 23:09:08 -0700412 // Number of threads used by Ceres to solve the Newton
413 // step. Currently only the SPARSE_SCHUR solver is capable of
414 // using this setting.
415 int num_linear_solver_threads;
416
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700417 // The order in which variables are eliminated in a linear solver
418 // can have a significant of impact on the efficiency and accuracy
419 // of the method. e.g., when doing sparse Cholesky factorization,
420 // there are matrices for which a good ordering will give a
421 // Cholesky factor with O(n) storage, where as a bad ordering will
422 // result in an completely dense factor.
423 //
424 // Ceres allows the user to provide varying amounts of hints to
425 // the solver about the variable elimination ordering to use. This
426 // can range from no hints, where the solver is free to decide the
427 // best possible ordering based on the user's choices like the
428 // linear solver being used, to an exact order in which the
429 // variables should be eliminated, and a variety of possibilities
430 // in between.
431 //
Sameer Agarwal2c94eed2012-10-01 16:34:37 -0700432 // Instances of the ParameterBlockOrdering class are used to
433 // communicate this information to Ceres.
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700434 //
Sameer Agarwal2c94eed2012-10-01 16:34:37 -0700435 // Formally an ordering is an ordered partitioning of the
436 // parameter blocks, i.e, each parameter block belongs to exactly
437 // one group, and each group has a unique non-negative integer
438 // associated with it, that determines its order in the set of
439 // groups.
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700440 //
441 // Given such an ordering, Ceres ensures that the parameter blocks in
442 // the lowest numbered group are eliminated first, and then the
443 // parmeter blocks in the next lowest numbered group and so on. Within
444 // each group, Ceres is free to order the parameter blocks as it
445 // chooses.
446 //
Sameer Agarwal65625f72012-09-17 12:06:57 -0700447 // If NULL, then all parameter blocks are assumed to be in the
448 // same group and the solver is free to decide the best
Sameer Agarwal2c94eed2012-10-01 16:34:37 -0700449 // ordering.
450 //
451 // e.g. Consider the linear system
452 //
453 // x + y = 3
454 // 2x + 3y = 7
455 //
456 // There are two ways in which it can be solved. First eliminating x
457 // from the two equations, solving for y and then back substituting
458 // for x, or first eliminating y, solving for x and back substituting
459 // for y. The user can construct three orderings here.
460 //
461 // {0: x}, {1: y} - eliminate x first.
462 // {0: y}, {1: x} - eliminate y first.
463 // {0: x, y} - Solver gets to decide the elimination order.
464 //
465 // Thus, to have Ceres determine the ordering automatically using
466 // heuristics, put all the variables in group 0 and to control the
467 // ordering for every variable, create groups 0..N-1, one per
468 // variable, in the desired order.
469 //
470 // Bundle Adjustment
471 // -----------------
472 //
473 // A particular case of interest is bundle adjustment, where the user
474 // has two options. The default is to not specify an ordering at all,
475 // the solver will see that the user wants to use a Schur type solver
476 // and figure out the right elimination ordering.
477 //
478 // But if the user already knows what parameter blocks are points and
479 // what are cameras, they can save preprocessing time by partitioning
480 // the parameter blocks into two groups, one for the points and one
481 // for the cameras, where the group containing the points has an id
482 // smaller than the group containing cameras.
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700483 //
484 // Once assigned, Solver::Options owns this pointer and will
485 // deallocate the memory when destroyed.
Sameer Agarwal68b32a92012-10-06 23:10:51 -0700486 ParameterBlockOrdering* linear_solver_ordering;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700487
Sameer Agarwal9189f4e2013-04-19 17:09:49 -0700488 // Sparse Cholesky factorization algorithms use a fill-reducing
489 // ordering to permute the columns of the Jacobian matrix. There
490 // are two ways of doing this.
491
492 // 1. Compute the Jacobian matrix in some order and then have the
493 // factorization algorithm permute the columns of the Jacobian.
494
495 // 2. Compute the Jacobian with its columns already permuted.
496
497 // The first option incurs a significant memory penalty. The
498 // factorization algorithm has to make a copy of the permuted
499 // Jacobian matrix, thus Ceres pre-permutes the columns of the
500 // Jacobian matrix and generally speaking, there is no performance
501 // penalty for doing so.
502
503 // In some rare cases, it is worth using a more complicated
504 // reordering algorithm which has slightly better runtime
505 // performance at the expense of an extra copy of the Jacobian
Sameer Agarwalcbdeb792013-04-22 10:18:18 -0700506 // matrix. Setting use_postordering to true enables this tradeoff.
Sameer Agarwal9189f4e2013-04-19 17:09:49 -0700507 bool use_postordering;
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700508
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700509 // Some non-linear least squares problems have additional
510 // structure in the way the parameter blocks interact that it is
511 // beneficial to modify the way the trust region step is computed.
512 //
513 // e.g., consider the following regression problem
514 //
515 // y = a_1 exp(b_1 x) + a_2 exp(b_3 x^2 + c_1)
516 //
517 // Given a set of pairs{(x_i, y_i)}, the user wishes to estimate
518 // a_1, a_2, b_1, b_2, and c_1.
519 //
520 // Notice here that the expression on the left is linear in a_1
521 // and a_2, and given any value for b_1, b_2 and c_1, it is
522 // possible to use linear regression to estimate the optimal
523 // values of a_1 and a_2. Indeed, its possible to analytically
524 // eliminate the variables a_1 and a_2 from the problem all
525 // together. Problems like these are known as separable least
526 // squares problem and the most famous algorithm for solving them
527 // is the Variable Projection algorithm invented by Golub &
528 // Pereyra.
529 //
530 // Similar structure can be found in the matrix factorization with
531 // missing data problem. There the corresponding algorithm is
532 // known as Wiberg's algorithm.
533 //
534 // Ruhe & Wedin (Algorithms for Separable Nonlinear Least Squares
535 // Problems, SIAM Reviews, 22(3), 1980) present an analyis of
536 // various algorithms for solving separable non-linear least
537 // squares problems and refer to "Variable Projection" as
538 // Algorithm I in their paper.
539 //
540 // Implementing Variable Projection is tedious and expensive, and
541 // they present a simpler algorithm, which they refer to as
542 // Algorithm II, where once the Newton/Trust Region step has been
543 // computed for the whole problem (a_1, a_2, b_1, b_2, c_1) and
544 // additional optimization step is performed to estimate a_1 and
545 // a_2 exactly.
546 //
547 // This idea can be generalized to cases where the residual is not
548 // linear in a_1 and a_2, i.e., Solve for the trust region step
549 // for the full problem, and then use it as the starting point to
550 // further optimize just a_1 and a_2. For the linear case, this
551 // amounts to doing a single linear least squares solve. For
552 // non-linear problems, any method for solving the a_1 and a_2
553 // optimization problems will do. The only constraint on a_1 and
554 // a_2 is that they do not co-occur in any residual block.
555 //
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700556 // This idea can be further generalized, by not just optimizing
557 // (a_1, a_2), but decomposing the graph corresponding to the
558 // Hessian matrix's sparsity structure in a collection of
559 // non-overlapping independent sets and optimizing each of them.
560 //
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700561 // Setting "use_inner_iterations" to true enables the use of this
562 // non-linear generalization of Ruhe & Wedin's Algorithm II. This
563 // version of Ceres has a higher iteration complexity, but also
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700564 // displays better convergence behaviour per iteration. Setting
565 // Solver::Options::num_threads to the maximum number possible is
566 // highly recommended.
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700567 bool use_inner_iterations;
568
569 // If inner_iterations is true, then the user has two choices.
570 //
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700571 // 1. Let the solver heuristically decide which parameter blocks
572 // to optimize in each inner iteration. To do this leave
573 // Solver::Options::inner_iteration_ordering untouched.
Sameer Agarwal9123e2f2012-09-18 21:49:06 -0700574 //
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700575 // 2. Specify a collection of of ordered independent sets. Where
576 // the lower numbered groups are optimized before the higher
Sameer Agarwal09396322013-05-28 22:29:36 -0700577 // number groups. Each group must be an independent set. Not
578 // all parameter blocks need to be present in the ordering.
Sameer Agarwalba8d9672012-10-02 00:48:57 -0700579 ParameterBlockOrdering* inner_iteration_ordering;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700580
Sameer Agarwald4cb94b2013-05-22 09:13:27 -0700581 // Generally speaking, inner iterations make significant progress
582 // in the early stages of the solve and then their contribution
583 // drops down sharply, at which point the time spent doing inner
584 // iterations is not worth it.
585 //
Sameer Agarwal09396322013-05-28 22:29:36 -0700586 // Once the relative decrease in the objective function due to
587 // inner iterations drops below inner_iteration_tolerance, the use
588 // of inner iterations in subsequent trust region minimizer
589 // iterations is disabled.
Sameer Agarwald4cb94b2013-05-22 09:13:27 -0700590 double inner_iteration_tolerance;
591
Keir Mierle8ebb0732012-04-30 23:09:08 -0700592 // Minimum number of iterations for which the linear solver should
593 // run, even if the convergence criterion is satisfied.
Sameer Agarwaleeedd2e2013-07-07 23:04:31 -0700594 int min_linear_solver_iterations;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700595
596 // Maximum number of iterations for which the linear solver should
597 // run. If the solver does not converge in less than
Sameer Agarwaleeedd2e2013-07-07 23:04:31 -0700598 // max_linear_solver_iterations, then it returns MAX_ITERATIONS,
599 // as its termination type.
600 int max_linear_solver_iterations;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700601
602 // Forcing sequence parameter. The truncated Newton solver uses
603 // this number to control the relative accuracy with which the
604 // Newton step is computed.
605 //
606 // This constant is passed to ConjugateGradientsSolver which uses
607 // it to terminate the iterations when
608 //
609 // (Q_i - Q_{i-1})/Q_i < eta/i
610 double eta;
611
612 // Normalize the jacobian using Jacobi scaling before calling
613 // the linear least squares solver.
614 bool jacobi_scaling;
615
616 // Logging options ---------------------------------------------------------
617
618 LoggingType logging_type;
619
620 // By default the Minimizer progress is logged to VLOG(1), which
621 // is sent to STDERR depending on the vlog level. If this flag is
622 // set to true, and logging_type is not SILENT, the logging output
623 // is sent to STDOUT.
624 bool minimizer_progress_to_stdout;
625
Sameer Agarwalc4a32912013-06-13 22:00:48 -0700626 // List of iterations at which the minimizer should dump the trust
627 // region problem. Useful for testing and benchmarking. If empty
628 // (default), no problems are dumped.
629 vector<int> trust_region_minimizer_iterations_to_dump;
630
631 // Directory to which the problems should be written to. Should be
632 // non-empty if trust_region_minimizer_iterations_to_dump is
633 // non-empty and trust_region_problem_dump_format_type is not
634 // CONSOLE.
635 string trust_region_problem_dump_directory;
636 DumpFormatType trust_region_problem_dump_format_type;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700637
Keir Mierle8ebb0732012-04-30 23:09:08 -0700638 // Finite differences options ----------------------------------------------
639
640 // Check all jacobians computed by each residual block with finite
641 // differences. This is expensive since it involves computing the
642 // derivative by normal means (e.g. user specified, autodiff,
643 // etc), then also computing it using finite differences. The
644 // results are compared, and if they differ substantially, details
645 // are printed to the log.
646 bool check_gradients;
647
648 // Relative precision to check for in the gradient checker. If the
649 // relative difference between an element in a jacobian exceeds
650 // this number, then the jacobian for that cost term is dumped.
651 double gradient_check_relative_precision;
652
653 // Relative shift used for taking numeric derivatives. For finite
654 // differencing, each dimension is evaluated at slightly shifted
655 // values; for the case of central difference, this is what gets
656 // evaluated:
657 //
658 // delta = numeric_derivative_relative_step_size;
659 // f_initial = f(x)
660 // f_forward = f((1 + delta) * x)
661 // f_backward = f((1 - delta) * x)
662 //
663 // The finite differencing is done along each dimension. The
664 // reason to use a relative (rather than absolute) step size is
665 // that this way, numeric differentation works for functions where
666 // the arguments are typically large (e.g. 1e9) and when the
667 // values are small (e.g. 1e-5). It is possible to construct
668 // "torture cases" which break this finite difference heuristic,
669 // but they do not come up often in practice.
670 //
671 // TODO(keir): Pick a smarter number than the default above! In
672 // theory a good choice is sqrt(eps) * x, which for doubles means
673 // about 1e-8 * x. However, I have found this number too
674 // optimistic. This number should be exposed for users to change.
675 double numeric_derivative_relative_step_size;
676
677 // If true, the user's parameter blocks are updated at the end of
678 // every Minimizer iteration, otherwise they are updated when the
679 // Minimizer terminates. This is useful if, for example, the user
680 // wishes to visualize the state of the optimization every
681 // iteration.
682 bool update_state_every_iteration;
683
684 // Callbacks that are executed at the end of each iteration of the
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700685 // Minimizer. An iteration may terminate midway, either due to
686 // numerical failures or because one of the convergence tests has
687 // been satisfied. In this case none of the callbacks are
688 // executed.
689
690 // Callbacks are executed in the order that they are specified in
691 // this vector. By default, parameter blocks are updated only at
692 // the end of the optimization, i.e when the Minimizer
693 // terminates. This behaviour is controlled by
Keir Mierle8ebb0732012-04-30 23:09:08 -0700694 // update_state_every_variable. If the user wishes to have access
695 // to the update parameter blocks when his/her callbacks are
696 // executed, then set update_state_every_iteration to true.
697 //
698 // The solver does NOT take ownership of these pointers.
699 vector<IterationCallback*> callbacks;
Sameer Agarwal1b7f3b52012-08-09 21:46:19 -0700700
701 // If non-empty, a summary of the execution of the solver is
702 // recorded to this file.
703 string solver_log;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700704 };
705
706 struct Summary {
707 Summary();
708
709 // A brief one line description of the state of the solver after
710 // termination.
711 string BriefReport() const;
712
713 // A full multiline description of the state of the solver after
714 // termination.
715 string FullReport() const;
716
717 // Minimizer summary -------------------------------------------------
Sameer Agarwald010de52013-02-15 14:26:56 -0800718 MinimizerType minimizer_type;
719
Keir Mierle8ebb0732012-04-30 23:09:08 -0700720 SolverTerminationType termination_type;
721
722 // If the solver did not run, or there was a failure, a
723 // description of the error.
724 string error;
725
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700726 // Cost of the problem (value of the objective function) before
727 // the optimization.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700728 double initial_cost;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700729
730 // Cost of the problem (value of the objective function) after the
731 // optimization.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700732 double final_cost;
733
734 // The part of the total cost that comes from residual blocks that
735 // were held fixed by the preprocessor because all the parameter
736 // blocks that they depend on were fixed.
737 double fixed_cost;
738
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700739 // IterationSummary for each minimizer iteration in order.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700740 vector<IterationSummary> iterations;
741
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700742 // Number of minimizer iterations in which the step was
743 // accepted. Unless use_non_monotonic_steps is true this is also
744 // the number of steps in which the objective function value/cost
745 // went down.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700746 int num_successful_steps;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700747
748 // Number of minimizer iterations in which the step was rejected
749 // either because it did not reduce the cost enough or the step
750 // was not numerically valid.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700751 int num_unsuccessful_steps;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700752
753 // Number of times inner iterations were performed.
Sameer Agarwal9f9488b2013-05-23 09:40:40 -0700754 int num_inner_iteration_steps;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700755
Sameer Agarwalf0b071b2013-05-31 13:22:51 -0700756 // All times reported below are wall times.
757
Sameer Agarwalfa015192012-06-11 14:21:42 -0700758 // When the user calls Solve, before the actual optimization
759 // occurs, Ceres performs a number of preprocessing steps. These
760 // include error checks, memory allocations, and reorderings. This
761 // time is accounted for as preprocessing time.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700762 double preprocessor_time_in_seconds;
Sameer Agarwalfa015192012-06-11 14:21:42 -0700763
764 // Time spent in the TrustRegionMinimizer.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700765 double minimizer_time_in_seconds;
Sameer Agarwalfa015192012-06-11 14:21:42 -0700766
767 // After the Minimizer is finished, some time is spent in
768 // re-evaluating residuals etc. This time is accounted for in the
769 // postprocessor time.
770 double postprocessor_time_in_seconds;
771
772 // Some total of all time spent inside Ceres when Solve is called.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700773 double total_time_in_seconds;
774
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700775 // Time (in seconds) spent in the linear solver computing the
776 // trust region step.
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800777 double linear_solver_time_in_seconds;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700778
779 // Time (in seconds) spent evaluating the residual vector.
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800780 double residual_evaluation_time_in_seconds;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700781
782 // Time (in seconds) spent evaluating the jacobian matrix.
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800783 double jacobian_evaluation_time_in_seconds;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700784
785 // Time (in seconds) spent doing inner iterations.
Sameer Agarwal9f9488b2013-05-23 09:40:40 -0700786 double inner_iteration_time_in_seconds;
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800787
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700788 // Number of parameter blocks in the problem.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700789 int num_parameter_blocks;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700790
791 // Number of parameters in the probem.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700792 int num_parameters;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700793
794 // Dimension of the tangent space of the problem (or the number of
795 // columns in the Jacobian for the problem). This is different
796 // from num_parameters if a parameter block is associated with a
797 // LocalParameterization
Keir Mierleba944212013-02-25 12:46:44 -0800798 int num_effective_parameters;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700799
800 // Number of residual blocks in the problem.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700801 int num_residual_blocks;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700802
803 // Number of residuals in the problem.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700804 int num_residuals;
805
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700806 // Number of parameter blocks in the problem after the inactive
807 // and constant parameter blocks have been removed. A parameter
808 // block is inactive if no residual block refers to it.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700809 int num_parameter_blocks_reduced;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700810
811 // Number of parameters in the reduced problem.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700812 int num_parameters_reduced;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700813
814 // Dimension of the tangent space of the reduced problem (or the
815 // number of columns in the Jacobian for the reduced
816 // problem). This is different from num_parameters_reduced if a
817 // parameter block in the reduced problem is associated with a
818 // LocalParameterization.
Keir Mierleba944212013-02-25 12:46:44 -0800819 int num_effective_parameters_reduced;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700820
821 // Number of residual blocks in the reduced problem.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700822 int num_residual_blocks_reduced;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700823
824 // Number of residuals in the reduced problem.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700825 int num_residuals_reduced;
826
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700827 // Number of threads specified by the user for Jacobian and
828 // residual evaluation.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700829 int num_threads_given;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700830
831 // Number of threads actually used by the solver for Jacobian and
832 // residual evaluation. This number is not equal to
833 // num_threads_given if OpenMP is not available.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700834 int num_threads_used;
835
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700836 // Number of threads specified by the user for solving the trust
837 // region problem.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700838 int num_linear_solver_threads_given;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700839
840 // Number of threads actually used by the solver for solving the
841 // trust region problem. This number is not equal to
842 // num_threads_given if OpenMP is not available.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700843 int num_linear_solver_threads_used;
844
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700845 // Type of the linear solver requested by the user.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700846 LinearSolverType linear_solver_type_given;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700847
848 // Type of the linear solver actually used. This may be different
849 // from linear_solver_type_given if Ceres determines that the
850 // problem structure is not compatible with the linear solver
851 // requested or if the linear solver requested by the user is not
852 // available, e.g. The user requested SPARSE_NORMAL_CHOLESKY but
853 // no sparse linear algebra library was available.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700854 LinearSolverType linear_solver_type_used;
855
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700856 // Size of the elimination groups given by the user as hints to
857 // the linear solver.
Sameer Agarwal977be7c2013-01-26 16:01:54 -0800858 vector<int> linear_solver_ordering_given;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700859
860 // Size of the parameter groups used by the solver when ordering
861 // the columns of the Jacobian. This maybe different from
862 // linear_solver_ordering_given if the user left
863 // linear_solver_ordering_given blank and asked for an automatic
864 // ordering, or if the problem contains some constant or inactive
865 // parameter blocks.
Sameer Agarwal977be7c2013-01-26 16:01:54 -0800866 vector<int> linear_solver_ordering_used;
867
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700868 // True if the user asked for inner iterations to be used as part
869 // of the optimization.
Sameer Agarwal9f9488b2013-05-23 09:40:40 -0700870 bool inner_iterations_given;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700871
872 // True if the user asked for inner iterations to be used as part
873 // of the optimization and the problem structure was such that
874 // they were actually performed. e.g., in a problem with just one
875 // parameter block, inner iterations are not performed.
Sameer Agarwal9f9488b2013-05-23 09:40:40 -0700876 bool inner_iterations_used;
877
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700878 // Size of the parameter groups given by the user for performing
879 // inner iterations.
Sameer Agarwal9f9488b2013-05-23 09:40:40 -0700880 vector<int> inner_iteration_ordering_given;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700881
882 // Size of the parameter groups given used by the solver for
883 // performing inner iterations. This maybe different from
884 // inner_iteration_ordering_given if the user left
885 // inner_iteration_ordering_given blank and asked for an automatic
886 // ordering, or if the problem contains some constant or inactive
887 // parameter blocks.
Sameer Agarwal9f9488b2013-05-23 09:40:40 -0700888 vector<int> inner_iteration_ordering_used;
889
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700890 // Type of preconditioner used for solving the trust region
891 // step. Only meaningful when an iterative linear solver is used.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700892 PreconditionerType preconditioner_type;
Sameer Agarwal97fb6d92012-06-17 10:08:19 -0700893
Sameer Agarwalf06b9fa2013-10-27 21:38:13 -0700894 // Type of clustering algorithm used for visibility based
895 // preconditioning. Only meaningful when the preconditioner_type
896 // is CLUSTER_JACOBI or CLUSTER_TRIDIAGONAL.
897 VisibilityClusteringType visibility_clustering_type;
898
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700899 // Type of trust region strategy.
Sameer Agarwal97fb6d92012-06-17 10:08:19 -0700900 TrustRegionStrategyType trust_region_strategy_type;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700901
902 // Type of dogleg strategy used for solving the trust region
903 // problem.
Sameer Agarwal1e289202012-08-21 18:00:54 -0700904 DoglegType dogleg_type;
Sameer Agarwald010de52013-02-15 14:26:56 -0800905
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700906 // Type of the dense linear algebra library used.
Sameer Agarwal367b65e2013-08-09 10:35:37 -0700907 DenseLinearAlgebraLibraryType dense_linear_algebra_library_type;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700908
909 // Type of the sparse linear algebra library used.
Sameer Agarwal367b65e2013-08-09 10:35:37 -0700910 SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type;
Sameer Agarwal977be7c2013-01-26 16:01:54 -0800911
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700912 // Type of line search direction used.
Sameer Agarwald010de52013-02-15 14:26:56 -0800913 LineSearchDirectionType line_search_direction_type;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700914
915 // Type of the line search algorithm used.
Sameer Agarwald010de52013-02-15 14:26:56 -0800916 LineSearchType line_search_type;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700917
918 // When performing line search, the degree of the polynomial used
919 // to approximate the objective function.
Sameer Agarwal4f010b22013-07-01 08:01:01 -0700920 LineSearchInterpolationType line_search_interpolation_type;
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700921
922 // If the line search direction is NONLINEAR_CONJUGATE_GRADIENT,
923 // then this indicates the particular variant of non-linear
924 // conjugate gradient used.
Sameer Agarwal67ccb732013-07-03 06:28:34 -0700925 NonlinearConjugateGradientType nonlinear_conjugate_gradient_type;
Sameer Agarwal09244012013-06-30 14:33:23 -0700926
Sameer Agarwal4ad80b72013-10-21 05:33:34 -0700927 // If the type of the line search direction is LBFGS, then this
928 // indicates the rank of the Hessian approximation.
Sameer Agarwald010de52013-02-15 14:26:56 -0800929 int max_lbfgs_rank;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700930 };
931
932 // Once a least squares problem has been built, this function takes
933 // the problem and optimizes it based on the values of the options
934 // parameters. Upon return, a detailed summary of the work performed
935 // by the preprocessor, the non-linear minmizer and the linear
936 // solver are reported in the summary object.
937 virtual void Solve(const Options& options,
938 Problem* problem,
939 Solver::Summary* summary);
940};
941
942// Helper function which avoids going through the interface.
943void Solve(const Solver::Options& options,
944 Problem* problem,
945 Solver::Summary* summary);
946
947} // namespace ceres
948
949#endif // CERES_PUBLIC_SOLVER_H_