| 1 | /* | 
|---|
| 2 | * MatrixContent.cpp | 
|---|
| 3 | * | 
|---|
| 4 | *  Created on: Nov 14, 2010 | 
|---|
| 5 | *      Author: heber | 
|---|
| 6 | */ | 
|---|
| 7 |  | 
|---|
| 8 |  | 
|---|
| 9 | // include config.h | 
|---|
| 10 | #ifdef HAVE_CONFIG_H | 
|---|
| 11 | #include <config.h> | 
|---|
| 12 | #endif | 
|---|
| 13 |  | 
|---|
| 14 | #include "CodePatterns/MemDebug.hpp" | 
|---|
| 15 |  | 
|---|
| 16 | #include "CodePatterns/Assert.hpp" | 
|---|
| 17 | #include "Exceptions/NotInvertibleException.hpp" | 
|---|
| 18 | #include "LinearAlgebra/defs.hpp" | 
|---|
| 19 | #include "LinearAlgebra/fast_functions.hpp" | 
|---|
| 20 | #include "LinearAlgebra/MatrixContent.hpp" | 
|---|
| 21 | #include "LinearAlgebra/RealSpaceMatrix.hpp" | 
|---|
| 22 | #include "LinearAlgebra/Vector.hpp" | 
|---|
| 23 | #include "LinearAlgebra/VectorContent.hpp" | 
|---|
| 24 |  | 
|---|
| 25 | #include <gsl/gsl_blas.h> | 
|---|
| 26 | #include <gsl/gsl_eigen.h> | 
|---|
| 27 | #include <gsl/gsl_linalg.h> | 
|---|
| 28 | #include <gsl/gsl_matrix.h> | 
|---|
| 29 | #include <gsl/gsl_multimin.h> | 
|---|
| 30 | #include <gsl/gsl_vector.h> | 
|---|
| 31 | #include <cmath> | 
|---|
| 32 | #include <cassert> | 
|---|
| 33 | #include <iostream> | 
|---|
| 34 | #include <limits> | 
|---|
| 35 | #include <set> | 
|---|
| 36 |  | 
|---|
| 37 | using namespace std; | 
|---|
| 38 |  | 
|---|
| 39 |  | 
|---|
| 40 | /** Constructor for class MatrixContent. | 
|---|
| 41 | * \param rows number of rows | 
|---|
| 42 | * \param columns number of columns | 
|---|
| 43 | */ | 
|---|
| 44 | MatrixContent::MatrixContent(size_t _rows, size_t _columns) : | 
|---|
| 45 | rows(_rows), | 
|---|
| 46 | columns(_columns), | 
|---|
| 47 | free_content_on_exit(true) | 
|---|
| 48 | { | 
|---|
| 49 | content = gsl_matrix_calloc(rows, columns); | 
|---|
| 50 | } | 
|---|
| 51 |  | 
|---|
| 52 | /** Constructor of class VectorContent. | 
|---|
| 53 | * We need this MatrixBaseCase for the VectorContentView class. | 
|---|
| 54 | * There no content should be allocated, as it is just a view with an internal | 
|---|
| 55 | * gsl_vector_view. Hence, MatrixBaseCase is just dummy class to give the | 
|---|
| 56 | * constructor a unique signature. | 
|---|
| 57 | * \param MatrixBaseCase | 
|---|
| 58 | */ | 
|---|
| 59 | MatrixContent::MatrixContent(size_t _rows, size_t _columns, MatrixBaseCase) : | 
|---|
| 60 | rows(_rows), | 
|---|
| 61 | columns(_columns), | 
|---|
| 62 | free_content_on_exit(true) | 
|---|
| 63 | {} | 
|---|
| 64 |  | 
|---|
| 65 | /** Constructor for class MatrixContent. | 
|---|
| 66 | * \param rows number of rows | 
|---|
| 67 | * \param columns number of columns | 
|---|
| 68 | * \param *src array with components to initialize matrix with | 
|---|
| 69 | */ | 
|---|
| 70 | MatrixContent::MatrixContent(size_t _rows, size_t _columns, const double *src) : | 
|---|
| 71 | rows(_rows), | 
|---|
| 72 | columns(_columns), | 
|---|
| 73 | free_content_on_exit(true) | 
|---|
| 74 | { | 
|---|
| 75 | content = gsl_matrix_calloc(rows, columns); | 
|---|
| 76 | set(0,0, src[0]); | 
|---|
| 77 | set(1,0, src[1]); | 
|---|
| 78 | set(2,0, src[2]); | 
|---|
| 79 |  | 
|---|
| 80 | set(0,1, src[3]); | 
|---|
| 81 | set(1,1, src[4]); | 
|---|
| 82 | set(2,1, src[5]); | 
|---|
| 83 |  | 
|---|
| 84 | set(0,2, src[6]); | 
|---|
| 85 | set(1,2, src[7]); | 
|---|
| 86 | set(2,2, src[8]); | 
|---|
| 87 | } | 
|---|
| 88 |  | 
|---|
| 89 | /** Constructor that parses square matrix from a stream. | 
|---|
| 90 | * | 
|---|
| 91 | * \note Matrix dimensions can be preparsed via | 
|---|
| 92 | * MatrixContent::preparseMatrixDimensions() without harming the stream. | 
|---|
| 93 | * | 
|---|
| 94 | * \param &inputstream stream to parse from | 
|---|
| 95 | */ | 
|---|
| 96 | MatrixContent::MatrixContent(size_t _row, size_t _column, std::istream &inputstream) : | 
|---|
| 97 | content(NULL), | 
|---|
| 98 | rows(_row), | 
|---|
| 99 | columns(_column), | 
|---|
| 100 | free_content_on_exit(true) | 
|---|
| 101 | { | 
|---|
| 102 | // allocate matrix and set contents | 
|---|
| 103 | content = gsl_matrix_alloc(rows, columns); | 
|---|
| 104 |  | 
|---|
| 105 | size_t row = 0; | 
|---|
| 106 | do { | 
|---|
| 107 | std::string line; | 
|---|
| 108 | getline(inputstream, line); | 
|---|
| 109 | //std::cout << line << std::endl; | 
|---|
| 110 | std::stringstream parseline(line); | 
|---|
| 111 | // skip comments | 
|---|
| 112 | if ((parseline.peek() == '#') || (parseline.str().empty())) | 
|---|
| 113 | continue; | 
|---|
| 114 | // break on empty lines | 
|---|
| 115 | if (parseline.peek() == '\n') | 
|---|
| 116 | break; | 
|---|
| 117 |  | 
|---|
| 118 | // parse line with values | 
|---|
| 119 | std::vector<double> line_of_values; | 
|---|
| 120 | do { | 
|---|
| 121 | double value; | 
|---|
| 122 | parseline >> value >> ws; | 
|---|
| 123 | line_of_values.push_back(value); | 
|---|
| 124 | } while (parseline.good()); | 
|---|
| 125 |  | 
|---|
| 126 | // check number of columns parsed | 
|---|
| 127 | ASSERT(line_of_values.size() == columns, | 
|---|
| 128 | "MatrixContent::MatrixContent() - row " | 
|---|
| 129 | +toString(row)+" has a different number of columns " | 
|---|
| 130 | +toString(line_of_values.size())+" than others before " | 
|---|
| 131 | +toString(columns)+"."); | 
|---|
| 132 |  | 
|---|
| 133 | for (size_t column = 0; column < columns; ++column) | 
|---|
| 134 | set(row, column, line_of_values[column]); | 
|---|
| 135 | ++row; | 
|---|
| 136 | } while (inputstream.good()); | 
|---|
| 137 | // check number of rows parsed | 
|---|
| 138 | ASSERT(row == rows, | 
|---|
| 139 | "MatrixContent::MatrixContent() - insufficent number of rows " | 
|---|
| 140 | +toString(row)+" < "+toString(rows)+" read."); | 
|---|
| 141 | } | 
|---|
| 142 |  | 
|---|
| 143 |  | 
|---|
| 144 | /** Constructor for class MatrixContent. | 
|---|
| 145 | * | 
|---|
| 146 | * \param *src source gsl_matrix vector to copy and put into this class | 
|---|
| 147 | */ | 
|---|
| 148 | MatrixContent::MatrixContent(gsl_matrix *&src) : | 
|---|
| 149 | rows(src->size1), | 
|---|
| 150 | columns(src->size2), | 
|---|
| 151 | free_content_on_exit(true) | 
|---|
| 152 | { | 
|---|
| 153 | content = gsl_matrix_alloc(src->size1, src->size2); | 
|---|
| 154 | gsl_matrix_memcpy(content,src); | 
|---|
| 155 | //  content = src; | 
|---|
| 156 | //  src = NULL; | 
|---|
| 157 | } | 
|---|
| 158 |  | 
|---|
| 159 | /** Copy constructor for class MatrixContent. | 
|---|
| 160 | * \param &src reference to source MatrixContent | 
|---|
| 161 | */ | 
|---|
| 162 | MatrixContent::MatrixContent(const MatrixContent &src) : | 
|---|
| 163 | rows(src.rows), | 
|---|
| 164 | columns(src.columns), | 
|---|
| 165 | free_content_on_exit(true) | 
|---|
| 166 | { | 
|---|
| 167 | content = gsl_matrix_alloc(src.rows, src.columns); | 
|---|
| 168 | gsl_matrix_memcpy(content,src.content); | 
|---|
| 169 | } | 
|---|
| 170 |  | 
|---|
| 171 | /** Copy constructor for class MatrixContent. | 
|---|
| 172 | * \param *src pointer to source MatrixContent | 
|---|
| 173 | */ | 
|---|
| 174 | MatrixContent::MatrixContent(const MatrixContent *src) : | 
|---|
| 175 | rows(src->rows), | 
|---|
| 176 | columns(src->columns), | 
|---|
| 177 | free_content_on_exit(true) | 
|---|
| 178 | { | 
|---|
| 179 | ASSERT(src != NULL, "MatrixContent::MatrixContent - pointer to source matrix is NULL!"); | 
|---|
| 180 | content = gsl_matrix_alloc(src->rows, src->columns); | 
|---|
| 181 | gsl_matrix_memcpy(content,src->content); | 
|---|
| 182 | } | 
|---|
| 183 |  | 
|---|
| 184 | /** Destructor for class MatrixContent. | 
|---|
| 185 | */ | 
|---|
| 186 | MatrixContent::~MatrixContent() | 
|---|
| 187 | { | 
|---|
| 188 | if (free_content_on_exit) | 
|---|
| 189 | gsl_matrix_free(content); | 
|---|
| 190 | } | 
|---|
| 191 |  | 
|---|
| 192 | /** Getter for MatrixContent::rows. | 
|---|
| 193 | * \return MatrixContent::rows | 
|---|
| 194 | */ | 
|---|
| 195 | const size_t MatrixContent::getRows() const | 
|---|
| 196 | { | 
|---|
| 197 | return rows; | 
|---|
| 198 | } | 
|---|
| 199 |  | 
|---|
| 200 | /** Getter for MatrixContent::columns. | 
|---|
| 201 | * \return MatrixContent::columns | 
|---|
| 202 | */ | 
|---|
| 203 | const size_t MatrixContent::getColumns() const | 
|---|
| 204 | { | 
|---|
| 205 | return columns; | 
|---|
| 206 | } | 
|---|
| 207 |  | 
|---|
| 208 | /** Return a VectorViewContent of the \a column -th column vector. | 
|---|
| 209 | * | 
|---|
| 210 | * @param column index of column | 
|---|
| 211 | * @return column of matrix as VectorContent | 
|---|
| 212 | */ | 
|---|
| 213 | VectorContent *MatrixContent::getColumnVector(size_t column) const | 
|---|
| 214 | { | 
|---|
| 215 | ASSERT(column < columns, | 
|---|
| 216 | "MatrixContent::getColumnVector() - requested column "+toString(column) | 
|---|
| 217 | +" greater than dimension "+toString(columns)); | 
|---|
| 218 | return (new VectorViewContent(gsl_matrix_column(content,column))); | 
|---|
| 219 | } | 
|---|
| 220 |  | 
|---|
| 221 | /** Returns a VectorViewContent of the \a row -th row vector. | 
|---|
| 222 | * @param row row index | 
|---|
| 223 | * @return VectorContent of row vector | 
|---|
| 224 | */ | 
|---|
| 225 | VectorContent *MatrixContent::getRowVector(size_t row) const | 
|---|
| 226 | { | 
|---|
| 227 | ASSERT(row < rows, | 
|---|
| 228 | "MatrixContent::getColumnVector() - requested row "+toString(row) | 
|---|
| 229 | +" greater than dimension "+toString(rows)); | 
|---|
| 230 | return (new VectorViewContent(gsl_matrix_row(content,row))); | 
|---|
| 231 | } | 
|---|
| 232 |  | 
|---|
| 233 | /** Returns the main diagonal of the matrix as VectorContent. | 
|---|
| 234 | * @return diagonal as VectorContent. | 
|---|
| 235 | */ | 
|---|
| 236 | VectorContent *MatrixContent::getDiagonalVector() const | 
|---|
| 237 | { | 
|---|
| 238 | return (new VectorViewContent(gsl_matrix_diagonal(content))); | 
|---|
| 239 | } | 
|---|
| 240 |  | 
|---|
| 241 | /** Set matrix to identity. | 
|---|
| 242 | */ | 
|---|
| 243 | void MatrixContent::setIdentity() | 
|---|
| 244 | { | 
|---|
| 245 | for(int i=rows;i--;){ | 
|---|
| 246 | for(int j=columns;j--;){ | 
|---|
| 247 | set(i,j,(double)(i==j)); | 
|---|
| 248 | } | 
|---|
| 249 | } | 
|---|
| 250 | } | 
|---|
| 251 |  | 
|---|
| 252 | /** Set all matrix components to zero. | 
|---|
| 253 | */ | 
|---|
| 254 | void MatrixContent::setZero() | 
|---|
| 255 | { | 
|---|
| 256 | for(int i=rows;i--;){ | 
|---|
| 257 | for(int j=columns;j--;){ | 
|---|
| 258 | set(i,j,0.); | 
|---|
| 259 | } | 
|---|
| 260 | } | 
|---|
| 261 | } | 
|---|
| 262 |  | 
|---|
| 263 | /** Set all matrix components to a given value. | 
|---|
| 264 | * \param _value value to set each component to | 
|---|
| 265 | */ | 
|---|
| 266 | void MatrixContent::setValue(double _value) | 
|---|
| 267 | { | 
|---|
| 268 | for(int i=rows;i--;){ | 
|---|
| 269 | for(int j=columns;j--;){ | 
|---|
| 270 | set(i,j,_value); | 
|---|
| 271 | } | 
|---|
| 272 | } | 
|---|
| 273 | } | 
|---|
| 274 |  | 
|---|
| 275 | /** Copy operator for MatrixContent with self-assignment check. | 
|---|
| 276 | * \param &src matrix to compare to | 
|---|
| 277 | * \return reference to this | 
|---|
| 278 | */ | 
|---|
| 279 | MatrixContent &MatrixContent::operator=(const MatrixContent &src) | 
|---|
| 280 | { | 
|---|
| 281 | if(&src!=this){ | 
|---|
| 282 | gsl_matrix_memcpy(content,src.content); | 
|---|
| 283 | } | 
|---|
| 284 | return *this; | 
|---|
| 285 | } | 
|---|
| 286 |  | 
|---|
| 287 | /** Addition operator. | 
|---|
| 288 | * \param &rhs matrix to add | 
|---|
| 289 | * \return reference to this | 
|---|
| 290 | */ | 
|---|
| 291 | const MatrixContent &MatrixContent::operator+=(const MatrixContent &rhs) | 
|---|
| 292 | { | 
|---|
| 293 | gsl_matrix_add(content, rhs.content); | 
|---|
| 294 | return *this; | 
|---|
| 295 | } | 
|---|
| 296 |  | 
|---|
| 297 | /** Subtraction operator. | 
|---|
| 298 | * \param &rhs matrix to subtract | 
|---|
| 299 | * \return reference to this | 
|---|
| 300 | */ | 
|---|
| 301 | const MatrixContent &MatrixContent::operator-=(const MatrixContent &rhs) | 
|---|
| 302 | { | 
|---|
| 303 | gsl_matrix_sub(content, rhs.content); | 
|---|
| 304 | return *this; | 
|---|
| 305 | } | 
|---|
| 306 |  | 
|---|
| 307 | /** Multiplication operator. | 
|---|
| 308 | * Note that here matrix have to have same dimensions. | 
|---|
| 309 | * \param &rhs matrix to multiply with | 
|---|
| 310 | * \return reference to this | 
|---|
| 311 | */ | 
|---|
| 312 | const MatrixContent &MatrixContent::operator*=(const MatrixContent &rhs) | 
|---|
| 313 | { | 
|---|
| 314 | ASSERT(rhs.columns == rhs.rows, | 
|---|
| 315 | "MatrixContent::operator*=() - rhs matrix is not square: "+toString(rhs.columns)+" != "+toString(rhs.rows)+"."); | 
|---|
| 316 | ASSERT(columns == rhs.rows, | 
|---|
| 317 | "MatrixContent::operator*=() - columns dimension differ: "+toString(columns)+" != "+toString(rhs.rows)+"."); | 
|---|
| 318 | (*this) = (*this)*rhs; | 
|---|
| 319 | return *this; | 
|---|
| 320 | } | 
|---|
| 321 |  | 
|---|
| 322 | /** Multiplication with copy operator. | 
|---|
| 323 | * \param &rhs matrix to multiply with | 
|---|
| 324 | * \return reference to newly allocated MatrixContent | 
|---|
| 325 | */ | 
|---|
| 326 | const MatrixContent MatrixContent::operator*(const MatrixContent &rhs) const | 
|---|
| 327 | { | 
|---|
| 328 | ASSERT (columns == rhs.rows, | 
|---|
| 329 | "MatrixContent::operator*() - dimensions not match for matrix product (a,b)*(b,c) = (a,c):" | 
|---|
| 330 | "("+toString(rows)+","+toString(columns)+")*("+toString(rhs.rows)+","+toString(rhs.columns)+")"); | 
|---|
| 331 | gsl_matrix *res = gsl_matrix_alloc(rows, rhs.columns); | 
|---|
| 332 | gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, content, rhs.content, 0.0, res); | 
|---|
| 333 | // gsl_matrix is taken over by constructor, hence no free | 
|---|
| 334 | MatrixContent tmp(res); | 
|---|
| 335 | gsl_matrix_free(res); | 
|---|
| 336 | return tmp; | 
|---|
| 337 | } | 
|---|
| 338 |  | 
|---|
| 339 | /** Hadamard multiplication with copy operator. | 
|---|
| 340 | * The Hadamard product is component-wise matrix product. | 
|---|
| 341 | * \param &rhs matrix to hadamard-multiply with | 
|---|
| 342 | * \return reference to newly allocated MatrixContent | 
|---|
| 343 | */ | 
|---|
| 344 | const MatrixContent MatrixContent::operator&(const MatrixContent &rhs) const | 
|---|
| 345 | { | 
|---|
| 346 | ASSERT ((rows == rhs.rows) && (columns == rhs.columns), | 
|---|
| 347 | "MatrixContent::operator&() - dimensions not match for matrix product (a,b) != (b,c):" | 
|---|
| 348 | "("+toString(rows)+","+toString(columns)+") != ("+toString(rhs.rows)+","+toString(rhs.columns)+")"); | 
|---|
| 349 | gsl_matrix *res = gsl_matrix_alloc(rows, rhs.columns); | 
|---|
| 350 | for (size_t i=0;i<rows;++i) | 
|---|
| 351 | for (size_t j=0;j<columns;++j) | 
|---|
| 352 | gsl_matrix_set(res, i,j, gsl_matrix_get(content, i,j)*gsl_matrix_get(rhs.content, i,j)); | 
|---|
| 353 | // gsl_matrix is taken over by constructor, hence no free | 
|---|
| 354 | MatrixContent tmp(res); | 
|---|
| 355 | gsl_matrix_free(res); | 
|---|
| 356 | return tmp; | 
|---|
| 357 | } | 
|---|
| 358 |  | 
|---|
| 359 | /** Hadamard multiplication with copy operator. | 
|---|
| 360 | * The Hadamard product is component-wise matrix product. | 
|---|
| 361 | * Note that Hadamard product can easily be done on top of \a *this matrix. | 
|---|
| 362 | * Hence, we don't need to use the multiply and copy operator as in the case of | 
|---|
| 363 | * MatrixContent::operator*=(). | 
|---|
| 364 | * \param &rhs matrix to hadamard-multiply with | 
|---|
| 365 | * \return reference to newly allocated MatrixContent | 
|---|
| 366 | */ | 
|---|
| 367 | const MatrixContent &MatrixContent::operator&=(const MatrixContent &rhs) | 
|---|
| 368 | { | 
|---|
| 369 | ASSERT ((rows == rhs.rows) && (columns == rhs.columns), | 
|---|
| 370 | "MatrixContent::operator&() - dimensions not match for matrix product (a,b) != (b,c):" | 
|---|
| 371 | "("+toString(rows)+","+toString(columns)+") != ("+toString(rhs.rows)+","+toString(rhs.columns)+")"); | 
|---|
| 372 | for (size_t i=0;i<rows;++i) | 
|---|
| 373 | for (size_t j=0;j<columns;++j) | 
|---|
| 374 | gsl_matrix_set(content, i,j, gsl_matrix_get(content, i,j)*gsl_matrix_get(rhs.content, i,j)); | 
|---|
| 375 | return *this; | 
|---|
| 376 | } | 
|---|
| 377 |  | 
|---|
| 378 | /* ========================== Accessing =============================== */ | 
|---|
| 379 |  | 
|---|
| 380 | /** Accessor for manipulating component (i,j). | 
|---|
| 381 | * \param i row number | 
|---|
| 382 | * \param j column number | 
|---|
| 383 | * \return reference to component (i,j) | 
|---|
| 384 | */ | 
|---|
| 385 | double &MatrixContent::at(size_t i, size_t j) | 
|---|
| 386 | { | 
|---|
| 387 | ASSERT((i>=0) && (i<rows), | 
|---|
| 388 | "MatrixContent::at() - Index i="+toString(i)+" for Matrix access out of range [0,"+toString(rows)+"]"); | 
|---|
| 389 | ASSERT((j>=0) && (j<columns), | 
|---|
| 390 | "MatrixContent::at() - Index j="+toString(j)+" for Matrix access out of range [0,"+toString(columns)+"]"); | 
|---|
| 391 | return *gsl_matrix_ptr (content, i, j); | 
|---|
| 392 | } | 
|---|
| 393 |  | 
|---|
| 394 | /** Constant accessor for (value of) component (i,j). | 
|---|
| 395 | * \param i row number | 
|---|
| 396 | * \param j column number | 
|---|
| 397 | * \return const component (i,j) | 
|---|
| 398 | */ | 
|---|
| 399 | const double MatrixContent::at(size_t i, size_t j) const | 
|---|
| 400 | { | 
|---|
| 401 | ASSERT((i>=0) && (i<rows), | 
|---|
| 402 | "MatrixContent::at() - Index i="+toString(i)+" for Matrix access out of range [0,"+toString(rows)+"]"); | 
|---|
| 403 | ASSERT((j>=0) && (j<columns), | 
|---|
| 404 | "MatrixContent::at() - Index j="+toString(j)+" for Matrix access out of range [0,"+toString(columns)+"]"); | 
|---|
| 405 | return gsl_matrix_get(content, i, j); | 
|---|
| 406 | } | 
|---|
| 407 |  | 
|---|
| 408 | /** These functions return a pointer to the \a m-th element of a matrix. | 
|---|
| 409 | *  If \a m or \a n lies outside the allowed range of 0 to MatrixContent::dimension-1 then the error handler is invoked and a null pointer is returned. | 
|---|
| 410 | * \param m index | 
|---|
| 411 | * \return pointer to \a m-th element | 
|---|
| 412 | */ | 
|---|
| 413 | double *MatrixContent::Pointer(size_t m, size_t n) | 
|---|
| 414 | { | 
|---|
| 415 | return gsl_matrix_ptr (content, m, n); | 
|---|
| 416 | }; | 
|---|
| 417 |  | 
|---|
| 418 | /** These functions return a constant pointer to the \a m-th element of a matrix. | 
|---|
| 419 | *  If \a m or \a n lies outside the allowed range of 0 to MatrixContent::dimension-1 then the error handler is invoked and a null pointer is returned. | 
|---|
| 420 | * \param m index | 
|---|
| 421 | * \return const pointer to \a m-th element | 
|---|
| 422 | */ | 
|---|
| 423 | const double *MatrixContent::const_Pointer(size_t m, size_t n) const | 
|---|
| 424 | { | 
|---|
| 425 | return gsl_matrix_const_ptr (content, m, n); | 
|---|
| 426 | }; | 
|---|
| 427 |  | 
|---|
| 428 | /* ========================== Initializing =============================== */ | 
|---|
| 429 |  | 
|---|
| 430 | /** Setter for component (i,j). | 
|---|
| 431 | * \param i row numbr | 
|---|
| 432 | * \param j column numnber | 
|---|
| 433 | * \param value value to set componnt (i,j) to | 
|---|
| 434 | */ | 
|---|
| 435 | void MatrixContent::set(size_t i, size_t j, const double value) | 
|---|
| 436 | { | 
|---|
| 437 | ASSERT((i>=0) && (i<rows), | 
|---|
| 438 | "MatrixContent::set() - Index i="+toString(i)+" for Matrix access out of range [0,"+toString(rows)+"]"); | 
|---|
| 439 | ASSERT((j>=0) && (j<columns), | 
|---|
| 440 | "MatrixContent::set() - Index j="+toString(j)+" for Matrix access out of range [0,"+toString(columns)+"]"); | 
|---|
| 441 | gsl_matrix_set(content,i,j,value); | 
|---|
| 442 | } | 
|---|
| 443 |  | 
|---|
| 444 | /** This function sets the matrix from a double array. | 
|---|
| 445 | * Creates a matrix view of the array and performs a memcopy. | 
|---|
| 446 | * \param *x array of values (no dimension check is performed) | 
|---|
| 447 | */ | 
|---|
| 448 | void MatrixContent::setFromDoubleArray(double * x) | 
|---|
| 449 | { | 
|---|
| 450 | gsl_matrix_view m = gsl_matrix_view_array (x, rows, columns); | 
|---|
| 451 | gsl_matrix_memcpy (content, &m.matrix); | 
|---|
| 452 | }; | 
|---|
| 453 |  | 
|---|
| 454 | /* ====================== Exchanging elements ============================ */ | 
|---|
| 455 | /** This function exchanges the \a i-th and \a j-th row of the matrix in-place. | 
|---|
| 456 | * \param i i-th row to swap with ... | 
|---|
| 457 | * \param j ... j-th row to swap against | 
|---|
| 458 | */ | 
|---|
| 459 | bool MatrixContent::SwapRows(size_t i, size_t j) | 
|---|
| 460 | { | 
|---|
| 461 | return (gsl_matrix_swap_rows (content, i, j) == GSL_SUCCESS); | 
|---|
| 462 | }; | 
|---|
| 463 |  | 
|---|
| 464 | /** This function exchanges the \a i-th and \a j-th column of the matrix in-place. | 
|---|
| 465 | * \param i i-th column to swap with ... | 
|---|
| 466 | * \param j ... j-th column to swap against | 
|---|
| 467 | */ | 
|---|
| 468 | bool MatrixContent::SwapColumns(size_t i, size_t j) | 
|---|
| 469 | { | 
|---|
| 470 | return (gsl_matrix_swap_columns (content, i, j) == GSL_SUCCESS); | 
|---|
| 471 | }; | 
|---|
| 472 |  | 
|---|
| 473 | /** This function exchanges the \a i-th row and \a j-th column of the matrix in-place. | 
|---|
| 474 | * The matrix must be square for this operation to be possible. | 
|---|
| 475 | * \param i i-th row to swap with ... | 
|---|
| 476 | * \param j ... j-th column to swap against | 
|---|
| 477 | */ | 
|---|
| 478 | bool MatrixContent::SwapRowColumn(size_t i, size_t j) | 
|---|
| 479 | { | 
|---|
| 480 | ASSERT (rows == columns, | 
|---|
| 481 | "MatrixContent::SwapRowColumn() - The matrix must be square for swapping row against column to be possible."); | 
|---|
| 482 | return (gsl_matrix_swap_rowcol (content, i, j) == GSL_SUCCESS); | 
|---|
| 483 | }; | 
|---|
| 484 |  | 
|---|
| 485 | /** Return transposed matrix. | 
|---|
| 486 | * \return new matrix that is transposed of this. | 
|---|
| 487 | */ | 
|---|
| 488 | MatrixContent MatrixContent::transpose() const | 
|---|
| 489 | { | 
|---|
| 490 | gsl_matrix *res = gsl_matrix_alloc(columns, rows);  // column and row dimensions exchanged! | 
|---|
| 491 | gsl_matrix_transpose_memcpy(res, content); | 
|---|
| 492 | MatrixContent newContent(res); | 
|---|
| 493 | gsl_matrix_free(res); | 
|---|
| 494 | return newContent; | 
|---|
| 495 | } | 
|---|
| 496 |  | 
|---|
| 497 | /** Turn this matrix into its transposed. | 
|---|
| 498 | * Note that this is only possible if rows == columns. | 
|---|
| 499 | */ | 
|---|
| 500 | MatrixContent &MatrixContent::transpose() | 
|---|
| 501 | { | 
|---|
| 502 | ASSERT( rows == columns, | 
|---|
| 503 | "MatrixContent::transpose() - cannot transpose onto itself as matrix not square: "+toString(rows)+"!="+toString(columns)+"!"); | 
|---|
| 504 | double tmp; | 
|---|
| 505 | for (size_t i=0;i<rows;i++) | 
|---|
| 506 | for (size_t j=i+1;j<rows;j++) { | 
|---|
| 507 | tmp = at(j,i); | 
|---|
| 508 | at(j,i) = at(i,j); | 
|---|
| 509 | at(i,j) = tmp; | 
|---|
| 510 | } | 
|---|
| 511 | return *this; | 
|---|
| 512 | } | 
|---|
| 513 |  | 
|---|
| 514 | /** Transform the matrix to its eigenbasis and return resulting eigenvalues. | 
|---|
| 515 | * Note that we only return real-space part in case of non-symmetric matrix. | 
|---|
| 516 | * \warn return vector has to be freed'd | 
|---|
| 517 | * TODO: encapsulate return value in boost::shared_ptr or in VectorContent. | 
|---|
| 518 | * \return gsl_vector pointer to vector of eigenvalues | 
|---|
| 519 | */ | 
|---|
| 520 | gsl_vector* MatrixContent::transformToEigenbasis() | 
|---|
| 521 | { | 
|---|
| 522 | if (rows == columns) { // symmetric | 
|---|
| 523 | gsl_eigen_symmv_workspace *T = gsl_eigen_symmv_alloc(rows); | 
|---|
| 524 | gsl_vector *eval = gsl_vector_alloc(rows); | 
|---|
| 525 | gsl_matrix *evec = gsl_matrix_alloc(rows, rows); | 
|---|
| 526 | gsl_eigen_symmv(content, eval, evec, T); | 
|---|
| 527 | gsl_eigen_symmv_free(T); | 
|---|
| 528 | gsl_matrix_memcpy(content, evec); | 
|---|
| 529 | gsl_matrix_free(evec); | 
|---|
| 530 | return eval; | 
|---|
| 531 | } else { // non-symmetric | 
|---|
| 532 | // blow up gsl_matrix in content to square matrix, fill other components with zero | 
|---|
| 533 | const size_t greaterDimension = rows > columns ? rows : columns; | 
|---|
| 534 | gsl_matrix *content_square = gsl_matrix_alloc(greaterDimension, greaterDimension); | 
|---|
| 535 | for (size_t i=0; i<greaterDimension; i++) { | 
|---|
| 536 | for (size_t j=0; j<greaterDimension; j++) { | 
|---|
| 537 | const double value = ((i < rows) && (j < columns)) ? gsl_matrix_get(content,i,j) : 0.; | 
|---|
| 538 | gsl_matrix_set(content_square, i,j, value); | 
|---|
| 539 | } | 
|---|
| 540 | } | 
|---|
| 541 |  | 
|---|
| 542 | // show squared matrix by putting it into a MatrixViewContent | 
|---|
| 543 | MatrixContent *ContentSquare = new MatrixViewContent(gsl_matrix_submatrix(content_square,0,0,content_square->size1, content_square->size2)); | 
|---|
| 544 | std::cout << "The squared matrix is " << *ContentSquare << std::endl; | 
|---|
| 545 |  | 
|---|
| 546 | // solve eigenvalue problem | 
|---|
| 547 | gsl_eigen_nonsymmv_workspace *T = gsl_eigen_nonsymmv_alloc(rows); | 
|---|
| 548 | gsl_vector_complex *eval = gsl_vector_complex_alloc(greaterDimension); | 
|---|
| 549 | gsl_matrix_complex *evec = gsl_matrix_complex_alloc(greaterDimension, greaterDimension); | 
|---|
| 550 | gsl_eigen_nonsymmv(content_square, eval, evec, T); | 
|---|
| 551 | gsl_eigen_nonsymmv_free(T); | 
|---|
| 552 |  | 
|---|
| 553 | // copy eigenvectors real-parts into content_square and ... | 
|---|
| 554 | for (size_t i=0; i<greaterDimension; i++) | 
|---|
| 555 | for (size_t j=0; j<greaterDimension; j++) | 
|---|
| 556 | gsl_matrix_set(content_square, i,j, GSL_REAL(gsl_matrix_complex_get(evec,i,j))); | 
|---|
| 557 |  | 
|---|
| 558 | // ... show complex-valued eigenvector matrix | 
|---|
| 559 | std::cout << "The real-value eigenvector matrix is " << *ContentSquare << std::endl; | 
|---|
| 560 | //    std::cout << "Resulting eigenvector matrix is ["; | 
|---|
| 561 | //    for (size_t i=0; i<greaterDimension; i++) { | 
|---|
| 562 | //      for (size_t j=0; j<greaterDimension; j++) { | 
|---|
| 563 | //        std::cout << "(" << GSL_REAL(gsl_matrix_complex_get(evec,i,j)) | 
|---|
| 564 | //            << "," << GSL_IMAG(gsl_matrix_complex_get(evec,i,j)) << ")"; | 
|---|
| 565 | //        if (j < greaterDimension-1) | 
|---|
| 566 | //          std::cout << " "; | 
|---|
| 567 | //      } | 
|---|
| 568 | //      if (i < greaterDimension-1) | 
|---|
| 569 | //        std::cout << "; "; | 
|---|
| 570 | //    } | 
|---|
| 571 | //    std::cout << "]" << std::endl; | 
|---|
| 572 |  | 
|---|
| 573 | // copy real-parts of complex eigenvalues and eigenvectors (column-wise orientation) | 
|---|
| 574 | gsl_vector *eval_real = gsl_vector_alloc(columns); | 
|---|
| 575 | size_t I=0; | 
|---|
| 576 | for (size_t i=0; i<greaterDimension; i++) { // only copy real space part | 
|---|
| 577 | if (fabs(GSL_REAL(gsl_vector_complex_get(eval,i))) > LINALG_MYEPSILON) { // only take eigenvectors with value > 0 | 
|---|
| 578 | std::cout << i << "th eigenvalue is (" << GSL_REAL(gsl_vector_complex_get(eval,i)) << "," << GSL_IMAG(gsl_vector_complex_get(eval,i)) << ")" << std::endl; | 
|---|
| 579 | for (size_t j=0; j<greaterDimension; j++) { | 
|---|
| 580 | if (fabs(GSL_IMAG(gsl_matrix_complex_get(evec,j,i))) > LINALG_MYEPSILON) | 
|---|
| 581 | std::cerr << "MatrixContent::transformToEigenbasis() - WARNING: eigenvectors are complex-valued!" << std::endl; | 
|---|
| 582 | gsl_matrix_set(content, j,I, GSL_REAL(gsl_matrix_complex_get(evec,j,i))); | 
|---|
| 583 | } | 
|---|
| 584 | if (fabs(GSL_IMAG(gsl_vector_complex_get(eval,I))) > LINALG_MYEPSILON) | 
|---|
| 585 | std::cerr << "MatrixContent::transformToEigenbasis() - WARNING: eigenvectors are complex-valued!" << std::endl; | 
|---|
| 586 | gsl_vector_set(eval_real, I, GSL_REAL(gsl_vector_complex_get(eval, i))); | 
|---|
| 587 | I++; | 
|---|
| 588 | } | 
|---|
| 589 | } | 
|---|
| 590 | gsl_matrix_complex_free(evec); | 
|---|
| 591 | gsl_vector_complex_free(eval); | 
|---|
| 592 | delete ContentSquare; | 
|---|
| 593 |  | 
|---|
| 594 | return eval_real; | 
|---|
| 595 | } | 
|---|
| 596 | } | 
|---|
| 597 |  | 
|---|
| 598 |  | 
|---|
| 599 | /** Sorts the eigenpairs in ascending order of the eigenvalues. | 
|---|
| 600 | * We assume that MatrixContent::transformToEigenbasis() has just been called. | 
|---|
| 601 | * @param eigenvalues vector of eigenvalue from | 
|---|
| 602 | *        MatrixContent::transformToEigenbasis() | 
|---|
| 603 | */ | 
|---|
| 604 | void MatrixContent::sortEigenbasis(gsl_vector *eigenvalues) | 
|---|
| 605 | { | 
|---|
| 606 | gsl_eigen_symmv_sort (eigenvalues, content, | 
|---|
| 607 | GSL_EIGEN_SORT_VAL_ASC); | 
|---|
| 608 | } | 
|---|
| 609 |  | 
|---|
| 610 | /* ============================ Properties ============================== */ | 
|---|
| 611 | /** Checks whether matrix' elements are strictly null. | 
|---|
| 612 | * \return true - is null, false - else | 
|---|
| 613 | */ | 
|---|
| 614 | bool MatrixContent::IsNull() const | 
|---|
| 615 | { | 
|---|
| 616 | return gsl_matrix_isnull (content); | 
|---|
| 617 | }; | 
|---|
| 618 |  | 
|---|
| 619 | /** Checks whether matrix' elements are strictly positive. | 
|---|
| 620 | * \return true - is positive, false - else | 
|---|
| 621 | */ | 
|---|
| 622 | bool MatrixContent::IsPositive() const | 
|---|
| 623 | { | 
|---|
| 624 | return gsl_matrix_ispos (content); | 
|---|
| 625 | }; | 
|---|
| 626 |  | 
|---|
| 627 | /** Checks whether matrix' elements are strictly negative. | 
|---|
| 628 | * \return true - is negative, false - else | 
|---|
| 629 | */ | 
|---|
| 630 | bool MatrixContent::IsNegative() const | 
|---|
| 631 | { | 
|---|
| 632 | return gsl_matrix_isneg (content); | 
|---|
| 633 | }; | 
|---|
| 634 |  | 
|---|
| 635 | /** Checks whether matrix' elements are strictly non-negative. | 
|---|
| 636 | * \return true - is non-negative, false - else | 
|---|
| 637 | */ | 
|---|
| 638 | bool MatrixContent::IsNonNegative() const | 
|---|
| 639 | { | 
|---|
| 640 | return gsl_matrix_isnonneg (content); | 
|---|
| 641 | }; | 
|---|
| 642 |  | 
|---|
| 643 | /** This function performs a Cholesky decomposition to determine whether matrix is positive definite. | 
|---|
| 644 | * We check whether GSL returns GSL_EDOM as error, indicating that decomposition failed due to matrix not being positive-definite. | 
|---|
| 645 | * \return true - matrix is positive-definite, false - else | 
|---|
| 646 | */ | 
|---|
| 647 | bool MatrixContent::IsPositiveDefinite() const | 
|---|
| 648 | { | 
|---|
| 649 | if (rows != columns)  // only possible for square matrices. | 
|---|
| 650 | return false; | 
|---|
| 651 | else | 
|---|
| 652 | return (gsl_linalg_cholesky_decomp (content) != GSL_EDOM); | 
|---|
| 653 | }; | 
|---|
| 654 |  | 
|---|
| 655 |  | 
|---|
| 656 | /** Calculates the determinant of the matrix. | 
|---|
| 657 | * if matrix is square, uses LU decomposition. | 
|---|
| 658 | */ | 
|---|
| 659 | double MatrixContent::Determinant() const | 
|---|
| 660 | { | 
|---|
| 661 | int signum = 0; | 
|---|
| 662 | ASSERT(rows == columns, | 
|---|
| 663 | "MatrixContent::Determinant() - determinant can only be calculated for square matrices."); | 
|---|
| 664 | gsl_permutation *p = gsl_permutation_alloc(rows); | 
|---|
| 665 | gsl_linalg_LU_decomp(content, p, &signum); | 
|---|
| 666 | gsl_permutation_free(p); | 
|---|
| 667 | return gsl_linalg_LU_det(content, signum); | 
|---|
| 668 | }; | 
|---|
| 669 |  | 
|---|
| 670 | /** Preparses a matrix in an input stream for row and column count. | 
|---|
| 671 | * | 
|---|
| 672 | *  \note This does not change the get position within the stream. | 
|---|
| 673 | * | 
|---|
| 674 | * @param inputstream to parse matrix from | 
|---|
| 675 | * @return pair with (row, column) | 
|---|
| 676 | */ | 
|---|
| 677 | std::pair<size_t, size_t> MatrixContent::preparseMatrixDimensions(std::istream &inputstream) | 
|---|
| 678 | { | 
|---|
| 679 | const std::streampos initial_position(inputstream.tellg()); | 
|---|
| 680 | //  std::cout << "We are at position " | 
|---|
| 681 | //      +toString((int)inputstream.tellg())+" from start of stream." << std::endl; | 
|---|
| 682 | size_t rows = 0; | 
|---|
| 683 | size_t columns = 0; | 
|---|
| 684 | do { | 
|---|
| 685 | std::string line; | 
|---|
| 686 | getline(inputstream, line); | 
|---|
| 687 | std::stringstream parseline(line); | 
|---|
| 688 | // skip comments | 
|---|
| 689 | if ((parseline.peek() == '#') || (line.empty())) | 
|---|
| 690 | continue; | 
|---|
| 691 | // break on empty lines | 
|---|
| 692 | if (parseline.peek() == '\n') | 
|---|
| 693 | break; | 
|---|
| 694 |  | 
|---|
| 695 | // parse line with values | 
|---|
| 696 | std::vector<double> line_of_values; | 
|---|
| 697 | do { | 
|---|
| 698 | double value; | 
|---|
| 699 | parseline >> value >> ws; | 
|---|
| 700 | line_of_values.push_back(value); | 
|---|
| 701 | } while (parseline.good()); | 
|---|
| 702 |  | 
|---|
| 703 | // set matrixdimension if first line | 
|---|
| 704 | if (columns == 0) { | 
|---|
| 705 | columns = line_of_values.size(); | 
|---|
| 706 | } else { // otherwise check | 
|---|
| 707 | ASSERT(columns == line_of_values.size(), | 
|---|
| 708 | "MatrixContent::MatrixContent() - row " | 
|---|
| 709 | +toString(rows)+" has a different number of columns " | 
|---|
| 710 | +toString(line_of_values.size())+" than this matrix has " | 
|---|
| 711 | +toString(columns)+"."); | 
|---|
| 712 | } | 
|---|
| 713 | ++rows; | 
|---|
| 714 | } while (inputstream.good()); | 
|---|
| 715 | // clear end of file flags | 
|---|
| 716 | inputstream.clear(); | 
|---|
| 717 |  | 
|---|
| 718 | // reset to initial position | 
|---|
| 719 | inputstream.seekg(initial_position); | 
|---|
| 720 | //  std::cout << "We are again at position " | 
|---|
| 721 | //      +toString((int)inputstream.tellg())+" from start of stream." << std::endl; | 
|---|
| 722 |  | 
|---|
| 723 | return make_pair(rows, columns); | 
|---|
| 724 | } | 
|---|
| 725 |  | 
|---|
| 726 | /** Writes matrix content to ostream in such a manner that it can be re-parsed | 
|---|
| 727 | * by the code in the constructor. | 
|---|
| 728 | * | 
|---|
| 729 | * @param ost stream to write to | 
|---|
| 730 | */ | 
|---|
| 731 | void MatrixContent::write(std::ostream &ost) const | 
|---|
| 732 | { | 
|---|
| 733 | for (size_t row = 0; row < rows; ++row) { | 
|---|
| 734 | for (size_t column = 0; column < columns; ++column) | 
|---|
| 735 | ost << at(row, column) << "\t"; | 
|---|
| 736 | ost << std::endl; | 
|---|
| 737 | } | 
|---|
| 738 | } | 
|---|
| 739 |  | 
|---|
| 740 |  | 
|---|
| 741 | /* ============================= Operators =============================== */ | 
|---|
| 742 |  | 
|---|
| 743 | /** Scalar multiplication operator. | 
|---|
| 744 | * \param factor factor to scale with | 
|---|
| 745 | */ | 
|---|
| 746 | const MatrixContent &MatrixContent::operator*=(const double factor) | 
|---|
| 747 | { | 
|---|
| 748 | gsl_matrix_scale(content, factor); | 
|---|
| 749 | return *this; | 
|---|
| 750 | } | 
|---|
| 751 |  | 
|---|
| 752 | /** Scalar multiplication and copy operator. | 
|---|
| 753 | * \param factor factor to scale with | 
|---|
| 754 | * \param &mat MatrixContent to scale | 
|---|
| 755 | * \return copied and scaled MatrixContent | 
|---|
| 756 | */ | 
|---|
| 757 | const MatrixContent operator*(const double factor,const MatrixContent& mat) | 
|---|
| 758 | { | 
|---|
| 759 | MatrixContent tmp = mat; | 
|---|
| 760 | tmp*=factor; | 
|---|
| 761 | return tmp; | 
|---|
| 762 | } | 
|---|
| 763 |  | 
|---|
| 764 | /** Scalar multiplication and copy operator (with operands exchanged). | 
|---|
| 765 | * \param &mat MatrixContent to scale | 
|---|
| 766 | * \param factor factor to scale with | 
|---|
| 767 | * \return copied and scaled MatrixContent | 
|---|
| 768 | */ | 
|---|
| 769 | const MatrixContent operator*(const MatrixContent &mat,const double factor) | 
|---|
| 770 | { | 
|---|
| 771 | return factor*mat; | 
|---|
| 772 | } | 
|---|
| 773 |  | 
|---|
| 774 | /** Sums two MatrixContents \a  and \b component-wise. | 
|---|
| 775 | * \param a first MatrixContent | 
|---|
| 776 | * \param b second MatrixContent | 
|---|
| 777 | * \return a + b | 
|---|
| 778 | */ | 
|---|
| 779 | MatrixContent const operator+(const MatrixContent& a, const MatrixContent& b) | 
|---|
| 780 | { | 
|---|
| 781 | ASSERT(a.rows == b.rows, "MatrixContent::operator+() - row counts have to match: " | 
|---|
| 782 | +toString(a.rows)+" != "+toString(b.rows)+"!"); | 
|---|
| 783 | ASSERT(a.columns == b.columns, "MatrixContent::operator+() - column counts have to match: " | 
|---|
| 784 | +toString(a.columns)+" != "+toString(b.columns)+"!"); | 
|---|
| 785 | MatrixContent x(a); | 
|---|
| 786 | x += b; | 
|---|
| 787 | return x; | 
|---|
| 788 | }; | 
|---|
| 789 |  | 
|---|
| 790 | /** Subtracts MatrixContent \a from \b component-wise. | 
|---|
| 791 | * \param a first MatrixContent | 
|---|
| 792 | * \param b second MatrixContent | 
|---|
| 793 | * \return a - b | 
|---|
| 794 | */ | 
|---|
| 795 | MatrixContent const operator-(const MatrixContent& a, const MatrixContent& b) | 
|---|
| 796 | { | 
|---|
| 797 | ASSERT(a.rows == b.rows, "MatrixContent::operator+() - row counts have to match: " | 
|---|
| 798 | +toString(a.rows)+" != "+toString(b.rows)+"!"); | 
|---|
| 799 | ASSERT(a.columns == b.columns, "MatrixContent::operator+() - column counts have to match: " | 
|---|
| 800 | +toString(a.columns)+" != "+toString(b.columns)+"!"); | 
|---|
| 801 | MatrixContent x(a); | 
|---|
| 802 | x -= b; | 
|---|
| 803 | return x; | 
|---|
| 804 | }; | 
|---|
| 805 |  | 
|---|
| 806 | /** Equality operator. | 
|---|
| 807 | * Note that we use numerical sensible checking, i.e. with threshold LINALG_MYEPSILON. | 
|---|
| 808 | * \param &rhs MatrixContent to checks against | 
|---|
| 809 | */ | 
|---|
| 810 | bool MatrixContent::operator==(const MatrixContent &rhs) const | 
|---|
| 811 | { | 
|---|
| 812 | if ((rows == rhs.rows) && (columns == rhs.columns)) { | 
|---|
| 813 | for(int i=rows;i--;){ | 
|---|
| 814 | for(int j=columns;j--;){ | 
|---|
| 815 | if(fabs(at(i,j)-rhs.at(i,j))>LINALG_MYEPSILON){ | 
|---|
| 816 | return false; | 
|---|
| 817 | } | 
|---|
| 818 | } | 
|---|
| 819 | } | 
|---|
| 820 | return true; | 
|---|
| 821 | } | 
|---|
| 822 | return false; | 
|---|
| 823 | } | 
|---|
| 824 |  | 
|---|
| 825 |  | 
|---|
| 826 | std::ostream & operator<<(std::ostream &ost, const MatrixContent &mat) | 
|---|
| 827 | { | 
|---|
| 828 | ost << "\\begin{pmatrix}"; | 
|---|
| 829 | for (size_t i=0;i<mat.rows;i++) { | 
|---|
| 830 | for (size_t j=0;j<mat.columns;j++) { | 
|---|
| 831 | ost << mat.at(i,j) << " "; | 
|---|
| 832 | if (j != mat.columns-1) | 
|---|
| 833 | ost << "& "; | 
|---|
| 834 | } | 
|---|
| 835 | if (i != mat.rows-1) | 
|---|
| 836 | ost << "\\\\ "; | 
|---|
| 837 | } | 
|---|
| 838 | ost << "\\end{pmatrix}"; | 
|---|
| 839 | return ost; | 
|---|
| 840 | } | 
|---|