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