| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * @file    Matrix.h | 
					
						
							|  |  |  |  * @brief   typedef and functions to augment Boost's ublas::matrix<double> | 
					
						
							|  |  |  |  * @author  Christian Potthast | 
					
						
							|  |  |  |  * @author  Kai Ni | 
					
						
							|  |  |  |  * @author  Frank Dellaert | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | // \callgraph
 | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | #pragma once
 | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-11-11 22:42:09 +08:00
										 |  |  | #include <list>
 | 
					
						
							| 
									
										
										
										
											2010-01-09 15:06:29 +08:00
										 |  |  | #include <boost/format.hpp>
 | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | #include <boost/numeric/ublas/matrix.hpp>
 | 
					
						
							| 
									
										
										
										
											2009-11-11 22:42:09 +08:00
										 |  |  | #include <boost/tuple/tuple.hpp>
 | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | #include "Vector.h"
 | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  * Vector is a *global* typedef | 
					
						
							|  |  |  |  * wrap-matlab does this typedef as well | 
					
						
							|  |  |  |  * we use the default < double,row_major,unbounded_array<double> > | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | #if ! defined (MEX_H)
 | 
					
						
							|  |  |  | typedef boost::numeric::ublas::matrix<double> Matrix; | 
					
						
							|  |  |  | #endif
 | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | namespace gtsam { | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  *  constructor with size and initial data, row order ! | 
					
						
							|  |  |  |  */ | 
					
						
							| 
									
										
										
										
											2009-12-11 12:56:29 +08:00
										 |  |  | Matrix Matrix_(size_t m, size_t n, const double* const data); | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  *  constructor with size and vector data, column order !!! | 
					
						
							|  |  |  |  */ | 
					
						
							| 
									
										
										
										
											2009-12-11 12:56:29 +08:00
										 |  |  | Matrix Matrix_(size_t m, size_t n, const Vector& v); | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  *  constructor from Vector yielding v.size()*1 vector | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | inline Matrix Matrix_(const Vector& v) { return Matrix_(v.size(),1,v);} | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  *  nice constructor, dangerous as number of arguments must be exactly right | 
					
						
							|  |  |  |  *  and you have to pass doubles !!! always use 0.0 never 0 | 
					
						
							|  |  |  | */ | 
					
						
							|  |  |  | Matrix Matrix_(size_t m, size_t n, ...); | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  * MATLAB like constructors | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | Matrix zeros(size_t m, size_t n); | 
					
						
							|  |  |  | Matrix eye(size_t m, size_t n); | 
					
						
							|  |  |  | inline Matrix eye( size_t m ) { return eye(m,m); } | 
					
						
							|  |  |  | Matrix diag(const Vector& v); | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  * equals with an tolerance | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | bool equal_with_abs_tol(const Matrix& A, const Matrix& B, double tol = 1e-9); | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  * equality is just equal_with_abs_tol 1e-9 | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | inline bool operator==(const Matrix& A, const Matrix& B) { | 
					
						
							|  |  |  |   return equal_with_abs_tol(A,B,1e-9); | 
					
						
							|  |  |  | } | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  * inequality  | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | inline bool operator!=(const Matrix& A, const Matrix& B) { | 
					
						
							|  |  |  |   return !(A==B); | 
					
						
							|  |  |  |  } | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  * equals with an tolerance, prints out message if unequal | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | bool assert_equal(const Matrix& A, const Matrix& B, double tol = 1e-9); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-02-28 08:48:42 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * equals with an tolerance, prints out message if unequal | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | bool assert_equal(const std::list<Matrix>& As, const std::list<Matrix>& Bs, double tol = 1e-9); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * overload * for matrix-vector multiplication (as BOOST does not) | 
					
						
							|  |  |  |  */ | 
					
						
							| 
									
										
										
										
											2010-01-31 10:53:03 +08:00
										 |  |  | inline Vector operator*(const Matrix& A, const Vector & v) { return prod(A,v);} | 
					
						
							| 
									
										
										
										
											2009-12-26 01:52:58 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-02-01 00:04:24 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * BLAS Level-2 style e <- e + alpha*A*x | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | void multiplyAdd(double alpha, const Matrix& A, const Vector& x, Vector& e); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-02-21 07:44:07 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * BLAS Level-2 style e <- e + A*x | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | void multiplyAdd(const Matrix& A, const Vector& x, Vector& e); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-12-26 01:52:58 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * overload ^ for trans(A)*v | 
					
						
							|  |  |  |  * We transpose the vectors for speed. | 
					
						
							|  |  |  |  */ | 
					
						
							| 
									
										
										
										
											2010-01-31 10:53:03 +08:00
										 |  |  | Vector operator^(const Matrix& A, const Vector & v); | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							| 
									
										
										
										
											2010-02-01 00:04:24 +08:00
										 |  |  |  * BLAS Level-2 style x <- x + alpha*A'*e | 
					
						
							| 
									
										
										
										
											2010-01-31 10:53:03 +08:00
										 |  |  |  */ | 
					
						
							| 
									
										
										
										
											2010-02-01 00:04:24 +08:00
										 |  |  | void transposeMultiplyAdd(double alpha, const Matrix& A, const Vector& e, Vector& x); | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-02-21 07:44:07 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * BLAS Level-2 style x <- x + A'*e | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | void transposeMultiplyAdd(const Matrix& A, const Vector& e, Vector& x); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-02-16 07:53:16 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * BLAS Level-2 style x <- x + alpha*A'*e | 
					
						
							|  |  |  |  */ | 
					
						
							| 
									
										
										
										
											2010-02-17 11:29:12 +08:00
										 |  |  | void transposeMultiplyAdd(double alpha, const Matrix& A, const Vector& e, SubVector x); | 
					
						
							| 
									
										
										
										
											2010-02-16 07:53:16 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-12-11 04:16:40 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * overload * for vector*matrix multiplication (as BOOST does not) | 
					
						
							|  |  |  |  */ | 
					
						
							| 
									
										
										
										
											2010-01-31 10:53:03 +08:00
										 |  |  | inline Vector operator*(const Vector & v, const Matrix& A) { return prod(v,A);} | 
					
						
							| 
									
										
										
										
											2009-12-11 04:16:40 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * overload * for matrix multiplication (as BOOST does not) | 
					
						
							|  |  |  |  */ | 
					
						
							| 
									
										
										
										
											2010-01-31 10:53:03 +08:00
										 |  |  | inline Matrix operator*(const Matrix& A, const Matrix& B) { return prod(A,B);} | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  * convert to column vector, column order !!! | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | Vector Vector_(const Matrix& A); | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  * print a matrix | 
					
						
							|  |  |  |  */ | 
					
						
							| 
									
										
										
										
											2010-02-24 08:12:10 +08:00
										 |  |  | void print(const Matrix& A, const std::string& s = "", std::ostream& stream = std::cout); | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  * save a matrix to file, which can be loaded by matlab | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | void save(const Matrix& A, const std::string &s, const std::string& filename); | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  * extract submatrix, slice semantics, i.e. range = [i1,i2[ excluding i2 | 
					
						
							|  |  |  |  * @param A matrix | 
					
						
							|  |  |  |  * @param i1 first row index | 
					
						
							|  |  |  |  * @param i2 last  row index + 1 | 
					
						
							|  |  |  |  * @param j1 first col index | 
					
						
							|  |  |  |  * @param j2 last  col index + 1 | 
					
						
							| 
									
										
										
										
											2009-11-05 04:59:16 +08:00
										 |  |  |  * @return submatrix A(i1:i2-1,j1:j2-1) | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  |  */ | 
					
						
							|  |  |  | Matrix sub(const Matrix& A, size_t i1, size_t i2, size_t j1, size_t j2); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-02-01 01:21:07 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * insert a submatrix IN PLACE at a specified location in a larger matrix | 
					
						
							|  |  |  |  * NOTE: there is no size checking | 
					
						
							|  |  |  |  * @param large matrix to be updated | 
					
						
							|  |  |  |  * @param small matrix to be inserted | 
					
						
							|  |  |  |  * @param i is the row of the upper left corner insert location | 
					
						
							|  |  |  |  * @param j is the column of the upper left corner insert location | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | void insertSub(Matrix& big, const Matrix& small, size_t i, size_t j); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-11-10 05:34:20 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * extracts a column from a matrix | 
					
						
							| 
									
										
										
										
											2010-01-21 00:08:14 +08:00
										 |  |  |  * NOTE: using this without the underscore is the ublas version! | 
					
						
							| 
									
										
										
										
											2009-11-10 05:34:20 +08:00
										 |  |  |  * @param matrix to extract column from | 
					
						
							|  |  |  |  * @param index of the column | 
					
						
							|  |  |  |  * @return the column in vector form | 
					
						
							|  |  |  |  */ | 
					
						
							| 
									
										
										
										
											2010-01-21 00:08:14 +08:00
										 |  |  | Vector column_(const Matrix& A, size_t j); | 
					
						
							| 
									
										
										
										
											2009-11-10 05:34:20 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-02-01 01:21:07 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * inserts a column into a matrix IN PLACE | 
					
						
							|  |  |  |  * NOTE: there is no size checking | 
					
						
							|  |  |  |  * Alternate form allows for vectors smaller than the whole column to be inserted | 
					
						
							|  |  |  |  * @param A matrix to be modified in place | 
					
						
							|  |  |  |  * @param col is the vector to be inserted | 
					
						
							|  |  |  |  * @param j is the index to insert the column | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | void insertColumn(Matrix& A, const Vector& col, size_t j); | 
					
						
							|  |  |  | void insertColumn(Matrix& A, const Vector& col, size_t i, size_t j); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-11-11 22:42:09 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * extracts a row from a matrix | 
					
						
							|  |  |  |  * @param matrix to extract row from | 
					
						
							|  |  |  |  * @param index of the row | 
					
						
							|  |  |  |  * @return the row in vector form | 
					
						
							|  |  |  |  */ | 
					
						
							| 
									
										
										
										
											2010-01-22 02:51:59 +08:00
										 |  |  | Vector row_(const Matrix& A, size_t j); | 
					
						
							| 
									
										
										
										
											2009-11-11 22:42:09 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | // left multiply with scalar
 | 
					
						
							|  |  |  | //inline Matrix operator*(double s, const Matrix& A) { return A*s;}
 | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  * solve AX=B via in-place Lu factorization and backsubstitution | 
					
						
							|  |  |  |  * After calling, A contains LU, B the solved RHS vectors | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | void solve(Matrix& A, Matrix& B); | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  * invert A | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | Matrix inverse(const Matrix& A); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-02-01 01:21:07 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * Perform updates of system matrices | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | static void updateAb(Matrix& A, Vector& b, int j, const Vector& a, | 
					
						
							|  |  |  | 		const Vector& r, double d); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * QR factorization, inefficient, best use imperative householder below | 
					
						
							|  |  |  |  * m*n matrix -> m*m Q, m*n R | 
					
						
							|  |  |  |  * @param A a matrix | 
					
						
							|  |  |  |  * @return <Q,R> rotation matrix Q, upper triangular R | 
					
						
							|  |  |  |  */  | 
					
						
							|  |  |  | std::pair<Matrix,Matrix> qr(const Matrix& A); | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  * Imperative version of Householder rank 1 update | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | void householder_update(Matrix &A, int j, double beta, const Vector& vjm); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-11-10 05:34:20 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * Imperative algorithm for in-place full elimination with | 
					
						
							|  |  |  |  * weights and constraint handling | 
					
						
							| 
									
										
										
										
											2009-11-11 22:42:09 +08:00
										 |  |  |  * @param A is a matrix to eliminate | 
					
						
							|  |  |  |  * @param b is the rhs | 
					
						
							| 
									
										
										
										
											2009-11-10 05:34:20 +08:00
										 |  |  |  * @param sigmas is a vector of the measurement standard deviation | 
					
						
							| 
									
										
										
										
											2009-11-11 22:42:09 +08:00
										 |  |  |  * @return list of r vectors, d  and sigma | 
					
						
							| 
									
										
										
										
											2009-10-27 22:21:22 +08:00
										 |  |  |  */ | 
					
						
							| 
									
										
										
										
											2009-11-11 22:42:09 +08:00
										 |  |  | std::list<boost::tuple<Vector, double, double> > | 
					
						
							|  |  |  | weighted_eliminate(Matrix& A, Vector& b, const Vector& sigmas); | 
					
						
							| 
									
										
										
										
											2009-10-27 22:21:22 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * Householder tranformation, Householder vectors below diagonal | 
					
						
							|  |  |  |  * @param k number of columns to zero out below diagonal | 
					
						
							|  |  |  |  * @param A matrix | 
					
						
							|  |  |  |  * @return nothing: in place !!! | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | void householder_(Matrix& A, size_t k); | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  * Householder tranformation, zeros below diagonal | 
					
						
							|  |  |  |  * @param k number of columns to zero out below diagonal | 
					
						
							|  |  |  |  * @param A matrix | 
					
						
							|  |  |  |  * @return nothing: in place !!! | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | void householder(Matrix& A, size_t k); | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							| 
									
										
										
										
											2009-12-29 21:59:34 +08:00
										 |  |  |  * backSubstitute U*x=b | 
					
						
							|  |  |  |  * @param U an upper triangular matrix | 
					
						
							| 
									
										
										
										
											2009-12-30 21:20:16 +08:00
										 |  |  |  * @param b an RHS vector | 
					
						
							|  |  |  |  * @param unit, set true if unit triangular | 
					
						
							| 
									
										
										
										
											2009-12-29 21:59:34 +08:00
										 |  |  |  * @return the solution x of U*x=b | 
					
						
							|  |  |  |  * TODO: use boost | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | Vector backSubstituteUpper(const Matrix& U, const Vector& b, bool unit=false); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-12-30 21:20:16 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * backSubstitute x'*U=b' | 
					
						
							|  |  |  |  * @param U an upper triangular matrix | 
					
						
							|  |  |  |  * @param b an RHS vector | 
					
						
							|  |  |  |  * @param unit, set true if unit triangular | 
					
						
							|  |  |  |  * @return the solution x of x'*U=b' | 
					
						
							|  |  |  |  * TODO: use boost | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | Vector backSubstituteUpper(const Vector& b, const Matrix& U, bool unit=false); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-12-29 21:59:34 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * backSubstitute L*x=b | 
					
						
							|  |  |  |  * @param L an lower triangular matrix | 
					
						
							| 
									
										
										
										
											2009-12-30 21:20:16 +08:00
										 |  |  |  * @param b an RHS vector | 
					
						
							|  |  |  |  * @param unit, set true if unit triangular | 
					
						
							| 
									
										
										
										
											2009-12-29 21:59:34 +08:00
										 |  |  |  * @return the solution x of L*x=b | 
					
						
							|  |  |  |  * TODO: use boost | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  |  */  | 
					
						
							| 
									
										
										
										
											2009-12-29 21:59:34 +08:00
										 |  |  | Vector backSubstituteLower(const Matrix& L, const Vector& d, bool unit=false); | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							| 
									
										
										
										
											2009-11-05 23:08:58 +08:00
										 |  |  |  * create a matrix by stacking other matrices | 
					
						
							|  |  |  |  * Given a set of matrices: A1, A2, A3... | 
					
						
							|  |  |  |  * @return combined matrix [A1; A2; A3] | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  |  */ | 
					
						
							|  |  |  | Matrix stack(size_t nrMatrices, ...); | 
					
						
							| 
									
										
										
										
											2009-11-05 23:08:58 +08:00
										 |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  * create a matrix by concatenating | 
					
						
							|  |  |  |  * Given a set of matrices: A1, A2, A3... | 
					
						
							| 
									
										
										
										
											2010-01-20 05:49:22 +08:00
										 |  |  |  * If all matrices have the same size, specifying single matrix dimensions | 
					
						
							|  |  |  |  * will avoid the lookup of dimensions | 
					
						
							|  |  |  |  * @param matrices is a vector of matrices in the order to be collected | 
					
						
							|  |  |  |  * @param m is the number of rows of a single matrix | 
					
						
							|  |  |  |  * @param n is the number of columns of a single matrix | 
					
						
							| 
									
										
										
										
											2009-11-05 23:08:58 +08:00
										 |  |  |  * @return combined matrix [A1 A2 A3] | 
					
						
							|  |  |  |  */ | 
					
						
							| 
									
										
										
										
											2010-01-20 05:49:22 +08:00
										 |  |  | Matrix collect(const std::vector<const Matrix *>& matrices, size_t m = 0, size_t n = 0); | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | Matrix collect(size_t nrMatrices, ...); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-11-06 21:43:39 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * scales a matrix row or column by the values in a vector | 
					
						
							| 
									
										
										
										
											2009-12-09 05:10:38 +08:00
										 |  |  |  * Arguments (Matrix, Vector) scales the columns, | 
					
						
							|  |  |  |  * (Vector, Matrix) scales the rows | 
					
						
							| 
									
										
										
										
											2009-11-06 21:43:39 +08:00
										 |  |  |  */ | 
					
						
							| 
									
										
										
										
											2010-02-24 14:12:56 +08:00
										 |  |  | void vector_scale_inplace(const Vector& v, Matrix& A); // row
 | 
					
						
							| 
									
										
										
										
											2009-12-09 05:10:38 +08:00
										 |  |  | Matrix vector_scale(const Vector& v, const Matrix& A); // row
 | 
					
						
							|  |  |  | Matrix vector_scale(const Matrix& A, const Vector& v); // column
 | 
					
						
							| 
									
										
										
										
											2009-11-06 21:43:39 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * skew symmetric matrix returns this: | 
					
						
							|  |  |  |  *   0  -wz   wy | 
					
						
							|  |  |  |  *  wz    0  -wx | 
					
						
							|  |  |  |  * -wy   wx    0 | 
					
						
							|  |  |  |  * @param wx 3 dimensional vector | 
					
						
							|  |  |  |  * @param wy | 
					
						
							|  |  |  |  * @param wz | 
					
						
							|  |  |  |  * @return a 3*3 skew symmetric matrix | 
					
						
							|  |  |  | */ | 
					
						
							|  |  |  | Matrix skewSymmetric(double wx, double wy, double wz); | 
					
						
							|  |  |  | inline Matrix skewSymmetric(const Vector& w) { return skewSymmetric(w(0),w(1),w(2));} | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /**
 | 
					
						
							|  |  |  |  * SVD computes economy SVD A=U*S*V' | 
					
						
							|  |  |  |  * @param A an m*n matrix | 
					
						
							| 
									
										
										
										
											2010-02-28 02:23:34 +08:00
										 |  |  |  * @param U output argument: rotation matrix | 
					
						
							|  |  |  |  * @param S output argument: vector of singular values, sorted by default, pass false as last argument to avoid sorting!!! | 
					
						
							|  |  |  |  * @param V output argument: rotation matrix | 
					
						
							| 
									
										
										
										
											2010-02-14 12:54:39 +08:00
										 |  |  |  * @param sort boolean flag to sort singular values and V | 
					
						
							| 
									
										
										
										
											2010-02-28 02:23:34 +08:00
										 |  |  |  * if m > n then U*S*V' = (m*n)*(n*n)*(n*n) (the m-n columns of U are of no use) | 
					
						
							|  |  |  |  * if m < n then U*S*V' = (m*m)*(m*m)*(m*n) (the n-m columns of V are of no use) | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  |  */  | 
					
						
							| 
									
										
										
										
											2010-02-14 12:54:39 +08:00
										 |  |  | void svd(const Matrix& A, Matrix& U, Vector& S, Matrix& V, bool sort=true); | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-02-28 02:23:34 +08:00
										 |  |  | /*
 | 
					
						
							|  |  |  |  * In place SVD, will throw an exception when m < n | 
					
						
							|  |  |  |  * @param A an m*n matrix is modified to contain U | 
					
						
							|  |  |  |  * @param S output argument: vector of singular values, sorted by default, pass false as last argument to avoid sorting!!! | 
					
						
							|  |  |  |  * @param V output argument: rotation matrix | 
					
						
							|  |  |  |  * @param sort boolean flag to sort singular values and V | 
					
						
							|  |  |  |  * if m > n then U*S*V' = (m*n)*(n*n)*(n*n) (the m-n columns of U are of no use) | 
					
						
							|  |  |  |  */ | 
					
						
							| 
									
										
										
										
											2010-02-14 12:54:39 +08:00
										 |  |  | void svd(Matrix& A, Vector& S, Matrix& V, bool sort=true); | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-12-09 04:48:13 +08:00
										 |  |  | /** Use SVD to calculate inverse square root of a matrix */ | 
					
						
							|  |  |  | Matrix inverse_square_root(const Matrix& A); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-01-17 02:39:39 +08:00
										 |  |  | /** Use SVD to calculate the positive square root of a matrix */ | 
					
						
							|  |  |  | Matrix square_root_positive(const Matrix& A); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-01-23 06:28:03 +08:00
										 |  |  | /** Calculate the LL^t decomposition of a S.P.D matrix */ | 
					
						
							| 
									
										
										
										
											2010-02-01 02:07:29 +08:00
										 |  |  | Matrix LLt(const Matrix& A); | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | /** Calculate the R^tR decomposition of a S.P.D matrix */ | 
					
						
							|  |  |  | Matrix RtR(const Matrix& A); | 
					
						
							| 
									
										
										
										
											2010-01-23 06:28:03 +08:00
										 |  |  | 
 | 
					
						
							|  |  |  | /** Return the inverse of a S.P.D. matrix.  Inversion is done via Cholesky decomposition. */ | 
					
						
							|  |  |  | Matrix cholesky_inverse(const Matrix &A); | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-03-02 10:24:38 +08:00
										 |  |  | /**
 | 
					
						
							|  |  |  |  * Numerical exponential map, naive approach, not industrial strength !!! | 
					
						
							|  |  |  |  * @param A matrix to exponentiate | 
					
						
							|  |  |  |  * @param K number of iterations | 
					
						
							|  |  |  |  */ | 
					
						
							|  |  |  | Matrix expm(const Matrix& A, int K=7); | 
					
						
							| 
									
										
										
										
											2010-01-23 06:28:03 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2009-08-22 06:23:24 +08:00
										 |  |  | // macro for unit tests
 | 
					
						
							|  |  |  | #define EQUALITY(expected,actual)\
 | 
					
						
							|  |  |  |   { if (!assert_equal(expected,actual)) \ | 
					
						
							|  |  |  |     result_.addFailure(Failure(name_, __FILE__, __LINE__, #expected, #actual)); } | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | } // namespace gtsam
 |