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


gtsam
Author(s):
autogenerated on Wed Jan 22 2025 04:02:10