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