base/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 #include <gtsam/config.h>
30 
31 #include <boost/format.hpp>
32 #include <boost/function.hpp>
33 #include <boost/tuple/tuple.hpp>
34 #include <boost/math/special_functions/fpclassify.hpp>
35 
41 namespace gtsam {
42 
43 typedef Eigen::MatrixXd Matrix;
45 
46 // Create handy typedefs and constants for square-size matrices
47 // MatrixMN, MatrixN = MatrixNN, I_NxN, and Z_NxN, for M,N=1..9
48 #define GTSAM_MAKE_MATRIX_DEFS(N) \
49 typedef Eigen::Matrix<double, N, N> Matrix##N; \
50 typedef Eigen::Matrix<double, 1, N> Matrix1##N; \
51 typedef Eigen::Matrix<double, 2, N> Matrix2##N; \
52 typedef Eigen::Matrix<double, 3, N> Matrix3##N; \
53 typedef Eigen::Matrix<double, 4, N> Matrix4##N; \
54 typedef Eigen::Matrix<double, 5, N> Matrix5##N; \
55 typedef Eigen::Matrix<double, 6, N> Matrix6##N; \
56 typedef Eigen::Matrix<double, 7, N> Matrix7##N; \
57 typedef Eigen::Matrix<double, 8, N> Matrix8##N; \
58 typedef Eigen::Matrix<double, 9, N> Matrix9##N; \
59 static const Eigen::MatrixBase<Matrix##N>::IdentityReturnType I_##N##x##N = Matrix##N::Identity(); \
60 static const Eigen::MatrixBase<Matrix##N>::ConstantReturnType Z_##N##x##N = Matrix##N::Zero();
61 
71 
72 // Matrix expressions for accessing parts of matrices
75 
76 // Matrix formatting arguments when printing.
77 // Akin to Matlab style.
79 
83 template <class MATRIX>
85 
86  const size_t n1 = A.cols(), m1 = A.rows();
87  const size_t n2 = B.cols(), m2 = B.rows();
88 
89  if(m1!=m2 || n1!=n2) return false;
90 
91  for(size_t i=0; i<m1; i++)
92  for(size_t j=0; j<n1; j++) {
93  if(!fpEqual(A(i,j), B(i,j), tol, false)) {
94  return false;
95  }
96  }
97  return true;
98 }
99 
103 inline bool operator==(const Matrix& A, const Matrix& B) {
104  return equal_with_abs_tol(A,B,1e-9);
105 }
106 
110 inline bool operator!=(const Matrix& A, const Matrix& B) {
111  return !(A==B);
112  }
113 
117 GTSAM_EXPORT bool assert_equal(const Matrix& A, const Matrix& B, double tol = 1e-9);
118 
122 GTSAM_EXPORT bool assert_inequal(const Matrix& A, const Matrix& B, double tol = 1e-9);
123 
127 GTSAM_EXPORT bool assert_equal(const std::list<Matrix>& As, const std::list<Matrix>& Bs, double tol = 1e-9);
128 
132 GTSAM_EXPORT bool linear_independent(const Matrix& A, const Matrix& B, double tol = 1e-9);
133 
137 GTSAM_EXPORT bool linear_dependent(const Matrix& A, const Matrix& B, double tol = 1e-9);
138 
143 GTSAM_EXPORT Vector operator^(const Matrix& A, const Vector & v);
144 
146 template<class MATRIX>
147 inline MATRIX prod(const MATRIX& A, const MATRIX&B) {
148  MATRIX result = A * B;
149  return result;
150 }
151 
155 GTSAM_EXPORT void print(const Matrix& A, const std::string& s, std::ostream& stream);
156 
160 GTSAM_EXPORT void print(const Matrix& A, const std::string& s = "");
161 
165 GTSAM_EXPORT void save(const Matrix& A, const std::string &s, const std::string& filename);
166 
172 GTSAM_EXPORT std::istream& operator>>(std::istream& inputStream, Matrix& destinationMatrix);
173 
183 template<class MATRIX>
184 Eigen::Block<const MATRIX> sub(const MATRIX& A, size_t i1, size_t i2, size_t j1, size_t j2) {
185  size_t m=i2-i1, n=j2-j1;
186  return A.block(i1,j1,m,n);
187 }
188 
197 template <typename Derived1, typename Derived2>
198 void insertSub(Eigen::MatrixBase<Derived1>& fullMatrix, const Eigen::MatrixBase<Derived2>& subMatrix, size_t i, size_t j) {
199  fullMatrix.block(i, j, subMatrix.rows(), subMatrix.cols()) = subMatrix;
200 }
201 
205 GTSAM_EXPORT Matrix diag(const std::vector<Matrix>& Hs);
206 
213 template<class MATRIX>
214 const typename MATRIX::ConstColXpr column(const MATRIX& A, size_t j) {
215  return A.col(j);
216 }
217 
224 template<class MATRIX>
225 const typename MATRIX::ConstRowXpr row(const MATRIX& A, size_t j) {
226  return A.row(j);
227 }
228 
234 template<class MATRIX>
235 void zeroBelowDiagonal(MATRIX& A, size_t cols=0) {
236  const size_t m = A.rows(), n = A.cols();
237  const size_t k = (cols) ? std::min(cols, std::min(m,n)) : std::min(m,n);
238  for (size_t j=0; j<k; ++j)
239  A.col(j).segment(j+1, m-(j+1)).setZero();
240 }
241 
245 inline Matrix trans(const Matrix& A) { return A.transpose(); }
246 
248 template <int OutM, int OutN, int OutOptions, int InM, int InN, int InOptions>
249 struct Reshape {
250  //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)
252  static inline ReshapedType reshape(const Eigen::Matrix<double, InM, InN, InOptions> & in) {
253  return in.data();
254  }
255 };
256 
258 template <int M, int InOptions>
259 struct Reshape<M, M, InOptions, M, M, InOptions> {
261  static inline ReshapedType reshape(const Eigen::Matrix<double, M, M, InOptions> & in) {
262  return in;
263  }
264 };
265 
267 template <int M, int N, int InOptions>
268 struct Reshape<M, N, InOptions, M, N, InOptions> {
270  static inline ReshapedType reshape(const Eigen::Matrix<double, M, N, InOptions> & in) {
271  return in;
272  }
273 };
274 
276 template <int M, int N, int InOptions>
277 struct Reshape<N, M, InOptions, M, N, InOptions> {
279  static inline ReshapedType reshape(const Eigen::Matrix<double, M, N, InOptions> & in) {
280  return in.transpose();
281  }
282 };
283 
284 template <int OutM, int OutN, int OutOptions, int InM, int InN, int InOptions>
286  BOOST_STATIC_ASSERT(InM * InN == OutM * OutN);
288 }
289 
296 GTSAM_EXPORT std::pair<Matrix,Matrix> qr(const Matrix& A);
297 
303 GTSAM_EXPORT void inplace_QR(Matrix& A);
304 
313 GTSAM_EXPORT std::list<boost::tuple<Vector, double, double> >
314 weighted_eliminate(Matrix& A, Vector& b, const Vector& sigmas);
315 
323 GTSAM_EXPORT void householder_(Matrix& A, size_t k, bool copy_vectors=true);
324 
331 GTSAM_EXPORT void householder(Matrix& A, size_t k);
332 
340 GTSAM_EXPORT Vector backSubstituteUpper(const Matrix& U, const Vector& b, bool unit=false);
341 
349 //TODO: is this function necessary? it isn't used
350 GTSAM_EXPORT Vector backSubstituteUpper(const Vector& b, const Matrix& U, bool unit=false);
351 
359 GTSAM_EXPORT Vector backSubstituteLower(const Matrix& L, const Vector& b, bool unit=false);
360 
367 GTSAM_EXPORT Matrix stack(size_t nrMatrices, ...);
368 GTSAM_EXPORT Matrix stack(const std::vector<Matrix>& blocks);
369 
380 GTSAM_EXPORT Matrix collect(const std::vector<const Matrix *>& matrices, size_t m = 0, size_t n = 0);
381 GTSAM_EXPORT Matrix collect(size_t nrMatrices, ...);
382 
389 GTSAM_EXPORT void vector_scale_inplace(const Vector& v, Matrix& A, bool inf_mask = false); // row
390 GTSAM_EXPORT Matrix vector_scale(const Vector& v, const Matrix& A, bool inf_mask = false); // row
391 GTSAM_EXPORT Matrix vector_scale(const Matrix& A, const Vector& v, bool inf_mask = false); // column
392 
404 inline Matrix3 skewSymmetric(double wx, double wy, double wz) {
405  return (Matrix3() << 0.0, -wz, +wy, +wz, 0.0, -wx, -wy, +wx, 0.0).finished();
406 }
407 
408 template <class Derived>
410  return skewSymmetric(w(0), w(1), w(2));
411 }
412 
414 GTSAM_EXPORT Matrix inverse_square_root(const Matrix& A);
415 
417 GTSAM_EXPORT Matrix cholesky_inverse(const Matrix &A);
418 
431 GTSAM_EXPORT void svd(const Matrix& A, Matrix& U, Vector& S, Matrix& V);
432 
440 GTSAM_EXPORT boost::tuple<int, double, Vector>
441 DLT(const Matrix& A, double rank_tol = 1e-9);
442 
448 GTSAM_EXPORT Matrix expm(const Matrix& A, size_t K=7);
449 
450 std::string formatMatrixIndented(const std::string& label, const Matrix& matrix, bool makeVectorHorizontal = false);
451 
458 template <int N>
462 
464  VectorN operator()(const MatrixN& A, const VectorN& b,
465  OptionalJacobian<N, N* N> H1 = boost::none,
466  OptionalJacobian<N, N> H2 = boost::none) const {
467  const MatrixN invA = A.inverse();
468  const VectorN c = invA * b;
469  // The derivative in A is just -[c[0]*invA c[1]*invA ... c[N-1]*invA]
470  if (H1)
471  for (size_t j = 0; j < N; j++)
472  H1->template middleCols<N>(N * j) = -c[j] * invA;
473  // The derivative in b is easy, as invA*b is just a linear map:
474  if (H2) *H2 = invA;
475  return c;
476  }
477 };
478 
484 template <typename T, int N>
489 
490  // The function phi should calculate f(a)*b, with derivatives in a and b.
491  // Naturally, the derivative in b is f(a).
492  typedef boost::function<VectorN(
493  const T&, const VectorN&, OptionalJacobian<N, M>, OptionalJacobian<N, N>)>
495 
497  MultiplyWithInverseFunction(const Operator& phi) : phi_(phi) {}
498 
500  VectorN operator()(const T& a, const VectorN& b,
501  OptionalJacobian<N, M> H1 = boost::none,
502  OptionalJacobian<N, N> H2 = boost::none) const {
503  MatrixN A;
504  phi_(a, b, boost::none, A); // get A = f(a) by calling f once
505  const MatrixN invA = A.inverse();
506  const VectorN c = invA * b;
507 
508  if (H1) {
510  phi_(a, c, H, boost::none); // get derivative H of forward mapping
511  *H1 = -invA* H;
512  }
513  if (H2) *H2 = invA;
514  return c;
515  }
516 
517  private:
518  const Operator phi_;
519 };
520 
521 GTSAM_EXPORT Matrix LLt(const Matrix& A);
522 
523 GTSAM_EXPORT Matrix RtR(const Matrix& A);
524 
525 GTSAM_EXPORT Vector columnNormSquare(const Matrix &A);
526 } // namespace gtsam
527 
528 #include <boost/serialization/nvp.hpp>
529 #include <boost/serialization/array.hpp>
530 #include <boost/serialization/split_free.hpp>
531 
532 namespace boost {
533  namespace serialization {
534 
549  // split version - sends sizes ahead
550  template<class Archive,
551  typename Scalar_,
552  int Rows_,
553  int Cols_,
554  int Ops_,
555  int MaxRows_,
556  int MaxCols_>
557  void save(Archive & ar,
559  const unsigned int /*version*/) {
560  const size_t rows = m.rows(), cols = m.cols();
561  ar << BOOST_SERIALIZATION_NVP(rows);
562  ar << BOOST_SERIALIZATION_NVP(cols);
563  ar << make_nvp("data", make_array(m.data(), m.size()));
564  }
565 
566  template<class Archive,
567  typename Scalar_,
568  int Rows_,
569  int Cols_,
570  int Ops_,
571  int MaxRows_,
572  int MaxCols_>
573  void load(Archive & ar,
575  const unsigned int /*version*/) {
576  size_t rows, cols;
577  ar >> BOOST_SERIALIZATION_NVP(rows);
578  ar >> BOOST_SERIALIZATION_NVP(cols);
579  m.resize(rows, cols);
580  ar >> make_nvp("data", make_array(m.data(), m.size()));
581  }
582 
583  // templated version of BOOST_SERIALIZATION_SPLIT_FREE(Eigen::Matrix);
584  template<class Archive,
585  typename Scalar_,
586  int Rows_,
587  int Cols_,
588  int Ops_,
589  int MaxRows_,
590  int MaxCols_>
591  void serialize(Archive & ar,
593  const unsigned int version) {
594  split_free(ar, m, version);
595  }
596 
597  // specialized to Matrix for MATLAB wrapper
598  template <class Archive>
599  void serialize(Archive& ar, gtsam::Matrix& m, const unsigned int version) {
600  split_free(ar, m, version);
601  }
602 
603  } // namespace serialization
604 } // namespace boost
void print(const Matrix &A, const string &s, ostream &stream)
Definition: Matrix.cpp:155
Matrix3f m
std::string formatMatrixIndented(const std::string &label, const Matrix &matrix, bool makeVectorHorizontal)
Definition: Matrix.cpp:598
void inplace_QR(Matrix &A)
Definition: Matrix.cpp:635
boost::tuple< int, double, Vector > DLT(const Matrix &A, double rank_tol)
Definition: Matrix.cpp:567
void zeroBelowDiagonal(MATRIX &A, size_t cols=0)
Definition: base/Matrix.h:235
Eigen::Matrix< double, M, N, InOptions >::ConstTransposeReturnType ReshapedType
Definition: base/Matrix.h:278
const MATRIX::ConstRowXpr row(const MATRIX &A, size_t j)
Definition: base/Matrix.h:225
Matrix< RealScalar, Dynamic, Dynamic > M
Definition: bench_gemm.cpp:38
Eigen::Matrix< double, N, N > MatrixN
Definition: base/Matrix.h:488
MATRIX prod(const MATRIX &A, const MATRIX &B)
Definition: base/Matrix.h:147
void save(const Matrix &A, const string &s, const string &filename)
Definition: Matrix.cpp:166
bool operator!=(const Matrix &A, const Matrix &B)
Definition: base/Matrix.h:110
boost::function< VectorN(const T &, const VectorN &, OptionalJacobian< N, M >, OptionalJacobian< N, N >)> Operator
Definition: base/Matrix.h:494
const Block< const Derived, internal::traits< Derived >::RowsAtCompileTime, 1,!IsRowMajor > ConstColXpr
Definition: BlockMethods.h:15
static ReshapedType reshape(const Eigen::Matrix< double, InM, InN, InOptions > &in)
Definition: base/Matrix.h:252
Matrix diag(const std::vector< Matrix > &Hs)
Definition: Matrix.cpp:206
std::string serialize(const T &input)
serializes to a string
void householder(Matrix &A, size_t k)
Definition: Matrix.cpp:353
const Eigen::Matrix< double, M, N, InOptions > & ReshapedType
Definition: base/Matrix.h:269
#define min(a, b)
Definition: datatypes.h:19
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:94
MatrixType m2(n_dims)
ArrayXcf v
Definition: Cwise_arg.cpp:1
GTSAM_MAKE_MATRIX_DEFS(1)
void vector_scale_inplace(const Vector &v, Matrix &A, bool inf_mask)
Definition: Matrix.cpp:481
int n
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:43
Scalar Scalar * c
Definition: benchVecAdd.cpp:17
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar * data() const
static ReshapedType reshape(const Eigen::Matrix< double, M, M, InOptions > &in)
Definition: base/Matrix.h:261
static Cal3_S2 K(500, 500, 0.1, 640/2, 480/2)
MatrixXd L
Definition: LLT_example.cpp:6
const Eigen::Matrix< double, M, M, InOptions > & ReshapedType
Definition: base/Matrix.h:260
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy y set format x g set format y g set format x2 g set format y2 g set format z g set angles radians set nogrid set key title set key left top Right noreverse box linetype linewidth samplen spacing width set nolabel set noarrow set nologscale set logscale x set set pointsize set encoding default set nopolar set noparametric set set set set surface set nocontour set clabel set mapping cartesian set nohidden3d set cntrparam order set cntrparam linear set cntrparam levels auto set cntrparam points set size set set xzeroaxis lt lw set x2zeroaxis lt lw set yzeroaxis lt lw set y2zeroaxis lt lw set tics in set ticslevel set tics set mxtics default set mytics default set mx2tics default set my2tics default set xtics border mirror norotate autofreq set ytics border mirror norotate autofreq set ztics border nomirror norotate autofreq set nox2tics set noy2tics set timestamp bottom norotate set rrange[*:*] noreverse nowriteback set trange[*:*] noreverse nowriteback set urange[*:*] noreverse nowriteback set vrange[*:*] noreverse nowriteback set xlabel matrix size set x2label set timefmt d m y n H
#define N
Definition: gksort.c:12
Base class for all dense matrices, vectors, and arrays.
Definition: DenseBase.h:41
Vector backSubstituteLower(const Matrix &L, const Vector &b, bool unit)
Definition: Matrix.cpp:366
T expm(const Vector &x, int K=7)
Definition: Lie.h:316
Matrix vector_scale(const Vector &v, const Matrix &A, bool inf_mask)
Definition: Matrix.cpp:497
Matrix trans(const Matrix &A)
Definition: base/Matrix.h:245
Matrix RtR(const Matrix &A)
Definition: Matrix.cpp:529
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: base/Matrix.h:464
bool operator==(const Matrix &A, const Matrix &B)
Definition: base/Matrix.h:103
Array33i a
Eigen::Block< const Matrix > ConstSubMatrix
Definition: base/Matrix.h:74
bool fpEqual(double a, double b, double tol, bool check_relative_also)
Definition: Vector.cpp:42
EIGEN_DOC_BLOCK_ADDONS_NOT_INNER_PANEL EIGEN_DEVICE_FUNC BlockXpr block(Index startRow, Index startCol, Index blockRows, Index blockCols)
Definition: DenseBase.h:65
Matrix LLt(const Matrix &A)
Definition: Matrix.cpp:522
const Eigen::IOFormat & matlabFormat()
Definition: Matrix.cpp:139
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > MatrixRowMajor
Definition: base/Matrix.h:44
Eigen::VectorXd Vector
Definition: Vector.h:38
istream & operator>>(istream &inputStream, Matrix &destinationMatrix)
Definition: Matrix.cpp:173
Reshape functor.
Definition: base/Matrix.h:249
Values result
MultiplyWithInverseFunction(const Operator &phi)
Construct with function as explained above.
Definition: base/Matrix.h:497
Matrix3d m1
Definition: IOFormat.cpp:2
Key S(std::uint64_t j)
Eigen::Matrix< double, N, 1 > VectorN
Definition: base/Matrix.h:487
static ReshapedType reshape(const Eigen::Matrix< double, M, N, InOptions > &in)
Definition: base/Matrix.h:270
pair< Matrix, Matrix > qr(const Matrix &A)
Definition: Matrix.cpp:234
void householder_(Matrix &A, size_t k, bool copy_vectors)
Definition: Matrix.cpp:326
list< boost::tuple< Vector, double, double > > weighted_eliminate(Matrix &A, Vector &b, const Vector &sigmas)
Definition: Matrix.cpp:272
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Eigen::Block< Matrix > SubMatrix
Definition: base/Matrix.h:73
RealScalar s
Eigen::Map< const Eigen::Matrix< double, OutM, OutN, OutOptions > > ReshapedType
Definition: base/Matrix.h:251
const Block< const Derived, 1, internal::traits< Derived >::ColsAtCompileTime, IsRowMajor > ConstRowXpr
Definition: BlockMethods.h:18
const G & b
Definition: Group.h:83
RowVector3d w
Vector operator^(const Matrix &A, const Vector &v)
Definition: Matrix.cpp:130
bool equal_with_abs_tol(const Eigen::DenseBase< MATRIX > &A, const Eigen::DenseBase< MATRIX > &B, double tol=1e-9)
Definition: base/Matrix.h:84
traits
Definition: chartTesting.h:28
typedef and functions to augment Eigen&#39;s VectorXd
bool linear_dependent(const Matrix &A, const Matrix &B, double tol)
Definition: Matrix.cpp:116
bool assert_equal(const Matrix &expected, const Matrix &actual, double tol)
Definition: Matrix.cpp:42
Matrix3 skewSymmetric(double wx, double wy, double wz)
Definition: base/Matrix.h:404
Matrix inverse_square_root(const Matrix &A)
Definition: Matrix.cpp:551
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:103
Vector columnNormSquare(const Matrix &A)
Definition: Matrix.cpp:223
void insertSub(Eigen::MatrixBase< Derived1 > &fullMatrix, const Eigen::MatrixBase< Derived2 > &subMatrix, size_t i, size_t j)
Definition: base/Matrix.h:198
void svd(const Matrix &A, Matrix &U, Vector &S, Matrix &V)
Definition: Matrix.cpp:559
Eigen::Matrix< double, N, 1 > VectorN
Definition: base/Matrix.h:460
Vector backSubstituteUpper(const Matrix &U, const Vector &b, bool unit)
Definition: Matrix.cpp:376
Matrix stack(size_t nrMatrices,...)
Definition: Matrix.cpp:396
bool assert_inequal(const Matrix &A, const Matrix &B, double tol)
Definition: Matrix.cpp:62
bool linear_independent(const Matrix &A, const Matrix &B, double tol)
Definition: Matrix.cpp:102
Eigen::Block< const MATRIX > sub(const MATRIX &A, size_t i1, size_t i2, size_t j1, size_t j2)
Definition: base/Matrix.h:184
static ReshapedType reshape(const Eigen::Matrix< double, M, N, InOptions > &in)
Definition: base/Matrix.h:279
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: base/Matrix.h:500
Matrix collect(const std::vector< const Matrix * > &matrices, size_t m, size_t n)
Definition: Matrix.cpp:442
Special class for optional Jacobian arguments.
const G double tol
Definition: Group.h:83
Eigen::Matrix< double, N, N > MatrixN
Definition: base/Matrix.h:461
void load(Archive &ar, Eigen::Matrix< Scalar_, Rows_, Cols_, Ops_, MaxRows_, MaxCols_ > &m, const unsigned int)
Definition: base/Matrix.h:573
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
The matrix class, also used for vectors and row-vectors.
Matrix cholesky_inverse(const Matrix &A)
Definition: Matrix.cpp:538
const MATRIX::ConstColXpr column(const MATRIX &A, size_t j)
Definition: base/Matrix.h:214
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Stores a set of parameters controlling the way matrices are printed.
Definition: IO.h:50
std::ptrdiff_t j


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:42:47