/* ---------------------------------------------------------------------------- * GTSAM Copyright 2010, Georgia Tech Research Corporation, * Atlanta, Georgia 30332-0415 * All Rights Reserved * Authors: Frank Dellaert, et al. (see THANKS for the full author list) * See LICENSE for the license information * -------------------------------------------------------------------------- */ /** * @file Matrix.cpp * @brief matrix class * @author Christian Potthast */ #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include using namespace std; namespace gtsam { /* ************************************************************************* */ Matrix zeros( size_t m, size_t n ) { return Matrix::Zero(m,n); } /* ************************************************************************* */ Matrix ones( size_t m, size_t n ) { return Matrix::Ones(m,n); } /* ************************************************************************* */ Matrix eye( size_t m, size_t n) { return Matrix::Identity(m, n); } /* ************************************************************************* */ Matrix diag(const Vector& v) { return v.asDiagonal(); } /* ************************************************************************* */ bool assert_equal(const Matrix& expected, const Matrix& actual, double tol) { if (equal_with_abs_tol(expected,actual,tol)) return true; size_t n1 = expected.cols(), m1 = expected.rows(); size_t n2 = actual.cols(), m2 = actual.rows(); cout << "not equal:" << endl; print(expected,"expected = "); print(actual,"actual = "); if(m1!=m2 || n1!=n2) cout << m1 << "," << n1 << " != " << m2 << "," << n2 << endl; else { Matrix diff = actual-expected; print(diff, "actual - expected = "); } return false; } /* ************************************************************************* */ bool assert_inequal(const Matrix& A, const Matrix& B, double tol) { if (!equal_with_abs_tol(A,B,tol)) return true; cout << "Erroneously equal:" << endl; print(A, "A = "); print(B, "B = "); return false; } /* ************************************************************************* */ bool assert_equal(const std::list& As, const std::list& Bs, double tol) { if (As.size() != Bs.size()) return false; list::const_iterator itA, itB; itA = As.begin(); itB = Bs.begin(); for (; itA != As.end(); itA++, itB++) if (!assert_equal(*itB, *itA, tol)) return false; return true; } /* ************************************************************************* */ static bool is_linear_dependent(const Matrix& A, const Matrix& B, double tol) { // This local static function is used by linear_independent and // linear_dependent just below. size_t n1 = A.cols(), m1 = A.rows(); size_t n2 = B.cols(), m2 = B.rows(); bool dependent = true; if(m1!=m2 || n1!=n2) dependent = false; for(size_t i=0; dependent && i& rowVec, coeffs) { destinationMatrix.row(row) = Eigen::Map(&rowVec[0], width); ++ row; } return inputStream; } /* ************************************************************************* */ Matrix diag(const std::vector& Hs) { size_t rows = 0, cols = 0; for (size_t i = 0; i qr(const Matrix& A) { const size_t m = A.rows(), n = A.cols(), kprime = min(m,n); Matrix Q=eye(m,m),R(A); Vector v(m); // loop over the kprime first columns for(size_t j=0; j < kprime; j++){ // we now work on the matrix (m-j)*(n-j) matrix A(j:end,j:end) const size_t mm=m-j; // copy column from matrix to xjm, i.e. x(j:m) = A(j:m,j) Vector xjm(mm); for(size_t k = 0 ; k < mm; k++) xjm(k) = R(j+k, j); // calculate the Householder vector v double beta; Vector vjm; boost::tie(beta,vjm) = house(xjm); // pad with zeros to get m-dimensional vector v for(size_t k = 0 ; k < m; k++) v(k) = k > weighted_eliminate(Matrix& A, Vector& b, const Vector& sigmas) { size_t m = A.rows(), n = A.cols(); // get size(A) size_t maxRank = min(m,n); // create list list > results; Vector pseudo(m); // allocate storage for pseudo-inverse Vector weights = sigmas.array().square().inverse(); // calculate weights once // We loop over all columns, because the columns that can be eliminated // are not necessarily contiguous. For each one, estimate the corresponding // scalar variable x as d-rS, with S the separator (remaining columns). // Then update A and b by substituting x with d-rS, zero-ing out x's column. for (size_t j=0; j=maxRank) break; // update A, b, expensive, using outer product // A' \define A_{S}-a*r and b'\define b-d*a A -= a * r.transpose(); b -= d * a; } return results; } /* ************************************************************************* */ /** Imperative version of Householder QR factorization, Golub & Van Loan p 224 * version with Householder vectors below diagonal, as in GVL */ /* ************************************************************************* */ void householder_(Matrix& A, size_t k, bool copy_vectors) { const size_t m = A.rows(), n = A.cols(), kprime = min(k,min(m,n)); // loop over the kprime first columns for(size_t j=0; j < kprime; j++) { // copy column from matrix to vjm, i.e. v(j:m) = A(j:m,j) Vector vjm = A.col(j).segment(j, m-j); // calculate the Householder vector, in place double beta = houseInPlace(vjm); // do outer product update A(j:m,:) = (I-beta vv')*A = A - v*(beta*A'*v)' = A - v*w' gttic(householder_update); // bottleneck for system // don't touch old columns Vector w = beta * A.block(j,j,m-j,n-j).transpose() * vjm; A.block(j,j,m-j,n-j) -= vjm * w.transpose(); gttoc(householder_update); // the Householder vector is copied in the zeroed out part if (copy_vectors) { gttic(householder_vector_copy); A.col(j).segment(j+1, m-(j+1)) = vjm.segment(1, m-(j+1)); gttoc(householder_vector_copy); } } // column j } /* ************************************************************************* */ void householder(Matrix& A, size_t k) { // version with zeros below diagonal gttic(householder_); householder_(A,k,false); gttoc(householder_); // gttic(householder_zero_fill); // const size_t m = A.rows(), n = A.cols(), kprime = min(k,min(m,n)); // for(size_t j=0; j < kprime; j++) // A.col(j).segment(j+1, m-(j+1)).setZero(); // gttoc(householder_zero_fill); } /* ************************************************************************* */ Vector backSubstituteLower(const Matrix& L, const Vector& b, bool unit) { // @return the solution x of L*x=b assert(L.rows() == L.cols()); if (unit) return L.triangularView().solve(b); else return L.triangularView().solve(b); } /* ************************************************************************* */ Vector backSubstituteUpper(const Matrix& U, const Vector& b, bool unit) { // @return the solution x of U*x=b assert(U.rows() == U.cols()); if (unit) return U.triangularView().solve(b); else return U.triangularView().solve(b); } /* ************************************************************************* */ Vector backSubstituteUpper(const Vector& b, const Matrix& U, bool unit) { // @return the solution x of x'*U=b' assert(U.rows() == U.cols()); if (unit) return U.triangularView().transpose().solve(b); else return U.triangularView().transpose().solve(b); } /* ************************************************************************* */ Matrix stack(size_t nrMatrices, ...) { size_t dimA1 = 0; size_t dimA2 = 0; va_list ap; va_start(ap, nrMatrices); for(size_t i = 0 ; i < nrMatrices ; i++) { Matrix *M = va_arg(ap, Matrix *); dimA1 += M->rows(); dimA2 = M->cols(); // TODO: should check if all the same ! } va_end(ap); va_start(ap, nrMatrices); Matrix A(dimA1, dimA2); size_t vindex = 0; for( size_t i = 0 ; i < nrMatrices ; i++) { Matrix *M = va_arg(ap, Matrix *); for(size_t d1 = 0; d1 < (size_t) M->rows(); d1++) for(size_t d2 = 0; d2 < (size_t) M->cols(); d2++) A(vindex+d1, d2) = (*M)(d1, d2); vindex += M->rows(); } return A; } /* ************************************************************************* */ Matrix stack(const std::vector& blocks) { if (blocks.size() == 1) return blocks.at(0); DenseIndex nrows = 0, ncols = blocks.at(0).cols(); BOOST_FOREACH(const Matrix& mat, blocks) { nrows += mat.rows(); if (ncols != mat.cols()) throw invalid_argument("Matrix::stack(): column size mismatch!"); } Matrix result(nrows, ncols); DenseIndex cur_row = 0; BOOST_FOREACH(const Matrix& mat, blocks) { result.middleRows(cur_row, mat.rows()) = mat; cur_row += mat.rows(); } return result; } /* ************************************************************************* */ Matrix collect(const std::vector& matrices, size_t m, size_t n) { // if we have known and constant dimensions, use them size_t dimA1 = m; size_t dimA2 = n*matrices.size(); if (!m && !n) { BOOST_FOREACH(const Matrix* M, matrices) { dimA1 = M->rows(); // TODO: should check if all the same ! dimA2 += M->cols(); } } // stl::copy version Matrix A(dimA1, dimA2); size_t hindex = 0; BOOST_FOREACH(const Matrix* M, matrices) { size_t row_len = M->cols(); A.block(0, hindex, dimA1, row_len) = *M; hindex += row_len; } return A; } /* ************************************************************************* */ Matrix collect(size_t nrMatrices, ...) { vector matrices; va_list ap; va_start(ap, nrMatrices); for( size_t i = 0 ; i < nrMatrices ; i++) { Matrix *M = va_arg(ap, Matrix *); matrices.push_back(M); } return collect(matrices); } /* ************************************************************************* */ // row scaling, in-place void vector_scale_inplace(const Vector& v, Matrix& A, bool inf_mask) { const DenseIndex m = A.rows(); if (inf_mask) { for (DenseIndex i=0; i llt(A); return llt.matrixL(); } /* ************************************************************************* */ Matrix RtR(const Matrix &A) { Eigen::LLT llt(A); return llt.matrixU(); } /* * This is not ultra efficient, but not terrible, either. */ Matrix cholesky_inverse(const Matrix &A) { Eigen::LLT llt(A); Matrix inv = eye(A.rows()); llt.matrixU().solveInPlace(inv); return inv*inv.transpose(); } /* ************************************************************************* */ // Semantics: // if B = inverse_square_root(A), then all of the following are true: // inv(B) * inv(B)' == A // inv(B' * B) == A Matrix inverse_square_root(const Matrix& A) { Eigen::LLT llt(A); Matrix inv = eye(A.rows()); llt.matrixU().solveInPlace(inv); return inv.transpose(); } /* ************************************************************************* */ void svd(const Matrix& A, Matrix& U, Vector& S, Matrix& V) { Eigen::JacobiSVD svd(A, Eigen::ComputeThinU | Eigen::ComputeThinV); U = svd.matrixU(); S = svd.singularValues(); V = svd.matrixV(); } /* ************************************************************************* */ boost::tuple DLT(const Matrix& A, double rank_tol) { // Check size of A size_t n = A.rows(), p = A.cols(), m = min(n,p); // Do SVD on A Eigen::JacobiSVD svd(A, Eigen::ComputeFullV); Vector s = svd.singularValues(); Matrix V = svd.matrixV(); // Find rank size_t rank = 0; for (size_t j = 0; j < m; j++) if (s(j) > rank_tol) rank++; // Return rank, error, and corresponding column of V double error = m

