blob: d477f31702594cea118031fb0143ad2ba18d1225 [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
Austin Schuh3de38b02024-06-25 18:25:10 -07002// Copyright 2023 Google Inc. All rights reserved.
Austin Schuh70cc9552019-01-21 19:46:48 -08003// http://ceres-solver.org/
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_COVARIANCE_H_
32#define CERES_PUBLIC_COVARIANCE_H_
33
34#include <memory>
35#include <utility>
36#include <vector>
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080037
Austin Schuh3de38b02024-06-25 18:25:10 -070038#include "ceres/internal/config.h"
Austin Schuh70cc9552019-01-21 19:46:48 -080039#include "ceres/internal/disable_warnings.h"
Austin Schuh3de38b02024-06-25 18:25:10 -070040#include "ceres/internal/export.h"
Austin Schuh70cc9552019-01-21 19:46:48 -080041#include "ceres/types.h"
42
43namespace ceres {
44
45class Problem;
46
47namespace internal {
48class CovarianceImpl;
49} // namespace internal
50
51// WARNING
52// =======
53// It is very easy to use this class incorrectly without understanding
54// the underlying mathematics. Please read and understand the
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080055// documentation completely before attempting to use it.
Austin Schuh70cc9552019-01-21 19:46:48 -080056//
57//
58// This class allows the user to evaluate the covariance for a
59// non-linear least squares problem and provides random access to its
60// blocks
61//
62// Background
63// ==========
64// One way to assess the quality of the solution returned by a
65// non-linear least squares solver is to analyze the covariance of the
66// solution.
67//
68// Let us consider the non-linear regression problem
69//
70// y = f(x) + N(0, I)
71//
72// i.e., the observation y is a random non-linear function of the
73// independent variable x with mean f(x) and identity covariance. Then
74// the maximum likelihood estimate of x given observations y is the
75// solution to the non-linear least squares problem:
76//
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080077// x* = arg min_x |f(x) - y|^2
Austin Schuh70cc9552019-01-21 19:46:48 -080078//
79// And the covariance of x* is given by
80//
81// C(x*) = inverse[J'(x*)J(x*)]
82//
83// Here J(x*) is the Jacobian of f at x*. The above formula assumes
84// that J(x*) has full column rank.
85//
86// If J(x*) is rank deficient, then the covariance matrix C(x*) is
87// also rank deficient and is given by
88//
89// C(x*) = pseudoinverse[J'(x*)J(x*)]
90//
91// Note that in the above, we assumed that the covariance
92// matrix for y was identity. This is an important assumption. If this
93// is not the case and we have
94//
95// y = f(x) + N(0, S)
96//
97// Where S is a positive semi-definite matrix denoting the covariance
98// of y, then the maximum likelihood problem to be solved is
99//
100// x* = arg min_x f'(x) inverse[S] f(x)
101//
102// and the corresponding covariance estimate of x* is given by
103//
104// C(x*) = inverse[J'(x*) inverse[S] J(x*)]
105//
106// So, if it is the case that the observations being fitted to have a
107// covariance matrix not equal to identity, then it is the user's
108// responsibility that the corresponding cost functions are correctly
109// scaled, e.g. in the above case the cost function for this problem
110// should evaluate S^{-1/2} f(x) instead of just f(x), where S^{-1/2}
111// is the inverse square root of the covariance matrix S.
112//
113// This class allows the user to evaluate the covariance for a
114// non-linear least squares problem and provides random access to its
115// blocks. The computation assumes that the CostFunctions compute
116// residuals such that their covariance is identity.
117//
118// Since the computation of the covariance matrix requires computing
119// the inverse of a potentially large matrix, this can involve a
120// rather large amount of time and memory. However, it is usually the
121// case that the user is only interested in a small part of the
122// covariance matrix. Quite often just the block diagonal. This class
123// allows the user to specify the parts of the covariance matrix that
124// she is interested in and then uses this information to only compute
125// and store those parts of the covariance matrix.
126//
127// Rank of the Jacobian
128// --------------------
129// As we noted above, if the jacobian is rank deficient, then the
130// inverse of J'J is not defined and instead a pseudo inverse needs to
131// be computed.
132//
133// The rank deficiency in J can be structural -- columns which are
134// always known to be zero or numerical -- depending on the exact
135// values in the Jacobian.
136//
137// Structural rank deficiency occurs when the problem contains
138// parameter blocks that are constant. This class correctly handles
139// structural rank deficiency like that.
140//
141// Numerical rank deficiency, where the rank of the matrix cannot be
142// predicted by its sparsity structure and requires looking at its
143// numerical values is more complicated. Here again there are two
144// cases.
145//
146// a. The rank deficiency arises from overparameterization. e.g., a
147// four dimensional quaternion used to parameterize SO(3), which is
148// a three dimensional manifold. In cases like this, the user should
Austin Schuh3de38b02024-06-25 18:25:10 -0700149// use an appropriate Manifold. Not only will this lead
Austin Schuh70cc9552019-01-21 19:46:48 -0800150// to better numerical behaviour of the Solver, it will also expose
151// the rank deficiency to the Covariance object so that it can
152// handle it correctly.
153//
154// b. More general numerical rank deficiency in the Jacobian
155// requires the computation of the so called Singular Value
156// Decomposition (SVD) of J'J. We do not know how to do this for
157// large sparse matrices efficiently. For small and moderate sized
158// problems this is done using dense linear algebra.
159//
160// Gauge Invariance
161// ----------------
162// In structure from motion (3D reconstruction) problems, the
163// reconstruction is ambiguous up to a similarity transform. This is
164// known as a Gauge Ambiguity. Handling Gauges correctly requires the
165// use of SVD or custom inversion algorithms. For small problems the
166// user can use the dense algorithm. For more details see
167//
168// Ken-ichi Kanatani, Daniel D. Morris: Gauges and gauge
169// transformations for uncertainty description of geometric structure
170// with indeterminacy. IEEE Transactions on Information Theory 47(5):
171// 2017-2028 (2001)
172//
173// Example Usage
174// =============
175//
176// double x[3];
177// double y[2];
178//
179// Problem problem;
180// problem.AddParameterBlock(x, 3);
181// problem.AddParameterBlock(y, 2);
182// <Build Problem>
183// <Solve Problem>
184//
185// Covariance::Options options;
186// Covariance covariance(options);
187//
188// std::vector<std::pair<const double*, const double*>> covariance_blocks;
189// covariance_blocks.push_back(make_pair(x, x));
190// covariance_blocks.push_back(make_pair(y, y));
191// covariance_blocks.push_back(make_pair(x, y));
192//
193// CHECK(covariance.Compute(covariance_blocks, &problem));
194//
195// double covariance_xx[3 * 3];
196// double covariance_yy[2 * 2];
197// double covariance_xy[3 * 2];
198// covariance.GetCovarianceBlock(x, x, covariance_xx)
199// covariance.GetCovarianceBlock(y, y, covariance_yy)
200// covariance.GetCovarianceBlock(x, y, covariance_xy)
201//
202class CERES_EXPORT Covariance {
203 public:
204 struct CERES_EXPORT Options {
205 // Sparse linear algebra library to use when a sparse matrix
206 // factorization is being used to compute the covariance matrix.
207 //
208 // Currently this only applies to SPARSE_QR.
209 SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type =
210#if !defined(CERES_NO_SUITESPARSE)
211 SUITE_SPARSE;
212#else
213 // Eigen's QR factorization is always available.
214 EIGEN_SPARSE;
215#endif
216
217 // Ceres supports two different algorithms for covariance
218 // estimation, which represent different tradeoffs in speed,
219 // accuracy and reliability.
220 //
221 // 1. DENSE_SVD uses Eigen's JacobiSVD to perform the
222 // computations. It computes the singular value decomposition
223 //
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800224 // U * D * V' = J
Austin Schuh70cc9552019-01-21 19:46:48 -0800225 //
226 // and then uses it to compute the pseudo inverse of J'J as
227 //
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800228 // pseudoinverse[J'J] = V * pseudoinverse[D^2] * V'
Austin Schuh70cc9552019-01-21 19:46:48 -0800229 //
230 // It is an accurate but slow method and should only be used
231 // for small to moderate sized problems. It can handle
232 // full-rank as well as rank deficient Jacobians.
233 //
234 // 2. SPARSE_QR uses the sparse QR factorization algorithm
235 // to compute the decomposition
236 //
237 // Q * R = J
238 //
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800239 // [J'J]^-1 = [R'*R]^-1
Austin Schuh70cc9552019-01-21 19:46:48 -0800240 //
241 // SPARSE_QR is not capable of computing the covariance if the
242 // Jacobian is rank deficient. Depending on the value of
243 // Covariance::Options::sparse_linear_algebra_library_type, either
244 // Eigen's Sparse QR factorization algorithm will be used or
245 // SuiteSparse's high performance SuiteSparseQR algorithm will be
246 // used.
247 CovarianceAlgorithmType algorithm_type = SPARSE_QR;
248
Austin Schuh3de38b02024-06-25 18:25:10 -0700249 // During QR factorization, if a column with Euclidean norm less
250 // than column_pivot_threshold is encountered it is treated as
251 // zero.
252 //
253 // If column_pivot_threshold < 0, then an automatic default value
254 // of 20*(m+n)*eps*sqrt(max(diag(J’*J))) is used. Here m and n are
255 // the number of rows and columns of the Jacobian (J)
256 // respectively.
257 //
258 // This is an advanced option meant for users who know enough
259 // about their Jacobian matrices that they can determine a value
260 // better than the default.
261 double column_pivot_threshold = -1;
262
Austin Schuh70cc9552019-01-21 19:46:48 -0800263 // If the Jacobian matrix is near singular, then inverting J'J
264 // will result in unreliable results, e.g, if
265 //
266 // J = [1.0 1.0 ]
267 // [1.0 1.0000001 ]
268 //
269 // which is essentially a rank deficient matrix, we have
270 //
271 // inv(J'J) = [ 2.0471e+14 -2.0471e+14]
272 // [-2.0471e+14 2.0471e+14]
273 //
274 // This is not a useful result. Therefore, by default
275 // Covariance::Compute will return false if a rank deficient
276 // Jacobian is encountered. How rank deficiency is detected
277 // depends on the algorithm being used.
278 //
279 // 1. DENSE_SVD
280 //
281 // min_sigma / max_sigma < sqrt(min_reciprocal_condition_number)
282 //
Austin Schuh3de38b02024-06-25 18:25:10 -0700283 // where min_sigma and max_sigma are the minimum and maximum
Austin Schuh70cc9552019-01-21 19:46:48 -0800284 // singular values of J respectively.
285 //
286 // 2. SPARSE_QR
287 //
288 // rank(J) < num_col(J)
289 //
290 // Here rank(J) is the estimate of the rank of J returned by the
291 // sparse QR factorization algorithm. It is a fairly reliable
292 // indication of rank deficiency.
293 //
294 double min_reciprocal_condition_number = 1e-14;
295
296 // When using DENSE_SVD, the user has more control in dealing with
297 // singular and near singular covariance matrices.
298 //
299 // As mentioned above, when the covariance matrix is near
300 // singular, instead of computing the inverse of J'J, the
301 // Moore-Penrose pseudoinverse of J'J should be computed.
302 //
303 // If J'J has the eigen decomposition (lambda_i, e_i), where
304 // lambda_i is the i^th eigenvalue and e_i is the corresponding
305 // eigenvector, then the inverse of J'J is
306 //
307 // inverse[J'J] = sum_i e_i e_i' / lambda_i
308 //
309 // and computing the pseudo inverse involves dropping terms from
310 // this sum that correspond to small eigenvalues.
311 //
312 // How terms are dropped is controlled by
313 // min_reciprocal_condition_number and null_space_rank.
314 //
315 // If null_space_rank is non-negative, then the smallest
316 // null_space_rank eigenvalue/eigenvectors are dropped
317 // irrespective of the magnitude of lambda_i. If the ratio of the
318 // smallest non-zero eigenvalue to the largest eigenvalue in the
319 // truncated matrix is still below
320 // min_reciprocal_condition_number, then the Covariance::Compute()
321 // will fail and return false.
322 //
323 // Setting null_space_rank = -1 drops all terms for which
324 //
325 // lambda_i / lambda_max < min_reciprocal_condition_number.
326 //
327 // This option has no effect on the SUITE_SPARSE_QR and
328 // EIGEN_SPARSE_QR algorithms.
329 int null_space_rank = 0;
330
331 int num_threads = 1;
332
333 // Even though the residual blocks in the problem may contain loss
334 // functions, setting apply_loss_function to false will turn off
335 // the application of the loss function to the output of the cost
336 // function and in turn its effect on the covariance.
337 //
338 // TODO(sameergaarwal): Expand this based on Jim's experiments.
339 bool apply_loss_function = true;
340 };
341
342 explicit Covariance(const Options& options);
343 ~Covariance();
344
345 // Compute a part of the covariance matrix.
346 //
347 // The vector covariance_blocks, indexes into the covariance matrix
348 // block-wise using pairs of parameter blocks. This allows the
349 // covariance estimation algorithm to only compute and store these
350 // blocks.
351 //
352 // Since the covariance matrix is symmetric, if the user passes
353 // (block1, block2), then GetCovarianceBlock can be called with
354 // block1, block2 as well as block2, block1.
355 //
356 // covariance_blocks cannot contain duplicates. Bad things will
357 // happen if they do.
358 //
359 // Note that the list of covariance_blocks is only used to determine
360 // what parts of the covariance matrix are computed. The full
361 // Jacobian is used to do the computation, i.e. they do not have an
362 // impact on what part of the Jacobian is used for computation.
363 //
364 // The return value indicates the success or failure of the
365 // covariance computation. Please see the documentation for
366 // Covariance::Options for more on the conditions under which this
367 // function returns false.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800368 bool Compute(const std::vector<std::pair<const double*, const double*>>&
369 covariance_blocks,
370 Problem* problem);
Austin Schuh70cc9552019-01-21 19:46:48 -0800371
372 // Compute a part of the covariance matrix.
373 //
374 // The vector parameter_blocks contains the parameter blocks that
375 // are used for computing the covariance matrix. From this vector
376 // all covariance pairs are generated. This allows the covariance
377 // estimation algorithm to only compute and store these blocks.
378 //
379 // parameter_blocks cannot contain duplicates. Bad things will
380 // happen if they do.
381 //
382 // Note that the list of covariance_blocks is only used to determine
383 // what parts of the covariance matrix are computed. The full
384 // Jacobian is used to do the computation, i.e. they do not have an
385 // impact on what part of the Jacobian is used for computation.
386 //
387 // The return value indicates the success or failure of the
388 // covariance computation. Please see the documentation for
389 // Covariance::Options for more on the conditions under which this
390 // function returns false.
391 bool Compute(const std::vector<const double*>& parameter_blocks,
392 Problem* problem);
393
394 // Return the block of the cross-covariance matrix corresponding to
395 // parameter_block1 and parameter_block2.
396 //
397 // Compute must be called before the first call to
398 // GetCovarianceBlock and the pair <parameter_block1,
399 // parameter_block2> OR the pair <parameter_block2,
400 // parameter_block1> must have been present in the vector
401 // covariance_blocks when Compute was called. Otherwise
402 // GetCovarianceBlock will return false.
403 //
404 // covariance_block must point to a memory location that can store a
405 // parameter_block1_size x parameter_block2_size matrix. The
406 // returned covariance will be a row-major matrix.
407 bool GetCovarianceBlock(const double* parameter_block1,
408 const double* parameter_block2,
409 double* covariance_block) const;
410
Austin Schuh3de38b02024-06-25 18:25:10 -0700411 // Returns the block of the cross-covariance in the tangent space if a
412 // manifold is associated with either parameter block; else returns
413 // cross-covariance in the ambient space.
Austin Schuh70cc9552019-01-21 19:46:48 -0800414 //
415 // Compute must be called before the first call to
416 // GetCovarianceBlock and the pair <parameter_block1,
417 // parameter_block2> OR the pair <parameter_block2,
418 // parameter_block1> must have been present in the vector
419 // covariance_blocks when Compute was called. Otherwise
420 // GetCovarianceBlock will return false.
421 //
422 // covariance_block must point to a memory location that can store a
423 // parameter_block1_local_size x parameter_block2_local_size matrix. The
424 // returned covariance will be a row-major matrix.
425 bool GetCovarianceBlockInTangentSpace(const double* parameter_block1,
426 const double* parameter_block2,
427 double* covariance_block) const;
428
429 // Return the covariance matrix corresponding to all parameter_blocks.
430 //
431 // Compute must be called before calling GetCovarianceMatrix and all
432 // parameter_blocks must have been present in the vector
433 // parameter_blocks when Compute was called. Otherwise
434 // GetCovarianceMatrix returns false.
435 //
436 // covariance_matrix must point to a memory location that can store
437 // the size of the covariance matrix. The covariance matrix will be
438 // a square matrix whose row and column count is equal to the sum of
439 // the sizes of the individual parameter blocks. The covariance
440 // matrix will be a row-major matrix.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800441 bool GetCovarianceMatrix(const std::vector<const double*>& parameter_blocks,
442 double* covariance_matrix) const;
Austin Schuh70cc9552019-01-21 19:46:48 -0800443
444 // Return the covariance matrix corresponding to parameter_blocks
Austin Schuh3de38b02024-06-25 18:25:10 -0700445 // in the tangent space if a manifold is associated with one of the parameter
446 // blocks else returns the covariance matrix in the ambient space.
Austin Schuh70cc9552019-01-21 19:46:48 -0800447 //
448 // Compute must be called before calling GetCovarianceMatrix and all
449 // parameter_blocks must have been present in the vector
450 // parameters_blocks when Compute was called. Otherwise
451 // GetCovarianceMatrix returns false.
452 //
453 // covariance_matrix must point to a memory location that can store
454 // the size of the covariance matrix. The covariance matrix will be
455 // a square matrix whose row and column count is equal to the sum of
456 // the sizes of the tangent spaces of the individual parameter
457 // blocks. The covariance matrix will be a row-major matrix.
458 bool GetCovarianceMatrixInTangentSpace(
459 const std::vector<const double*>& parameter_blocks,
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800460 double* covariance_matrix) const;
Austin Schuh70cc9552019-01-21 19:46:48 -0800461
462 private:
463 std::unique_ptr<internal::CovarianceImpl> impl_;
464};
465
466} // namespace ceres
467
468#include "ceres/internal/reenable_warnings.h"
469
470#endif // CERES_PUBLIC_COVARIANCE_H_