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 /* ************************************************************************* */
129 Vector operator^(const Matrix& A, const Vector & v) {
130  if (A.rows()!=v.size()) {
131  throw std::invalid_argument("Matrix operator^ : A.m(" + std::to_string(A.rows()) + ")!=v.size(" +
132  std::to_string(v.size()) + ")");
133  }
134 // Vector vt = v.transpose();
135 // Vector vtA = vt * A;
136 // return vtA.transpose();
137  return A.transpose() * v;
138 }
139 
141  static const Eigen::IOFormat matlab(
142  Eigen::StreamPrecision, // precision
143  Eigen::DontAlignCols, // flags set such that rowSpacers are not added
144  ", ", // coeffSeparator
145  ";\n", // rowSeparator
146  "\t", // rowPrefix
147  "", // rowSuffix
148  "[\n", // matPrefix
149  "\n]" // matSuffix
150  );
151  return matlab;
152 }
153 
154 /* ************************************************************************* */
155 //3 argument call
156 void print(const Matrix& A, const string &s, ostream& stream) {
157  stream << s << A.format(matlabFormat()) << endl;
158 }
159 
160 /* ************************************************************************* */
161 //1 or 2 argument call
162 void print(const Matrix& A, const string &s){
163  print(A, s, cout);
164 }
165 
166 /* ************************************************************************* */
167 void save(const Matrix& A, const string &s, const string& filename) {
168  fstream stream(filename.c_str(), fstream::out | fstream::app);
169  print(A, s + "=", stream);
170  stream.close();
171 }
172 
173 /* ************************************************************************* */
174 istream& operator>>(istream& inputStream, Matrix& destinationMatrix) {
175  string line;
176  FastList<vector<double> > coeffs;
177  bool first = true;
178  size_t width = 0;
179  size_t height = 0;
180  while(getline(inputStream, line)) {
181  // Read coefficients from file
182  coeffs.push_back(vector<double>());
183  if(!first)
184  coeffs.back().reserve(width);
185  stringstream lineStream(line);
186  std::copy(istream_iterator<double>(lineStream), istream_iterator<double>(),
187  back_insert_iterator<vector<double> >(coeffs.back()));
188  if(first)
189  width = coeffs.back().size();
190  if(coeffs.back().size() != width)
191  throw runtime_error("Error reading matrix from input stream, inconsistent numbers of elements in rows");
192  ++ height;
193  }
194 
195  // Copy coefficients to matrix
196  destinationMatrix.resize(height, width);
197  int row = 0;
198  for(const vector<double>& rowVec: coeffs) {
199  destinationMatrix.row(row) = Eigen::Map<const Eigen::RowVectorXd>(&rowVec[0], width);
200  ++ row;
201  }
202 
203  return inputStream;
204 }
205 
206 /* ************************************************************************* */
207 Matrix diag(const std::vector<Matrix>& Hs) {
208  size_t rows = 0, cols = 0;
209  for (size_t i = 0; i<Hs.size(); ++i) {
210  rows+= Hs[i].rows();
211  cols+= Hs[i].cols();
212  }
213  Matrix results = Matrix::Zero(rows,cols);
214  size_t r = 0, c = 0;
215  for (size_t i = 0; i<Hs.size(); ++i) {
216  insertSub(results, Hs[i], r, c);
217  r+=Hs[i].rows();
218  c+=Hs[i].cols();
219  }
220  return results;
221 }
222 
223 /* ************************************************************************* */
225  Vector v (A.cols()) ;
226  for ( size_t i = 0 ; i < (size_t) A.cols() ; ++i ) {
227  v[i] = A.col(i).dot(A.col(i));
228  }
229  return v ;
230 }
231 
232 /* ************************************************************************* */
234 /* ************************************************************************* */
235 pair<Matrix,Matrix> qr(const Matrix& A) {
236  const size_t m = A.rows(), n = A.cols(), kprime = min(m,n);
237 
238  Matrix Q=Matrix::Identity(m,m),R(A);
239  Vector v(m);
240 
241  // loop over the kprime first columns
242  for(size_t j=0; j < kprime; j++){
243 
244  // we now work on the matrix (m-j)*(n-j) matrix A(j:end,j:end)
245  const size_t mm=m-j;
246 
247  // copy column from matrix to xjm, i.e. x(j:m) = A(j:m,j)
248  Vector xjm(mm);
249  for(size_t k = 0 ; k < mm; k++)
250  xjm(k) = R(j+k, j);
251 
252  // calculate the Householder vector v
253  const auto [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<std::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<std::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(std::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 /* ************************************************************************* */
522 Matrix LLt(const Matrix& A)
523 {
525  return llt.matrixL();
526 }
527 
528 /* ************************************************************************* */
529 Matrix RtR(const Matrix &A)
530 {
532  return llt.matrixU();
533 }
534 
535 /*
536  * This is not ultra efficient, but not terrible, either.
537  */
539 {
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
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 std::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 std::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 
616  // Split the matrix string into lines and indent them
617  std::string line;
618  std::istringstream iss(matrixStr);
619  DenseIndex row = 0;
620  while (std::getline(iss, line)) {
621  assert(row < effectiveRows);
622  if(row > 0)
623  ss << padding;
624  ss << "[ " << line << " ]";
625  if(row < effectiveRows - 1)
626  ss << "\n";
627  ++ row;
628  }
629 
630  } else {
631  ss << "Empty (" << matrix.rows() << "x" << matrix.cols() << ")";
632  }
633  return ss.str();
634 }
635 
636 /* ************************************************************************* */
638  size_t rows = A.rows();
639  size_t cols = A.cols();
640  size_t size = std::min(rows,cols);
641 
643  typedef Eigen::internal::plain_row_type<Matrix>::type RowVectorType;
644  HCoeffsType hCoeffs(size);
645  RowVectorType temp(cols);
646 
647 #if !EIGEN_VERSION_AT_LEAST(3,2,5)
649 #else
651 #endif
652 
654 }
655 
656 } // namespace gtsam
gttoc
#define gttoc(label)
Definition: timing.h:296
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:272
gtsam::DLT
std::tuple< int, double, Vector > DLT(const Matrix &A, double rank_tol)
Definition: Matrix.cpp:567
gtsam::operator>>
istream & operator>>(istream &inputStream, Matrix &destinationMatrix)
Definition: Matrix.cpp:174
operator^
const EIGEN_DEVICE_FUNC CwiseBinaryOp< internal::scalar_boolean_xor_op, const Derived, const OtherDerived > operator^(const EIGEN_CURRENT_STORAGE_BASE_CLASS< OtherDerived > &other) const
Definition: ArrayCwiseBinaryOps.h:310
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:210
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:231
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:551
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:505
gtsam::diag
Matrix diag(const std::vector< Matrix > &Hs)
Definition: Matrix.cpp:207
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:326
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:423
gtsam::matlabFormat
const Eigen::IOFormat & matlabFormat()
Definition: Matrix.cpp:140
gtsam::qr
pair< Matrix, Matrix > qr(const Matrix &A)
Definition: Matrix.cpp:235
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:481
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:529
gtsam::cholesky_inverse
Matrix cholesky_inverse(const Matrix &A)
Definition: Matrix.cpp:538
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:588
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:224
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:167
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:637
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:366
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:194
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:386
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:162
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:295
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:467
gtsam::formatMatrixIndented
std::string formatMatrixIndented(const std::string &label, const Matrix &matrix, bool makeVectorHorizontal)
Definition: Matrix.cpp:598
gtsam::LLt
Matrix LLt(const Matrix &A)
Definition: Matrix.cpp:522
M
Matrix< RealScalar, Dynamic, Dynamic > M
Definition: bench_gemm.cpp:51


gtsam
Author(s):
autogenerated on Sun Dec 22 2024 04:12:12