blob: e3dac3c8541d2786d2989b653d7c0048a0f890bc [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.
Keir Mierle8ebb0732012-04-30 23:09:08 -070048typedef int int32;
49
50// Argument type used in interfaces that can optionally take ownership
51// of a passed in argument. If TAKE_OWNERSHIP is passed, the called
52// object takes ownership of the pointer argument, and will call
53// delete on it upon completion.
54enum Ownership {
55 DO_NOT_TAKE_OWNERSHIP,
56 TAKE_OWNERSHIP
57};
58
59// TODO(keir): Considerably expand the explanations of each solver type.
60enum LinearSolverType {
61 // These solvers are for general rectangular systems formed from the
62 // normal equations A'A x = A'b. They are direct solvers and do not
63 // assume any special problem structure.
64
Sameer Agarwalb9f15a52012-08-18 13:06:19 -070065 // Solve the normal equations using a dense Cholesky solver; based
66 // on Eigen.
67 DENSE_NORMAL_CHOLESKY,
Keir Mierle8ebb0732012-04-30 23:09:08 -070068
69 // Solve the normal equations using a dense QR solver; based on
70 // Eigen.
71 DENSE_QR,
72
Sameer Agarwalb9f15a52012-08-18 13:06:19 -070073 // Solve the normal equations using a sparse cholesky solver; requires
74 // SuiteSparse or CXSparse.
75 SPARSE_NORMAL_CHOLESKY,
76
Keir Mierle8ebb0732012-04-30 23:09:08 -070077 // Specialized solvers, specific to problems with a generalized
78 // bi-partitite structure.
79
80 // Solves the reduced linear system using a dense Cholesky solver;
81 // based on Eigen.
82 DENSE_SCHUR,
83
84 // Solves the reduced linear system using a sparse Cholesky solver;
85 // based on CHOLMOD.
86 SPARSE_SCHUR,
87
88 // Solves the reduced linear system using Conjugate Gradients, based
89 // on a new Ceres implementation. Suitable for large scale
90 // problems.
91 ITERATIVE_SCHUR,
92
Keir Mierlee2a6cdc2012-05-07 06:39:56 -070093 // Conjugate gradients on the normal equations.
94 CGNR
Keir Mierle8ebb0732012-04-30 23:09:08 -070095};
96
97enum PreconditionerType {
98 // Trivial preconditioner - the identity matrix.
99 IDENTITY,
100
101 // Block diagonal of the Gauss-Newton Hessian.
102 JACOBI,
103
Sameer Agarwalf06b9fa2013-10-27 21:38:13 -0700104 // Note: The following three preconditioners can only be used with
105 // the ITERATIVE_SCHUR solver. They are well suited for Structure
106 // from Motion problems.
107
Keir Mierle8ebb0732012-04-30 23:09:08 -0700108 // Block diagonal of the Schur complement. This preconditioner may
Sameer Agarwal290b9752013-02-17 16:50:37 -0800109 // only be used with the ITERATIVE_SCHUR solver.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700110 SCHUR_JACOBI,
111
112 // Visibility clustering based preconditioners.
113 //
Sameer Agarwalf06b9fa2013-10-27 21:38:13 -0700114 // The following two preconditioners use the visibility structure of
115 // the scene to determine the sparsity structure of the
116 // preconditioner. This is done using a clustering algorithm. The
117 // available visibility clustering algorithms are described below.
118 //
119 // Note: Requires SuiteSparse.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700120 CLUSTER_JACOBI,
121 CLUSTER_TRIDIAGONAL
122};
123
Sameer Agarwalf06b9fa2013-10-27 21:38:13 -0700124enum VisibilityClusteringType {
125 // Canonical views algorithm as described in
126 //
127 // "Scene Summarization for Online Image Collections", Ian Simon, Noah
128 // Snavely, Steven M. Seitz, ICCV 2007.
129 //
130 // This clustering algorithm can be quite slow, but gives high
131 // quality clusters. The original visibility based clustering paper
132 // used this algorithm.
133 CANONICAL_VIEWS,
134
135 // The classic single linkage algorithm. It is extremely fast as
Sameer Agarwal9ba0b352013-11-05 13:04:56 -0800136 // compared to CANONICAL_VIEWS, but can give slightly poorer
Sameer Agarwalf06b9fa2013-10-27 21:38:13 -0700137 // results. For problems with large number of cameras though, this
138 // is generally a pretty good option.
139 //
140 // If you are using SCHUR_JACOBI preconditioner and have SuiteSparse
141 // available, CLUSTER_JACOBI and CLUSTER_TRIDIAGONAL in combination
142 // with the SINGLE_LINKAGE algorithm will generally give better
143 // results.
144 SINGLE_LINKAGE
145};
146
Sameer Agarwalb0518732012-05-29 00:27:57 -0700147enum SparseLinearAlgebraLibraryType {
148 // High performance sparse Cholesky factorization and approximate
149 // minimum degree ordering.
150 SUITE_SPARSE,
151
152 // A lightweight replacment for SuiteSparse.
153 CX_SPARSE
154};
155
Sameer Agarwal367b65e2013-08-09 10:35:37 -0700156enum DenseLinearAlgebraLibraryType {
157 EIGEN,
158 LAPACK
159};
160
Keir Mierle8ebb0732012-04-30 23:09:08 -0700161// Logging options
162// The options get progressively noisier.
163enum LoggingType {
164 SILENT,
Sameer Agarwal64472192012-05-03 21:53:07 -0700165 PER_MINIMIZER_ITERATION
Keir Mierle8ebb0732012-04-30 23:09:08 -0700166};
167
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800168enum MinimizerType {
169 LINE_SEARCH,
170 TRUST_REGION
171};
172
173enum LineSearchDirectionType {
174 // Negative of the gradient.
175 STEEPEST_DESCENT,
176
177 // A generalization of the Conjugate Gradient method to non-linear
178 // functions. The generalization can be performed in a number of
179 // different ways, resulting in a variety of search directions. The
180 // precise choice of the non-linear conjugate gradient algorithm
Sameer Agarwalc89ea4b2013-01-09 16:09:35 -0800181 // used is determined by NonlinerConjuateGradientType.
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800182 NONLINEAR_CONJUGATE_GRADIENT,
Sameer Agarwal3e8d1922012-11-28 17:20:22 -0800183
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100184 // BFGS, and it's limited memory approximation L-BFGS, are quasi-Newton
185 // algorithms that approximate the Hessian matrix by iteratively refining
186 // an initial estimate with rank-one updates using the gradient at each
187 // iteration. They are a generalisation of the Secant method and satisfy
188 // the Secant equation. The Secant equation has an infinium of solutions
189 // in multiple dimensions, as there are N*(N+1)/2 degrees of freedom in a
190 // symmetric matrix but only N conditions are specified by the Secant
191 // equation. The requirement that the Hessian approximation be positive
192 // definite imposes another N additional constraints, but that still leaves
193 // remaining degrees-of-freedom. (L)BFGS methods uniquely deteremine the
194 // approximate Hessian by imposing the additional constraints that the
195 // approximation at the next iteration must be the 'closest' to the current
196 // approximation (the nature of how this proximity is measured is actually
197 // the defining difference between a family of quasi-Newton methods including
198 // (L)BFGS & DFP). (L)BFGS is currently regarded as being the best known
199 // general quasi-Newton method.
Sameer Agarwal3e8d1922012-11-28 17:20:22 -0800200 //
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100201 // The principal difference between BFGS and L-BFGS is that whilst BFGS
202 // maintains a full, dense approximation to the (inverse) Hessian, L-BFGS
203 // maintains only a window of the last M observations of the parameters and
204 // gradients. Using this observation history, the calculation of the next
205 // search direction can be computed without requiring the construction of the
206 // full dense inverse Hessian approximation. This is particularly important
207 // for problems with a large number of parameters, where storage of an N-by-N
208 // matrix in memory would be prohibitive.
209 //
210 // For more details on BFGS see:
211 //
212 // Broyden, C.G., "The Convergence of a Class of Double-rank Minimization
213 // Algorithms,"; J. Inst. Maths. Applics., Vol. 6, pp 76–90, 1970.
214 //
215 // Fletcher, R., "A New Approach to Variable Metric Algorithms,"
216 // Computer Journal, Vol. 13, pp 317–322, 1970.
217 //
218 // Goldfarb, D., "A Family of Variable Metric Updates Derived by Variational
219 // Means," Mathematics of Computing, Vol. 24, pp 23–26, 1970.
220 //
221 // Shanno, D.F., "Conditioning of Quasi-Newton Methods for Function
222 // Minimization," Mathematics of Computing, Vol. 24, pp 647–656, 1970.
223 //
224 // For more details on L-BFGS see:
Sameer Agarwal3e8d1922012-11-28 17:20:22 -0800225 //
226 // Nocedal, J. (1980). "Updating Quasi-Newton Matrices with Limited
227 // Storage". Mathematics of Computation 35 (151): 773–782.
228 //
229 // Byrd, R. H.; Nocedal, J.; Schnabel, R. B. (1994).
230 // "Representations of Quasi-Newton Matrices and their use in
231 // Limited Memory Methods". Mathematical Programming 63 (4):
232 // 129–156.
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100233 //
234 // A general reference for both methods:
235 //
236 // Nocedal J., Wright S., Numerical Optimization, 2nd Ed. Springer, 1999.
Sameer Agarwal3e8d1922012-11-28 17:20:22 -0800237 LBFGS,
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100238 BFGS,
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800239};
240
241// Nonliner conjugate gradient methods are a generalization of the
242// method of Conjugate Gradients for linear systems. The
243// generalization can be carried out in a number of different ways
244// leading to number of different rules for computing the search
245// direction. Ceres provides a number of different variants. For more
246// details see Numerical Optimization by Nocedal & Wright.
247enum NonlinearConjugateGradientType {
248 FLETCHER_REEVES,
249 POLAK_RIBIRERE,
Sameer Agarwal1afd4982012-11-29 10:28:11 -0800250 HESTENES_STIEFEL,
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800251};
252
253enum LineSearchType {
254 // Backtracking line search with polynomial interpolation or
255 // bisection.
Sameer Agarwal1afd4982012-11-29 10:28:11 -0800256 ARMIJO,
Alex Stewart9aa0e3c2013-07-05 20:22:37 +0100257 WOLFE,
Sameer Agarwalf4d01642012-11-26 12:55:58 -0800258};
259
Sameer Agarwalfa015192012-06-11 14:21:42 -0700260// Ceres supports different strategies for computing the trust region
261// step.
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700262enum TrustRegionStrategyType {
Sameer Agarwalfa015192012-06-11 14:21:42 -0700263 // The default trust region strategy is to use the step computation
264 // used in the Levenberg-Marquardt algorithm. For more details see
265 // levenberg_marquardt_strategy.h
Sameer Agarwalaa9a83c2012-05-29 17:40:17 -0700266 LEVENBERG_MARQUARDT,
Sameer Agarwalfa015192012-06-11 14:21:42 -0700267
268 // Powell's dogleg algorithm interpolates between the Cauchy point
269 // and the Gauss-Newton step. It is particularly useful if the
270 // LEVENBERG_MARQUARDT algorithm is making a large number of
271 // unsuccessful steps. For more details see dogleg_strategy.h.
272 //
273 // NOTES:
274 //
275 // 1. This strategy has not been experimented with or tested as
276 // extensively as LEVENBERG_MARQUARDT, and therefore it should be
277 // considered EXPERIMENTAL for now.
278 //
279 // 2. For now this strategy should only be used with exact
280 // factorization based linear solvers, i.e., SPARSE_SCHUR,
281 // DENSE_SCHUR, DENSE_QR and SPARSE_NORMAL_CHOLESKY.
282 DOGLEG
Keir Mierle8ebb0732012-04-30 23:09:08 -0700283};
284
Markus Moll51cf7cb2012-08-20 20:10:20 +0200285// Ceres supports two different dogleg strategies.
286// The "traditional" dogleg method by Powell and the
287// "subspace" method described in
288// R. H. Byrd, R. B. Schnabel, and G. A. Shultz,
289// "Approximate solution of the trust region problem by minimization
290// over two-dimensional subspaces", Mathematical Programming,
291// 40 (1988), pp. 247--263
292enum DoglegType {
293 // The traditional approach constructs a dogleg path
294 // consisting of two line segments and finds the furthest
295 // point on that path that is still inside the trust region.
296 TRADITIONAL_DOGLEG,
297
298 // The subspace approach finds the exact minimum of the model
299 // constrained to the subspace spanned by the dogleg path.
300 SUBSPACE_DOGLEG
301};
302
Sameer Agarwaldcee1202013-12-07 21:48:56 -0800303enum TerminationType {
304 // Minimizer terminated because one of the convergence criterion set
305 // by the user was satisfied.
306 //
307 // 1. (new_cost - old_cost) < function_tolerance * old_cost;
308 // 2. max_i |gradient_i| < gradient_tolerance * max_i|initial_gradient_i|
309 // 3. |step|_2 <= parameter_tolerance * ( |x|_2 + parameter_tolerance)
310 //
311 // The user's parameter blocks will be updated with the solution.
312 CONVERGENCE,
Keir Mierle8ebb0732012-04-30 23:09:08 -0700313
Sameer Agarwaldcee1202013-12-07 21:48:56 -0800314 // The solver ran for maximum number of iterations or maximum amount
315 // of time specified by the user, but none of the convergence
316 // criterion specified by the user were met. The user's parameter
317 // blocks will be updated with the solution found so far.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700318 NO_CONVERGENCE,
319
Sameer Agarwaldcee1202013-12-07 21:48:56 -0800320 // The minimizer terminated because of an error. The user's
321 // parameter blocks will not be updated.
322 FAILURE,
Keir Mierle8ebb0732012-04-30 23:09:08 -0700323
324 // Using an IterationCallback object, user code can control the
325 // minimizer. The following enums indicate that the user code was
326 // responsible for termination.
Sameer Agarwaldcee1202013-12-07 21:48:56 -0800327 //
328 // Minimizer terminated successfully because a user
329 // IterationCallback returned SOLVER_TERMINATE_SUCCESSFULLY.
330 //
331 // The user's parameter blocks will be updated with the solution.
332 USER_SUCCESS,
Keir Mierle8ebb0732012-04-30 23:09:08 -0700333
Sameer Agarwaldcee1202013-12-07 21:48:56 -0800334 // Minimizer terminated because because a user IterationCallback
335 // returned SOLVER_ABORT.
336 //
337 // The user's parameter blocks will not be updated.
338 USER_FAILURE
Keir Mierle8ebb0732012-04-30 23:09:08 -0700339};
340
341// Enums used by the IterationCallback instances to indicate to the
342// solver whether it should continue solving, the user detected an
343// error or the solution is good enough and the solver should
344// terminate.
345enum CallbackReturnType {
346 // Continue solving to next iteration.
347 SOLVER_CONTINUE,
348
349 // Terminate solver, and do not update the parameter blocks upon
350 // return. Unless the user has set
351 // Solver:Options:::update_state_every_iteration, in which case the
352 // state would have been updated every iteration
353 // anyways. Solver::Summary::termination_type is set to USER_ABORT.
354 SOLVER_ABORT,
355
356 // Terminate solver, update state and
357 // return. Solver::Summary::termination_type is set to USER_SUCCESS.
358 SOLVER_TERMINATE_SUCCESSFULLY
359};
360
Sameer Agarwal82f4b882012-05-06 21:05:28 -0700361// The format in which linear least squares problems should be logged
362// when Solver::Options::lsqp_iterations_to_dump is non-empty.
363enum DumpFormatType {
364 // Print the linear least squares problem in a human readable format
365 // to stderr. The Jacobian is printed as a dense matrix. The vectors
366 // D, x and f are printed as dense vectors. This should only be used
367 // for small problems.
368 CONSOLE,
369
370 // Write out the linear least squares problem to the directory
Sameer Agarwal82f4b882012-05-06 21:05:28 -0700371 // pointed to by Solver::Options::lsqp_dump_directory as text files
372 // which can be read into MATLAB/Octave. The Jacobian is dumped as a
373 // text file containing (i,j,s) triplets, the vectors D, x and f are
374 // dumped as text files containing a list of their values.
375 //
376 // A MATLAB/octave script called lm_iteration_???.m is also output,
377 // which can be used to parse and load the problem into memory.
378 TEXTFILE
379};
380
Keir Mierlefdeb5772012-05-09 07:38:07 -0700381// For SizedCostFunction and AutoDiffCostFunction, DYNAMIC can be specified for
382// the number of residuals. If specified, then the number of residuas for that
383// cost function can vary at runtime.
384enum DimensionType {
385 DYNAMIC = -1
386};
387
Sameer Agarwal2fc0ed62013-01-15 11:34:10 -0800388enum NumericDiffMethod {
389 CENTRAL,
390 FORWARD
391};
392
Sameer Agarwal09244012013-06-30 14:33:23 -0700393enum LineSearchInterpolationType {
394 BISECTION,
395 QUADRATIC,
396 CUBIC
397};
398
Sameer Agarwal5a974712013-06-07 22:38:30 -0700399enum CovarianceAlgorithmType {
400 DENSE_SVD,
401 SPARSE_CHOLESKY,
402 SPARSE_QR
403};
404
Keir Mierle8ebb0732012-04-30 23:09:08 -0700405const char* LinearSolverTypeToString(LinearSolverType type);
Sameer Agarwalcbae8562012-09-02 13:50:43 -0700406bool StringToLinearSolverType(string value, LinearSolverType* type);
407
Keir Mierle8ebb0732012-04-30 23:09:08 -0700408const char* PreconditionerTypeToString(PreconditionerType type);
Sameer Agarwalcbae8562012-09-02 13:50:43 -0700409bool StringToPreconditionerType(string value, PreconditionerType* type);
410
Sameer Agarwalf06b9fa2013-10-27 21:38:13 -0700411const char* VisibilityClusteringTypeToString(VisibilityClusteringType type);
412bool StringToVisibilityClusteringType(string value,
413 VisibilityClusteringType* type);
414
Sameer Agarwalb0518732012-05-29 00:27:57 -0700415const char* SparseLinearAlgebraLibraryTypeToString(
416 SparseLinearAlgebraLibraryType type);
Sameer Agarwalcbae8562012-09-02 13:50:43 -0700417bool StringToSparseLinearAlgebraLibraryType(
418 string value,
419 SparseLinearAlgebraLibraryType* type);
420
Sameer Agarwal367b65e2013-08-09 10:35:37 -0700421const char* DenseLinearAlgebraLibraryTypeToString(
422 DenseLinearAlgebraLibraryType type);
423bool StringToDenseLinearAlgebraLibraryType(
424 string value,
425 DenseLinearAlgebraLibraryType* type);
426
Sameer Agarwalcbae8562012-09-02 13:50:43 -0700427const char* TrustRegionStrategyTypeToString(TrustRegionStrategyType type);
428bool StringToTrustRegionStrategyType(string value,
429 TrustRegionStrategyType* type);
430
431const char* DoglegTypeToString(DoglegType type);
432bool StringToDoglegType(string value, DoglegType* type);
433
Sameer Agarwal1afd4982012-11-29 10:28:11 -0800434const char* MinimizerTypeToString(MinimizerType type);
435bool StringToMinimizerType(string value, MinimizerType* type);
436
437const char* LineSearchDirectionTypeToString(LineSearchDirectionType type);
438bool StringToLineSearchDirectionType(string value,
439 LineSearchDirectionType* type);
440
441const char* LineSearchTypeToString(LineSearchType type);
442bool StringToLineSearchType(string value, LineSearchType* type);
443
444const char* NonlinearConjugateGradientTypeToString(
445 NonlinearConjugateGradientType type);
446bool StringToNonlinearConjugateGradientType(
Sameer Agarwal09244012013-06-30 14:33:23 -0700447 string value,
448 NonlinearConjugateGradientType* type);
449
450const char* LineSearchInterpolationTypeToString(
451 LineSearchInterpolationType type);
452bool StringToLineSearchInterpolationType(
453 string value,
454 LineSearchInterpolationType* type);
Sameer Agarwal1afd4982012-11-29 10:28:11 -0800455
Sameer Agarwal5a974712013-06-07 22:38:30 -0700456const char* CovarianceAlgorithmTypeToString(
457 CovarianceAlgorithmType type);
458bool StringToCovarianceAlgorithmType(
459 string value,
460 CovarianceAlgorithmType* type);
461
Sameer Agarwaldcee1202013-12-07 21:48:56 -0800462const char* TerminationTypeToString(TerminationType type);
Sameer Agarwal14ee7952012-09-06 11:05:32 -0700463
Keir Mierle8ebb0732012-04-30 23:09:08 -0700464bool IsSchurType(LinearSolverType type);
Sameer Agarwal14ee7952012-09-06 11:05:32 -0700465bool IsSparseLinearAlgebraLibraryTypeAvailable(
466 SparseLinearAlgebraLibraryType type);
Sameer Agarwal367b65e2013-08-09 10:35:37 -0700467bool IsDenseLinearAlgebraLibraryTypeAvailable(
468 DenseLinearAlgebraLibraryType type);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700469
470} // namespace ceres
471
472#endif // CERES_PUBLIC_TYPES_H_