Matrix.cpp
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 
18 #include <gtsam/global_includes.h>
19 #include <gtsam/base/Matrix.h>
20 #include <gtsam/base/timing.h>
21 #include <gtsam/base/Vector.h>
22 #include <gtsam/base/FastList.h>
23 #include <Eigen/SVD>
24 #include <Eigen/LU>
25 
26 #include <boost/tuple/tuple.hpp>
27 #include <boost/tokenizer.hpp>
28 
29 #include <cstdarg>
30 #include <cstring>
31 #include <iomanip>
32 #include <list>
33 #include <fstream>
34 #include <limits>
35 #include <iostream>
36 
37 using namespace std;
38 
39 namespace gtsam {
40 
41 /* ************************************************************************* */
42 bool assert_equal(const Matrix& expected, const Matrix& actual, double tol) {
43 
44  if (equal_with_abs_tol(expected,actual,tol)) return true;
45 
46  size_t n1 = expected.cols(), m1 = expected.rows();
47  size_t n2 = actual.cols(), m2 = actual.rows();
48 
49  cout << "not equal:" << endl;
50  print(expected,"expected = ");
51  print(actual,"actual = ");
52  if(m1!=m2 || n1!=n2)
53  cout << m1 << "," << n1 << " != " << m2 << "," << n2 << endl;
54  else {
55  Matrix diff = actual-expected;
56  print(diff, "actual - expected = ");
57  }
58  return false;
59 }
60 
61 /* ************************************************************************* */
62 bool assert_inequal(const Matrix& A, const Matrix& B, double tol) {
63  if (!equal_with_abs_tol(A,B,tol)) return true;
64  cout << "Erroneously equal:" << endl;
65  print(A, "A = ");
66  print(B, "B = ");
67  return false;
68 }
69 
70 /* ************************************************************************* */
71 bool assert_equal(const std::list<Matrix>& As, const std::list<Matrix>& Bs, double tol) {
72  if (As.size() != Bs.size()) return false;
73 
74  list<Matrix>::const_iterator itA, itB;
75  itA = As.begin(); itB = Bs.begin();
76  for (; itA != As.end(); itA++, itB++)
77  if (!assert_equal(*itB, *itA, tol))
78  return false;
79 
80  return true;
81 }
82 
83 /* ************************************************************************* */
84 static bool is_linear_dependent(const Matrix& A, const Matrix& B, double tol) {
85  // This local static function is used by linear_independent and
86  // linear_dependent just below.
87  size_t n1 = A.cols(), m1 = A.rows();
88  size_t n2 = B.cols(), m2 = B.rows();
89 
90  bool dependent = true;
91  if(m1!=m2 || n1!=n2) dependent = false;
92 
93  for(size_t i=0; dependent && i<m1; i++) {
94  if (!gtsam::linear_dependent(Vector(row(A,i)), Vector(row(B,i)), tol))
95  dependent = false;
96  }
97 
98  return dependent;
99 }
100 
101 /* ************************************************************************* */
102 bool linear_independent(const Matrix& A, const Matrix& B, double tol) {
103  if(!is_linear_dependent(A, B, tol))
104  return true;
105  else {
106  cout << "not linearly dependent:" << endl;
107  print(A,"A = ");
108  print(B,"B = ");
109  if(A.rows()!=B.rows() || A.cols()!=B.cols())
110  cout << A.rows() << "x" << A.cols() << " != " << B.rows() << "x" << B.cols() << endl;
111  return false;
112  }
113 }
114 
115 /* ************************************************************************* */
116 bool linear_dependent(const Matrix& A, const Matrix& B, double tol) {
117  if(is_linear_dependent(A, B, tol))
118  return true;
119  else {
120  cout << "not linearly dependent:" << endl;
121  print(A,"A = ");
122  print(B,"B = ");
123  if(A.rows()!=B.rows() || A.cols()!=B.cols())
124  cout << A.rows() << "x" << A.cols() << " != " << B.rows() << "x" << B.cols() << endl;
125  return false;
126  }
127 }
128 
129 /* ************************************************************************* */
130 Vector operator^(const Matrix& A, const Vector & v) {
131  if (A.rows()!=v.size()) throw std::invalid_argument(
132  boost::str(boost::format("Matrix operator^ : A.m(%d)!=v.size(%d)") % A.rows() % v.size()));
133 // Vector vt = v.transpose();
134 // Vector vtA = vt * A;
135 // return vtA.transpose();
136  return A.transpose() * v;
137 }
138 
140  static const Eigen::IOFormat matlab(
141  Eigen::StreamPrecision, // precision
142  Eigen::DontAlignCols, // flags set such that rowSpacers are not added
143  ", ", // coeffSeparator
144  ";\n", // rowSeparator
145  "\t", // rowPrefix
146  "", // rowSuffix
147  "[\n", // matPrefix
148  "\n]" // matSuffix
149  );
150  return matlab;
151 }
152 
153 /* ************************************************************************* */
154 //3 argument call
155 void print(const Matrix& A, const string &s, ostream& stream) {
156  stream << s << A.format(matlabFormat()) << endl;
157 }
158 
159 /* ************************************************************************* */
160 //1 or 2 argument call
161 void print(const Matrix& A, const string &s){
162  print(A, s, cout);
163 }
164 
165 /* ************************************************************************* */
166 void save(const Matrix& A, const string &s, const string& filename) {
167  fstream stream(filename.c_str(), fstream::out | fstream::app);
168  print(A, s + "=", stream);
169  stream.close();
170 }
171 
172 /* ************************************************************************* */
173 istream& operator>>(istream& inputStream, Matrix& destinationMatrix) {
174  string line;
175  FastList<vector<double> > coeffs;
176  bool first = true;
177  size_t width = 0;
178  size_t height = 0;
179  while(getline(inputStream, line)) {
180  // Read coefficients from file
181  coeffs.push_back(vector<double>());
182  if(!first)
183  coeffs.back().reserve(width);
184  stringstream lineStream(line);
185  std::copy(istream_iterator<double>(lineStream), istream_iterator<double>(),
186  back_insert_iterator<vector<double> >(coeffs.back()));
187  if(first)
188  width = coeffs.back().size();
189  if(coeffs.back().size() != width)
190  throw runtime_error("Error reading matrix from input stream, inconsistent numbers of elements in rows");
191  ++ height;
192  }
193 
194  // Copy coefficients to matrix
195  destinationMatrix.resize(height, width);
196  int row = 0;
197  for(const vector<double>& rowVec: coeffs) {
198  destinationMatrix.row(row) = Eigen::Map<const Eigen::RowVectorXd>(&rowVec[0], width);
199  ++ row;
200  }
201 
202  return inputStream;
203 }
204 
205 /* ************************************************************************* */
206 Matrix diag(const std::vector<Matrix>& Hs) {
207  size_t rows = 0, cols = 0;
208  for (size_t i = 0; i<Hs.size(); ++i) {
209  rows+= Hs[i].rows();
210  cols+= Hs[i].cols();
211  }
212  Matrix results = Matrix::Zero(rows,cols);
213  size_t r = 0, c = 0;
214  for (size_t i = 0; i<Hs.size(); ++i) {
215  insertSub(results, Hs[i], r, c);
216  r+=Hs[i].rows();
217  c+=Hs[i].cols();
218  }
219  return results;
220 }
221 
222 /* ************************************************************************* */
224  Vector v (A.cols()) ;
225  for ( size_t i = 0 ; i < (size_t) A.cols() ; ++i ) {
226  v[i] = A.col(i).dot(A.col(i));
227  }
228  return v ;
229 }
230 
231 /* ************************************************************************* */
233 /* ************************************************************************* */
234 pair<Matrix,Matrix> qr(const Matrix& A) {
235  const size_t m = A.rows(), n = A.cols(), kprime = min(m,n);
236 
237  Matrix Q=Matrix::Identity(m,m),R(A);
238  Vector v(m);
239 
240  // loop over the kprime first columns
241  for(size_t j=0; j < kprime; j++){
242 
243  // we now work on the matrix (m-j)*(n-j) matrix A(j:end,j:end)
244  const size_t mm=m-j;
245 
246  // copy column from matrix to xjm, i.e. x(j:m) = A(j:m,j)
247  Vector xjm(mm);
248  for(size_t k = 0 ; k < mm; k++)
249  xjm(k) = R(j+k, j);
250 
251  // calculate the Householder vector v
252  double beta; Vector vjm;
253  boost::tie(beta,vjm) = house(xjm);
254 
255  // pad with zeros to get m-dimensional vector v
256  for(size_t k = 0 ; k < m; k++)
257  v(k) = k<j ? 0.0 : vjm(k-j);
258 
259  // create Householder reflection matrix Qj = I-beta*v*v'
260  Matrix Qj = Matrix::Identity(m,m) - beta * v * v.transpose();
261 
262  R = Qj * R; // update R
263  Q = Q * Qj; // update Q
264 
265  } // column j
266 
267  return make_pair(Q,R);
268 }
269 
270 /* ************************************************************************* */
271 list<boost::tuple<Vector, double, double> >
273  size_t m = A.rows(), n = A.cols(); // get size(A)
274  size_t maxRank = min(m,n);
275 
276  // create list
277  list<boost::tuple<Vector, double, double> > results;
278 
279  Vector pseudo(m); // allocate storage for pseudo-inverse
280  Vector weights = sigmas.array().square().inverse(); // calculate weights once
281 
282  // We loop over all columns, because the columns that can be eliminated
283  // are not necessarily contiguous. For each one, estimate the corresponding
284  // scalar variable x as d-rS, with S the separator (remaining columns).
285  // Then update A and b by substituting x with d-rS, zero-ing out x's column.
286  for (size_t j=0; j<n; ++j) {
287  // extract the first column of A
288  Vector a(column(A, j));
289 
290  // Calculate weighted pseudo-inverse and corresponding precision
291  double precision = weightedPseudoinverse(a, weights, pseudo);
292 
293  // if precision is zero, no information on this column
294  if (precision < 1e-8) continue;
295 
296  // create solution and copy into r
297  Vector r(Vector::Unit(n,j));
298  for (size_t j2=j+1; j2<n; ++j2)
299  r(j2) = pseudo.dot(A.col(j2));
300 
301  // create the rhs
302  double d = pseudo.dot(b);
303 
304  // construct solution (r, d, sigma)
305  // TODO: avoid sqrt, store precision or at least variance
306  results.push_back(boost::make_tuple(r, d, 1./sqrt(precision)));
307 
308  // exit after rank exhausted
309  if (results.size()>=maxRank) break;
310 
311  // update A, b, expensive, using outer product
312  // A' \define A_{S}-a*r and b'\define b-d*a
313  A -= a * r.transpose();
314  b -= d * a;
315  }
316 
317  return results;
318 }
319 
320 /* ************************************************************************* */
325 /* ************************************************************************* */
326 void householder_(Matrix& A, size_t k, bool copy_vectors) {
327  const size_t m = A.rows(), n = A.cols(), kprime = min(k,min(m,n));
328  // loop over the kprime first columns
329  for(size_t j=0; j < kprime; j++) {
330  // copy column from matrix to vjm, i.e. v(j:m) = A(j:m,j)
331  Vector vjm = A.col(j).segment(j, m-j);
332 
333  // calculate the Householder vector, in place
334  double beta = houseInPlace(vjm);
335 
336  // do outer product update A(j:m,:) = (I-beta vv')*A = A - v*(beta*A'*v)' = A - v*w'
337  gttic(householder_update); // bottleneck for system
338  // don't touch old columns
339  Vector w = beta * A.block(j,j,m-j,n-j).transpose() * vjm;
340  A.block(j,j,m-j,n-j) -= vjm * w.transpose();
341  gttoc(householder_update);
342 
343  // the Householder vector is copied in the zeroed out part
344  if (copy_vectors) {
345  gttic(householder_vector_copy);
346  A.col(j).segment(j+1, m-(j+1)) = vjm.segment(1, m-(j+1));
347  gttoc(householder_vector_copy);
348  }
349  } // column j
350 }
351 
352 /* ************************************************************************* */
353 void householder(Matrix& A, size_t k) {
354  // version with zeros below diagonal
356  householder_(A,k,false);
358 // gttic(householder_zero_fill);
359 // const size_t m = A.rows(), n = A.cols(), kprime = min(k,min(m,n));
360 // for(size_t j=0; j < kprime; j++)
361 // A.col(j).segment(j+1, m-(j+1)).setZero();
362 // gttoc(householder_zero_fill);
363 }
364 
365 /* ************************************************************************* */
366 Vector backSubstituteLower(const Matrix& L, const Vector& b, bool unit) {
367  // @return the solution x of L*x=b
368  assert(L.rows() == L.cols());
369  if (unit)
370  return L.triangularView<Eigen::UnitLower>().solve(b);
371  else
372  return L.triangularView<Eigen::Lower>().solve(b);
373 }
374 
375 /* ************************************************************************* */
376 Vector backSubstituteUpper(const Matrix& U, const Vector& b, bool unit) {
377  // @return the solution x of U*x=b
378  assert(U.rows() == U.cols());
379  if (unit)
380  return U.triangularView<Eigen::UnitUpper>().solve(b);
381  else
382  return U.triangularView<Eigen::Upper>().solve(b);
383 }
384 
385 /* ************************************************************************* */
386 Vector backSubstituteUpper(const Vector& b, const Matrix& U, bool unit) {
387  // @return the solution x of x'*U=b'
388  assert(U.rows() == U.cols());
389  if (unit)
390  return U.triangularView<Eigen::UnitUpper>().transpose().solve<Eigen::OnTheLeft>(b);
391  else
392  return U.triangularView<Eigen::Upper>().transpose().solve<Eigen::OnTheLeft>(b);
393 }
394 
395 /* ************************************************************************* */
396 Matrix stack(size_t nrMatrices, ...)
397 {
398  size_t dimA1 = 0;
399  size_t dimA2 = 0;
400  va_list ap;
401  va_start(ap, nrMatrices);
402  for(size_t i = 0 ; i < nrMatrices ; i++) {
403  Matrix *M = va_arg(ap, Matrix *);
404  dimA1 += M->rows();
405  dimA2 = M->cols(); // TODO: should check if all the same !
406  }
407  va_end(ap);
408  va_start(ap, nrMatrices);
409  Matrix A(dimA1, dimA2);
410  size_t vindex = 0;
411  for( size_t i = 0 ; i < nrMatrices ; i++) {
412  Matrix *M = va_arg(ap, Matrix *);
413  for(size_t d1 = 0; d1 < (size_t) M->rows(); d1++)
414  for(size_t d2 = 0; d2 < (size_t) M->cols(); d2++)
415  A(vindex+d1, d2) = (*M)(d1, d2);
416  vindex += M->rows();
417  }
418 
419  return A;
420 }
421 
422 /* ************************************************************************* */
423 Matrix stack(const std::vector<Matrix>& blocks) {
424  if (blocks.size() == 1) return blocks.at(0);
425  DenseIndex nrows = 0, ncols = blocks.at(0).cols();
426  for(const Matrix& mat: blocks) {
427  nrows += mat.rows();
428  if (ncols != mat.cols())
429  throw invalid_argument("Matrix::stack(): column size mismatch!");
430  }
431  Matrix result(nrows, ncols);
432 
433  DenseIndex cur_row = 0;
434  for(const Matrix& mat: blocks) {
435  result.middleRows(cur_row, mat.rows()) = mat;
436  cur_row += mat.rows();
437  }
438  return result;
439 }
440 
441 /* ************************************************************************* */
442 Matrix collect(const std::vector<const Matrix *>& matrices, size_t m, size_t n)
443 {
444  // if we have known and constant dimensions, use them
445  size_t dimA1 = m;
446  size_t dimA2 = n*matrices.size();
447  if (!m && !n) {
448  for(const Matrix* M: matrices) {
449  dimA1 = M->rows(); // TODO: should check if all the same !
450  dimA2 += M->cols();
451  }
452  }
453 
454  // stl::copy version
455  Matrix A(dimA1, dimA2);
456  size_t hindex = 0;
457  for(const Matrix* M: matrices) {
458  size_t row_len = M->cols();
459  A.block(0, hindex, dimA1, row_len) = *M;
460  hindex += row_len;
461  }
462 
463  return A;
464 }
465 
466 /* ************************************************************************* */
467 Matrix collect(size_t nrMatrices, ...)
468 {
469  vector<const Matrix *> matrices;
470  va_list ap;
471  va_start(ap, nrMatrices);
472  for( size_t i = 0 ; i < nrMatrices ; i++) {
473  Matrix *M = va_arg(ap, Matrix *);
474  matrices.push_back(M);
475  }
476  return collect(matrices);
477 }
478 
479 /* ************************************************************************* */
480 // row scaling, in-place
481 void vector_scale_inplace(const Vector& v, Matrix& A, bool inf_mask) {
482  const DenseIndex m = A.rows();
483  if (inf_mask) {
484  for (DenseIndex i=0; i<m; ++i) {
485  const double& vi = v(i);
486  if (std::isfinite(vi))
487  A.row(i) *= vi;
488  }
489  } else {
490  for (DenseIndex i=0; i<m; ++i)
491  A.row(i) *= v(i);
492  }
493 }
494 
495 /* ************************************************************************* */
496 // row scaling
497 Matrix vector_scale(const Vector& v, const Matrix& A, bool inf_mask) {
498  Matrix M(A);
499  vector_scale_inplace(v, M, inf_mask);
500  return M;
501 }
502 
503 /* ************************************************************************* */
504 // column scaling
505 Matrix vector_scale(const Matrix& A, const Vector& v, bool inf_mask) {
506  Matrix M(A);
507  const size_t n = A.cols();
508  if (inf_mask) {
509  for (size_t j=0; j<n; ++j) {
510  const double& vj = v(j);
511  if (std::isfinite(vj))
512  M.col(j) *= vj;
513  }
514  } else {
515  for (size_t j=0; j<n; ++j)
516  M.col(j) *= v(j);
517  }
518  return M;
519 }
520 
521 /* ************************************************************************* */
523 {
524  Eigen::LLT<Matrix> llt(A);
525  return llt.matrixL();
526 }
527 
528 /* ************************************************************************* */
530 {
531  Eigen::LLT<Matrix> llt(A);
532  return llt.matrixU();
533 }
534 
535 /*
536  * This is not ultra efficient, but not terrible, either.
537  */
539 {
540  Eigen::LLT<Matrix> llt(A);
541  Matrix inv = Matrix::Identity(A.rows(),A.rows());
542  llt.matrixU().solveInPlace<Eigen::OnTheRight>(inv);
543  return inv*inv.transpose();
544 }
545 
546 /* ************************************************************************* */
547 // Semantics:
548 // if B = inverse_square_root(A), then all of the following are true:
549 // inv(B) * inv(B)' == A
550 // inv(B' * B) == A
552  Eigen::LLT<Matrix> llt(A);
553  Matrix inv = Matrix::Identity(A.rows(),A.rows());
554  llt.matrixU().solveInPlace<Eigen::OnTheRight>(inv);
555  return inv.transpose();
556 }
557 
558 /* ************************************************************************* */
559 void svd(const Matrix& A, Matrix& U, Vector& S, Matrix& V) {
561  U = svd.matrixU();
562  S = svd.singularValues();
563  V = svd.matrixV();
564 }
565 
566 /* ************************************************************************* */
567 boost::tuple<int, double, Vector> DLT(const Matrix& A, double rank_tol) {
568 
569  // Check size of A
570  size_t n = A.rows(), p = A.cols(), m = min(n,p);
571 
572  // Do SVD on A
574  Vector s = svd.singularValues();
575  Matrix V = svd.matrixV();
576 
577  // Find rank
578  size_t rank = 0;
579  for (size_t j = 0; j < m; j++)
580  if (s(j) > rank_tol) rank++;
581 
582  // Return rank, error, and corresponding column of V
583  double error = m<p ? 0 : s(m-1);
584  return boost::tuple<int, double, Vector>((int)rank, error, Vector(column(V, p-1)));
585 }
586 
587 /* ************************************************************************* */
588 Matrix expm(const Matrix& A, size_t K) {
589  Matrix E = Matrix::Identity(A.rows(),A.rows()), A_k = Matrix::Identity(A.rows(),A.rows());
590  for(size_t k=1;k<=K;k++) {
591  A_k = A_k*A/double(k);
592  E = E + A_k;
593  }
594  return E;
595 }
596 
597 /* ************************************************************************* */
598 std::string formatMatrixIndented(const std::string& label, const Matrix& matrix, bool makeVectorHorizontal)
599 {
600  stringstream ss;
601  const string firstline = label;
602  ss << firstline;
603  const string padding(firstline.size(), ' ');
604  const bool transposeMatrix = makeVectorHorizontal && matrix.cols() == 1 && matrix.rows() > 1;
605  const DenseIndex effectiveRows = transposeMatrix ? matrix.cols() : matrix.rows();
606 
607  if(matrix.rows() > 0 && matrix.cols() > 0)
608  {
609  stringstream matrixPrinted;
610  if(transposeMatrix)
611  matrixPrinted << matrix.transpose();
612  else
613  matrixPrinted << matrix;
614  const std::string matrixStr = matrixPrinted.str();
615  boost::tokenizer<boost::char_separator<char> > tok(matrixStr, boost::char_separator<char>("\n"));
616 
617  DenseIndex row = 0;
618  for(const std::string& line: tok)
619  {
620  assert(row < effectiveRows);
621  if(row > 0)
622  ss << padding;
623  ss << "[ " << line << " ]";
624  if(row < effectiveRows - 1)
625  ss << "\n";
626  ++ row;
627  }
628  } else {
629  ss << "Empty (" << matrix.rows() << "x" << matrix.cols() << ")";
630  }
631  return ss.str();
632 }
633 
634 /* ************************************************************************* */
636  size_t rows = A.rows();
637  size_t cols = A.cols();
638  size_t size = std::min(rows,cols);
639 
641  typedef Eigen::internal::plain_row_type<Matrix>::type RowVectorType;
642  HCoeffsType hCoeffs(size);
643  RowVectorType temp(cols);
644 
645 #if !EIGEN_VERSION_AT_LEAST(3,2,5)
647 #else
649 #endif
650 
652 }
653 
654 } // namespace gtsam
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
cout<< "Here is the matrix m:"<< endl<< m<< endl;JacobiSVD< MatrixXf > svd(m, ComputeThinU|ComputeThinV)
Matrix< RealScalar, Dynamic, Dynamic > M
Definition: bench_gemm.cpp:38
static void run(MatrixQR &mat, HCoeffs &hCoeffs, Index maxBlockSize=32, typename MatrixQR::Scalar *tempData=0)
void save(const Matrix &A, const string &s, const string &filename)
Definition: Matrix.cpp:166
Key E(std::uint64_t j)
double weightedPseudoinverse(const Vector &a, const Vector &weights, Vector &pseudo)
Definition: Vector.cpp:246
A thin wrapper around std::list that uses boost&#39;s fast_pool_allocator.
Scalar * b
Definition: benchVecAdd.cpp:17
Matrix vector_scale(const Matrix &A, const Vector &v, bool inf_mask)
Definition: Matrix.cpp:505
return int(ret)+1
EIGEN_DEVICE_FUNC const CwiseBinaryOp< internal::scalar_boolean_xor_op, const Derived, const OtherDerived > operator^(const EIGEN_CURRENT_STORAGE_BASE_CLASS< OtherDerived > &other) const
Matrix diag(const std::vector< Matrix > &Hs)
Definition: Matrix.cpp:206
const SingularValuesType & singularValues() const
Definition: SVDBase.h:111
#define min(a, b)
Definition: datatypes.h:19
Matrix expected
Definition: testMatrix.cpp:974
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:94
bool assert_equal(const std::list< Matrix > &As, const std::list< Matrix > &Bs, double tol)
Definition: Matrix.cpp:71
MatrixType m2(n_dims)
ArrayXcf v
Definition: Cwise_arg.cpp:1
Traits::MatrixU matrixU() const
Definition: LLT.h:119
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
Rot2 R(Rot2::fromAngle(0.1))
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
static Cal3_S2 K(500, 500, 0.1, 640/2, 480/2)
MatrixXd L
Definition: LLT_example.cpp:6
Definition: Half.h:150
#define isfinite(X)
Definition: main.h:74
Vector backSubstituteLower(const Matrix &L, const Vector &b, bool unit)
Definition: Matrix.cpp:366
void print(const Matrix &A, const string &s)
Definition: Matrix.cpp:161
ptrdiff_t DenseIndex
The index type for Eigen objects.
Definition: types.h:67
Matrix RtR(const Matrix &A)
Definition: Matrix.cpp:529
Included from all GTSAM files.
Array33i a
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
Matrix expm(const Matrix &A, size_t K)
Definition: Matrix.cpp:588
Matrix LLt(const Matrix &A)
Definition: Matrix.cpp:522
#define gttic(label)
Definition: timing.h:280
const Eigen::IOFormat & matlabFormat()
Definition: Matrix.cpp:139
void householder(const MatrixType &m)
Definition: householder.cpp:13
Vector backSubstituteUpper(const Vector &b, const Matrix &U, bool unit)
Definition: Matrix.cpp:386
constexpr int first(int i)
Implementation details for constexpr functions.
Eigen::VectorXd Vector
Definition: Vector.h:38
Tuple< Args... > make_tuple(Args...args)
Creates a tuple object, deducing the target type from the types of arguments.
istream & operator>>(istream &inputStream, Matrix &destinationMatrix)
Definition: Matrix.cpp:173
Values result
Matrix3d m1
Definition: IOFormat.cpp:2
Key S(std::uint64_t j)
m row(1)
pair< Matrix, Matrix > qr(const Matrix &A)
Definition: Matrix.cpp:234
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition: LLT.h:56
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
static bool is_linear_dependent(const Matrix &A, const Matrix &B, double tol)
Definition: Matrix.cpp:84
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Matrix collect(size_t nrMatrices,...)
Definition: Matrix.cpp:467
RealScalar s
RowVector3d w
static std::stringstream ss
Definition: testBTree.cpp:33
Matrix stack(const std::vector< Matrix > &blocks)
Definition: Matrix.cpp:423
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
Matrix inverse_square_root(const Matrix &A)
Definition: Matrix.cpp:551
Vector columnNormSquare(const Matrix &A)
Definition: Matrix.cpp:223
#define gttoc(label)
Definition: timing.h:281
void insertSub(Eigen::MatrixBase< Derived1 > &fullMatrix, const Eigen::MatrixBase< Derived2 > &subMatrix, size_t i, size_t j)
Definition: base/Matrix.h:198
The quaternion class used to represent 3D orientations and rotations.
std::map< std::string, Array< float, 1, 8, DontAlign|RowMajor > > results
cout precision(2)
Two-sided Jacobi SVD decomposition of a rectangular matrix.
float * p
bool assert_inequal(const Matrix &A, const Matrix &B, double tol)
Definition: Matrix.cpp:62
pair< double, Vector > house(const Vector &x)
Definition: Vector.cpp:237
bool linear_independent(const Matrix &A, const Matrix &B, double tol)
Definition: Matrix.cpp:102
static double error
Definition: testRot3.cpp:39
const G double tol
Definition: Group.h:83
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Traits::MatrixL matrixL() const
Definition: LLT.h:126
Matrix cholesky_inverse(const Matrix &A)
Definition: Matrix.cpp:538
Eigen::Matrix< double, Eigen::Dynamic, 1 > Vector
const MATRIX::ConstColXpr column(const MATRIX &A, size_t j)
Definition: base/Matrix.h:214
Stores a set of parameters controlling the way matrices are printed.
Definition: IO.h:50
double houseInPlace(Vector &v)
Definition: Vector.cpp:212
std::ptrdiff_t j
Timing utilities.


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