source: src/LinearAlgebra/Matrix.cpp@ 6e5084

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 6e5084 was a439e5, checked in by Frederik Heber <heber@…>, 15 years ago

New functions vor class Vector and Matrix.

Matrix:

  • transpose() (returning copy and on itself)
  • zero() (initialize to zero)
  • transformToEigenbasis() (columns are eigenvectors, return value has eigenvalues)

Vector:

  • Property mode set to 100644
File size: 8.0 KB
Line 
1/*
2 * Matrix.cpp
3 *
4 * Created on: Jun 25, 2010
5 * Author: crueger
6 */
7
8#include "LinearAlgebra/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 "LinearAlgebra/Vector.hpp"
14#include "VectorContent.hpp"
15#include "MatrixContent.hpp"
16
17#include <gsl/gsl_blas.h>
18#include <gsl/gsl_eigen.h>
19#include <gsl/gsl_matrix.h>
20#include <gsl/gsl_multimin.h>
21#include <gsl/gsl_vector.h>
22#include <cmath>
23#include <iostream>
24
25using namespace std;
26
27Matrix::Matrix()
28{
29 content = new MatrixContent();
30 content->content = gsl_matrix_calloc(NDIM, NDIM);
31 createViews();
32}
33
34Matrix::Matrix(const double* src){
35 content = new MatrixContent();
36 content->content = gsl_matrix_alloc(NDIM, NDIM);
37 set(0,0, src[0]);
38 set(1,0, src[1]);
39 set(2,0, src[2]);
40
41 set(0,1, src[3]);
42 set(1,1, src[4]);
43 set(2,1, src[5]);
44
45 set(0,2, src[6]);
46 set(1,2, src[7]);
47 set(2,2, src[8]);
48 createViews();
49}
50
51Matrix::Matrix(const Matrix &src){
52 content = new MatrixContent();
53 content->content = gsl_matrix_alloc(NDIM, NDIM);
54 gsl_matrix_memcpy(content->content,src.content->content);
55 createViews();
56}
57
58Matrix::Matrix(MatrixContent* _content) :
59 content(_content)
60{
61 createViews();
62}
63
64Matrix::~Matrix()
65{
66 // delete all views
67 for(int i=NDIM;i--;){
68 delete rows_ptr[i];
69 }
70 for(int i=NDIM;i--;){
71 delete columns_ptr[i];
72 }
73 delete diagonal_ptr;
74 gsl_matrix_free(content->content);
75 delete content;
76}
77
78void Matrix::createViews(){
79 // create row views
80 for(int i=NDIM;i--;){
81 VectorContent *rowContent = new VectorViewContent(gsl_matrix_row(content->content,i));
82 rows_ptr[i] = new Vector(rowContent);
83 }
84 // create column views
85 for(int i=NDIM;i--;){
86 VectorContent *columnContent = new VectorViewContent(gsl_matrix_column(content->content,i));
87 columns_ptr[i] = new Vector(columnContent);
88 }
89 // create diagonal view
90 VectorContent *diagonalContent = new VectorViewContent(gsl_matrix_diagonal(content->content));
91 diagonal_ptr = new Vector(diagonalContent);
92}
93
94void Matrix::one(){
95 for(int i=NDIM;i--;){
96 for(int j=NDIM;j--;){
97 set(i,j,i==j);
98 }
99 }
100}
101
102void Matrix::zero(){
103 for(int i=NDIM;i--;){
104 for(int j=NDIM;j--;){
105 set(i,j,0.);
106 }
107 }
108}
109
110Matrix &Matrix::operator=(const Matrix &src){
111 if(&src!=this){
112 gsl_matrix_memcpy(content->content,src.content->content);
113 }
114 return *this;
115}
116
117Matrix &Matrix::operator+=(const Matrix &rhs){
118 gsl_matrix_add(content->content, rhs.content->content);
119 return *this;
120}
121
122Matrix &Matrix::operator-=(const Matrix &rhs){
123 gsl_matrix_sub(content->content, rhs.content->content);
124 return *this;
125}
126
127Matrix &Matrix::operator*=(const Matrix &rhs){
128 (*this) = (*this)*rhs;
129 return *this;
130}
131
132Matrix Matrix::operator+(const Matrix &rhs) const{
133 Matrix tmp = *this;
134 tmp+=rhs;
135 return tmp;
136}
137
138Matrix Matrix::operator-(const Matrix &rhs) const{
139 Matrix tmp = *this;
140 tmp-=rhs;
141 return tmp;
142}
143
144Matrix Matrix::operator*(const Matrix &rhs) const{
145 gsl_matrix *res = gsl_matrix_alloc(NDIM, NDIM);
146 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, content->content, rhs.content->content, 0.0, res);
147 MatrixContent *content= new MatrixContent();
148 content->content = res;
149 return Matrix(content);
150}
151
152double &Matrix::at(size_t i, size_t j){
153 ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
154 ASSERT(j>=0&&j<NDIM,"Index j for Matrix access out of range");
155 return *gsl_matrix_ptr (content->content, i, j);
156}
157
158const double Matrix::at(size_t i, size_t j) const{
159 ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
160 ASSERT(j>=0&&j<NDIM,"Index j for Matrix access out of range");
161 return gsl_matrix_get(content->content, i, j);
162}
163
164Vector &Matrix::row(size_t i){
165 ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
166 return *rows_ptr[i];
167}
168
169const Vector &Matrix::row(size_t i) const{
170 ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
171 return *rows_ptr[i];
172}
173
174Vector &Matrix::column(size_t i){
175 ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
176 return *columns_ptr[i];
177}
178
179const Vector &Matrix::column(size_t i) const{
180 ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
181 return *columns_ptr[i];
182}
183
184Vector &Matrix::diagonal(){
185 return *diagonal_ptr;
186}
187
188const Vector &Matrix::diagonal() const{
189 return *diagonal_ptr;
190}
191
192void Matrix::set(size_t i, size_t j, const double value){
193 ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
194 ASSERT(j>=0&&j<NDIM,"Index j for Matrix access out of range");
195 gsl_matrix_set(content->content,i,j,value);
196}
197
198double Matrix::determinant() const{
199 return at(0,0)*at(1,1)*at(2,2)
200 + at(0,1)*at(1,2)*at(2,0)
201 + at(0,2)*at(1,0)*at(2,1)
202 - at(2,0)*at(1,1)*at(0,2)
203 - at(2,1)*at(1,2)*at(0,0)
204 - at(2,2)*at(1,0)*at(0,1);
205}
206
207Matrix Matrix::transpose() const
208{
209 Matrix copy(*this);
210 copy.transpose();
211 return copy;
212}
213
214
215void Matrix::transpose()
216{
217 double tmp;
218 for (int i=0;i<NDIM;i++)
219 for (int j=i+1;j<NDIM;j++) {
220 tmp = at(j,i);
221 at(i,j) = tmp;
222 at(j,i) = tmp;
223 }
224}
225
226
227Matrix Matrix::invert() const{
228 double det = determinant();
229 if(fabs(det)<MYEPSILON){
230 throw NotInvertibleException(__FILE__,__LINE__);
231 }
232
233 double detReci = 1./det;
234 Matrix res;
235 res.set(0,0, detReci*RDET2(at(1,1),at(2,1),at(1,2),at(2,2))); // A_11
236 res.set(1,0, -detReci*RDET2(at(1,0),at(2,0),at(1,2),at(2,2))); // A_21
237 res.set(2,0, detReci*RDET2(at(1,0),at(2,0),at(1,1),at(2,1))); // A_31
238 res.set(0,1, -detReci*RDET2(at(0,1),at(2,1),at(0,2),at(2,2))); // A_12
239 res.set(1,1, detReci*RDET2(at(0,0),at(2,0),at(0,2),at(2,2))); // A_22
240 res.set(2,1, -detReci*RDET2(at(0,0),at(2,0),at(0,1),at(2,1))); // A_32
241 res.set(0,2, detReci*RDET2(at(0,1),at(1,1),at(0,2),at(1,2))); // A_13
242 res.set(1,2, -detReci*RDET2(at(0,0),at(1,0),at(0,2),at(1,2))); // A_23
243 res.set(2,2, detReci*RDET2(at(0,0),at(1,0),at(0,1),at(1,1))); // A_33
244 return res;
245}
246
247Vector Matrix::transformToEigenbasis()
248{
249 gsl_eigen_symmv_workspace *T = gsl_eigen_symmv_alloc(NDIM);
250 gsl_vector *eval = gsl_vector_alloc(NDIM);
251 gsl_matrix *evec = gsl_matrix_alloc(NDIM, NDIM);
252 gsl_eigen_symmv(content->content, eval, evec, T);
253 gsl_eigen_symmv_free(T);
254 gsl_matrix_memcpy(content->content, evec);
255 Vector evalues(gsl_vector_get(eval,0), gsl_vector_get(eval,1), gsl_vector_get(eval,2));
256 return evalues;
257}
258
259Matrix &Matrix::operator*=(const double factor){
260 gsl_matrix_scale(content->content, factor);
261 return *this;
262}
263
264Matrix operator*(const double factor,const Matrix& mat){
265 Matrix tmp = mat;
266 tmp*=factor;
267 return tmp;
268}
269
270Matrix operator*(const Matrix &mat,const double factor){
271 return factor*mat;
272}
273
274bool Matrix::operator==(const Matrix &rhs) const{
275 for(int i=NDIM;i--;){
276 for(int j=NDIM;j--;){
277 if(fabs(at(i,j)-rhs.at(i,j))>MYEPSILON){
278 return false;
279 }
280 }
281 }
282 return true;
283}
284
285/** Blows the 6-dimensional \a cell_size array up to a full NDIM by NDIM matrix.
286 * \param *symm 6-dim array of unique symmetric matrix components
287 * \return allocated NDIM*NDIM array with the symmetric matrix
288 */
289Matrix ReturnFullMatrixforSymmetric(const double * const symm)
290{
291 Matrix matrix;
292 matrix.set(0,0, symm[0]);
293 matrix.set(1,0, symm[1]);
294 matrix.set(2,0, symm[3]);
295 matrix.set(0,1, symm[1]);
296 matrix.set(1,1, symm[2]);
297 matrix.set(2,1, symm[4]);
298 matrix.set(0,2, symm[3]);
299 matrix.set(1,2, symm[4]);
300 matrix.set(2,2, symm[5]);
301 return matrix;
302};
303
304ostream &operator<<(ostream &ost,const Matrix &mat){
305 for(int i = 0;i<NDIM;++i){
306 ost << "\n";
307 for(int j = 0; j<NDIM;++j){
308 ost << mat.at(i,j);
309 if(j!=NDIM-1)
310 ost << "; ";
311 }
312 }
313 return ost;
314}
315
316Vector operator*(const Matrix &mat,const Vector &vec){
317 gsl_vector *res = gsl_vector_calloc(NDIM);
318 gsl_blas_dgemv( CblasNoTrans, 1.0, mat.content->content, vec.content->content, 0.0, res);
319 VectorContent *content = new VectorContent();
320 content->content = res;
321 return Vector(content);
322}
323
324Vector &operator*=(Vector& lhs,const Matrix &mat){
325 lhs = mat*lhs;
326 return lhs;
327}
328
Note: See TracBrowser for help on using the repository browser.