covariance.h 19 KB

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