1 | /*
|
---|
2 | * Project: MoleCuilder
|
---|
3 | * Description: creates and alters molecular systems
|
---|
4 | * Copyright (C) 2010 University of Bonn. All rights reserved.
|
---|
5 | * Please see the LICENSE file or "Copyright notice" in builder.cpp for details.
|
---|
6 | */
|
---|
7 |
|
---|
8 | /*
|
---|
9 | * RealSpaceMatrix.cpp
|
---|
10 | *
|
---|
11 | * Created on: Jun 25, 2010
|
---|
12 | * Author: crueger
|
---|
13 | */
|
---|
14 |
|
---|
15 | // include config.h
|
---|
16 | #ifdef HAVE_CONFIG_H
|
---|
17 | #include <config.h>
|
---|
18 | #endif
|
---|
19 |
|
---|
20 | #include "CodePatterns/MemDebug.hpp"
|
---|
21 |
|
---|
22 | #include "Exceptions/NotInvertibleException.hpp"
|
---|
23 | #include "CodePatterns/Assert.hpp"
|
---|
24 | #include "Helpers/defs.hpp"
|
---|
25 | #include "Helpers/fast_functions.hpp"
|
---|
26 | #include "LinearAlgebra/defs.hpp"
|
---|
27 | #include "LinearAlgebra/MatrixContent.hpp"
|
---|
28 | #include "LinearAlgebra/RealSpaceMatrix.hpp"
|
---|
29 | #include "LinearAlgebra/Vector.hpp"
|
---|
30 | #include "LinearAlgebra/VectorContent.hpp"
|
---|
31 | #include "RandomNumbers/RandomNumberGeneratorFactory.hpp"
|
---|
32 | #include "RandomNumbers/RandomNumberGenerator.hpp"
|
---|
33 |
|
---|
34 | #include <gsl/gsl_blas.h>
|
---|
35 | #include <gsl/gsl_eigen.h>
|
---|
36 | #include <gsl/gsl_matrix.h>
|
---|
37 | #include <gsl/gsl_multimin.h>
|
---|
38 | #include <gsl/gsl_vector.h>
|
---|
39 | #include <cmath>
|
---|
40 | #include <iostream>
|
---|
41 | #include <limits>
|
---|
42 |
|
---|
43 | using namespace std;
|
---|
44 |
|
---|
45 | RealSpaceMatrix::RealSpaceMatrix()
|
---|
46 | {
|
---|
47 | content = new MatrixContent(NDIM, NDIM);
|
---|
48 | createViews();
|
---|
49 | }
|
---|
50 |
|
---|
51 | RealSpaceMatrix::RealSpaceMatrix(const double* src)
|
---|
52 | {
|
---|
53 | content = new MatrixContent(NDIM, NDIM, src);
|
---|
54 | createViews();
|
---|
55 | }
|
---|
56 |
|
---|
57 | RealSpaceMatrix::RealSpaceMatrix(const RealSpaceMatrix &src)
|
---|
58 | {
|
---|
59 | content = new MatrixContent(src.content);
|
---|
60 | createViews();
|
---|
61 | }
|
---|
62 |
|
---|
63 | RealSpaceMatrix::RealSpaceMatrix(const MatrixContent &src)
|
---|
64 | {
|
---|
65 | content = new MatrixContent(src);
|
---|
66 | createViews();
|
---|
67 | }
|
---|
68 |
|
---|
69 | RealSpaceMatrix::RealSpaceMatrix(MatrixContent* _content)
|
---|
70 | {
|
---|
71 | content = new MatrixContent(_content);
|
---|
72 | createViews();
|
---|
73 | }
|
---|
74 |
|
---|
75 | RealSpaceMatrix::~RealSpaceMatrix()
|
---|
76 | {
|
---|
77 | // delete all views
|
---|
78 | for(int i=NDIM;i--;){
|
---|
79 | delete rows_ptr[i];
|
---|
80 | }
|
---|
81 | for(int i=NDIM;i--;){
|
---|
82 | delete columns_ptr[i];
|
---|
83 | }
|
---|
84 | delete diagonal_ptr;
|
---|
85 | delete content;
|
---|
86 | }
|
---|
87 |
|
---|
88 | void RealSpaceMatrix::createViews(){
|
---|
89 | // create row views
|
---|
90 | for(int i=NDIM;i--;){
|
---|
91 | VectorContent *rowContent = content->getRowVector(i);
|
---|
92 | rows_ptr[i] = new Vector(rowContent);
|
---|
93 | ASSERT(rowContent == NULL, "RealSpaceMatrix::createViews() - rowContent was not taken over as supposed to happen!");
|
---|
94 | }
|
---|
95 | // create column views
|
---|
96 | for(int i=NDIM;i--;){
|
---|
97 | VectorContent *columnContent = content->getColumnVector(i);
|
---|
98 | columns_ptr[i] = new Vector(columnContent);
|
---|
99 | ASSERT(columnContent == NULL, "RealSpaceMatrix::createViews() - columnContent was not taken over as supposed to happen!");
|
---|
100 | }
|
---|
101 | // create diagonal view
|
---|
102 | VectorContent *diagonalContent = content->getDiagonalVector();
|
---|
103 | diagonal_ptr = new Vector(diagonalContent);
|
---|
104 | ASSERT(diagonalContent == NULL, "RealSpaceMatrix::createViews() - diagonalContent was not taken over as supposed to happen!");
|
---|
105 | }
|
---|
106 |
|
---|
107 | void RealSpaceMatrix::setIdentity(){
|
---|
108 | content->setIdentity();
|
---|
109 | }
|
---|
110 |
|
---|
111 | void RealSpaceMatrix::setZero(){
|
---|
112 | content->setZero();
|
---|
113 | }
|
---|
114 |
|
---|
115 | void RealSpaceMatrix::setRotation(const double x, const double y, const double z)
|
---|
116 | {
|
---|
117 | set(0,0, cos(y)*cos(z));
|
---|
118 | set(0,1, cos(z)*sin(x)*sin(y) - cos(x)*sin(z));
|
---|
119 | set(0,2, cos(x)*cos(z)*sin(y) + sin(x) * sin(z));
|
---|
120 | set(1,0, cos(y)*sin(z));
|
---|
121 | set(1,1, cos(x)*cos(z) + sin(x)*sin(y)*sin(z));
|
---|
122 | set(1,2, -cos(z)*sin(x) + cos(x)*sin(y)*sin(z));
|
---|
123 | set(2,0, -sin(y));
|
---|
124 | set(2,1, cos(y)*sin(x));
|
---|
125 | set(2,2, cos(x)*cos(y));
|
---|
126 | }
|
---|
127 |
|
---|
128 | void RealSpaceMatrix::setRandomRotation()
|
---|
129 | {
|
---|
130 | double phi[NDIM];
|
---|
131 | RandomNumberGenerator &random = RandomNumberGeneratorFactory::getInstance().makeRandomNumberGenerator();
|
---|
132 | const double rng_min = random.min();
|
---|
133 | const double rng_max = random.max();
|
---|
134 |
|
---|
135 |
|
---|
136 | for (int i=0;i<NDIM;i++) {
|
---|
137 | phi[i] = (random()/(rng_max-rng_min))*(2.*M_PI);
|
---|
138 | std::cout << "Random angle is " << phi[i] << std::endl;
|
---|
139 | }
|
---|
140 |
|
---|
141 | set(0,0, cos(phi[0]) *cos(phi[2]) + (sin(phi[0])*sin(phi[1])*sin(phi[2])));
|
---|
142 | set(0,1, sin(phi[0]) *cos(phi[2]) - (cos(phi[0])*sin(phi[1])*sin(phi[2])));
|
---|
143 | set(0,2, cos(phi[1])*sin(phi[2]) );
|
---|
144 | set(1,0, -sin(phi[0])*cos(phi[1]) );
|
---|
145 | set(1,1, cos(phi[0])*cos(phi[1]) );
|
---|
146 | set(1,2, sin(phi[1]) );
|
---|
147 | set(2,0, -cos(phi[0]) *sin(phi[2]) + (sin(phi[0])*sin(phi[1])*cos(phi[2])));
|
---|
148 | set(2,1, -sin(phi[0]) *sin(phi[2]) - (cos(phi[0])*sin(phi[1])*cos(phi[2])));
|
---|
149 | set(2,2, cos(phi[1])*cos(phi[2]) );
|
---|
150 | }
|
---|
151 |
|
---|
152 |
|
---|
153 | RealSpaceMatrix &RealSpaceMatrix::operator=(const RealSpaceMatrix &src)
|
---|
154 | {
|
---|
155 | // MatrixContent checks for self-assignment
|
---|
156 | *content = *(src.content);
|
---|
157 | return *this;
|
---|
158 | }
|
---|
159 |
|
---|
160 | const RealSpaceMatrix &RealSpaceMatrix::operator+=(const RealSpaceMatrix &rhs)
|
---|
161 | {
|
---|
162 | *content += *(rhs.content);
|
---|
163 | return *this;
|
---|
164 | }
|
---|
165 |
|
---|
166 | const RealSpaceMatrix &RealSpaceMatrix::operator-=(const RealSpaceMatrix &rhs)
|
---|
167 | {
|
---|
168 | *content -= *(rhs.content);
|
---|
169 | return *this;
|
---|
170 | }
|
---|
171 |
|
---|
172 | const RealSpaceMatrix &RealSpaceMatrix::operator*=(const RealSpaceMatrix &rhs)
|
---|
173 | {
|
---|
174 | (*content) *= (*rhs.content);
|
---|
175 | return *this;
|
---|
176 | }
|
---|
177 |
|
---|
178 | const RealSpaceMatrix RealSpaceMatrix::operator+(const RealSpaceMatrix &rhs) const
|
---|
179 | {
|
---|
180 | RealSpaceMatrix tmp = *this;
|
---|
181 | tmp+=rhs;
|
---|
182 | return tmp;
|
---|
183 | }
|
---|
184 |
|
---|
185 | const RealSpaceMatrix RealSpaceMatrix::operator-(const RealSpaceMatrix &rhs) const
|
---|
186 | {
|
---|
187 | RealSpaceMatrix tmp = *this;
|
---|
188 | tmp-=rhs;
|
---|
189 | return tmp;
|
---|
190 | }
|
---|
191 |
|
---|
192 | const RealSpaceMatrix RealSpaceMatrix::operator*(const RealSpaceMatrix &rhs) const
|
---|
193 | {
|
---|
194 | RealSpaceMatrix tmp(content);
|
---|
195 | tmp *= rhs;
|
---|
196 | return tmp;
|
---|
197 | }
|
---|
198 |
|
---|
199 | double &RealSpaceMatrix::at(size_t i, size_t j)
|
---|
200 | {
|
---|
201 | return content->at(i,j);
|
---|
202 | }
|
---|
203 |
|
---|
204 | const double RealSpaceMatrix::at(size_t i, size_t j) const
|
---|
205 | {
|
---|
206 | return content->at(i,j);
|
---|
207 | }
|
---|
208 |
|
---|
209 | Vector &RealSpaceMatrix::row(size_t i)
|
---|
210 | {
|
---|
211 | ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
|
---|
212 | return *rows_ptr[i];
|
---|
213 | }
|
---|
214 |
|
---|
215 | const Vector &RealSpaceMatrix::row(size_t i) const
|
---|
216 | {
|
---|
217 | ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
|
---|
218 | return *rows_ptr[i];
|
---|
219 | }
|
---|
220 |
|
---|
221 | Vector &RealSpaceMatrix::column(size_t i)
|
---|
222 | {
|
---|
223 | ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
|
---|
224 | return *columns_ptr[i];
|
---|
225 | }
|
---|
226 |
|
---|
227 | const Vector &RealSpaceMatrix::column(size_t i) const
|
---|
228 | {
|
---|
229 | ASSERT(i>=0&&i<NDIM,"Index i for Matrix access out of range");
|
---|
230 | return *columns_ptr[i];
|
---|
231 | }
|
---|
232 |
|
---|
233 | Vector &RealSpaceMatrix::diagonal()
|
---|
234 | {
|
---|
235 | return *diagonal_ptr;
|
---|
236 | }
|
---|
237 |
|
---|
238 | const Vector &RealSpaceMatrix::diagonal() const
|
---|
239 | {
|
---|
240 | return *diagonal_ptr;
|
---|
241 | }
|
---|
242 |
|
---|
243 | void RealSpaceMatrix::set(size_t i, size_t j, const double value)
|
---|
244 | {
|
---|
245 | content->set(i,j, value);
|
---|
246 | }
|
---|
247 |
|
---|
248 | double RealSpaceMatrix::determinant() const{
|
---|
249 | return at(0,0)*at(1,1)*at(2,2)
|
---|
250 | + at(0,1)*at(1,2)*at(2,0)
|
---|
251 | + at(0,2)*at(1,0)*at(2,1)
|
---|
252 | - at(2,0)*at(1,1)*at(0,2)
|
---|
253 | - at(2,1)*at(1,2)*at(0,0)
|
---|
254 | - at(2,2)*at(1,0)*at(0,1);
|
---|
255 | }
|
---|
256 |
|
---|
257 |
|
---|
258 | RealSpaceMatrix RealSpaceMatrix::invert() const{
|
---|
259 | double det = determinant();
|
---|
260 | if(fabs(det) <= LINALG_MYEPSILON){
|
---|
261 | throw NotInvertibleException(__FILE__,__LINE__);
|
---|
262 | }
|
---|
263 |
|
---|
264 | double detReci = 1./det;
|
---|
265 | RealSpaceMatrix res;
|
---|
266 | res.set(0,0, detReci*RDET2(at(1,1),at(2,1),at(1,2),at(2,2))); // A_11
|
---|
267 | res.set(1,0, -detReci*RDET2(at(1,0),at(2,0),at(1,2),at(2,2))); // A_21
|
---|
268 | res.set(2,0, detReci*RDET2(at(1,0),at(2,0),at(1,1),at(2,1))); // A_31
|
---|
269 | res.set(0,1, -detReci*RDET2(at(0,1),at(2,1),at(0,2),at(2,2))); // A_12
|
---|
270 | res.set(1,1, detReci*RDET2(at(0,0),at(2,0),at(0,2),at(2,2))); // A_22
|
---|
271 | res.set(2,1, -detReci*RDET2(at(0,0),at(2,0),at(0,1),at(2,1))); // A_32
|
---|
272 | res.set(0,2, detReci*RDET2(at(0,1),at(1,1),at(0,2),at(1,2))); // A_13
|
---|
273 | res.set(1,2, -detReci*RDET2(at(0,0),at(1,0),at(0,2),at(1,2))); // A_23
|
---|
274 | res.set(2,2, detReci*RDET2(at(0,0),at(1,0),at(0,1),at(1,1))); // A_33
|
---|
275 | return res;
|
---|
276 | }
|
---|
277 |
|
---|
278 | RealSpaceMatrix RealSpaceMatrix::transpose() const
|
---|
279 | {
|
---|
280 | RealSpaceMatrix res = RealSpaceMatrix(const_cast<const MatrixContent *>(content)->transpose());
|
---|
281 | return res;
|
---|
282 | }
|
---|
283 |
|
---|
284 | RealSpaceMatrix &RealSpaceMatrix::transpose()
|
---|
285 | {
|
---|
286 | content->transpose();
|
---|
287 | return *this;
|
---|
288 | }
|
---|
289 |
|
---|
290 | Vector RealSpaceMatrix::transformToEigenbasis()
|
---|
291 | {
|
---|
292 | gsl_vector *eval = content->transformToEigenbasis();
|
---|
293 | Vector evalues(gsl_vector_get(eval,0), gsl_vector_get(eval,1), gsl_vector_get(eval,2));
|
---|
294 | gsl_vector_free(eval);
|
---|
295 | return evalues;
|
---|
296 | }
|
---|
297 |
|
---|
298 | const RealSpaceMatrix &RealSpaceMatrix::operator*=(const double factor)
|
---|
299 | {
|
---|
300 | *content *= factor;
|
---|
301 | return *this;
|
---|
302 | }
|
---|
303 |
|
---|
304 | const RealSpaceMatrix operator*(const double factor,const RealSpaceMatrix& mat)
|
---|
305 | {
|
---|
306 | RealSpaceMatrix tmp = mat;
|
---|
307 | return tmp *= factor;
|
---|
308 | }
|
---|
309 |
|
---|
310 | const RealSpaceMatrix operator*(const RealSpaceMatrix &mat,const double factor)
|
---|
311 | {
|
---|
312 | return factor*mat;
|
---|
313 | }
|
---|
314 |
|
---|
315 | bool RealSpaceMatrix::operator==(const RealSpaceMatrix &rhs) const
|
---|
316 | {
|
---|
317 | return (*content == *(rhs.content));
|
---|
318 | }
|
---|
319 |
|
---|
320 | /** Blows the 6-dimensional \a cell_size array up to a full NDIM by NDIM matrix.
|
---|
321 | * \param *symm 6-dim array of unique symmetric matrix components
|
---|
322 | * \return allocated NDIM*NDIM array with the symmetric matrix
|
---|
323 | */
|
---|
324 | RealSpaceMatrix ReturnFullMatrixforSymmetric(const double * const symm)
|
---|
325 | {
|
---|
326 | RealSpaceMatrix matrix;
|
---|
327 | matrix.set(0,0, symm[0]);
|
---|
328 | matrix.set(1,0, symm[1]);
|
---|
329 | matrix.set(2,0, symm[3]);
|
---|
330 | matrix.set(0,1, symm[1]);
|
---|
331 | matrix.set(1,1, symm[2]);
|
---|
332 | matrix.set(2,1, symm[4]);
|
---|
333 | matrix.set(0,2, symm[3]);
|
---|
334 | matrix.set(1,2, symm[4]);
|
---|
335 | matrix.set(2,2, symm[5]);
|
---|
336 | return matrix;
|
---|
337 | };
|
---|
338 |
|
---|
339 | ostream &operator<<(ostream &ost,const RealSpaceMatrix &mat)
|
---|
340 | {
|
---|
341 | for(int i = 0;i<NDIM;++i){
|
---|
342 | ost << "\n";
|
---|
343 | for(int j = 0; j<NDIM;++j){
|
---|
344 | ost << mat.at(i,j);
|
---|
345 | if(j!=NDIM-1)
|
---|
346 | ost << "; ";
|
---|
347 | }
|
---|
348 | }
|
---|
349 | return ost;
|
---|
350 | }
|
---|
351 |
|
---|
352 | Vector operator*(const RealSpaceMatrix &mat,const Vector &vec)
|
---|
353 | {
|
---|
354 | return (*mat.content) * vec;
|
---|
355 | }
|
---|
356 |
|
---|
357 | Vector &operator*=(Vector& lhs,const RealSpaceMatrix &mat)
|
---|
358 | {
|
---|
359 | lhs = mat*lhs;
|
---|
360 | return lhs;
|
---|
361 | }
|
---|
362 |
|
---|