((int)rank, error, Vector(column(V, p-1))); } /* ************************************************************************* */ Matrix expm(const Matrix& A, size_t K) { Matrix E = eye(A.rows()), A_k = eye(A.rows()); for(size_t k=1;k<=K;k++) { A_k = A_k*A/double(k); E = E + A_k; } return E; } /* ************************************************************************* */ std::string formatMatrixIndented(const std::string& label, const Matrix& matrix, bool makeVectorHorizontal) { stringstream ss; const string firstline = label; ss << firstline; const string padding(firstline.size(), ' '); const bool transposeMatrix = makeVectorHorizontal && matrix.cols() == 1 && matrix.rows() > 1; const DenseIndex effectiveRows = transposeMatrix ? matrix.cols() : matrix.rows(); if(matrix.rows() > 0 && matrix.cols() > 0) { stringstream matrixPrinted; if(transposeMatrix) matrixPrinted << matrix.transpose(); else matrixPrinted << matrix; const std::string matrixStr = matrixPrinted.str(); boost::tokenizer > tok(matrixStr, boost::char_separator("\n")); DenseIndex row = 0; BOOST_FOREACH(const std::string& line, tok) { assert(row < effectiveRows); if(row > 0) ss << padding; ss << "[ " << line << " ]"; if(row < effectiveRows - 1) ss << "\n"; ++ row; } } else { ss << "Empty (" << matrix.rows() << "x" << matrix.cols() << ")"; } return ss.str(); } /* ************************************************************************* */ void inplace_QR(Matrix& A){ size_t rows = A.rows(); size_t cols = A.cols(); size_t size = std::min(rows,cols); typedef Eigen::internal::plain_diag_type::type HCoeffsType; typedef Eigen::internal::plain_row_type::type RowVectorType; HCoeffsType hCoeffs(size); RowVectorType temp(cols); #if !EIGEN_VERSION_AT_LEAST(3,2,5) Eigen::internal::householder_qr_inplace_blocked(A, hCoeffs, 48, temp.data()); #else Eigen::internal::householder_qr_inplace_blocked::run(A, hCoeffs, 48, temp.data()); #endif zeroBelowDiagonal(A); } } // namespace gtsam