blob: 600d226ab162a0a4476b78b1412d620897dd35b6 [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// Enums and other top level class definitions.
32//
33// Note: internal/types.cc defines stringification routines for some
34// of these enums. Please update those routines if you extend or
35// remove enums from here.
36
37#ifndef CERES_PUBLIC_TYPES_H_
38#define CERES_PUBLIC_TYPES_H_
39
Alex Stewart9aa0e3c2013-07-05 20:22:37 +010040#include <string>
41
Sameer Agarwalcbae8562012-09-02 13:50:43 -070042#include "ceres/internal/port.h"
43
Keir Mierle8ebb0732012-04-30 23:09:08 -070044namespace ceres {
45
46// Basic integer types. These typedefs are in the Ceres namespace to avoid
47// conflicts with other packages having similar typedefs.
48typedef short int16;
49typedef int int32;
50
51// Argument type used in interfaces that can optionally take ownership
52// of a passed in argument. If TAKE_OWNERSHIP is passed, the called
53// object takes ownership of the pointer argument, and will call
54// delete on it upon completion.
55enum Ownership {
56 DO_NOT_TAKE_OWNERSHIP,
57 TAKE_OWNERSHIP
58};
59
60// TODO(keir): Considerably expand the explanations of each solver type.
61enum LinearSolverType {
62 // These solvers are for general rectangular systems formed from the
63 // normal equations A'A x = A'b. They are direct solvers and do not
64 // assume any special problem structure.
65
Sameer Agarwalb9f15a52012-08-18 13:06:19 -070066 // Solve the normal equations using a dense Cholesky solver; based
67 // on Eigen.
68 DENSE_NORMAL_CHOLESKY,
Keir Mierle8ebb0732012-04-30 23:09:08 -070069
70 // Solve the normal equations using a dense QR solver; based on
71 // Eigen.
72 DENSE_QR,
73
Sameer Agarwalb9f15a52012-08-18 13:06:19 -070074 // Solve the normal equations using a sparse cholesky solver; requires
75 // SuiteSparse or CXSparse.
76 SPARSE_NORMAL_CHOLESKY,
77
Keir Mierle8ebb0732012-04-30 23:09:08 -070078 // Specialized solvers, specific to problems with a generalized
79 // bi-partitite structure.
80
81 // Solves the reduced linear system using a dense Cholesky solver;
82 // based on Eigen.
83 DENSE_SCHUR,
84
85 // Solves the reduced linear system using a sparse Cholesky solver;
86 // based on CHOLMOD.
87 SPARSE_SCHUR,
88
89 // Solves the reduced linear system using Conjugate Gradients, based
90 // on a new Ceres implementation. Suitable for large scale
91 // problems.
92 ITERATIVE_SCHUR,
93
Keir Mierlee2a6cdc2012-05-07 06:39:56 -070094 // Conjugate gradients on the normal equations.
95 CGNR
Keir Mierle8ebb0732012-04-30 23:09:08 -070096};
97
98enum PreconditionerType {
99 // Trivial preconditioner - the identity matrix.
100 IDENTITY,
101
102 // Block diagonal of the Gauss-Newton Hessian.
103 JACOBI,
104
Sameer Agarwalf06b9fa2013-10-27 21:38:13 -0700105 // Note: The following three preconditioners can only be used with
106 // the ITERATIVE_SCHUR solver. They are well suited for Structure
107 // from Motion problems.
108
Keir Mierle8ebb0732012-04-30 23:09:08 -0700109 // Block diagonal of the Schur complement. This preconditioner may
Sameer Agarwal290b9752013-02-17 16:50:37 -0800110 // only be used with the ITERATIVE_SCHUR solver.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700111 SCHUR_JACOBI,
112
113 // Visibility clustering based preconditioners.
114 //
Sameer Agarwalf06b9fa2013-10-27 21:38:13 -0700115 // The following two preconditioners use the visibility structure of
116 // the scene to determine the sparsity structure of the
117 // preconditioner. This is done using a clustering algorithm. The
118 // available visibility clustering algorithms are described below.
119 //
120 // Note: Requires SuiteSparse.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700121 CLUSTER_JACOBI,
122 CLUSTER_TRIDIAGONAL
123};
124
Sameer Agarwalf06b9fa2013-10-27 21:38:13 -0700125enum VisibilityClusteringType {
126 // Canonical views algorithm as described in
127 //
128 // "Scene Summarization for Online Image Collections", Ian Simon, Noah
129 // Snavely, Steven M. Seitz, ICCV 2007.
130 //
131 // This clustering algorithm can be quite slow, but gives high
132 // quality clusters. The original visibility based clustering paper
133 // used this algorithm.
134 CANONICAL_VIEWS,
135
136 // The classic single linkage algorithm. It is extremely fast as
Sameer Agarwal9ba0b352013-11-05 13:04:56 -0800137 // compared to CANONICAL_VIEWS, but can give slightly poorer
Sameer Agarwalf06b9fa2013-10-27 21:38:13 -0700138 // results. For problems with large number of cameras though, this
139 // is generally a pretty good option.
140 //
141 // If you are using SCHUR_JACOBI preconditioner and have SuiteSparse
142 // available, CLUSTER_JACOBI and CLUSTER_TRIDIAGONAL in combination
143 // with the SINGLE_LINKAGE algorithm will generally give better
144 // results.
145 SINGLE_LINKAGE
146};
147
Sameer Agarwalb0518732012-05-29 00:27:57 -0700148enum SparseLinearAlgebraLibraryType {
149 // High performance sparse Cholesky factorization and approximate
150 // minimum degree ordering.
151 SUITE_SPARSE,
152
153 // A lightweight replacment for SuiteSparse.
154 CX_SPARSE
155};
156
Sameer Agarwal367b65e2013-08-09 10:35:37 -0700157enum DenseLinearAlgebraLibraryType {
158 EIGEN,
159 LAPACK
160};
161
Keir Mierle8ebb0732012-04-30 23:09:08 -0700162// Logging options
163// The options get progressively noisier.
164enum LoggingType {
165 SILENT,
Sameer Agarwal64472192012-05-03 21:53:07 -0700166 PER_MINIMIZER_ITERATION
Keir Mierle8ebb0732012-04-30 23:09:08 -0700167};
168
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800169enum MinimizerType {
170 LINE_SEARCH,
171 TRUST_REGION
172};
173
174enum LineSearchDirectionType {
175 // Negative of the gradient.
176 STEEPEST_DESCENT,
177
178 // A generalization of the Conjugate Gradient method to non-linear
179 // functions. The generalization can be performed in a number of
180 // different ways, resulting in a variety of search directions. The
181 // precise choice of the non-linear conjugate gradient algorithm
Sameer Agarwalc89ea4b2013-01-09 16:09:35 -0800182 // used is determined by NonlinerConjuateGradientType.
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800183 NONLINEAR_CONJUGATE_GRADIENT,
Sameer Agarwal3e8d1922012-11-28 17:20:22 -0800184
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100185 // BFGS, and it's limited memory approximation L-BFGS, are quasi-Newton
186 // algorithms that approximate the Hessian matrix by iteratively refining
187 // an initial estimate with rank-one updates using the gradient at each
188 // iteration. They are a generalisation of the Secant method and satisfy
189 // the Secant equation. The Secant equation has an infinium of solutions
190 // in multiple dimensions, as there are N*(N+1)/2 degrees of freedom in a
191 // symmetric matrix but only N conditions are specified by the Secant
192 // equation. The requirement that the Hessian approximation be positive
193 // definite imposes another N additional constraints, but that still leaves
194 // remaining degrees-of-freedom. (L)BFGS methods uniquely deteremine the
195 // approximate Hessian by imposing the additional constraints that the
196 // approximation at the next iteration must be the 'closest' to the current
197 // approximation (the nature of how this proximity is measured is actually
198 // the defining difference between a family of quasi-Newton methods including
199 // (L)BFGS & DFP). (L)BFGS is currently regarded as being the best known
200 // general quasi-Newton method.
Sameer Agarwal3e8d1922012-11-28 17:20:22 -0800201 //
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100202 // The principal difference between BFGS and L-BFGS is that whilst BFGS
203 // maintains a full, dense approximation to the (inverse) Hessian, L-BFGS
204 // maintains only a window of the last M observations of the parameters and
205 // gradients. Using this observation history, the calculation of the next
206 // search direction can be computed without requiring the construction of the
207 // full dense inverse Hessian approximation. This is particularly important
208 // for problems with a large number of parameters, where storage of an N-by-N
209 // matrix in memory would be prohibitive.
210 //
211 // For more details on BFGS see:
212 //
213 // Broyden, C.G., "The Convergence of a Class of Double-rank Minimization
214 // Algorithms,"; J. Inst. Maths. Applics., Vol. 6, pp 76–90, 1970.
215 //
216 // Fletcher, R., "A New Approach to Variable Metric Algorithms,"
217 // Computer Journal, Vol. 13, pp 317–322, 1970.
218 //
219 // Goldfarb, D., "A Family of Variable Metric Updates Derived by Variational
220 // Means," Mathematics of Computing, Vol. 24, pp 23–26, 1970.
221 //
222 // Shanno, D.F., "Conditioning of Quasi-Newton Methods for Function
223 // Minimization," Mathematics of Computing, Vol. 24, pp 647–656, 1970.
224 //
225 // For more details on L-BFGS see:
Sameer Agarwal3e8d1922012-11-28 17:20:22 -0800226 //
227 // Nocedal, J. (1980). "Updating Quasi-Newton Matrices with Limited
228 // Storage". Mathematics of Computation 35 (151): 773–782.
229 //
230 // Byrd, R. H.; Nocedal, J.; Schnabel, R. B. (1994).
231 // "Representations of Quasi-Newton Matrices and their use in
232 // Limited Memory Methods". Mathematical Programming 63 (4):
233 // 129–156.
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100234 //
235 // A general reference for both methods:
236 //
237 // Nocedal J., Wright S., Numerical Optimization, 2nd Ed. Springer, 1999.
Sameer Agarwal3e8d1922012-11-28 17:20:22 -0800238 LBFGS,
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100239 BFGS,
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800240};
241
242// Nonliner conjugate gradient methods are a generalization of the
243// method of Conjugate Gradients for linear systems. The
244// generalization can be carried out in a number of different ways
245// leading to number of different rules for computing the search
246// direction. Ceres provides a number of different variants. For more
247// details see Numerical Optimization by Nocedal & Wright.
248enum NonlinearConjugateGradientType {
249 FLETCHER_REEVES,
250 POLAK_RIBIRERE,
Sameer Agarwal1afd4982012-11-29 10:28:11 -0800251 HESTENES_STIEFEL,
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800252};
253
254enum LineSearchType {
255 // Backtracking line search with polynomial interpolation or
256 // bisection.
Sameer Agarwal1afd4982012-11-29 10:28:11 -0800257 ARMIJO,
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100258 WOLFE,
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800259};
260
Sameer Agarwalfa015192012-06-11 14:21:42 -0700261// Ceres supports different strategies for computing the trust region
262// step.
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700263enum TrustRegionStrategyType {
Sameer Agarwalfa015192012-06-11 14:21:42 -0700264 // The default trust region strategy is to use the step computation
265 // used in the Levenberg-Marquardt algorithm. For more details see
266 // levenberg_marquardt_strategy.h
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700267 LEVENBERG_MARQUARDT,
Sameer Agarwalfa015192012-06-11 14:21:42 -0700268
269 // Powell's dogleg algorithm interpolates between the Cauchy point
270 // and the Gauss-Newton step. It is particularly useful if the
271 // LEVENBERG_MARQUARDT algorithm is making a large number of
272 // unsuccessful steps. For more details see dogleg_strategy.h.
273 //
274 // NOTES:
275 //
276 // 1. This strategy has not been experimented with or tested as
277 // extensively as LEVENBERG_MARQUARDT, and therefore it should be
278 // considered EXPERIMENTAL for now.
279 //
280 // 2. For now this strategy should only be used with exact
281 // factorization based linear solvers, i.e., SPARSE_SCHUR,
282 // DENSE_SCHUR, DENSE_QR and SPARSE_NORMAL_CHOLESKY.
283 DOGLEG
Keir Mierle8ebb0732012-04-30 23:09:08 -0700284};
285
Markus Moll51cf7cb2012-08-20 20:10:20 +0200286// Ceres supports two different dogleg strategies.
287// The "traditional" dogleg method by Powell and the
288// "subspace" method described in
289// R. H. Byrd, R. B. Schnabel, and G. A. Shultz,
290// "Approximate solution of the trust region problem by minimization
291// over two-dimensional subspaces", Mathematical Programming,
292// 40 (1988), pp. 247--263
293enum DoglegType {
294 // The traditional approach constructs a dogleg path
295 // consisting of two line segments and finds the furthest
296 // point on that path that is still inside the trust region.
297 TRADITIONAL_DOGLEG,
298
299 // The subspace approach finds the exact minimum of the model
300 // constrained to the subspace spanned by the dogleg path.
301 SUBSPACE_DOGLEG
302};
303
Sameer Agarwaldcee1202013-12-07 21:48:56 -0800304enum TerminationType {
305 // Minimizer terminated because one of the convergence criterion set
306 // by the user was satisfied.
307 //
308 // 1. (new_cost - old_cost) < function_tolerance * old_cost;
309 // 2. max_i |gradient_i| < gradient_tolerance * max_i|initial_gradient_i|
310 // 3. |step|_2 <= parameter_tolerance * ( |x|_2 + parameter_tolerance)
311 //
312 // The user's parameter blocks will be updated with the solution.
313 CONVERGENCE,
Keir Mierle8ebb0732012-04-30 23:09:08 -0700314
Sameer Agarwaldcee1202013-12-07 21:48:56 -0800315 // The solver ran for maximum number of iterations or maximum amount
316 // of time specified by the user, but none of the convergence
317 // criterion specified by the user were met. The user's parameter
318 // blocks will be updated with the solution found so far.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700319 NO_CONVERGENCE,
320
Sameer Agarwaldcee1202013-12-07 21:48:56 -0800321 // The minimizer terminated because of an error. The user's
322 // parameter blocks will not be updated.
323 FAILURE,
Keir Mierle8ebb0732012-04-30 23:09:08 -0700324
325 // Using an IterationCallback object, user code can control the
326 // minimizer. The following enums indicate that the user code was
327 // responsible for termination.
Sameer Agarwaldcee1202013-12-07 21:48:56 -0800328 //
329 // Minimizer terminated successfully because a user
330 // IterationCallback returned SOLVER_TERMINATE_SUCCESSFULLY.
331 //
332 // The user's parameter blocks will be updated with the solution.
333 USER_SUCCESS,
Keir Mierle8ebb0732012-04-30 23:09:08 -0700334
Sameer Agarwaldcee1202013-12-07 21:48:56 -0800335 // Minimizer terminated because because a user IterationCallback
336 // returned SOLVER_ABORT.
337 //
338 // The user's parameter blocks will not be updated.
339 USER_FAILURE
Keir Mierle8ebb0732012-04-30 23:09:08 -0700340};
341
342// Enums used by the IterationCallback instances to indicate to the
343// solver whether it should continue solving, the user detected an
344// error or the solution is good enough and the solver should
345// terminate.
346enum CallbackReturnType {
347 // Continue solving to next iteration.
348 SOLVER_CONTINUE,
349
350 // Terminate solver, and do not update the parameter blocks upon
351 // return. Unless the user has set
352 // Solver:Options:::update_state_every_iteration, in which case the
353 // state would have been updated every iteration
354 // anyways. Solver::Summary::termination_type is set to USER_ABORT.
355 SOLVER_ABORT,
356
357 // Terminate solver, update state and
358 // return. Solver::Summary::termination_type is set to USER_SUCCESS.
359 SOLVER_TERMINATE_SUCCESSFULLY
360};
361
Sameer Agarwal82f4b882012-05-06 21:05:28 -0700362// The format in which linear least squares problems should be logged
363// when Solver::Options::lsqp_iterations_to_dump is non-empty.
364enum DumpFormatType {
365 // Print the linear least squares problem in a human readable format
366 // to stderr. The Jacobian is printed as a dense matrix. The vectors
367 // D, x and f are printed as dense vectors. This should only be used
368 // for small problems.
369 CONSOLE,
370
371 // Write out the linear least squares problem to the directory
Sameer Agarwal82f4b882012-05-06 21:05:28 -0700372 // pointed to by Solver::Options::lsqp_dump_directory as text files
373 // which can be read into MATLAB/Octave. The Jacobian is dumped as a
374 // text file containing (i,j,s) triplets, the vectors D, x and f are
375 // dumped as text files containing a list of their values.
376 //
377 // A MATLAB/octave script called lm_iteration_???.m is also output,
378 // which can be used to parse and load the problem into memory.
379 TEXTFILE
380};
381
Keir Mierlefdeb5772012-05-09 07:38:07 -0700382// For SizedCostFunction and AutoDiffCostFunction, DYNAMIC can be specified for
383// the number of residuals. If specified, then the number of residuas for that
384// cost function can vary at runtime.
385enum DimensionType {
386 DYNAMIC = -1
387};
388
Sameer Agarwal2fc0ed62013-01-15 11:34:10 -0800389enum NumericDiffMethod {
390 CENTRAL,
391 FORWARD
392};
393
Sameer Agarwal09244012013-06-30 14:33:23 -0700394enum LineSearchInterpolationType {
395 BISECTION,
396 QUADRATIC,
397 CUBIC
398};
399
Sameer Agarwal5a974712013-06-07 22:38:30 -0700400enum CovarianceAlgorithmType {
401 DENSE_SVD,
402 SPARSE_CHOLESKY,
403 SPARSE_QR
404};
405
Keir Mierle8ebb0732012-04-30 23:09:08 -0700406const char* LinearSolverTypeToString(LinearSolverType type);
Sameer Agarwalcbae8562012-09-02 13:50:43 -0700407bool StringToLinearSolverType(string value, LinearSolverType* type);
408
Keir Mierle8ebb0732012-04-30 23:09:08 -0700409const char* PreconditionerTypeToString(PreconditionerType type);
Sameer Agarwalcbae8562012-09-02 13:50:43 -0700410bool StringToPreconditionerType(string value, PreconditionerType* type);
411
Sameer Agarwalf06b9fa2013-10-27 21:38:13 -0700412const char* VisibilityClusteringTypeToString(VisibilityClusteringType type);
413bool StringToVisibilityClusteringType(string value,
414 VisibilityClusteringType* type);
415
Sameer Agarwalb0518732012-05-29 00:27:57 -0700416const char* SparseLinearAlgebraLibraryTypeToString(
417 SparseLinearAlgebraLibraryType type);
Sameer Agarwalcbae8562012-09-02 13:50:43 -0700418bool StringToSparseLinearAlgebraLibraryType(
419 string value,
420 SparseLinearAlgebraLibraryType* type);
421
Sameer Agarwal367b65e2013-08-09 10:35:37 -0700422const char* DenseLinearAlgebraLibraryTypeToString(
423 DenseLinearAlgebraLibraryType type);
424bool StringToDenseLinearAlgebraLibraryType(
425 string value,
426 DenseLinearAlgebraLibraryType* type);
427
Sameer Agarwalcbae8562012-09-02 13:50:43 -0700428const char* TrustRegionStrategyTypeToString(TrustRegionStrategyType type);
429bool StringToTrustRegionStrategyType(string value,
430 TrustRegionStrategyType* type);
431
432const char* DoglegTypeToString(DoglegType type);
433bool StringToDoglegType(string value, DoglegType* type);
434
Sameer Agarwal1afd4982012-11-29 10:28:11 -0800435const char* MinimizerTypeToString(MinimizerType type);
436bool StringToMinimizerType(string value, MinimizerType* type);
437
438const char* LineSearchDirectionTypeToString(LineSearchDirectionType type);
439bool StringToLineSearchDirectionType(string value,
440 LineSearchDirectionType* type);
441
442const char* LineSearchTypeToString(LineSearchType type);
443bool StringToLineSearchType(string value, LineSearchType* type);
444
445const char* NonlinearConjugateGradientTypeToString(
446 NonlinearConjugateGradientType type);
447bool StringToNonlinearConjugateGradientType(
Sameer Agarwal09244012013-06-30 14:33:23 -0700448 string value,
449 NonlinearConjugateGradientType* type);
450
451const char* LineSearchInterpolationTypeToString(
452 LineSearchInterpolationType type);
453bool StringToLineSearchInterpolationType(
454 string value,
455 LineSearchInterpolationType* type);
Sameer Agarwal1afd4982012-11-29 10:28:11 -0800456
Sameer Agarwal5a974712013-06-07 22:38:30 -0700457const char* CovarianceAlgorithmTypeToString(
458 CovarianceAlgorithmType type);
459bool StringToCovarianceAlgorithmType(
460 string value,
461 CovarianceAlgorithmType* type);
462
Sameer Agarwaldcee1202013-12-07 21:48:56 -0800463const char* TerminationTypeToString(TerminationType type);
Sameer Agarwal14ee7952012-09-06 11:05:32 -0700464
Keir Mierle8ebb0732012-04-30 23:09:08 -0700465bool IsSchurType(LinearSolverType type);
Sameer Agarwal14ee7952012-09-06 11:05:32 -0700466bool IsSparseLinearAlgebraLibraryTypeAvailable(
467 SparseLinearAlgebraLibraryType type);
Sameer Agarwal367b65e2013-08-09 10:35:37 -0700468bool IsDenseLinearAlgebraLibraryTypeAvailable(
469 DenseLinearAlgebraLibraryType type);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700470
471} // namespace ceres
472
473#endif // CERES_PUBLIC_TYPES_H_