blob: 219939edcc3c589eb38ae929c70b8ee207db9134 [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// 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>
Sameer Agarwal0c52f1e2012-09-17 11:30:14 -070038#include <vector>
Keir Mierle8ebb0732012-04-30 23:09:08 -070039
40#include <glog/logging.h>
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"
Sameer Agarwalbdd87c02013-01-29 16:24:31 -080045#include "ceres/execution_summary.h"
Sameer Agarwal42a84b82013-02-01 12:22:53 -080046#include "ceres/execution_summary.h"
Keir Mierle8ebb0732012-04-30 23:09:08 -070047#include "ceres/triplet_sparse_matrix.h"
48#include "ceres/types.h"
49
50namespace ceres {
51namespace internal {
52
53class LinearOperator;
54
55// Abstract base class for objects that implement algorithms for
56// solving linear systems
57//
58// Ax = b
59//
60// It is expected that a single instance of a LinearSolver object
Sameer Agarwala9d8ef82012-05-14 02:28:05 -070061// maybe used multiple times for solving multiple linear systems with
62// the same sparsity structure. This allows them to cache and reuse
63// information across solves. This means that calling Solve on the
64// same LinearSolver instance with two different linear systems will
65// result in undefined behaviour.
Keir Mierle8ebb0732012-04-30 23:09:08 -070066//
67// Subclasses of LinearSolver use two structs to configure themselves.
68// The Options struct configures the LinearSolver object for its
69// lifetime. The PerSolveOptions struct is used to specify options for
70// a particular Solve call.
71class LinearSolver {
72 public:
73 struct Options {
74 Options()
75 : type(SPARSE_NORMAL_CHOLESKY),
76 preconditioner_type(JACOBI),
Sameer Agarwalb0518732012-05-29 00:27:57 -070077 sparse_linear_algebra_library(SUITE_SPARSE),
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -070078 use_block_amd(true),
Keir Mierle8ebb0732012-04-30 23:09:08 -070079 min_num_iterations(1),
80 max_num_iterations(1),
81 num_threads(1),
Keir Mierle8ebb0732012-04-30 23:09:08 -070082 residual_reset_period(10),
83 row_block_size(Dynamic),
84 e_block_size(Dynamic),
85 f_block_size(Dynamic) {
86 }
87
88 LinearSolverType type;
89
90 PreconditionerType preconditioner_type;
91
Sameer Agarwalb0518732012-05-29 00:27:57 -070092 SparseLinearAlgebraLibraryType sparse_linear_algebra_library;
93
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -070094 // See solver.h for explanation of this option.
95 bool use_block_amd;
96
Keir Mierle8ebb0732012-04-30 23:09:08 -070097 // Number of internal iterations that the solver uses. This
98 // parameter only makes sense for iterative solvers like CG.
99 int min_num_iterations;
100 int max_num_iterations;
101
102 // If possible, how many threads can the solver use.
103 int num_threads;
104
Sameer Agarwal0c52f1e2012-09-17 11:30:14 -0700105 // Hints about the order in which the parameter blocks should be
106 // eliminated by the linear solver.
107 //
108 // For example if elimination_groups is a vector of size k, then
109 // the linear solver is informed that it should eliminate the
110 // parameter blocks 0 - elimination_groups[0] - 1 first, and then
111 // elimination_groups[0] - elimination_groups[1] and so on. Within
112 // each elimination group, the linear solver is free to choose how
113 // the parameter blocks are ordered. Different linear solvers have
114 // differing requirements on elimination_groups.
115 //
116 // The most common use is for Schur type solvers, where there
117 // should be at least two elimination groups and the first
118 // elimination group must form an independent set in the normal
119 // equations. The first elimination group corresponds to the
120 // num_eliminate_blocks in the Schur type solvers.
121 vector<int> elimination_groups;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700122
123 // Iterative solvers, e.g. Preconditioned Conjugate Gradients
124 // maintain a cheap estimate of the residual which may become
125 // inaccurate over time. Thus for non-zero values of this
126 // parameter, the solver can be told to recalculate the value of
127 // the residual using a |b - Ax| evaluation.
128 int residual_reset_period;
129
130 // If the block sizes in a BlockSparseMatrix are fixed, then in
131 // some cases the Schur complement based solvers can detect and
132 // specialize on them.
133 //
134 // It is expected that these parameters are set programmatically
135 // rather than manually.
136 //
Sameer Agarwala9d8ef82012-05-14 02:28:05 -0700137 // Please see schur_complement_solver.h and schur_eliminator.h for
138 // more details.
Keir Mierle8ebb0732012-04-30 23:09:08 -0700139 int row_block_size;
140 int e_block_size;
141 int f_block_size;
142 };
143
144 // Options for the Solve method.
145 struct PerSolveOptions {
146 PerSolveOptions()
147 : D(NULL),
148 preconditioner(NULL),
149 r_tolerance(0.0),
150 q_tolerance(0.0) {
151 }
152
153 // This option only makes sense for unsymmetric linear solvers
154 // that can solve rectangular linear systems.
155 //
156 // Given a matrix A, an optional diagonal matrix D as a vector,
157 // and a vector b, the linear solver will solve for
158 //
159 // | A | x = | b |
160 // | D | | 0 |
161 //
162 // If D is null, then it is treated as zero, and the solver returns
163 // the solution to
164 //
165 // A x = b
166 //
167 // In either case, x is the vector that solves the following
168 // optimization problem.
169 //
Keir Mierlef7898fb2012-05-05 20:55:08 -0700170 // arg min_x ||Ax - b||^2 + ||Dx||^2
Keir Mierle8ebb0732012-04-30 23:09:08 -0700171 //
172 // Here A is a matrix of size m x n, with full column rank. If A
173 // does not have full column rank, the results returned by the
174 // solver cannot be relied on. D, if it is not null is an array of
175 // size n. b is an array of size m and x is an array of size n.
176 double * D;
177
178 // This option only makes sense for iterative solvers.
179 //
180 // In general the performance of an iterative linear solver
181 // depends on the condition number of the matrix A. For example
182 // the convergence rate of the conjugate gradients algorithm
183 // is proportional to the square root of the condition number.
184 //
185 // One particularly useful technique for improving the
186 // conditioning of a linear system is to precondition it. In its
187 // simplest form a preconditioner is a matrix M such that instead
188 // of solving Ax = b, we solve the linear system AM^{-1} y = b
189 // instead, where M is such that the condition number k(AM^{-1})
190 // is smaller than the conditioner k(A). Given the solution to
191 // this system, x = M^{-1} y. The iterative solver takes care of
192 // the mechanics of solving the preconditioned system and
193 // returning the corrected solution x. The user only needs to
194 // supply a linear operator.
195 //
196 // A null preconditioner is equivalent to an identity matrix being
197 // used a preconditioner.
198 LinearOperator* preconditioner;
199
200
201 // The following tolerance related options only makes sense for
202 // iterative solvers. Direct solvers ignore them.
203
204 // Solver terminates when
205 //
206 // |Ax - b| <= r_tolerance * |b|.
207 //
208 // This is the most commonly used termination criterion for
209 // iterative solvers.
210 double r_tolerance;
211
212 // For PSD matrices A, let
213 //
214 // Q(x) = x'Ax - 2b'x
215 //
216 // be the cost of the quadratic function defined by A and b. Then,
217 // the solver terminates at iteration i if
218 //
219 // i * (Q(x_i) - Q(x_i-1)) / Q(x_i) < q_tolerance.
220 //
221 // This termination criterion is more useful when using CG to
222 // solve the Newton step. This particular convergence test comes
223 // from Stephen Nash's work on truncated Newton
224 // methods. References:
225 //
226 // 1. Stephen G. Nash & Ariela Sofer, Assessing A Search
227 // Direction Within A Truncated Newton Method, Operation
228 // Research Letters 9(1990) 219-221.
229 //
230 // 2. Stephen G. Nash, A Survey of Truncated Newton Methods,
231 // Journal of Computational and Applied Mathematics,
232 // 124(1-2), 45-59, 2000.
233 //
234 double q_tolerance;
235 };
236
237 // Summary of a call to the Solve method. We should move away from
238 // the true/false method for determining solver success. We should
239 // let the summary object do the talking.
240 struct Summary {
241 Summary()
242 : residual_norm(0.0),
243 num_iterations(-1),
244 termination_type(FAILURE) {
245 }
246
247 double residual_norm;
248 int num_iterations;
249 LinearSolverTerminationType termination_type;
250 };
251
252 virtual ~LinearSolver();
253
254 // Solve Ax = b.
255 virtual Summary Solve(LinearOperator* A,
256 const double* b,
257 const PerSolveOptions& per_solve_options,
258 double* x) = 0;
259
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800260 // The following two methods return copies instead of references so
261 // that the base class implementation does not have to worry about
262 // life time issues. Further, these calls are not expected to be
263 // frequent or performance sensitive.
264 virtual map<string, int> CallStatistics() const {
265 return map<string, int>();
266 }
267
268 virtual map<string, double> TimeStatistics() const {
269 return map<string, double>();
270 }
271
Sameer Agarwala9d8ef82012-05-14 02:28:05 -0700272 // Factory
Keir Mierle8ebb0732012-04-30 23:09:08 -0700273 static LinearSolver* Create(const Options& options);
274};
275
276// This templated subclass of LinearSolver serves as a base class for
277// other linear solvers that depend on the particular matrix layout of
278// the underlying linear operator. For example some linear solvers
279// need low level access to the TripletSparseMatrix implementing the
280// LinearOperator interface. This class hides those implementation
281// details behind a private virtual method, and has the Solve method
282// perform the necessary upcasting.
283template <typename MatrixType>
284class TypedLinearSolver : public LinearSolver {
285 public:
286 virtual ~TypedLinearSolver() {}
287 virtual LinearSolver::Summary Solve(
288 LinearOperator* A,
289 const double* b,
290 const LinearSolver::PerSolveOptions& per_solve_options,
291 double* x) {
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800292 ScopedExecutionTimer total_time("LinearSolver::Solve", &execution_summary_);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700293 CHECK_NOTNULL(A);
294 CHECK_NOTNULL(b);
295 CHECK_NOTNULL(x);
296 return SolveImpl(down_cast<MatrixType*>(A), b, per_solve_options, x);
297 }
298
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800299 virtual map<string, int> CallStatistics() const {
300 return execution_summary_.calls();
301 }
302
303 virtual map<string, double> TimeStatistics() const {
304 return execution_summary_.times();
305 }
306
Keir Mierle8ebb0732012-04-30 23:09:08 -0700307 private:
308 virtual LinearSolver::Summary SolveImpl(
309 MatrixType* A,
310 const double* b,
311 const LinearSolver::PerSolveOptions& per_solve_options,
312 double* x) = 0;
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800313
314 ExecutionSummary execution_summary_;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700315};
316
317// Linear solvers that depend on acccess to the low level structure of
318// a SparseMatrix.
319typedef TypedLinearSolver<BlockSparseMatrix> BlockSparseMatrixSolver; // NOLINT
320typedef TypedLinearSolver<BlockSparseMatrixBase> BlockSparseMatrixBaseSolver; // NOLINT
321typedef TypedLinearSolver<CompressedRowSparseMatrix> CompressedRowSparseMatrixSolver; // NOLINT
322typedef TypedLinearSolver<DenseSparseMatrix> DenseSparseMatrixSolver; // NOLINT
323typedef TypedLinearSolver<TripletSparseMatrix> TripletSparseMatrixSolver; // NOLINT
324
325} // namespace internal
326} // namespace ceres
327
328#endif // CERES_INTERNAL_LINEAR_SOLVER_H_