gtsam  4.0.0
gtsam
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 
22 // \callgraph
23 
24 #pragma once
26 #include <gtsam/base/Vector.h>
27 #include <gtsam/config.h>
28 #ifdef GTSAM_ALLOW_DEPRECATED_SINCE_V4
29 #include <Eigen/Core>
30 #include <Eigen/Cholesky>
31 #include <Eigen/LU>
32 #endif
33 #include <boost/format.hpp>
34 #include <boost/function.hpp>
35 #include <boost/tuple/tuple.hpp>
36 #include <boost/math/special_functions/fpclassify.hpp>
37 
43 namespace gtsam {
44 
45 typedef Eigen::MatrixXd Matrix;
46 typedef Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> MatrixRowMajor;
47 
48 // Create handy typedefs and constants for square-size matrices
49 // MatrixMN, MatrixN = MatrixNN, I_NxN, and Z_NxN, for M,N=1..9
50 #define GTSAM_MAKE_MATRIX_DEFS(N) \
51 typedef Eigen::Matrix<double, N, N> Matrix##N; \
52 typedef Eigen::Matrix<double, 1, N> Matrix1##N; \
53 typedef Eigen::Matrix<double, 2, N> Matrix2##N; \
54 typedef Eigen::Matrix<double, 3, N> Matrix3##N; \
55 typedef Eigen::Matrix<double, 4, N> Matrix4##N; \
56 typedef Eigen::Matrix<double, 5, N> Matrix5##N; \
57 typedef Eigen::Matrix<double, 6, N> Matrix6##N; \
58 typedef Eigen::Matrix<double, 7, N> Matrix7##N; \
59 typedef Eigen::Matrix<double, 8, N> Matrix8##N; \
60 typedef Eigen::Matrix<double, 9, N> Matrix9##N; \
61 static const Eigen::MatrixBase<Matrix##N>::IdentityReturnType I_##N##x##N = Matrix##N::Identity(); \
62 static const Eigen::MatrixBase<Matrix##N>::ConstantReturnType Z_##N##x##N = Matrix##N::Zero();
63 
64 GTSAM_MAKE_MATRIX_DEFS(1);
65 GTSAM_MAKE_MATRIX_DEFS(2);
66 GTSAM_MAKE_MATRIX_DEFS(3);
67 GTSAM_MAKE_MATRIX_DEFS(4);
68 GTSAM_MAKE_MATRIX_DEFS(5);
69 GTSAM_MAKE_MATRIX_DEFS(6);
70 GTSAM_MAKE_MATRIX_DEFS(7);
71 GTSAM_MAKE_MATRIX_DEFS(8);
72 GTSAM_MAKE_MATRIX_DEFS(9);
73 
74 // Matrix expressions for accessing parts of matrices
75 typedef Eigen::Block<Matrix> SubMatrix;
76 typedef Eigen::Block<const Matrix> ConstSubMatrix;
77 
81 template <class MATRIX>
82 bool equal_with_abs_tol(const Eigen::DenseBase<MATRIX>& A, const Eigen::DenseBase<MATRIX>& B, double tol = 1e-9) {
83 
84  const size_t n1 = A.cols(), m1 = A.rows();
85  const size_t n2 = B.cols(), m2 = B.rows();
86 
87  if(m1!=m2 || n1!=n2) return false;
88 
89  for(size_t i=0; i<m1; i++)
90  for(size_t j=0; j<n1; j++) {
91  if(boost::math::isnan(A(i,j)) ^ boost::math::isnan(B(i,j)))
92  return false;
93  else if(fabs(A(i,j) - B(i,j)) > tol)
94  return false;
95  }
96  return true;
97 }
98 
102 inline bool operator==(const Matrix& A, const Matrix& B) {
103  return equal_with_abs_tol(A,B,1e-9);
104 }
105 
109 inline bool operator!=(const Matrix& A, const Matrix& B) {
110  return !(A==B);
111  }
112 
116 GTSAM_EXPORT bool assert_equal(const Matrix& A, const Matrix& B, double tol = 1e-9);
117 
121 GTSAM_EXPORT bool assert_inequal(const Matrix& A, const Matrix& B, double tol = 1e-9);
122 
126 GTSAM_EXPORT bool assert_equal(const std::list<Matrix>& As, const std::list<Matrix>& Bs, double tol = 1e-9);
127 
131 GTSAM_EXPORT bool linear_independent(const Matrix& A, const Matrix& B, double tol = 1e-9);
132 
136 GTSAM_EXPORT bool linear_dependent(const Matrix& A, const Matrix& B, double tol = 1e-9);
137 
142 GTSAM_EXPORT Vector operator^(const Matrix& A, const Vector & v);
143 
145 template<class MATRIX>
146 inline MATRIX prod(const MATRIX& A, const MATRIX&B) {
147  MATRIX result = A * B;
148  return result;
149 }
150 
154 GTSAM_EXPORT void print(const Matrix& A, const std::string& s, std::ostream& stream);
155 
159 GTSAM_EXPORT void print(const Matrix& A, const std::string& s = "");
160 
164 GTSAM_EXPORT void save(const Matrix& A, const std::string &s, const std::string& filename);
165 
171 GTSAM_EXPORT std::istream& operator>>(std::istream& inputStream, Matrix& destinationMatrix);
172 
182 template<class MATRIX>
183 Eigen::Block<const MATRIX> sub(const MATRIX& A, size_t i1, size_t i2, size_t j1, size_t j2) {
184  size_t m=i2-i1, n=j2-j1;
185  return A.block(i1,j1,m,n);
186 }
187 
196 template <typename Derived1, typename Derived2>
197 void insertSub(Eigen::MatrixBase<Derived1>& fullMatrix, const Eigen::MatrixBase<Derived2>& subMatrix, size_t i, size_t j) {
198  fullMatrix.block(i, j, subMatrix.rows(), subMatrix.cols()) = subMatrix;
199 }
200 
204 GTSAM_EXPORT Matrix diag(const std::vector<Matrix>& Hs);
205 
212 template<class MATRIX>
213 const typename MATRIX::ConstColXpr column(const MATRIX& A, size_t j) {
214  return A.col(j);
215 }
216 
223 template<class MATRIX>
224 const typename MATRIX::ConstRowXpr row(const MATRIX& A, size_t j) {
225  return A.row(j);
226 }
227 
233 template<class MATRIX>
234 void zeroBelowDiagonal(MATRIX& A, size_t cols=0) {
235  const size_t m = A.rows(), n = A.cols();
236  const size_t k = (cols) ? std::min(cols, std::min(m,n)) : std::min(m,n);
237  for (size_t j=0; j<k; ++j)
238  A.col(j).segment(j+1, m-(j+1)).setZero();
239 }
240 
244 inline Matrix trans(const Matrix& A) { return A.transpose(); }
245 
247 template <int OutM, int OutN, int OutOptions, int InM, int InN, int InOptions>
248 struct Reshape {
249  //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)
250  typedef Eigen::Map<const Eigen::Matrix<double, OutM, OutN, OutOptions> > ReshapedType;
251  static inline ReshapedType reshape(const Eigen::Matrix<double, InM, InN, InOptions> & in) {
252  return in.data();
253  }
254 };
255 
257 template <int M, int InOptions>
258 struct Reshape<M, M, InOptions, M, M, InOptions> {
259  typedef const Eigen::Matrix<double, M, M, InOptions> & ReshapedType;
260  static inline ReshapedType reshape(const Eigen::Matrix<double, M, M, InOptions> & in) {
261  return in;
262  }
263 };
264 
266 template <int M, int N, int InOptions>
267 struct Reshape<M, N, InOptions, M, N, InOptions> {
268  typedef const Eigen::Matrix<double, M, N, InOptions> & ReshapedType;
269  static inline ReshapedType reshape(const Eigen::Matrix<double, M, N, InOptions> & in) {
270  return in;
271  }
272 };
273 
275 template <int M, int N, int InOptions>
276 struct Reshape<N, M, InOptions, M, N, InOptions> {
277  typedef typename Eigen::Matrix<double, M, N, InOptions>::ConstTransposeReturnType ReshapedType;
278  static inline ReshapedType reshape(const Eigen::Matrix<double, M, N, InOptions> & in) {
279  return in.transpose();
280  }
281 };
282 
283 template <int OutM, int OutN, int OutOptions, int InM, int InN, int InOptions>
284 inline typename Reshape<OutM, OutN, OutOptions, InM, InN, InOptions>::ReshapedType reshape(const Eigen::Matrix<double, InM, InN, InOptions> & m){
285  BOOST_STATIC_ASSERT(InM * InN == OutM * OutN);
287 }
288 
295 GTSAM_EXPORT std::pair<Matrix,Matrix> qr(const Matrix& A);
296 
302 void inplace_QR(Matrix& A);
303 
312 GTSAM_EXPORT std::list<boost::tuple<Vector, double, double> >
313 weighted_eliminate(Matrix& A, Vector& b, const Vector& sigmas);
314 
322 GTSAM_EXPORT void householder_(Matrix& A, size_t k, bool copy_vectors=true);
323 
330 GTSAM_EXPORT void householder(Matrix& A, size_t k);
331 
339 GTSAM_EXPORT Vector backSubstituteUpper(const Matrix& U, const Vector& b, bool unit=false);
340 
348 //TODO: is this function necessary? it isn't used
349 GTSAM_EXPORT Vector backSubstituteUpper(const Vector& b, const Matrix& U, bool unit=false);
350 
358 GTSAM_EXPORT Vector backSubstituteLower(const Matrix& L, const Vector& b, bool unit=false);
359 
366 GTSAM_EXPORT Matrix stack(size_t nrMatrices, ...);
367 GTSAM_EXPORT Matrix stack(const std::vector<Matrix>& blocks);
368 
379 GTSAM_EXPORT Matrix collect(const std::vector<const Matrix *>& matrices, size_t m = 0, size_t n = 0);
380 GTSAM_EXPORT Matrix collect(size_t nrMatrices, ...);
381 
388 GTSAM_EXPORT void vector_scale_inplace(const Vector& v, Matrix& A, bool inf_mask = false); // row
389 GTSAM_EXPORT Matrix vector_scale(const Vector& v, const Matrix& A, bool inf_mask = false); // row
390 GTSAM_EXPORT Matrix vector_scale(const Matrix& A, const Vector& v, bool inf_mask = false); // column
391 
403 inline Matrix3 skewSymmetric(double wx, double wy, double wz) {
404  return (Matrix3() << 0.0, -wz, +wy, +wz, 0.0, -wx, -wy, +wx, 0.0).finished();
405 }
406 
407 template <class Derived>
408 inline Matrix3 skewSymmetric(const Eigen::MatrixBase<Derived>& w) {
409  return skewSymmetric(w(0), w(1), w(2));
410 }
411 
413 GTSAM_EXPORT Matrix inverse_square_root(const Matrix& A);
414 
416 GTSAM_EXPORT Matrix cholesky_inverse(const Matrix &A);
417 
430 GTSAM_EXPORT void svd(const Matrix& A, Matrix& U, Vector& S, Matrix& V);
431 
439 GTSAM_EXPORT boost::tuple<int, double, Vector>
440 DLT(const Matrix& A, double rank_tol = 1e-9);
441 
447 GTSAM_EXPORT Matrix expm(const Matrix& A, size_t K=7);
448 
449 std::string formatMatrixIndented(const std::string& label, const Matrix& matrix, bool makeVectorHorizontal = false);
450 
457 template <int N>
459  typedef Eigen::Matrix<double, N, 1> VectorN;
460  typedef Eigen::Matrix<double, N, N> MatrixN;
461 
463  VectorN operator()(const MatrixN& A, const VectorN& b,
464  OptionalJacobian<N, N* N> H1 = boost::none,
465  OptionalJacobian<N, N> H2 = boost::none) const {
466  const MatrixN invA = A.inverse();
467  const VectorN c = invA * b;
468  // The derivative in A is just -[c[0]*invA c[1]*invA ... c[N-1]*invA]
469  if (H1)
470  for (size_t j = 0; j < N; j++)
471  H1->template middleCols<N>(N * j) = -c[j] * invA;
472  // The derivative in b is easy, as invA*b is just a linear map:
473  if (H2) *H2 = invA;
474  return c;
475  }
476 };
477 
483 template <typename T, int N>
485  enum { M = traits<T>::dimension };
486  typedef Eigen::Matrix<double, N, 1> VectorN;
487  typedef Eigen::Matrix<double, N, N> MatrixN;
488 
489  // The function phi should calculate f(a)*b, with derivatives in a and b.
490  // Naturally, the derivative in b is f(a).
491  typedef boost::function<VectorN(
492  const T&, const VectorN&, OptionalJacobian<N, M>, OptionalJacobian<N, N>)>
493  Operator;
494 
496  MultiplyWithInverseFunction(const Operator& phi) : phi_(phi) {}
497 
499  VectorN operator()(const T& a, const VectorN& b,
500  OptionalJacobian<N, M> H1 = boost::none,
501  OptionalJacobian<N, N> H2 = boost::none) const {
502  MatrixN A;
503  phi_(a, b, boost::none, A); // get A = f(a) by calling f once
504  const MatrixN invA = A.inverse();
505  const VectorN c = invA * b;
506 
507  if (H1) {
508  Eigen::Matrix<double, N, M> H;
509  phi_(a, c, H, boost::none); // get derivative H of forward mapping
510  *H1 = -invA* H;
511  }
512  if (H2) *H2 = invA;
513  return c;
514  }
515 
516  private:
517  const Operator phi_;
518 };
519 
520 #ifdef GTSAM_ALLOW_DEPRECATED_SINCE_V4
521 inline Matrix zeros( size_t m, size_t n ) { return Matrix::Zero(m,n); }
522 inline Matrix ones( size_t m, size_t n ) { return Matrix::Ones(m,n); }
523 inline Matrix eye( size_t m, size_t n) { return Matrix::Identity(m, n); }
524 inline Matrix eye( size_t m ) { return eye(m,m); }
525 inline Matrix diag(const Vector& v) { return v.asDiagonal(); }
526 inline void multiplyAdd(double alpha, const Matrix& A, const Vector& x, Vector& e) { e += alpha * A * x; }
527 inline void multiplyAdd(const Matrix& A, const Vector& x, Vector& e) { e += A * x; }
528 inline void transposeMultiplyAdd(double alpha, const Matrix& A, const Vector& e, Vector& x) { x += alpha * A.transpose() * e; }
529 inline void transposeMultiplyAdd(const Matrix& A, const Vector& e, Vector& x) { x += A.transpose() * e; }
530 inline void transposeMultiplyAdd(double alpha, const Matrix& A, const Vector& e, SubVector x) { x += alpha * A.transpose() * e; }
531 inline void insertColumn(Matrix& A, const Vector& col, size_t j) { A.col(j) = col; }
532 inline void insertColumn(Matrix& A, const Vector& col, size_t i, size_t j) { A.col(j).segment(i, col.size()) = col; }
533 inline void solve(Matrix& A, Matrix& B) { B = A.fullPivLu().solve(B); }
534 inline Matrix inverse(const Matrix& A) { return A.inverse(); }
535 #endif
536 
537 GTSAM_EXPORT Matrix LLt(const Matrix& A);
538 
539 GTSAM_EXPORT Matrix RtR(const Matrix& A);
540 
541 GTSAM_EXPORT Vector columnNormSquare(const Matrix &A);
542 } // namespace gtsam
543 
544 #include <boost/serialization/nvp.hpp>
545 #include <boost/serialization/array.hpp>
546 #include <boost/serialization/split_free.hpp>
547 
548 namespace boost {
549  namespace serialization {
550 
551  // split version - sends sizes ahead
552  template<class Archive>
553  void save(Archive & ar, const gtsam::Matrix & m, unsigned int /*version*/) {
554  const size_t rows = m.rows(), cols = m.cols();
555  ar << BOOST_SERIALIZATION_NVP(rows);
556  ar << BOOST_SERIALIZATION_NVP(cols);
557  ar << make_nvp("data", make_array(m.data(), m.size()));
558  }
559 
560  template<class Archive>
561  void load(Archive & ar, gtsam::Matrix & m, unsigned int /*version*/) {
562  size_t rows, cols;
563  ar >> BOOST_SERIALIZATION_NVP(rows);
564  ar >> BOOST_SERIALIZATION_NVP(cols);
565  m.resize(rows, cols);
566  ar >> make_nvp("data", make_array(m.data(), m.size()));
567  }
568 
569  } // namespace serialization
570 } // namespace boost
571 
572 BOOST_SERIALIZATION_SPLIT_FREE(gtsam::Matrix);
573 
Matrix3 skewSymmetric(double wx, double wy, double wz)
skew symmetric matrix returns this: 0 -wz wy wz 0 -wx -wy wx 0
Definition: Matrix.h:403
const MATRIX::ConstColXpr column(const MATRIX &A, size_t j)
Extracts a column view from a matrix that avoids a copy.
Definition: Matrix.h:213
Vector operator^(const Matrix &A, const Vector &v)
overload ^ for trans(A)*v We transpose the vectors for speed.
Definition: Matrix.cpp:130
Vector backSubstituteUpper(const Matrix &U, const Vector &b, bool unit)
backSubstitute U*x=b
Definition: Matrix.cpp:372
Matrix diag(const std::vector< Matrix > &Hs)
Create a matrix with submatrices along its diagonal.
Definition: Matrix.cpp:202
void vector_scale_inplace(const Vector &v, Matrix &A, bool inf_mask)
scales a matrix row or column by the values in a vector Arguments (Matrix, Vector) scales the columns...
Definition: Matrix.cpp:477
void print(const Matrix &A, const string &s, ostream &stream)
print without optional string, must specify cout yourself
Definition: Matrix.cpp:141
bool linear_dependent(const Matrix &A, const Matrix &B, double tol)
check whether the rows of two matrices are linear dependent
Definition: Matrix.cpp:116
bool equal_with_abs_tol(const Eigen::DenseBase< MATRIX > &A, const Eigen::DenseBase< MATRIX > &B, double tol=1e-9)
equals with a tolerance
Definition: Matrix.h:82
bool assert_equal(const Matrix &expected, const Matrix &actual, double tol)
equals with an tolerance, prints out message if unequal
Definition: Matrix.cpp:42
void householder_(Matrix &A, size_t k, bool copy_vectors)
Imperative version of Householder QR factorization, Golub & Van Loan p 224 version with Householder v...
Definition: Matrix.cpp:322
void insertSub(Eigen::MatrixBase< Derived1 > &fullMatrix, const Eigen::MatrixBase< Derived2 > &subMatrix, size_t i, size_t j)
insert a submatrix IN PLACE at a specified location in a larger matrix NOTE: there is no size checkin...
Definition: Matrix.h:197
Matrix cholesky_inverse(const Matrix &A)
Return the inverse of a S.P.D.
Definition: Matrix.cpp:534
Reshape functor.
Definition: Matrix.h:248
void svd(const Matrix &A, Matrix &U, Vector &S, Matrix &V)
SVD computes economy SVD A=U*S*V'.
Definition: Matrix.cpp:555
pair< Matrix, Matrix > qr(const Matrix &A)
Householder QR factorization, Golub & Van Loan p 224, explicit version.
Definition: Matrix.cpp:230
Vector backSubstituteLower(const Matrix &L, const Vector &b, bool unit)
backSubstitute L*x=b
Definition: Matrix.cpp:362
OptionalJacobian is an Eigen::Ref like class that can take be constructed using either a fixed size o...
Definition: OptionalJacobian.h:39
void save(const Matrix &A, const string &s, const string &filename)
save a matrix to file, which can be loaded by matlab
Definition: Matrix.cpp:162
Functor that implements multiplication with the inverse of a matrix, itself the result of a function ...
Definition: Matrix.h:484
boost::tuple< int, double, Vector > DLT(const Matrix &A, double rank_tol)
Direct linear transform algorithm that calls svd to find a vector v that minimizes the algebraic erro...
Definition: Matrix.cpp:563
T inverse(const T &t)
unary functions
Definition: lieProxies.h:43
MultiplyWithInverseFunction(const Operator &phi)
Construct with function as explained above.
Definition: Matrix.h:496
Functor that implements multiplication of a vector b with the inverse of a matrix A.
Definition: Matrix.h:458
bool operator!=(const Matrix &A, const Matrix &B)
inequality
Definition: Matrix.h:109
MATRIX prod(const MATRIX &A, const MATRIX &B)
products using old-style format to improve compatibility
Definition: Matrix.h:146
A manifold defines a space in which there is a notion of a linear tangent space that can be centered ...
Definition: concepts.h:30
void inplace_QR(Matrix &A)
QR factorization using Eigen's internal block QR algorithm.
Definition: Matrix.cpp:631
VectorN operator()(const MatrixN &A, const VectorN &b, OptionalJacobian< N, N *N > H1=boost::none, OptionalJacobian< N, N > H2=boost::none) const
A.inverse() * b, with optional derivatives.
Definition: Matrix.h:463
T expm(const Vector &x, int K=7)
Exponential map given exponential coordinates class T needs a wedge<> function and a constructor from...
Definition: Lie.h:325
Matrix collect(const std::vector< const Matrix * > &matrices, size_t m, size_t n)
create a matrix by concatenating Given a set of matrices: A1, A2, A3...
Definition: Matrix.cpp:438
list< boost::tuple< Vector, double, double > > weighted_eliminate(Matrix &A, Vector &b, const Vector &sigmas)
Imperative algorithm for in-place full elimination with weights and constraint handling.
Definition: Matrix.cpp:268
Matrix stack(size_t nrMatrices,...)
create a matrix by stacking other matrices Given a set of matrices: A1, A2, A3...
Definition: Matrix.cpp:392
Eigen::Block< const MATRIX > sub(const MATRIX &A, size_t i1, size_t i2, size_t j1, size_t j2)
extract submatrix, slice semantics, i.e.
Definition: Matrix.h:183
Global functions in a separate testing namespace.
Definition: chartTesting.h:28
void zeroBelowDiagonal(MATRIX &A, size_t cols=0)
Zeros all of the elements below the diagonal of a matrix, in place.
Definition: Matrix.h:234
bool assert_inequal(const Matrix &A, const Matrix &B, double tol)
inequals with an tolerance, prints out message if within tolerance
Definition: Matrix.cpp:62
const MATRIX::ConstRowXpr row(const MATRIX &A, size_t j)
Extracts a row view from a matrix that avoids a copy.
Definition: Matrix.h:224
VectorN operator()(const T &a, const VectorN &b, OptionalJacobian< N, M > H1=boost::none, OptionalJacobian< N, N > H2=boost::none) const
f(a).inverse() * b, with optional derivatives
Definition: Matrix.h:499
Matrix trans(const Matrix &A)
static transpose function, just calls Eigen transpose member function
Definition: Matrix.h:244
void householder(Matrix &A, size_t k)
Householder tranformation, zeros below diagonal.
Definition: Matrix.cpp:349
Matrix inverse_square_root(const Matrix &A)
Use Cholesky to calculate inverse square root of a matrix.
Definition: Matrix.cpp:547
bool linear_independent(const Matrix &A, const Matrix &B, double tol)
check whether the rows of two matrices are linear independent
Definition: Matrix.cpp:102
bool operator==(const Matrix &A, const Matrix &B)
equality is just equal_with_abs_tol 1e-9
Definition: Matrix.h:102
typedef and functions to augment Eigen's VectorXd
Special class for optional Jacobian arguments.
istream & operator>>(istream &inputStream, Matrix &destinationMatrix)
Read a matrix from an input stream, such as a file.
Definition: Matrix.cpp:169