GTSAM  4.0.2
C++ library for smoothing and mapping (SAM)
Matrix.h
Go to the documentation of this file.
1 /* ----------------------------------------------------------------------------
2 
3  * GTSAM Copyright 2010, Georgia Tech Research Corporation,
4  * Atlanta, Georgia 30332-0415
5  * All Rights Reserved
6  * Authors: Frank Dellaert, et al. (see THANKS for the full author list)
7 
8  * See LICENSE for the license information
9 
10  * -------------------------------------------------------------------------- */
11 
23 // \callgraph
24 
25 #pragma once
26 
28 #include <gtsam/base/Vector.h>
29 
30 #include <vector>
31 
37 namespace gtsam {
38 
39 typedef Eigen::MatrixXd Matrix;
40 typedef Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> MatrixRowMajor;
41 
42 // Create handy typedefs and constants for square-size matrices
43 // MatrixMN, MatrixN = MatrixNN, I_NxN, and Z_NxN, for M,N=1..9
44 #define GTSAM_MAKE_MATRIX_DEFS(N) \
45 using Matrix##N = Eigen::Matrix<double, N, N>; \
46 using Matrix1##N = Eigen::Matrix<double, 1, N>; \
47 using Matrix2##N = Eigen::Matrix<double, 2, N>; \
48 using Matrix3##N = Eigen::Matrix<double, 3, N>; \
49 using Matrix4##N = Eigen::Matrix<double, 4, N>; \
50 using Matrix5##N = Eigen::Matrix<double, 5, N>; \
51 using Matrix6##N = Eigen::Matrix<double, 6, N>; \
52 using Matrix7##N = Eigen::Matrix<double, 7, N>; \
53 using Matrix8##N = Eigen::Matrix<double, 8, N>; \
54 using Matrix9##N = Eigen::Matrix<double, 9, N>; \
55 static const Eigen::MatrixBase<Matrix##N>::IdentityReturnType I_##N##x##N = Matrix##N::Identity(); \
56 static const Eigen::MatrixBase<Matrix##N>::ConstantReturnType Z_##N##x##N = Matrix##N::Zero();
57 
58 GTSAM_MAKE_MATRIX_DEFS(1)
59 GTSAM_MAKE_MATRIX_DEFS(2)
60 GTSAM_MAKE_MATRIX_DEFS(3)
61 GTSAM_MAKE_MATRIX_DEFS(4)
62 GTSAM_MAKE_MATRIX_DEFS(5)
63 GTSAM_MAKE_MATRIX_DEFS(6)
64 GTSAM_MAKE_MATRIX_DEFS(7)
65 GTSAM_MAKE_MATRIX_DEFS(8)
66 GTSAM_MAKE_MATRIX_DEFS(9)
67 
68 // Matrix expressions for accessing parts of matrices
69 typedef Eigen::Block<Matrix> SubMatrix;
70 typedef Eigen::Block<const Matrix> ConstSubMatrix;
71 
72 // Matrix formatting arguments when printing.
73 // Akin to Matlab style.
74 const Eigen::IOFormat& matlabFormat();
75 
79 template <class MATRIX>
80 bool equal_with_abs_tol(const Eigen::DenseBase<MATRIX>& A, const Eigen::DenseBase<MATRIX>& B, double tol = 1e-9) {
81 
82  const size_t n1 = A.cols(), m1 = A.rows();
83  const size_t n2 = B.cols(), m2 = B.rows();
84 
85  if(m1!=m2 || n1!=n2) return false;
86 
87  for(size_t i=0; i<m1; i++)
88  for(size_t j=0; j<n1; j++) {
89  if(!fpEqual(A(i,j), B(i,j), tol, false)) {
90  return false;
91  }
92  }
93  return true;
94 }
95 
99 inline bool operator==(const Matrix& A, const Matrix& B) {
100  return equal_with_abs_tol(A,B,1e-9);
101 }
102 
106 inline bool operator!=(const Matrix& A, const Matrix& B) {
107  return !(A==B);
108  }
109 
113 GTSAM_EXPORT bool assert_equal(const Matrix& A, const Matrix& B, double tol = 1e-9);
114 
118 GTSAM_EXPORT bool assert_inequal(const Matrix& A, const Matrix& B, double tol = 1e-9);
119 
123 GTSAM_EXPORT bool assert_equal(const std::list<Matrix>& As, const std::list<Matrix>& Bs, double tol = 1e-9);
124 
128 GTSAM_EXPORT bool linear_independent(const Matrix& A, const Matrix& B, double tol = 1e-9);
129 
133 GTSAM_EXPORT bool linear_dependent(const Matrix& A, const Matrix& B, double tol = 1e-9);
134 
139 GTSAM_EXPORT Vector operator^(const Matrix& A, const Vector & v);
140 
142 template<class MATRIX>
143 inline MATRIX prod(const MATRIX& A, const MATRIX&B) {
144  MATRIX result = A * B;
145  return result;
146 }
147 
151 GTSAM_EXPORT void print(const Matrix& A, const std::string& s, std::ostream& stream);
152 
156 GTSAM_EXPORT void print(const Matrix& A, const std::string& s = "");
157 
161 GTSAM_EXPORT void save(const Matrix& A, const std::string &s, const std::string& filename);
162 
168 GTSAM_EXPORT std::istream& operator>>(std::istream& inputStream, Matrix& destinationMatrix);
169 
179 template<class MATRIX>
180 Eigen::Block<const MATRIX> sub(const MATRIX& A, size_t i1, size_t i2, size_t j1, size_t j2) {
181  size_t m=i2-i1, n=j2-j1;
182  return A.block(i1,j1,m,n);
183 }
184 
193 template <typename Derived1, typename Derived2>
194 void insertSub(Eigen::MatrixBase<Derived1>& fullMatrix, const Eigen::MatrixBase<Derived2>& subMatrix, size_t i, size_t j) {
195  fullMatrix.block(i, j, subMatrix.rows(), subMatrix.cols()) = subMatrix;
196 }
197 
201 GTSAM_EXPORT Matrix diag(const std::vector<Matrix>& Hs);
202 
209 template<class MATRIX>
210 const typename MATRIX::ConstColXpr column(const MATRIX& A, size_t j) {
211  return A.col(j);
212 }
213 
220 template<class MATRIX>
221 const typename MATRIX::ConstRowXpr row(const MATRIX& A, size_t j) {
222  return A.row(j);
223 }
224 
230 template<class MATRIX>
231 void zeroBelowDiagonal(MATRIX& A, size_t cols=0) {
232  const size_t m = A.rows(), n = A.cols();
233  const size_t k = (cols) ? std::min(cols, std::min(m,n)) : std::min(m,n);
234  for (size_t j=0; j<k; ++j)
235  A.col(j).segment(j+1, m-(j+1)).setZero();
236 }
237 
241 inline Matrix trans(const Matrix& A) { return A.transpose(); }
242 
244 template <int OutM, int OutN, int OutOptions, int InM, int InN, int InOptions>
245 struct Reshape {
246  //TODO replace this with Eigen's reshape function as soon as available. (There is a PR already pending : https://bitbucket.org/eigen/eigen/pull-request/41/reshape/diff)
247  typedef Eigen::Map<const Eigen::Matrix<double, OutM, OutN, OutOptions> > ReshapedType;
248  static inline ReshapedType reshape(const Eigen::Matrix<double, InM, InN, InOptions> & in) {
249  return in.data();
250  }
251 };
252 
254 template <int M, int InOptions>
255 struct Reshape<M, M, InOptions, M, M, InOptions> {
256  typedef const Eigen::Matrix<double, M, M, InOptions> & ReshapedType;
257  static inline ReshapedType reshape(const Eigen::Matrix<double, M, M, InOptions> & in) {
258  return in;
259  }
260 };
261 
263 template <int M, int N, int InOptions>
264 struct Reshape<M, N, InOptions, M, N, InOptions> {
265  typedef const Eigen::Matrix<double, M, N, InOptions> & ReshapedType;
266  static inline ReshapedType reshape(const Eigen::Matrix<double, M, N, InOptions> & in) {
267  return in;
268  }
269 };
270 
272 template <int M, int N, int InOptions>
273 struct Reshape<N, M, InOptions, M, N, InOptions> {
274  typedef typename Eigen::Matrix<double, M, N, InOptions>::ConstTransposeReturnType ReshapedType;
275  static inline ReshapedType reshape(const Eigen::Matrix<double, M, N, InOptions> & in) {
276  return in.transpose();
277  }
278 };
279 
280 template <int OutM, int OutN, int OutOptions, int InM, int InN, int InOptions>
281 inline typename Reshape<OutM, OutN, OutOptions, InM, InN, InOptions>::ReshapedType reshape(const Eigen::Matrix<double, InM, InN, InOptions> & m){
282  static_assert(InM * InN == OutM * OutN);
284 }
285 
292 GTSAM_EXPORT std::pair<Matrix,Matrix> qr(const Matrix& A);
293 
299 GTSAM_EXPORT void inplace_QR(Matrix& A);
300 
309 GTSAM_EXPORT std::list<std::tuple<Vector, double, double> >
310 weighted_eliminate(Matrix& A, Vector& b, const Vector& sigmas);
311 
319 GTSAM_EXPORT void householder_(Matrix& A, size_t k, bool copy_vectors=true);
320 
327 GTSAM_EXPORT void householder(Matrix& A, size_t k);
328 
336 GTSAM_EXPORT Vector backSubstituteUpper(const Matrix& U, const Vector& b, bool unit=false);
337 
345 //TODO: is this function necessary? it isn't used
346 GTSAM_EXPORT Vector backSubstituteUpper(const Vector& b, const Matrix& U, bool unit=false);
347 
355 GTSAM_EXPORT Vector backSubstituteLower(const Matrix& L, const Vector& b, bool unit=false);
356 
363 GTSAM_EXPORT Matrix stack(size_t nrMatrices, ...);
364 GTSAM_EXPORT Matrix stack(const std::vector<Matrix>& blocks);
365 
376 GTSAM_EXPORT Matrix collect(const std::vector<const Matrix *>& matrices, size_t m = 0, size_t n = 0);
377 GTSAM_EXPORT Matrix collect(size_t nrMatrices, ...);
378 
385 GTSAM_EXPORT void vector_scale_inplace(const Vector& v, Matrix& A, bool inf_mask = false); // row
386 GTSAM_EXPORT Matrix vector_scale(const Vector& v, const Matrix& A, bool inf_mask = false); // row
387 GTSAM_EXPORT Matrix vector_scale(const Matrix& A, const Vector& v, bool inf_mask = false); // column
388 
400 inline Matrix3 skewSymmetric(double wx, double wy, double wz) {
401  return (Matrix3() << 0.0, -wz, +wy, +wz, 0.0, -wx, -wy, +wx, 0.0).finished();
402 }
403 
404 template <class Derived>
405 inline Matrix3 skewSymmetric(const Eigen::MatrixBase<Derived>& w) {
406  return skewSymmetric(w(0), w(1), w(2));
407 }
408 
410 GTSAM_EXPORT Matrix inverse_square_root(const Matrix& A);
411 
413 GTSAM_EXPORT Matrix cholesky_inverse(const Matrix &A);
414 
427 GTSAM_EXPORT void svd(const Matrix& A, Matrix& U, Vector& S, Matrix& V);
428 
436 GTSAM_EXPORT std::tuple<int, double, Vector>
437 DLT(const Matrix& A, double rank_tol = 1e-9);
438 
444 GTSAM_EXPORT Matrix expm(const Matrix& A, size_t K=7);
445 
446 std::string formatMatrixIndented(const std::string& label, const Matrix& matrix, bool makeVectorHorizontal = false);
447 
454 template <int N>
456  typedef Eigen::Matrix<double, N, 1> VectorN;
457  typedef Eigen::Matrix<double, N, N> MatrixN;
458 
460  VectorN operator()(const MatrixN& A, const VectorN& b,
462  OptionalJacobian<N, N> H2 = {}) const {
463  const MatrixN invA = A.inverse();
464  const VectorN c = invA * b;
465  // The derivative in A is just -[c[0]*invA c[1]*invA ... c[N-1]*invA]
466  if (H1)
467  for (size_t j = 0; j < N; j++)
468  H1->template middleCols<N>(N * j) = -c[j] * invA;
469  // The derivative in b is easy, as invA*b is just a linear map:
470  if (H2) *H2 = invA;
471  return c;
472  }
473 };
474 
480 template <typename T, int N>
482  enum { M = traits<T>::dimension };
483  typedef Eigen::Matrix<double, N, 1> VectorN;
484  typedef Eigen::Matrix<double, N, N> MatrixN;
485 
486  // The function phi should calculate f(a)*b, with derivatives in a and b.
487  // Naturally, the derivative in b is f(a).
488  typedef std::function<VectorN(
489  const T&, const VectorN&, OptionalJacobian<N, M>, OptionalJacobian<N, N>)>
490  Operator;
491 
493  MultiplyWithInverseFunction(const Operator& phi) : phi_(phi) {}
494 
496  VectorN operator()(const T& a, const VectorN& b,
497  OptionalJacobian<N, M> H1 = {},
498  OptionalJacobian<N, N> H2 = {}) const {
499  MatrixN A;
500  phi_(a, b, {}, A); // get A = f(a) by calling f once
501  const MatrixN invA = A.inverse();
502  const VectorN c = invA * b;
503 
504  if (H1) {
505  Eigen::Matrix<double, N, M> H;
506  phi_(a, c, H, {}); // get derivative H of forward mapping
507  *H1 = -invA* H;
508  }
509  if (H2) *H2 = invA;
510  return c;
511  }
512 
513  private:
514  const Operator phi_;
515 };
516 
517 GTSAM_EXPORT Matrix LLt(const Matrix& A);
518 
519 GTSAM_EXPORT Matrix RtR(const Matrix& A);
520 
521 GTSAM_EXPORT Vector columnNormSquare(const Matrix &A);
522 } // namespace gtsam
GTSAM_EXPORT void inplace_QR(Matrix &A)
void zeroBelowDiagonal(MATRIX &A, size_t cols=0)
Definition: Matrix.h:231
VectorN operator()(const T &a, const VectorN &b, OptionalJacobian< N, M > H1={}, OptionalJacobian< N, N > H2={}) const
f(a).inverse() * b, with optional derivatives
Definition: Matrix.h:496
VectorN operator()(const MatrixN &A, const VectorN &b, OptionalJacobian< N, N *N > H1={}, OptionalJacobian< N, N > H2={}) const
A.inverse() * b, with optional derivatives.
Definition: Matrix.h:460
const MATRIX::ConstRowXpr row(const MATRIX &A, size_t j)
Definition: Matrix.h:221
MATRIX prod(const MATRIX &A, const MATRIX &B)
Definition: Matrix.h:143
GTSAM_EXPORT bool assert_inequal(const Matrix &A, const Matrix &B, double tol=1e-9)
bool operator!=(const Matrix &A, const Matrix &B)
Definition: Matrix.h:106
GTSAM_EXPORT Vector backSubstituteLower(const Matrix &L, const Vector &b, bool unit=false)
GTSAM_EXPORT std::pair< Matrix, Matrix > qr(const Matrix &A)
GTSAM_EXPORT Matrix diag(const std::vector< Matrix > &Hs)
GTSAM_EXPORT Matrix collect(const std::vector< const Matrix *> &matrices, size_t m=0, size_t n=0)
GTSAM_EXPORT void householder(Matrix &A, size_t k)
GTSAM_EXPORT bool linear_dependent(const Matrix &A, const Matrix &B, double tol=1e-9)
Definition: Group.h:43
GTSAM_EXPORT bool fpEqual(double a, double b, double tol, bool check_relative_also=true)
bool operator==(const Matrix &A, const Matrix &B)
Definition: Matrix.h:99
Definition: Matrix.h:455
GTSAM_EXPORT std::tuple< int, double, Vector > DLT(const Matrix &A, double rank_tol=1e-9)
T expm(const Vector &x, int K=7)
Definition: Lie.h:317
GTSAM_EXPORT void householder_(Matrix &A, size_t k, bool copy_vectors=true)
Matrix trans(const Matrix &A)
Definition: Matrix.h:241
Reshape functor.
Definition: Matrix.h:245
MultiplyWithInverseFunction(const Operator &phi)
Construct with function as explained above.
Definition: Matrix.h:493
GTSAM_EXPORT void print(const Matrix &A, const std::string &s, std::ostream &stream)
GTSAM_EXPORT void save(const Matrix &A, const std::string &s, const std::string &filename)
GTSAM_EXPORT Vector operator^(const Matrix &A, const Vector &v)
GTSAM_EXPORT std::list< std::tuple< Vector, double, double > > weighted_eliminate(Matrix &A, Vector &b, const Vector &sigmas)
bool equal_with_abs_tol(const Eigen::DenseBase< MATRIX > &A, const Eigen::DenseBase< MATRIX > &B, double tol=1e-9)
Definition: Matrix.h:80
Definition: chartTesting.h:28
typedef and functions to augment Eigen&#39;s VectorXd
GTSAM_EXPORT bool assert_equal(const Matrix &A, const Matrix &B, double tol=1e-9)
Matrix3 skewSymmetric(double wx, double wy, double wz)
Definition: Matrix.h:400
GTSAM_EXPORT Matrix inverse_square_root(const Matrix &A)
Definition: OptionalJacobian.h:38
void insertSub(Eigen::MatrixBase< Derived1 > &fullMatrix, const Eigen::MatrixBase< Derived2 > &subMatrix, size_t i, size_t j)
Definition: Matrix.h:194
GTSAM_EXPORT std::istream & operator>>(std::istream &inputStream, Matrix &destinationMatrix)
GTSAM_EXPORT void svd(const Matrix &A, Matrix &U, Vector &S, Matrix &V)
GTSAM_EXPORT Matrix stack(size_t nrMatrices,...)
GTSAM_EXPORT Vector backSubstituteUpper(const Matrix &U, const Vector &b, bool unit=false)
Eigen::Block< const MATRIX > sub(const MATRIX &A, size_t i1, size_t i2, size_t j1, size_t j2)
Definition: Matrix.h:180
Special class for optional Jacobian arguments.
Definition: Matrix.h:481
GTSAM_EXPORT Matrix cholesky_inverse(const Matrix &A)
const MATRIX::ConstColXpr column(const MATRIX &A, size_t j)
Definition: Matrix.h:210
GTSAM_EXPORT bool linear_independent(const Matrix &A, const Matrix &B, double tol=1e-9)
GTSAM_EXPORT void vector_scale_inplace(const Vector &v, Matrix &A, bool inf_mask=false)