source: src/LinearAlgebra/MatrixContent.cpp@ 0e4aa4

Action_Thermostats Add_AtomRandomPerturbation Add_FitFragmentPartialChargesAction Add_RotateAroundBondAction Add_SelectAtomByNameAction Added_ParseSaveFragmentResults AddingActions_SaveParseParticleParameters Adding_Graph_to_ChangeBondActions Adding_MD_integration_tests Adding_ParticleName_to_Atom Adding_StructOpt_integration_tests AtomFragments Automaking_mpqc_open AutomationFragmentation_failures Candidate_v1.5.4 Candidate_v1.6.0 Candidate_v1.6.1 ChangeBugEmailaddress ChangingTestPorts ChemicalSpaceEvaluator CombiningParticlePotentialParsing Combining_Subpackages Debian_Package_split Debian_package_split_molecuildergui_only Disabling_MemDebug Docu_Python_wait EmpiricalPotential_contain_HomologyGraph EmpiricalPotential_contain_HomologyGraph_documentation Enable_parallel_make_install Enhance_userguide Enhanced_StructuralOptimization Enhanced_StructuralOptimization_continued Example_ManyWaysToTranslateAtom Exclude_Hydrogens_annealWithBondGraph FitPartialCharges_GlobalError Fix_BoundInBox_CenterInBox_MoleculeActions Fix_ChargeSampling_PBC Fix_ChronosMutex Fix_FitPartialCharges Fix_FitPotential_needs_atomicnumbers Fix_ForceAnnealing Fix_IndependentFragmentGrids Fix_ParseParticles Fix_ParseParticles_split_forward_backward_Actions Fix_PopActions Fix_QtFragmentList_sorted_selection Fix_Restrictedkeyset_FragmentMolecule Fix_StatusMsg Fix_StepWorldTime_single_argument Fix_Verbose_Codepatterns Fix_fitting_potentials Fixes ForceAnnealing_goodresults ForceAnnealing_oldresults ForceAnnealing_tocheck ForceAnnealing_with_BondGraph ForceAnnealing_with_BondGraph_continued ForceAnnealing_with_BondGraph_continued_betteresults ForceAnnealing_with_BondGraph_contraction-expansion FragmentAction_writes_AtomFragments FragmentMolecule_checks_bonddegrees GeometryObjects Gui_Fixes Gui_displays_atomic_force_velocity ImplicitCharges IndependentFragmentGrids IndependentFragmentGrids_IndividualZeroInstances IndependentFragmentGrids_IntegrationTest IndependentFragmentGrids_Sole_NN_Calculation JobMarket_RobustOnKillsSegFaults JobMarket_StableWorkerPool JobMarket_unresolvable_hostname_fix MoreRobust_FragmentAutomation ODR_violation_mpqc_open PartialCharges_OrthogonalSummation PdbParser_setsAtomName PythonUI_with_named_parameters QtGui_reactivate_TimeChanged_changes Recreated_GuiChecks Rewrite_FitPartialCharges RotateToPrincipalAxisSystem_UndoRedo SaturateAtoms_findBestMatching SaturateAtoms_singleDegree StoppableMakroAction Subpackage_CodePatterns Subpackage_JobMarket Subpackage_LinearAlgebra Subpackage_levmar Subpackage_mpqc_open Subpackage_vmg Switchable_LogView ThirdParty_MPQC_rebuilt_buildsystem TrajectoryDependenant_MaxOrder TremoloParser_IncreasedPrecision TremoloParser_MultipleTimesteps TremoloParser_setsAtomName Ubuntu_1604_changes stable
Last change on this file since 0e4aa4 was 0e4aa4, checked in by Frederik Heber <heber@…>, 14 years ago

FIX: MatrixContent::sortEigenbasis() sorted by magnitude not value.

  • Property mode set to 100644
