/* * MatrixContent.cpp * * Created on: Nov 14, 2010 * Author: heber */ // include config.h #ifdef HAVE_CONFIG_H #include #endif #include "CodePatterns/MemDebug.hpp" #include "CodePatterns/Assert.hpp" #include "Exceptions/NotInvertibleException.hpp" #include "LinearAlgebra/defs.hpp" #include "LinearAlgebra/fast_functions.hpp" #include "LinearAlgebra/MatrixContent.hpp" #include "LinearAlgebra/RealSpaceMatrix.hpp" #include "LinearAlgebra/Vector.hpp" #include "LinearAlgebra/VectorContent.hpp" #include #include #include #include #include #include #include #include #include #include #include using namespace std; /** Constructor for class MatrixContent. * \param rows number of rows * \param columns number of columns */ MatrixContent::MatrixContent(size_t _rows, size_t _columns) : rows(_rows), columns(_columns), free_content_on_exit(true) { content = gsl_matrix_calloc(rows, columns); } /** Constructor of class VectorContent. * We need this MatrixBaseCase for the VectorContentView class. * There no content should be allocated, as it is just a view with an internal * gsl_vector_view. Hence, MatrixBaseCase is just dummy class to give the * constructor a unique signature. * \param MatrixBaseCase */ MatrixContent::MatrixContent(size_t _rows, size_t _columns, MatrixBaseCase) : rows(_rows), columns(_columns), free_content_on_exit(true) {} /** Constructor for class MatrixContent. * \param rows number of rows * \param columns number of columns * \param *src array with components to initialize matrix with */ MatrixContent::MatrixContent(size_t _rows, size_t _columns, const double *src) : rows(_rows), columns(_columns), free_content_on_exit(true) { content = gsl_matrix_calloc(rows, columns); set(0,0, src[0]); set(1,0, src[1]); set(2,0, src[2]); set(0,1, src[3]); set(1,1, src[4]); set(2,1, src[5]); set(0,2, src[6]); set(1,2, src[7]); set(2,2, src[8]); } /** Constructor that parses square matrix from a stream. * * \note Matrix dimensions can be preparsed via * MatrixContent::preparseMatrixDimensions() without harming the stream. * * \param &inputstream stream to parse from */ MatrixContent::MatrixContent(size_t _row, size_t _column, std::istream &inputstream) : content(NULL), rows(_row), columns(_column), free_content_on_exit(true) { // allocate matrix and set contents content = gsl_matrix_alloc(rows, columns); size_t row = 0; do { std::string line; getline(inputstream, line); //std::cout << line << std::endl; std::stringstream parseline(line); // skip comments if ((parseline.peek() == '#') || (parseline.str().empty())) continue; // break on empty lines if (parseline.peek() == '\n') break; // parse line with values std::vector line_of_values; do { double value; parseline >> value >> ws; line_of_values.push_back(value); } while (parseline.good()); // check number of columns parsed ASSERT(line_of_values.size() == columns, "MatrixContent::MatrixContent() - row " +toString(row)+" has a different number of columns " +toString(line_of_values.size())+" than others before " +toString(columns)+"."); for (size_t column = 0; column < columns; ++column) set(row, column, line_of_values[column]); ++row; } while (inputstream.good()); // check number of rows parsed ASSERT(row == rows, "MatrixContent::MatrixContent() - insufficent number of rows " +toString(row)+" < "+toString(rows)+" read."); } /** Constructor for class MatrixContent. * * \param *src source gsl_matrix vector to copy and put into this class */ MatrixContent::MatrixContent(gsl_matrix *&src) : rows(src->size1), columns(src->size2), free_content_on_exit(true) { content = gsl_matrix_alloc(src->size1, src->size2); gsl_matrix_memcpy(content,src); // content = src; // src = NULL; } /** Copy constructor for class MatrixContent. * \param &src reference to source MatrixContent */ MatrixContent::MatrixContent(const MatrixContent &src) : rows(src.rows), columns(src.columns), free_content_on_exit(true) { content = gsl_matrix_alloc(src.rows, src.columns); gsl_matrix_memcpy(content,src.content); } /** Copy constructor for class MatrixContent. * \param *src pointer to source MatrixContent */ MatrixContent::MatrixContent(const MatrixContent *src) : rows(src->rows), columns(src->columns), free_content_on_exit(true) { ASSERT(src != NULL, "MatrixContent::MatrixContent - pointer to source matrix is NULL!"); content = gsl_matrix_alloc(src->rows, src->columns); gsl_matrix_memcpy(content,src->content); } /** Destructor for class MatrixContent. */ MatrixContent::~MatrixContent() { if (free_content_on_exit) gsl_matrix_free(content); } /** Getter for MatrixContent::rows. * \return MatrixContent::rows */ const size_t MatrixContent::getRows() const { return rows; } /** Getter for MatrixContent::columns. * \return MatrixContent::columns */ const size_t MatrixContent::getColumns() const { return columns; } /** Return a VectorViewContent of the \a column -th column vector. * * @param column index of column * @return column of matrix as VectorContent */ VectorContent *MatrixContent::getColumnVector(size_t column) const { ASSERT(column < columns, "MatrixContent::getColumnVector() - requested column "+toString(column) +" greater than dimension "+toString(columns)); return (new VectorViewContent(gsl_matrix_column(content,column))); } /** Returns a VectorViewContent of the \a row -th row vector. * @param row row index * @return VectorContent of row vector */ VectorContent *MatrixContent::getRowVector(size_t row) const { ASSERT(row < rows, "MatrixContent::getColumnVector() - requested row "+toString(row) +" greater than dimension "+toString(rows)); return (new VectorViewContent(gsl_matrix_row(content,row))); } /** Returns the main diagonal of the matrix as VectorContent. * @return diagonal as VectorContent. */ VectorContent *MatrixContent::getDiagonalVector() const { return (new VectorViewContent(gsl_matrix_diagonal(content))); } /** Set matrix to identity. */ void MatrixContent::setIdentity() { for(int i=rows;i--;){ for(int j=columns;j--;){ set(i,j,(double)(i==j)); } } } /** Set all matrix components to zero. */ void MatrixContent::setZero() { for(int i=rows;i--;){ for(int j=columns;j--;){ set(i,j,0.); } } } /** Set all matrix components to a given value. * \param _value value to set each component to */ void MatrixContent::setValue(double _value) { for(int i=rows;i--;){ for(int j=columns;j--;){ set(i,j,_value); } } } /** Copy operator for MatrixContent with self-assignment check. * \param &src matrix to compare to * \return reference to this */ MatrixContent &MatrixContent::operator=(const MatrixContent &src) { if(&src!=this){ gsl_matrix_memcpy(content,src.content); } return *this; } /** Addition operator. * \param &rhs matrix to add * \return reference to this */ const MatrixContent &MatrixContent::operator+=(const MatrixContent &rhs) { gsl_matrix_add(content, rhs.content); return *this; } /** Subtraction operator. * \param &rhs matrix to subtract * \return reference to this */ const MatrixContent &MatrixContent::operator-=(const MatrixContent &rhs) { gsl_matrix_sub(content, rhs.content); return *this; } /** Multiplication operator. * Note that here matrix have to have same dimensions. * \param &rhs matrix to multiply with * \return reference to this */ const MatrixContent &MatrixContent::operator*=(const MatrixContent &rhs) { ASSERT(rhs.columns == rhs.rows, "MatrixContent::operator*=() - rhs matrix is not square: "+toString(rhs.columns)+" != "+toString(rhs.rows)+"."); ASSERT(columns == rhs.rows, "MatrixContent::operator*=() - columns dimension differ: "+toString(columns)+" != "+toString(rhs.rows)+"."); (*this) = (*this)*rhs; return *this; } /** Multiplication with copy operator. * \param &rhs matrix to multiply with * \return reference to newly allocated MatrixContent */ const MatrixContent MatrixContent::operator*(const MatrixContent &rhs) const { ASSERT (columns == rhs.rows, "MatrixContent::operator*() - dimensions not match for matrix product (a,b)*(b,c) = (a,c):" "("+toString(rows)+","+toString(columns)+")*("+toString(rhs.rows)+","+toString(rhs.columns)+")"); gsl_matrix *res = gsl_matrix_alloc(rows, rhs.columns); gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, content, rhs.content, 0.0, res); // gsl_matrix is taken over by constructor, hence no free MatrixContent tmp(res); gsl_matrix_free(res); return tmp; } /** Hadamard multiplication with copy operator. * The Hadamard product is component-wise matrix product. * \param &rhs matrix to hadamard-multiply with * \return reference to newly allocated MatrixContent */ const MatrixContent MatrixContent::operator&(const MatrixContent &rhs) const { ASSERT ((rows == rhs.rows) && (columns == rhs.columns), "MatrixContent::operator&() - dimensions not match for matrix product (a,b) != (b,c):" "("+toString(rows)+","+toString(columns)+") != ("+toString(rhs.rows)+","+toString(rhs.columns)+")"); gsl_matrix *res = gsl_matrix_alloc(rows, rhs.columns); for (size_t i=0;i=0) && (i=0) && (j=0) && (i=0) && (j=0) && (i=0) && (j columns ? rows : columns; gsl_matrix *content_square = gsl_matrix_alloc(greaterDimension, greaterDimension); for (size_t i=0; isize1, content_square->size2)); std::cout << "The squared matrix is " << *ContentSquare << std::endl; // solve eigenvalue problem gsl_eigen_nonsymmv_workspace *T = gsl_eigen_nonsymmv_alloc(rows); gsl_vector_complex *eval = gsl_vector_complex_alloc(greaterDimension); gsl_matrix_complex *evec = gsl_matrix_complex_alloc(greaterDimension, greaterDimension); gsl_eigen_nonsymmv(content_square, eval, evec, T); gsl_eigen_nonsymmv_free(T); // copy eigenvectors real-parts into content_square and ... for (size_t i=0; i LINALG_MYEPSILON()) { // only take eigenvectors with value > 0 std::cout << i << "th eigenvalue is (" << GSL_REAL(gsl_vector_complex_get(eval,i)) << "," << GSL_IMAG(gsl_vector_complex_get(eval,i)) << ")" << std::endl; for (size_t j=0; j LINALG_MYEPSILON()) std::cerr << "MatrixContent::transformToEigenbasis() - WARNING: eigenvectors are complex-valued!" << std::endl; gsl_matrix_set(content, j,I, GSL_REAL(gsl_matrix_complex_get(evec,j,i))); } if (fabs(GSL_IMAG(gsl_vector_complex_get(eval,I))) > LINALG_MYEPSILON()) std::cerr << "MatrixContent::transformToEigenbasis() - WARNING: eigenvectors are complex-valued!" << std::endl; gsl_vector_set(eval_real, I, GSL_REAL(gsl_vector_complex_get(eval, i))); I++; } } gsl_matrix_complex_free(evec); gsl_vector_complex_free(eval); delete ContentSquare; return eval_real; } } /** Sorts the eigenpairs in ascending order of the eigenvalues. * We assume that MatrixContent::transformToEigenbasis() has just been called. * @param eigenvalues vector of eigenvalue from * MatrixContent::transformToEigenbasis() */ void MatrixContent::sortEigenbasis(gsl_vector *eigenvalues) { gsl_eigen_symmv_sort (eigenvalues, content, GSL_EIGEN_SORT_VAL_ASC); } /* ============================ Properties ============================== */ /** Checks whether matrix' elements are strictly null. * \return true - is null, false - else */ bool MatrixContent::IsNull() const { return gsl_matrix_isnull (content); }; /** Checks whether matrix' elements are strictly positive. * \return true - is positive, false - else */ bool MatrixContent::IsPositive() const { return gsl_matrix_ispos (content); }; /** Checks whether matrix' elements are strictly negative. * \return true - is negative, false - else */ bool MatrixContent::IsNegative() const { return gsl_matrix_isneg (content); }; /** Checks whether matrix' elements are strictly non-negative. * \return true - is non-negative, false - else */ bool MatrixContent::IsNonNegative() const { return gsl_matrix_isnonneg (content); }; /** This function performs a Cholesky decomposition to determine whether matrix is positive definite. * We check whether GSL returns GSL_EDOM as error, indicating that decomposition failed due to matrix not being positive-definite. * \return true - matrix is positive-definite, false - else */ bool MatrixContent::IsPositiveDefinite() const { if (rows != columns) // only possible for square matrices. return false; else return (gsl_linalg_cholesky_decomp (content) != GSL_EDOM); }; /** Calculates the determinant of the matrix. * if matrix is square, uses LU decomposition. */ double MatrixContent::Determinant() const { int signum = 0; ASSERT(rows == columns, "MatrixContent::Determinant() - determinant can only be calculated for square matrices."); gsl_permutation *p = gsl_permutation_alloc(rows); gsl_linalg_LU_decomp(content, p, &signum); gsl_permutation_free(p); return gsl_linalg_LU_det(content, signum); }; /** Preparses a matrix in an input stream for row and column count. * * \note This does not change the get position within the stream. * * @param inputstream to parse matrix from * @return pair with (row, column) */ std::pair MatrixContent::preparseMatrixDimensions(std::istream &inputstream) { const std::streampos initial_position(inputstream.tellg()); // std::cout << "We are at position " // +toString((int)inputstream.tellg())+" from start of stream." << std::endl; size_t rows = 0; size_t columns = 0; do { std::string line; getline(inputstream, line); std::stringstream parseline(line); // skip comments if ((parseline.peek() == '#') || (line.empty())) continue; // break on empty lines if (parseline.peek() == '\n') break; // parse line with values std::vector line_of_values; do { double value; parseline >> value >> ws; line_of_values.push_back(value); } while (parseline.good()); // set matrixdimension if first line if (columns == 0) { columns = line_of_values.size(); } else { // otherwise check ASSERT(columns == line_of_values.size(), "MatrixContent::MatrixContent() - row " +toString(rows)+" has a different number of columns " +toString(line_of_values.size())+" than this matrix has " +toString(columns)+"."); } ++rows; } while (inputstream.good()); // clear end of file flags inputstream.clear(); // reset to initial position inputstream.seekg(initial_position); // std::cout << "We are again at position " // +toString((int)inputstream.tellg())+" from start of stream." << std::endl; return make_pair(rows, columns); } /** Writes matrix content to ostream in such a manner that it can be re-parsed * by the code in the constructor. * * @param ost stream to write to */ void MatrixContent::write(std::ostream &ost) const { for (size_t row = 0; row < rows; ++row) { for (size_t column = 0; column < columns; ++column) ost << at(row, column) << "\t"; ost << std::endl; } } /* ============================= Operators =============================== */ /** Scalar multiplication operator. * \param factor factor to scale with */ const MatrixContent &MatrixContent::operator*=(const double factor) { gsl_matrix_scale(content, factor); return *this; } /** Scalar multiplication and copy operator. * \param factor factor to scale with * \param &mat MatrixContent to scale * \return copied and scaled MatrixContent */ const MatrixContent operator*(const double factor,const MatrixContent& mat) { MatrixContent tmp = mat; tmp*=factor; return tmp; } /** Scalar multiplication and copy operator (with operands exchanged). * \param &mat MatrixContent to scale * \param factor factor to scale with * \return copied and scaled MatrixContent */ const MatrixContent operator*(const MatrixContent &mat,const double factor) { return factor*mat; } /** Sums two MatrixContents \a and \b component-wise. * \param a first MatrixContent * \param b second MatrixContent * \return a + b */ MatrixContent const operator+(const MatrixContent& a, const MatrixContent& b) { ASSERT(a.rows == b.rows, "MatrixContent::operator+() - row counts have to match: " +toString(a.rows)+" != "+toString(b.rows)+"!"); ASSERT(a.columns == b.columns, "MatrixContent::operator+() - column counts have to match: " +toString(a.columns)+" != "+toString(b.columns)+"!"); MatrixContent x(a); x += b; return x; }; /** Subtracts MatrixContent \a from \b component-wise. * \param a first MatrixContent * \param b second MatrixContent * \return a - b */ MatrixContent const operator-(const MatrixContent& a, const MatrixContent& b) { ASSERT(a.rows == b.rows, "MatrixContent::operator+() - row counts have to match: " +toString(a.rows)+" != "+toString(b.rows)+"!"); ASSERT(a.columns == b.columns, "MatrixContent::operator+() - column counts have to match: " +toString(a.columns)+" != "+toString(b.columns)+"!"); MatrixContent x(a); x -= b; return x; }; /** Equality operator. * Note that we use numerical sensible checking, i.e. with threshold LINALG_MYEPSILON(). * \param &rhs MatrixContent to checks against */ bool MatrixContent::operator==(const MatrixContent &rhs) const { if ((rows == rhs.rows) && (columns == rhs.columns)) { for(int i=rows;i--;){ for(int j=columns;j--;){ if(fabs(at(i,j)-rhs.at(i,j))>LINALG_MYEPSILON()){ return false; } } } return true; } return false; } std::ostream & operator<<(std::ostream &ost, const MatrixContent &mat) { ost << "\\begin{pmatrix}"; for (size_t i=0;i