source: src/Matrix.cpp@ 5e8e02

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 5e8e02 was fee079, checked in by Tillmann Crueger <crueger@…>, 15 years ago

Added a "forward declaration" of the gsl_matrix struct in file Matrix.hpp

  • Property mode set to 100644
File size: 5.5 KB
Line 
1/*
2 * Matrix.cpp
3 *
4 * Created on: Jun 25, 2010
5 * Author: crueger
6 */
7
8#include "Matrix.hpp"
9#include "Helpers/Assert.hpp"
10#include "Exceptions/NotInvertibleException.hpp"
11#include "Helpers/fast_functions.hpp"
12#include "Helpers/Assert.hpp"
13#include "vector.hpp"
14#include "VectorContent.hpp"
15#include "MatrixContent.hpp"
16
17#include <gsl/gsl_blas.h>
18#include <cmath>
19#include <iostream>
20
21using namespace std;
22
23Matrix::Matrix()
24{
25 content = new MatrixContent();
26 content->content = gsl_matrix_calloc(NDIM, NDIM);
27}
28
29Matrix::Matrix(const double* src){
30 content = new MatrixContent();
31 content->content = gsl_matrix_alloc(NDIM, NDIM);
32 set(0,0, src[0]);
33 set(1,0, src[1]);
34 set(2,0, src[2]);
35
36 set(0,1, src[3]);
37 set(1,1, src[4]);
38 set(2,1, src[5]);
39
40 set(0,2, src[6]);
41 set(1,2, src[7]);
42 set(2,2, src[8]);
43}
44
45Matrix::Matrix(const Matrix &src){
46 content = new MatrixContent();
47 content->content = gsl_matrix_alloc(NDIM, NDIM);
48 gsl_matrix_memcpy(content->content,src.content->content);
49}
50
51Matrix::Matrix(MatrixContent* _content) :
52 content(_content)
53{}
54
55Matrix::~Matrix()
56{
57 gsl_matrix_free(content->content);
58 delete content;
59}
60
61void Matrix::one(){
62 gsl_matrix_free(content->content);
63 content->content = gsl_matrix_calloc(NDIM, NDIM);
64 for(int i = NDIM;i--;){
65 set(i,i,1.);
66 }
67}
68
69Matrix &Matrix::operator=(const Matrix &src){
70 if(&src!=this){
71 gsl_matrix_memcpy(content->content,src.content->content);
72 }
73 return *this;
74}
75
76Matrix &Matrix::operator+=(const Matrix &rhs){
77 gsl_matrix_add(content->content, rhs.content->content);
78 return *this;
79}
80
81Matrix &Matrix::operator-=(const Matrix &rhs){
82 gsl_matrix_sub(content->content, rhs.content->content);
83 return *this;
84}
85
86Matrix &Matrix::operator*=(const Matrix &rhs){
87 (*this) = (*this)*rhs;
88 return *this;
89}
90
91Matrix Matrix::operator+(const Matrix &rhs) const{
92 Matrix tmp = *this;
93 tmp+=rhs;
94 return tmp;
95}
96
97Matrix Matrix::operator-(const Matrix &rhs) const{
98 Matrix tmp = *this;
99 tmp-=rhs;
100 return tmp;
101}
102
103Matrix Matrix::operator*(const Matrix &rhs) const{
104 gsl_matrix *res = gsl_matrix_alloc(NDIM, NDIM);
105 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, content->content, rhs.content->content, 0.0, res);
106 MatrixContent *content= new MatrixContent();
107 content->content = res;
108 return Matrix(content);
109}
110
111double &Matrix::at(size_t i, size_t j){
112 ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
113 ASSERT(j>=0&&j<NDIM,"Index j for Matrix access out of range");
114 return *gsl_matrix_ptr (content->content, i, j);
115}
116
117const double Matrix::at(size_t i, size_t j) const{
118 ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
119 ASSERT(j>=0&&j<NDIM,"Index j for Matrix access out of range");
120 return gsl_matrix_get(content->content, i, j);
121}
122
123void Matrix::set(size_t i, size_t j, const double value){
124 ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
125 ASSERT(j>=0&&j<NDIM,"Index j for Matrix access out of range");
126 gsl_matrix_set(content->content,i,j,value);
127}
128
129double Matrix::determinant() const{
130 return at(0,0)*at(1,1)*at(2,2)
131 + at(0,1)*at(1,2)*at(2,0)
132 + at(0,2)*at(1,0)*at(2,1)
133 - at(2,0)*at(1,1)*at(0,2)
134 - at(2,1)*at(1,2)*at(0,0)
135 - at(2,2)*at(1,0)*at(0,1);
136}
137
138Matrix Matrix::invert() const{
139 double det = determinant();
140 if(fabs(det)<MYEPSILON){
141 throw NotInvertibleException(__FILE__,__LINE__);
142 }
143
144 double detReci = 1./det;
145 Matrix res;
146 res.set(0,0, detReci*RDET2(at(1,1),at(2,1),at(1,2),at(2,2))); // A_11
147 res.set(1,0, -detReci*RDET2(at(1,0),at(2,0),at(1,2),at(2,2))); // A_21
148 res.set(2,0, detReci*RDET2(at(1,0),at(2,0),at(1,1),at(2,1))); // A_31
149 res.set(0,1, -detReci*RDET2(at(0,1),at(2,1),at(0,2),at(2,2))); // A_12
150 res.set(1,1, detReci*RDET2(at(0,0),at(2,0),at(0,2),at(2,2))); // A_22
151 res.set(2,1, -detReci*RDET2(at(0,0),at(2,0),at(0,1),at(2,1))); // A_32
152 res.set(0,2, detReci*RDET2(at(0,1),at(1,1),at(0,2),at(1,2))); // A_13
153 res.set(1,2, -detReci*RDET2(at(0,0),at(1,0),at(0,2),at(1,2))); // A_23
154 res.set(2,2, detReci*RDET2(at(0,0),at(1,0),at(0,1),at(1,1))); // A_33
155 return res;
156}
157
158Matrix &Matrix::operator*=(const double factor){
159 gsl_matrix_scale(content->content, factor);
160 return *this;
161}
162
163Matrix operator*(const double factor,const Matrix& mat){
164 Matrix tmp = mat;
165 tmp*=factor;
166 return tmp;
167}
168
169Matrix operator*(const Matrix &mat,const double factor){
170 return factor*mat;
171}
172
173/** Blows the 6-dimensional \a cell_size array up to a full NDIM by NDIM matrix.
174 * \param *symm 6-dim array of unique symmetric matrix components
175 * \return allocated NDIM*NDIM array with the symmetric matrix
176 */
177Matrix ReturnFullMatrixforSymmetric(const double * const symm)
178{
179 Matrix matrix;
180 matrix.set(0,0, symm[0]);
181 matrix.set(1,0, symm[1]);
182 matrix.set(2,0, symm[3]);
183 matrix.set(0,1, symm[1]);
184 matrix.set(1,1, symm[2]);
185 matrix.set(2,1, symm[4]);
186 matrix.set(0,2, symm[3]);
187 matrix.set(1,2, symm[4]);
188 matrix.set(2,2, symm[5]);
189 return matrix;
190};
191
192ostream &operator<<(ostream &ost,const Matrix &mat){
193 for(int i = 0;i<NDIM;++i){
194 ost << "\n";
195 for(int j = 0; j<NDIM;++j){
196 ost << mat.at(i,j);
197 if(j!=NDIM-1)
198 ost << "; ";
199 }
200 }
201 return ost;
202}
203
204Vector operator*(const Matrix &mat,const Vector &vec){
205 gsl_vector *res = gsl_vector_calloc(NDIM);
206 gsl_blas_dgemv( CblasNoTrans, 1.0, mat.content->content, vec.content->content, 0.0, res);
207 VectorContent *content = new VectorContent();
208 content->content = res;
209 return Vector(content);
210}
211
212Vector &operator*=(Vector& lhs,const Matrix &mat){
213 lhs = mat*lhs;
214 return lhs;
215}
216
Note: See TracBrowser for help on using the repository browser.