File size: 22.5 KB
Line 
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
37using namespace std;
38
39
40/** Constructor for class MatrixContent.
41 * \param rows number of rows
42 * \param columns number of columns
43 */
44MatrixContent::MatrixContent(size_t _rows, size_t _columns) :
45 rows(_rows),
46 columns(_columns)
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 */
58MatrixContent::MatrixContent(size_t _rows, size_t _columns, MatrixBaseCase) :
59 rows(_rows),
60 columns(_columns)
61{}
62
63/** Constructor for class MatrixContent.
64 * \param rows number of rows
65 * \param columns number of columns
66 * \param *src array with components to initialize matrix with
67 */
68MatrixContent::MatrixContent(size_t _rows, size_t _columns, const double *src) :
69 rows(_rows),
70 columns(_columns)
71{
72 content = gsl_matrix_calloc(rows, columns);
73 set(0,0, src[0]);
74 set(1,0, src[1]);
75 set(2,0, src[2]);
76
77 set(0,1, src[3]);
78 set(1,1, src[4]);
79 set(2,1, src[5]);
80
81 set(0,2, src[6]);
82 set(1,2, src[7]);
83 set(2,2, src[8]);
84}
85
86/** Constructor for class MatrixContent.
87 * We embed the given gls_matrix pointer within this class and set it to NULL
88 * afterwards.
89 * \param *src source gsl_matrix vector to embed within this class
90 */
91MatrixContent::MatrixContent(gsl_matrix *&src) :
92 rows(src->size1),
93 columns(src->size2)
94{
95 content = gsl_matrix_alloc(src->size1, src->size2);
96 gsl_matrix_memcpy(content,src);
97// content = src;
98// src = NULL;
99}
100
101/** Copy constructor for class MatrixContent.
102 * \param &src reference to source MatrixContent
103 */
104MatrixContent::MatrixContent(const MatrixContent &src) :
105 rows(src.rows),
106 columns(src.columns)
107{
108 content = gsl_matrix_alloc(src.rows, src.columns);
109 gsl_matrix_memcpy(content,src.content);
110}
111
112/** Copy constructor for class MatrixContent.
113 * \param *src pointer to source MatrixContent
114 */
115MatrixContent::MatrixContent(const MatrixContent *src) :
116 rows(src->rows),
117 columns(src->columns)
118{
119 ASSERT(src != NULL, "MatrixContent::MatrixContent - pointer to source matrix is NULL!");
120 content = gsl_matrix_alloc(src->rows, src->columns);
121 gsl_matrix_memcpy(content,src->content);
122}
123
124/** Destructor for class MatrixContent.
125 */
126MatrixContent::~MatrixContent()
127{
128 gsl_matrix_free(content);
129}
130
131/** Getter for MatrixContent::rows.
132 * \return MatrixContent::rows
133 */
134const size_t MatrixContent::getRows() const
135{
136 return rows;
137}
138
139/** Getter for MatrixContent::columns.
140 * \return MatrixContent::columns
141 */
142const size_t MatrixContent::getColumns() const
143{
144 return columns;
145}
146
147/** Return a VectorViewContent of the \a column -th column vector.
148 *
149 * @param column index of column
150 * @return column of matrix as VectorContent
151 */
152VectorContent *MatrixContent::getColumnVector(size_t column) const
153{
154 ASSERT(column < columns,
155 "MatrixContent::getColumnVector() - requested column "+toString(column)
156 +" greater than dimension "+toString(columns));
157 return (new VectorViewContent(gsl_matrix_column(content,column)));
158}
159
160/** Returns a VectorViewContent of the \a row -th row vector.
161 * @param row row index
162 * @return VectorContent of row vector
163 */
164VectorContent *MatrixContent::getRowVector(size_t row) const
165{
166 ASSERT(row < rows,
167 "MatrixContent::getColumnVector() - requested row "+toString(row)
168 +" greater than dimension "+toString(rows));
169 return (new VectorViewContent(gsl_matrix_row(content,row)));
170}
171
172/** Returns the main diagonal of the matrix as VectorContent.
173 * @return diagonal as VectorContent.
174 */
175VectorContent *MatrixContent::getDiagonalVector() const
176{
177 return (new VectorViewContent(gsl_matrix_diagonal(content)));
178}
179
180/** Set matrix to identity.
181 */
182void MatrixContent::setIdentity()
183{
184 for(int i=rows;i--;){
185 for(int j=columns;j--;){
186 set(i,j,(double)(i==j));
187 }
188 }
189}
190
191/** Set all matrix components to zero.
192 */
193void MatrixContent::setZero()
194{
195 for(int i=rows;i--;){
196 for(int j=columns;j--;){
197 set(i,j,0.);
198 }
199 }
200}
201
202/** Set all matrix components to a given value.
203 * \param _value value to set each component to
204 */
205void MatrixContent::setValue(double _value)
206{
207 for(int i=rows;i--;){
208 for(int j=columns;j--;){
209 set(i,j,_value);
210 }
211 }
212}
213
214/** Copy operator for MatrixContent with self-assignment check.
215 * \param &src matrix to compare to
216 * \return reference to this
217 */
218MatrixContent &MatrixContent::operator=(const MatrixContent &src)
219{
220 if(&src!=this){
221 gsl_matrix_memcpy(content,src.content);
222 }
223 return *this;
224}
225
226/** Addition operator.
227 * \param &rhs matrix to add
228 * \return reference to this
229 */
230const MatrixContent &MatrixContent::operator+=(const MatrixContent &rhs)
231{
232 gsl_matrix_add(content, rhs.content);
233 return *this;
234}
235
236/** Subtraction operator.
237 * \param &rhs matrix to subtract
238 * \return reference to this
239 */
240const MatrixContent &MatrixContent::operator-=(const MatrixContent &rhs)
241 {
242 gsl_matrix_sub(content, rhs.content);
243 return *this;
244}
245
246/** Multiplication operator.
247 * Note that here matrix have to have same dimensions.
248 * \param &rhs matrix to multiply with
249 * \return reference to this
250 */
251const MatrixContent &MatrixContent::operator*=(const MatrixContent &rhs)
252{
253 ASSERT(rhs.columns == rhs.rows,
254 "MatrixContent::operator*=() - rhs matrix is not square: "+toString(rhs.columns)+" != "+toString(rhs.rows)+".");
255 ASSERT(columns == rhs.rows,
256 "MatrixContent::operator*=() - columns dimension differ: "+toString(columns)+" != "+toString(rhs.rows)+".");
257 (*this) = (*this)*rhs;
258 return *this;
259}
260
261/** Multiplication with copy operator.
262 * \param &rhs matrix to multiply with
263 * \return reference to newly allocated MatrixContent
264 */
265const MatrixContent MatrixContent::operator*(const MatrixContent &rhs) const
266{
267 ASSERT (columns == rhs.rows,
268 "MatrixContent::operator*() - dimensions not match for matrix product (a,b)*(b,c) = (a,c):"
269 "("+toString(rows)+","+toString(columns)+")*("+toString(rhs.rows)+","+toString(rhs.columns)+")");
270 gsl_matrix *res = gsl_matrix_alloc(rows, rhs.columns);
271 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, content, rhs.content, 0.0, res);
272 // gsl_matrix is taken over by constructor, hence no free
273 MatrixContent tmp(res);
274 gsl_matrix_free(res);
275 return tmp;
276}
277
278/** Hadamard multiplication with copy operator.
279 * The Hadamard product is component-wise matrix product.
280 * \param &rhs matrix to hadamard-multiply with
281 * \return reference to newly allocated MatrixContent
282 */
283const MatrixContent MatrixContent::operator&(const MatrixContent &rhs) const
284{
285 ASSERT ((rows == rhs.rows) && (columns == rhs.columns),
286 "MatrixContent::operator&() - dimensions not match for matrix product (a,b) != (b,c):"
287 "("+toString(rows)+","+toString(columns)+") != ("+toString(rhs.rows)+","+toString(rhs.columns)+")");
288 gsl_matrix *res = gsl_matrix_alloc(rows, rhs.columns);
289 for (size_t i=0;i<rows;++i)
290 for (size_t j=0;j<columns;++j)
291 gsl_matrix_set(res, i,j, gsl_matrix_get(content, i,j)*gsl_matrix_get(rhs.content, i,j));
292 // gsl_matrix is taken over by constructor, hence no free
293 MatrixContent tmp(res);
294 gsl_matrix_free(res);
295 return tmp;
296}
297
298/** Hadamard multiplication with copy operator.
299 * The Hadamard product is component-wise matrix product.
300 * Note that Hadamard product can easily be done on top of \a *this matrix.
301 * Hence, we don't need to use the multiply and copy operator as in the case of
302 * MatrixContent::operator*=().
303 * \param &rhs matrix to hadamard-multiply with
304 * \return reference to newly allocated MatrixContent
305 */
306const MatrixContent &MatrixContent::operator&=(const MatrixContent &rhs)
307{
308 ASSERT ((rows == rhs.rows) && (columns == rhs.columns),
309 "MatrixContent::operator&() - dimensions not match for matrix product (a,b) != (b,c):"
310 "("+toString(rows)+","+toString(columns)+") != ("+toString(rhs.rows)+","+toString(rhs.columns)+")");
311 for (size_t i=0;i<rows;++i)
312 for (size_t j=0;j<columns;++j)
313 gsl_matrix_set(content, i,j, gsl_matrix_get(content, i,j)*gsl_matrix_get(rhs.content, i,j));
314 return *this;
315}
316
317/* ========================== Accessing =============================== */
318
319/** Accessor for manipulating component (i,j).
320 * \param i row number
321 * \param j column number
322 * \return reference to component (i,j)
323 */
324double &MatrixContent::at(size_t i, size_t j)
325{
326 ASSERT((i>=0) && (i<rows),
327 "MatrixContent::at() - Index i="+toString(i)+" for Matrix access out of range [0,"+toString(rows)+"]");
328 ASSERT((j>=0) && (j<columns),
329 "MatrixContent::at() - Index j="+toString(j)+" for Matrix access out of range [0,"+toString(columns)+"]");
330 return *gsl_matrix_ptr (content, i, j);
331}
332
333/** Constant accessor for (value of) component (i,j).
334 * \param i row number
335 * \param j column number
336 * \return const component (i,j)
337 */
338const double MatrixContent::at(size_t i, size_t j) const
339{
340 ASSERT((i>=0) && (i<rows),
341 "MatrixContent::at() - Index i="+toString(i)+" for Matrix access out of range [0,"+toString(rows)+"]");
342 ASSERT((j>=0) && (j<columns),
343 "MatrixContent::at() - Index j="+toString(j)+" for Matrix access out of range [0,"+toString(columns)+"]");
344 return gsl_matrix_get(content, i, j);
345}
346
347/** These functions return a pointer to the \a m-th element of a matrix.
348 * 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.
349 * \param m index
350 * \return pointer to \a m-th element
351 */
352double *MatrixContent::Pointer(size_t m, size_t n)
353{
354 return gsl_matrix_ptr (content, m, n);
355};
356
357/** These functions return a constant pointer to the \a m-th element of a matrix.
358 * 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.
359 * \param m index
360 * \return const pointer to \a m-th element
361 */
362const double *MatrixContent::const_Pointer(size_t m, size_t n) const
363{
364 return gsl_matrix_const_ptr (content, m, n);
365};
366
367/* ========================== Initializing =============================== */
368
369/** Setter for component (i,j).
370 * \param i row numbr
371 * \param j column numnber
372 * \param value value to set componnt (i,j) to
373 */
374void MatrixContent::set(size_t i, size_t j, const double value)
375{
376 ASSERT((i>=0) && (i<rows),
377 "MatrixContent::set() - Index i="+toString(i)+" for Matrix access out of range [0,"+toString(rows)+"]");
378 ASSERT((j>=0) && (j<columns),
379 "MatrixContent::set() - Index j="+toString(j)+" for Matrix access out of range [0,"+toString(columns)+"]");
380 gsl_matrix_set(content,i,j,value);
381}
382
383/** This function sets the matrix from a double array.
384 * Creates a matrix view of the array and performs a memcopy.
385 * \param *x array of values (no dimension check is performed)
386 */
387void MatrixContent::setFromDoubleArray(double * x)
388{
389 gsl_matrix_view m = gsl_matrix_view_array (x, rows, columns);
390 gsl_matrix_memcpy (content, &m.matrix);
391};
392
393/* ====================== Exchanging elements ============================ */
394/** This function exchanges the \a i-th and \a j-th row of the matrix in-place.
395 * \param i i-th row to swap with ...
396 * \param j ... j-th row to swap against
397 */
398bool MatrixContent::SwapRows(size_t i, size_t j)
399{
400 return (gsl_matrix_swap_rows (content, i, j) == GSL_SUCCESS);
401};
402
403/** This function exchanges the \a i-th and \a j-th column of the matrix in-place.
404 * \param i i-th column to swap with ...
405 * \param j ... j-th column to swap against
406 */
407bool MatrixContent::SwapColumns(size_t i, size_t j)
408{
409 return (gsl_matrix_swap_columns (content, i, j) == GSL_SUCCESS);
410};
411
412/** This function exchanges the \a i-th row and \a j-th column of the matrix in-place.
413 * The matrix must be square for this operation to be possible.
414 * \param i i-th row to swap with ...
415 * \param j ... j-th column to swap against
416 */
417bool MatrixContent::SwapRowColumn(size_t i, size_t j)
418{
419 ASSERT (rows == columns,
420 "MatrixContent::SwapRowColumn() - The matrix must be square for swapping row against column to be possible.");
421 return (gsl_matrix_swap_rowcol (content, i, j) == GSL_SUCCESS);
422};
423
424/** Return transposed matrix.
425 * \return new matrix that is transposed of this.
426 */
427MatrixContent MatrixContent::transpose() const
428{
429 gsl_matrix *res = gsl_matrix_alloc(columns, rows); // column and row dimensions exchanged!
430 gsl_matrix_transpose_memcpy(res, content);
431 MatrixContent newContent(res);
432 gsl_matrix_free(res);
433 return newContent;
434}
435
436/** Turn this matrix into its transposed.
437 * Note that this is only possible if rows == columns.
438 */
439MatrixContent &MatrixContent::transpose()
440{
441 ASSERT( rows == columns,
442 "MatrixContent::transpose() - cannot transpose onto itself as matrix not square: "+toString(rows)+"!="+toString(columns)+"!");
443 double tmp;
444 for (size_t i=0;i<rows;i++)
445 for (size_t j=i+1;j<rows;j++) {
446 tmp = at(j,i);
447 at(j,i) = at(i,j);
448 at(i,j) = tmp;
449 }
450 return *this;
451}
452
453/** Transform the matrix to its eigenbasis and return resulting eigenvalues.
454 * Note that we only return real-space part in case of non-symmetric matrix.
455 * \warn return vector has to be freed'd
456 * TODO: encapsulate return value in boost::shared_ptr or in VectorContent.
457 * \return gsl_vector pointer to vector of eigenvalues
458 */
459gsl_vector* MatrixContent::transformToEigenbasis()
460{
461 if (rows == columns) { // symmetric
462 gsl_eigen_symmv_workspace *T = gsl_eigen_symmv_alloc(rows);
463 gsl_vector *eval = gsl_vector_alloc(rows);
464 gsl_matrix *evec = gsl_matrix_alloc(rows, rows);
465 gsl_eigen_symmv(content, eval, evec, T);
466 gsl_eigen_symmv_free(T);
467 gsl_matrix_memcpy(content, evec);
468 gsl_matrix_free(evec);
469 return eval;
470 } else { // non-symmetric
471 // blow up gsl_matrix in content to square matrix, fill other components with zero
472 const size_t greaterDimension = rows > columns ? rows : columns;
473 gsl_matrix *content_square = gsl_matrix_alloc(greaterDimension, greaterDimension);
474 for (size_t i=0; i<greaterDimension; i++) {
475 for (size_t j=0; j<greaterDimension; j++) {
476 const double value = ((i < rows) && (j < columns)) ? gsl_matrix_get(content,i,j) : 0.;
477 gsl_matrix_set(content_square, i,j, value);
478 }
479 }
480
481 // show squared matrix by putting it into a MatrixViewContent
482 MatrixContent *ContentSquare = new MatrixViewContent(gsl_matrix_submatrix(content_square,0,0,content_square->size1, content_square->size2));
483 std::cout << "The squared matrix is " << *ContentSquare << std::endl;
484
485 // solve eigenvalue problem
486 gsl_eigen_nonsymmv_workspace *T = gsl_eigen_nonsymmv_alloc(rows);
487 gsl_vector_complex *eval = gsl_vector_complex_alloc(greaterDimension);
488 gsl_matrix_complex *evec = gsl_matrix_complex_alloc(greaterDimension, greaterDimension);
489 gsl_eigen_nonsymmv(content_square, eval, evec, T);
490 gsl_eigen_nonsymmv_free(T);
491
492 // copy eigenvectors real-parts into content_square and ...
493 for (size_t i=0; i<greaterDimension; i++)
494 for (size_t j=0; j<greaterDimension; j++)
495 gsl_matrix_set(content_square, i,j, GSL_REAL(gsl_matrix_complex_get(evec,i,j)));
496
497 // ... show complex-valued eigenvector matrix
498 std::cout << "The real-value eigenvector matrix is " << *ContentSquare << std::endl;
499// std::cout << "Resulting eigenvector matrix is [";
500// for (size_t i=0; i<greaterDimension; i++) {
501// for (size_t j=0; j<greaterDimension; j++) {
502// std::cout << "(" << GSL_REAL(gsl_matrix_complex_get(evec,i,j))
503// << "," << GSL_IMAG(gsl_matrix_complex_get(evec,i,j)) << ")";
504// if (j < greaterDimension-1)
505// std::cout << " ";
506// }
507// if (i < greaterDimension-1)
508// std::cout << "; ";
509// }
510// std::cout << "]" << std::endl;
511
512 // copy real-parts of complex eigenvalues and eigenvectors (column-wise orientation)
513 gsl_vector *eval_real = gsl_vector_alloc(columns);
514 size_t I=0;
515 for (size_t i=0; i<greaterDimension; i++) { // only copy real space part
516 if (fabs(GSL_REAL(gsl_vector_complex_get(eval,i))) > LINALG_MYEPSILON) { // only take eigenvectors with value > 0
517 std::cout << i << "th eigenvalue is (" << GSL_REAL(gsl_vector_complex_get(eval,i)) << "," << GSL_IMAG(gsl_vector_complex_get(eval,i)) << ")" << std::endl;
518 for (size_t j=0; j<greaterDimension; j++) {
519 if (fabs(GSL_IMAG(gsl_matrix_complex_get(evec,j,i))) > LINALG_MYEPSILON)
520 std::cerr << "MatrixContent::transformToEigenbasis() - WARNING: eigenvectors are complex-valued!" << std::endl;
521 gsl_matrix_set(content, j,I, GSL_REAL(gsl_matrix_complex_get(evec,j,i)));
522 }
523 if (fabs(GSL_IMAG(gsl_vector_complex_get(eval,I))) > LINALG_MYEPSILON)
524 std::cerr << "MatrixContent::transformToEigenbasis() - WARNING: eigenvectors are complex-valued!" << std::endl;
525 gsl_vector_set(eval_real, I, GSL_REAL(gsl_vector_complex_get(eval, i)));
526 I++;
527 }
528 }
529 gsl_matrix_complex_free(evec);
530 gsl_vector_complex_free(eval);
531 delete ContentSquare;
532
533 return eval_real;
534 }
535}
536
537
538/** Sorts the eigenpairs in ascending order of the eigenvalues.
539 * We assume that MatrixContent::transformToEigenbasis() has just been called.
540 * @param eigenvalues vector of eigenvalue from
541 * MatrixContent::transformToEigenbasis()
542 */
543void MatrixContent::sortEigenbasis(gsl_vector *eigenvalues)
544{
545 gsl_eigen_symmv_sort (eigenvalues, content,
546 GSL_EIGEN_SORT_VAL_ASC);
547}
548
549/* ============================ Properties ============================== */
550/** Checks whether matrix' elements are strictly null.
551 * \return true - is null, false - else
552 */
553bool MatrixContent::IsNull() const
554{
555 return gsl_matrix_isnull (content);
556};
557
558/** Checks whether matrix' elements are strictly positive.
559 * \return true - is positive, false - else
560 */
561bool MatrixContent::IsPositive() const
562{
563 return gsl_matrix_ispos (content);
564};
565
566/** Checks whether matrix' elements are strictly negative.
567 * \return true - is negative, false - else
568 */
569bool MatrixContent::IsNegative() const
570{
571 return gsl_matrix_isneg (content);
572};
573
574/** Checks whether matrix' elements are strictly non-negative.
575 * \return true - is non-negative, false - else
576 */
577bool MatrixContent::IsNonNegative() const
578{
579 return gsl_matrix_isnonneg (content);
580};
581
582/** This function performs a Cholesky decomposition to determine whether matrix is positive definite.
583 * We check whether GSL returns GSL_EDOM as error, indicating that decomposition failed due to matrix not being positive-definite.
584 * \return true - matrix is positive-definite, false - else
585 */
586bool MatrixContent::IsPositiveDefinite() const
587{
588 if (rows != columns) // only possible for square matrices.
589 return false;
590 else
591 return (gsl_linalg_cholesky_decomp (content) != GSL_EDOM);
592};
593
594
595/** Calculates the determinant of the matrix.
596 * if matrix is square, uses LU decomposition.
597 */
598double MatrixContent::Determinant() const
599{
600 int signum = 0;
601 ASSERT(rows == columns,
602 "MatrixContent::Determinant() - determinant can only be calculated for square matrices.");
603 gsl_permutation *p = gsl_permutation_alloc(rows);
604 gsl_linalg_LU_decomp(content, p, &signum);
605 gsl_permutation_free(p);
606 return gsl_linalg_LU_det(content, signum);
607};
608
609/* ============================= Operators =============================== */
610
611/** Scalar multiplication operator.
612 * \param factor factor to scale with
613 */
614const MatrixContent &MatrixContent::operator*=(const double factor)
615{
616 gsl_matrix_scale(content, factor);
617 return *this;
618}
619
620/** Scalar multiplication and copy operator.
621 * \param factor factor to scale with
622 * \param &mat MatrixContent to scale
623 * \return copied and scaled MatrixContent
624 */
625const MatrixContent operator*(const double factor,const MatrixContent& mat)
626{
627 MatrixContent tmp = mat;
628 tmp*=factor;
629 return tmp;
630}
631
632/** Scalar multiplication and copy operator (with operands exchanged).
633 * \param &mat MatrixContent to scale
634 * \param factor factor to scale with
635 * \return copied and scaled MatrixContent
636 */
637const MatrixContent operator*(const MatrixContent &mat,const double factor)
638{
639 return factor*mat;
640}
641
642/** Sums two MatrixContents \a and \b component-wise.
643 * \param a first MatrixContent
644 * \param b second MatrixContent
645 * \return a + b
646 */
647MatrixContent const operator+(const MatrixContent& a, const MatrixContent& b)
648{
649 ASSERT(a.rows == b.rows, "MatrixContent::operator+() - row counts have to match: "
650 +toString(a.rows)+" != "+toString(b.rows)+"!");
651 ASSERT(a.columns == b.columns, "MatrixContent::operator+() - column counts have to match: "
652 +toString(a.columns)+" != "+toString(b.columns)+"!");
653 MatrixContent x(a);
654 x += b;
655 return x;
656};
657
658/** Subtracts MatrixContent \a from \b component-wise.
659 * \param a first MatrixContent
660 * \param b second MatrixContent
661 * \return a - b
662 */
663MatrixContent const operator-(const MatrixContent& a, const MatrixContent& b)
664{
665 ASSERT(a.rows == b.rows, "MatrixContent::operator+() - row counts have to match: "
666 +toString(a.rows)+" != "+toString(b.rows)+"!");
667 ASSERT(a.columns == b.columns, "MatrixContent::operator+() - column counts have to match: "
668 +toString(a.columns)+" != "+toString(b.columns)+"!");
669 MatrixContent x(a);
670 x -= b;
671 return x;
672};
673
674/** Equality operator.
675 * Note that we use numerical sensible checking, i.e. with threshold LINALG_MYEPSILON.
676 * \param &rhs MatrixContent to checks against
677 */
678bool MatrixContent::operator==(const MatrixContent &rhs) const
679 {
680 if ((rows == rhs.rows) && (columns == rhs.columns)) {
681 for(int i=rows;i--;){
682 for(int j=columns;j--;){
683 if(fabs(at(i,j)-rhs.at(i,j))>LINALG_MYEPSILON){
684 return false;
685 }
686 }
687 }
688 return true;
689 }
690 return false;
691}
692
693
694std::ostream & operator<<(std::ostream &ost, const MatrixContent &mat)
695{
696 ost << "\\begin{pmatrix}";
697 for (size_t i=0;i<mat.rows;i++) {
698 for (size_t j=0;j<mat.columns;j++) {
699 ost << mat.at(i,j) << " ";
700 if (j != mat.columns-1)
701 ost << "& ";
702 }
703 if (i != mat.rows-1)
704 ost << "\\\\ ";
705 }
706 ost << "\\end{pmatrix}";
707 return ost;
708}
Note: See TracBrowser for help on using the repository browser.