source: LinearAlgebra/src/LinearAlgebra/MatrixContent.hpp@ c8b17b

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 c8b17b was 26108c, checked in by Frederik Heber <heber@…>, 13 years ago

FIX: Changed ambigious interface for MatrixContent::transpose().

  • now, if we transpose in place, the function is transpose() (active voice).
  • if we transpose to other MatrixContent, it is transposed() (passive voice).
  • this also changed the interface for RealSpaceMatrix::transpose().
  • fixes in several unit tests where before we always had to explicitly cast the matrix to const in order to make use of the passive voice form.
  • Property mode set to 100644
File size: 5.8 KB
Line 
1/*
2 * MatrixContent.hpp
3 *
4 * Created on: Jul 2, 2010
5 * Author: crueger
6 */
7
8#ifndef MATRIXCONTENT_HPP_
9#define MATRIXCONTENT_HPP_
10
11// include config.h
12#ifdef HAVE_CONFIG_H
13#include <config.h>
14#endif
15
16
17/** MatrixContent is a wrapper for gsl_matrix.
18 *
19 * The GNU Scientific Library contaisn some very well written routines for
20 * linear algebra problems. However, it's syntax is C and hence it does not
21 * lend itself to readable written code, i.e. C = A * B, where A, B, and C
22 * are all matrices. Writing code this way is very convenient, concise and
23 * also in the same manner as one would type in MatLab.
24 * However, with C++ and its feature to overload functions and its operator
25 * functions such syntax is easy to create.
26 *
27 * Hence, this class is a C++ wrapper to gsl_matrix. There already some other
28 * libraries, however they fall short for the following reasons:
29 * -# gslwrap: not very well written and uses floats instead of doubles
30 * -# gslcpp: last change is from 2007 and only very few commits
31 * -# o2scl: basically, the same functionality as gsl only written in C++,
32 * however it uses GPLv3 license which is inconvenient for MoleCuilder.
33 *
34 * <h1>Howto use MatrixContent</h1>
35 *
36 * One should not use MatrixContent directly but either have it as a member
37 * variable in a specialized class or inherit from it.
38 *
39 */
40
41#include <gsl/gsl_matrix.h>
42
43#include <iosfwd>
44
45#include "MatrixVector_ops.hpp"
46
47class VectorContent;
48
49/** Dummy structure to create a unique signature.
50 *
51 */
52struct MatrixBaseCase{};
53
54class MatrixContent
55{
56 friend std::ostream & operator<<(std::ostream &ost, const MatrixContent &mat);
57 friend class RealSpaceMatrix;
58 friend class LinearSystemOfEquations;
59
60 // matrix vector products
61 friend VectorContent const operator*(const VectorContent& vec, const MatrixContent& mat);
62 friend VectorContent const operator*(const MatrixContent& mat, const VectorContent& vec);
63 friend MatrixContent const operator+(const MatrixContent& a, const MatrixContent& b);
64 friend MatrixContent const operator-(const MatrixContent& a, const MatrixContent& b);
65
66 // unit tests
67 friend class MatrixContentTest;
68 friend class MatrixContentSymmetricTest;
69public:
70 MatrixContent(size_t rows, size_t columns);
71 MatrixContent(size_t _rows, size_t _columns, MatrixBaseCase);
72 MatrixContent(size_t rows, size_t columns, const double *src);
73 MatrixContent(size_t rows, size_t columns, std::istream &inputstream);
74 MatrixContent(gsl_matrix *&src);
75 MatrixContent(const MatrixContent &src);
76 MatrixContent(const MatrixContent *src);
77 virtual ~MatrixContent();
78
79 // Accessing
80 double &at(size_t i, size_t j);
81 const double at(size_t i, size_t j) const;
82 void set(size_t i, size_t j, const double value);
83
84 // Initializing
85 void setFromDoubleArray(double * x);
86 void setIdentity();
87 void setValue(double _value);
88 void setZero();
89
90 // Exchanging elements
91 bool SwapRows(size_t i, size_t j);
92 bool SwapColumns(size_t i, size_t j);
93 bool SwapRowColumn(size_t i, size_t j);
94
95 // Transformations
96 MatrixContent transposed() const;
97 MatrixContent &transpose();
98 gsl_vector* transformToEigenbasis();
99 void sortEigenbasis(gsl_vector *eigenvalues);
100 void performSingularValueDecomposition(MatrixContent &V, VectorContent &S);
101 VectorContent solveOverdeterminedLinearEquation(const VectorContent &b);
102 VectorContent solveOverdeterminedLinearEquation(MatrixContent &V, VectorContent &S, const VectorContent &b);
103
104 // Properties
105 bool IsNull() const;
106 bool IsPositive() const;
107 bool IsNegative() const;
108 bool IsNonNegative() const;
109 bool IsPositiveDefinite() const;
110 double Determinant() const;
111
112 // Operators
113 MatrixContent &operator=(const MatrixContent &src);
114 const MatrixContent &operator+=(const MatrixContent &rhs);
115 const MatrixContent &operator-=(const MatrixContent &rhs);
116 const MatrixContent &operator*=(const MatrixContent &rhs);
117 const MatrixContent operator*(const MatrixContent &rhs) const;
118 const MatrixContent &operator&=(const MatrixContent &rhs);
119 const MatrixContent operator&(const MatrixContent &rhs) const;
120 const MatrixContent &operator*=(const double factor);
121 bool operator==(const MatrixContent &rhs) const;
122
123 const size_t getRows() const;
124 const size_t getColumns() const;
125
126 VectorContent *getColumnVector(size_t column) const;
127 VectorContent *getRowVector(size_t row) const;
128 VectorContent *getDiagonalVector() const;
129
130 static std::pair<size_t, size_t> preparseMatrixDimensions(std::istream &inputstream);
131 void write(std::ostream &ost) const;
132
133protected:
134 double *Pointer(size_t m, size_t n);
135 const double *const_Pointer(size_t m, size_t n) const;
136
137 gsl_matrix * content;
138 const size_t rows; // store for internal purposes
139 const size_t columns; // store for internal purposes
140
141private:
142 bool free_content_on_exit;
143};
144
145std::ostream & operator<<(std::ostream &ost, const MatrixContent &mat);
146const MatrixContent operator*(const double factor,const MatrixContent& mat);
147const MatrixContent operator*(const MatrixContent &mat,const double factor);
148MatrixContent const operator+(const MatrixContent& a, const MatrixContent& b);
149MatrixContent const operator-(const MatrixContent& a, const MatrixContent& b);
150
151/** Matrix view.
152 * Extension of MatrixContent to contain not a gsl_matrix but only a view on a
153 * gsl_matrix
154 *
155 * We need the above MatrixBaseCase here:
156 * content, i.e. the gsl_matrix, must not be allocated, as it is just a view
157 * with an internal gsl_matrix_view. Hence, MatrixBaseCase is just dummy class
158 * to give the constructor a unique signature.
159 */
160struct MatrixViewContent : public MatrixContent
161{
162 MatrixViewContent(gsl_matrix_view _view) :
163 MatrixContent(_view.matrix.size1, _view.matrix.size2, MatrixBaseCase()),
164 view(_view)
165 {
166 content = &view.matrix;
167 }
168 ~MatrixViewContent(){
169 content=0;
170 }
171 gsl_matrix_view view;
172};
173
174#endif /* MATRIXCONTENT_HPP_ */
Note: See TracBrowser for help on using the repository browser.