• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1 // Ceres Solver - A fast non-linear least squares minimizer
2 // Copyright 2010, 2011, 2012 Google Inc. All rights reserved.
3 // http://code.google.com/p/ceres-solver/
4 //
5 // Redistribution and use in source and binary forms, with or without
6 // modification, are permitted provided that the following conditions are met:
7 //
8 // * Redistributions of source code must retain the above copyright notice,
9 //   this list of conditions and the following disclaimer.
10 // * Redistributions in binary form must reproduce the above copyright notice,
11 //   this list of conditions and the following disclaimer in the documentation
12 //   and/or other materials provided with the distribution.
13 // * Neither the name of Google Inc. nor the names of its contributors may be
14 //   used to endorse or promote products derived from this software without
15 //   specific prior written permission.
16 //
17 // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
18 // AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
19 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
20 // ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
21 // LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
22 // CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
23 // SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
24 // INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
25 // CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
26 // ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
27 // POSSIBILITY OF SUCH DAMAGE.
28 //
29 // Author: sameeragarwal@google.com (Sameer Agarwal)
30 //
31 // Abstract interface for objects solving linear systems of various
32 // kinds.
33 
34 #ifndef CERES_INTERNAL_LINEAR_SOLVER_H_
35 #define CERES_INTERNAL_LINEAR_SOLVER_H_
36 
37 #include <cstddef>
38 #include <map>
39 #include <string>
40 #include <vector>
41 #include "ceres/block_sparse_matrix.h"
42 #include "ceres/casts.h"
43 #include "ceres/compressed_row_sparse_matrix.h"
44 #include "ceres/dense_sparse_matrix.h"
45 #include "ceres/execution_summary.h"
46 #include "ceres/triplet_sparse_matrix.h"
47 #include "ceres/types.h"
48 #include "glog/logging.h"
49 
50 namespace ceres {
51 namespace internal {
52 
53 enum LinearSolverTerminationType {
54   // Termination criterion was met.
55   LINEAR_SOLVER_SUCCESS,
56 
57   // Solver ran for max_num_iterations and terminated before the
58   // termination tolerance could be satisfied.
59   LINEAR_SOLVER_NO_CONVERGENCE,
60 
61   // Solver was terminated due to numerical problems, generally due to
62   // the linear system being poorly conditioned.
63   LINEAR_SOLVER_FAILURE,
64 
65   // Solver failed with a fatal error that cannot be recovered from,
66   // e.g. CHOLMOD ran out of memory when computing the symbolic or
67   // numeric factorization or an underlying library was called with
68   // the wrong arguments.
69   LINEAR_SOLVER_FATAL_ERROR
70 };
71 
72 
73 class LinearOperator;
74 
75 // Abstract base class for objects that implement algorithms for
76 // solving linear systems
77 //
78 //   Ax = b
79 //
80 // It is expected that a single instance of a LinearSolver object
81 // maybe used multiple times for solving multiple linear systems with
82 // the same sparsity structure. This allows them to cache and reuse
83 // information across solves. This means that calling Solve on the
84 // same LinearSolver instance with two different linear systems will
85 // result in undefined behaviour.
86 //
87 // Subclasses of LinearSolver use two structs to configure themselves.
88 // The Options struct configures the LinearSolver object for its
89 // lifetime. The PerSolveOptions struct is used to specify options for
90 // a particular Solve call.
91 class LinearSolver {
92  public:
93   struct Options {
OptionsOptions94     Options()
95         : type(SPARSE_NORMAL_CHOLESKY),
96           preconditioner_type(JACOBI),
97           visibility_clustering_type(CANONICAL_VIEWS),
98           dense_linear_algebra_library_type(EIGEN),
99           sparse_linear_algebra_library_type(SUITE_SPARSE),
100           use_postordering(false),
101           dynamic_sparsity(false),
102           min_num_iterations(1),
103           max_num_iterations(1),
104           num_threads(1),
105           residual_reset_period(10),
106           row_block_size(Eigen::Dynamic),
107           e_block_size(Eigen::Dynamic),
108           f_block_size(Eigen::Dynamic) {
109     }
110 
111     LinearSolverType type;
112     PreconditionerType preconditioner_type;
113     VisibilityClusteringType visibility_clustering_type;
114     DenseLinearAlgebraLibraryType dense_linear_algebra_library_type;
115     SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type;
116 
117     // See solver.h for information about this flag.
118     bool use_postordering;
119     bool dynamic_sparsity;
120 
121     // Number of internal iterations that the solver uses. This
122     // parameter only makes sense for iterative solvers like CG.
123     int min_num_iterations;
124     int max_num_iterations;
125 
126     // If possible, how many threads can the solver use.
127     int num_threads;
128 
129     // Hints about the order in which the parameter blocks should be
130     // eliminated by the linear solver.
131     //
132     // For example if elimination_groups is a vector of size k, then
133     // the linear solver is informed that it should eliminate the
134     // parameter blocks 0 ... elimination_groups[0] - 1 first, and
135     // then elimination_groups[0] ... elimination_groups[1] - 1 and so
136     // on. Within each elimination group, the linear solver is free to
137     // choose how the parameter blocks are ordered. Different linear
138     // solvers have differing requirements on elimination_groups.
139     //
140     // The most common use is for Schur type solvers, where there
141     // should be at least two elimination groups and the first
142     // elimination group must form an independent set in the normal
143     // equations. The first elimination group corresponds to the
144     // num_eliminate_blocks in the Schur type solvers.
145     vector<int> elimination_groups;
146 
147     // Iterative solvers, e.g. Preconditioned Conjugate Gradients
148     // maintain a cheap estimate of the residual which may become
149     // inaccurate over time. Thus for non-zero values of this
150     // parameter, the solver can be told to recalculate the value of
151     // the residual using a |b - Ax| evaluation.
152     int residual_reset_period;
153 
154     // If the block sizes in a BlockSparseMatrix are fixed, then in
155     // some cases the Schur complement based solvers can detect and
156     // specialize on them.
157     //
158     // It is expected that these parameters are set programmatically
159     // rather than manually.
160     //
161     // Please see schur_complement_solver.h and schur_eliminator.h for
162     // more details.
163     int row_block_size;
164     int e_block_size;
165     int f_block_size;
166   };
167 
168   // Options for the Solve method.
169   struct PerSolveOptions {
PerSolveOptionsPerSolveOptions170     PerSolveOptions()
171         : D(NULL),
172           preconditioner(NULL),
173           r_tolerance(0.0),
174           q_tolerance(0.0) {
175     }
176 
177     // This option only makes sense for unsymmetric linear solvers
178     // that can solve rectangular linear systems.
179     //
180     // Given a matrix A, an optional diagonal matrix D as a vector,
181     // and a vector b, the linear solver will solve for
182     //
183     //   | A | x = | b |
184     //   | D |     | 0 |
185     //
186     // If D is null, then it is treated as zero, and the solver returns
187     // the solution to
188     //
189     //   A x = b
190     //
191     // In either case, x is the vector that solves the following
192     // optimization problem.
193     //
194     //   arg min_x ||Ax - b||^2 + ||Dx||^2
195     //
196     // Here A is a matrix of size m x n, with full column rank. If A
197     // does not have full column rank, the results returned by the
198     // solver cannot be relied on. D, if it is not null is an array of
199     // size n.  b is an array of size m and x is an array of size n.
200     double * D;
201 
202     // This option only makes sense for iterative solvers.
203     //
204     // In general the performance of an iterative linear solver
205     // depends on the condition number of the matrix A. For example
206     // the convergence rate of the conjugate gradients algorithm
207     // is proportional to the square root of the condition number.
208     //
209     // One particularly useful technique for improving the
210     // conditioning of a linear system is to precondition it. In its
211     // simplest form a preconditioner is a matrix M such that instead
212     // of solving Ax = b, we solve the linear system AM^{-1} y = b
213     // instead, where M is such that the condition number k(AM^{-1})
214     // is smaller than the conditioner k(A). Given the solution to
215     // this system, x = M^{-1} y. The iterative solver takes care of
216     // the mechanics of solving the preconditioned system and
217     // returning the corrected solution x. The user only needs to
218     // supply a linear operator.
219     //
220     // A null preconditioner is equivalent to an identity matrix being
221     // used a preconditioner.
222     LinearOperator* preconditioner;
223 
224 
225     // The following tolerance related options only makes sense for
226     // iterative solvers. Direct solvers ignore them.
227 
228     // Solver terminates when
229     //
230     //   |Ax - b| <= r_tolerance * |b|.
231     //
232     // This is the most commonly used termination criterion for
233     // iterative solvers.
234     double r_tolerance;
235 
236     // For PSD matrices A, let
237     //
238     //   Q(x) = x'Ax - 2b'x
239     //
240     // be the cost of the quadratic function defined by A and b. Then,
241     // the solver terminates at iteration i if
242     //
243     //   i * (Q(x_i) - Q(x_i-1)) / Q(x_i) < q_tolerance.
244     //
245     // This termination criterion is more useful when using CG to
246     // solve the Newton step. This particular convergence test comes
247     // from Stephen Nash's work on truncated Newton
248     // methods. References:
249     //
250     //   1. Stephen G. Nash & Ariela Sofer, Assessing A Search
251     //      Direction Within A Truncated Newton Method, Operation
252     //      Research Letters 9(1990) 219-221.
253     //
254     //   2. Stephen G. Nash, A Survey of Truncated Newton Methods,
255     //      Journal of Computational and Applied Mathematics,
256     //      124(1-2), 45-59, 2000.
257     //
258     double q_tolerance;
259   };
260 
261   // Summary of a call to the Solve method. We should move away from
262   // the true/false method for determining solver success. We should
263   // let the summary object do the talking.
264   struct Summary {
SummarySummary265     Summary()
266         : residual_norm(0.0),
267           num_iterations(-1),
268           termination_type(LINEAR_SOLVER_FAILURE) {
269     }
270 
271     double residual_norm;
272     int num_iterations;
273     LinearSolverTerminationType termination_type;
274     string message;
275   };
276 
277   // If the optimization problem is such that there are no remaining
278   // e-blocks, a Schur type linear solver cannot be used. If the
279   // linear solver is of Schur type, this function implements a policy
280   // to select an alternate nearest linear solver to the one selected
281   // by the user. The input linear_solver_type is returned otherwise.
282   static LinearSolverType LinearSolverForZeroEBlocks(
283       LinearSolverType linear_solver_type);
284 
285   virtual ~LinearSolver();
286 
287   // Solve Ax = b.
288   virtual Summary Solve(LinearOperator* A,
289                         const double* b,
290                         const PerSolveOptions& per_solve_options,
291                         double* x) = 0;
292 
293   // The following two methods return copies instead of references so
294   // that the base class implementation does not have to worry about
295   // life time issues. Further, these calls are not expected to be
296   // frequent or performance sensitive.
CallStatistics()297   virtual map<string, int> CallStatistics() const {
298     return map<string, int>();
299   }
300 
TimeStatistics()301   virtual map<string, double> TimeStatistics() const {
302     return map<string, double>();
303   }
304 
305   // Factory
306   static LinearSolver* Create(const Options& options);
307 };
308 
309 // This templated subclass of LinearSolver serves as a base class for
310 // other linear solvers that depend on the particular matrix layout of
311 // the underlying linear operator. For example some linear solvers
312 // need low level access to the TripletSparseMatrix implementing the
313 // LinearOperator interface. This class hides those implementation
314 // details behind a private virtual method, and has the Solve method
315 // perform the necessary upcasting.
316 template <typename MatrixType>
317 class TypedLinearSolver : public LinearSolver {
318  public:
~TypedLinearSolver()319   virtual ~TypedLinearSolver() {}
Solve(LinearOperator * A,const double * b,const LinearSolver::PerSolveOptions & per_solve_options,double * x)320   virtual LinearSolver::Summary Solve(
321       LinearOperator* A,
322       const double* b,
323       const LinearSolver::PerSolveOptions& per_solve_options,
324       double* x) {
325     ScopedExecutionTimer total_time("LinearSolver::Solve", &execution_summary_);
326     CHECK_NOTNULL(A);
327     CHECK_NOTNULL(b);
328     CHECK_NOTNULL(x);
329     return SolveImpl(down_cast<MatrixType*>(A), b, per_solve_options, x);
330   }
331 
CallStatistics()332   virtual map<string, int> CallStatistics() const {
333     return execution_summary_.calls();
334   }
335 
TimeStatistics()336   virtual map<string, double> TimeStatistics() const {
337     return execution_summary_.times();
338   }
339 
340  private:
341   virtual LinearSolver::Summary SolveImpl(
342       MatrixType* A,
343       const double* b,
344       const LinearSolver::PerSolveOptions& per_solve_options,
345       double* x) = 0;
346 
347   ExecutionSummary execution_summary_;
348 };
349 
350 // Linear solvers that depend on acccess to the low level structure of
351 // a SparseMatrix.
352 typedef TypedLinearSolver<BlockSparseMatrix>         BlockSparseMatrixSolver;          // NOLINT
353 typedef TypedLinearSolver<CompressedRowSparseMatrix> CompressedRowSparseMatrixSolver;  // NOLINT
354 typedef TypedLinearSolver<DenseSparseMatrix>         DenseSparseMatrixSolver;          // NOLINT
355 typedef TypedLinearSolver<TripletSparseMatrix>       TripletSparseMatrixSolver;        // NOLINT
356 
357 }  // namespace internal
358 }  // namespace ceres
359 
360 #endif  // CERES_INTERNAL_LINEAR_SOLVER_H_
361