source: src/LinearAlgebra/MatrixContent.hpp@ ef0c6b

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 ef0c6b was 8e9ce1, checked in by Frederik Heber <heber@…>, 14 years ago

Added MatrixViewContent to allow views on blocks of a gsl_matrix.

  • Property mode set to 100644
File size: 4.4 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/** MatrixContent is a wrapper for gsl_matrix.
12 *
13 * The GNU Scientific Library contaisn some very well written routines for
14 * linear algebra problems. However, it's syntax is C and hence it does not
15 * lend itself to readable written code, i.e. C = A * B, where A, B, and C
16 * are all matrices. Writing code this way is very convenient, concise and
17 * also in the same manner as one would type in MatLab.
18 * However, with C++ and its feature to overload functions and its operator
19 * functions such syntax is easy to create.
20 *
21 * Hence, this class is a C++ wrapper to gsl_matrix. There already some other
22 * libraries, however they fall short for the following reasons:
23 * -# gslwrap: not very well written and uses floats instead of doubles
24 * -# gslcpp: last change is from 2007 and only very few commits
25 * -# o2scl: basically, the same functionality as gsl only written in C++,
26 * however it uses GPLv3 license which is inconvenient for MoleCuilder.
27 *
28 * <h1>Howto use MatrixContent</h1>
29 *
30 * One should not use MatrixContent directly but either have it as a member
31 * variable in a specialized class or inherit from it.
32 *
33 */
34
35#include <gsl/gsl_matrix.h>
36
37#include <iosfwd>
38
39class Vector;
40
41/** Dummy structure to create a unique signature.
42 *
43 */
44struct MatrixBaseCase{};
45
46class MatrixContent
47{
48 friend Vector operator*(const MatrixContent &mat,const Vector &vec);
49 friend std::ostream & operator<<(std::ostream &ost, const MatrixContent &mat);
50 friend class RealSpaceMatrix;
51 friend class LinearSystemOfEquations;
52
53 friend class MatrixContentTest;
54 friend class MatrixContentSymmetricTest;
55public:
56 MatrixContent(size_t rows, size_t columns);
57 MatrixContent(size_t _rows, size_t _columns, MatrixBaseCase);
58 MatrixContent(size_t rows, size_t columns, const double *src);
59 MatrixContent(gsl_matrix *&src);
60 MatrixContent(const MatrixContent &src);
61 MatrixContent(const MatrixContent *src);
62 virtual ~MatrixContent();
63
64 // Accessing
65 double &at(size_t i, size_t j);
66 const double at(size_t i, size_t j) const;
67 void set(size_t i, size_t j, const double value);
68
69 // Initializing
70 void setFromDoubleArray(double * x);
71 void setIdentity();
72 void setValue(double _value);
73 void setZero();
74
75 // Exchanging elements
76 bool SwapRows(size_t i, size_t j);
77 bool SwapColumns(size_t i, size_t j);
78 bool SwapRowColumn(size_t i, size_t j);
79
80 // Transformations
81 MatrixContent transpose() const;
82 MatrixContent &transpose();
83 gsl_vector* transformToEigenbasis();
84
85 // Properties
86 bool IsNull() const;
87 bool IsPositive() const;
88 bool IsNegative() const;
89 bool IsNonNegative() const;
90 bool IsPositiveDefinite() const;
91 double Determinant() const;
92
93 // Operators
94 MatrixContent &operator=(const MatrixContent &src);
95 const MatrixContent &operator+=(const MatrixContent &rhs);
96 const MatrixContent &operator-=(const MatrixContent &rhs);
97 const MatrixContent &operator*=(const MatrixContent &rhs);
98 const MatrixContent operator*(const MatrixContent &rhs) const;
99 const MatrixContent &operator*=(const double factor);
100 bool operator==(const MatrixContent &rhs) const;
101
102protected:
103 double *Pointer(size_t m, size_t n);
104 const double *const_Pointer(size_t m, size_t n) const;
105
106 gsl_matrix * content;
107 const size_t rows; // store for internal purposes
108 const size_t columns; // store for internal purposes
109
110private:
111};
112
113const MatrixContent operator*(const double factor,const MatrixContent& mat);
114const MatrixContent operator*(const MatrixContent &mat,const double factor);
115Vector operator*(const MatrixContent &mat,const Vector &vec);
116std::ostream & operator<<(std::ostream &ost, const MatrixContent &mat);
117
118/** Matrix view.
119 * Extension of MatrixContent to contain not a gsl_matrix but only a view on a
120 * gsl_matrix
121 *
122 * We need the above MatrixBaseCase here:
123 * content, i.e. the gsl_matrix, must not be allocated, as it is just a view
124 * with an internal gsl_matrix_view. Hence, MatrixBaseCase is just dummy class
125 * to give the constructor a unique signature.
126 */
127struct MatrixViewContent : public MatrixContent
128{
129 MatrixViewContent(gsl_matrix_view _view) :
130 MatrixContent(_view.matrix.size1, _view.matrix.size2, MatrixBaseCase()),
131 view(_view)
132 {
133 content = &view.matrix;
134 }
135 ~MatrixViewContent(){
136 content=0;
137 }
138 gsl_matrix_view view;
139};
140
141#endif /* MATRIXCONTENT_HPP_ */
Note: See TracBrowser for help on using the repository browser.