| 1 | /*
 | 
|---|
| 2 |  * molecule_geometry.cpp
 | 
|---|
| 3 |  *
 | 
|---|
| 4 |  *  Created on: Oct 5, 2009
 | 
|---|
| 5 |  *      Author: heber
 | 
|---|
| 6 |  */
 | 
|---|
| 7 | 
 | 
|---|
| 8 | #include "atom.hpp"
 | 
|---|
| 9 | #include "bond.hpp"
 | 
|---|
| 10 | #include "config.hpp"
 | 
|---|
| 11 | #include "element.hpp"
 | 
|---|
| 12 | #include "helpers.hpp"
 | 
|---|
| 13 | #include "leastsquaremin.hpp"
 | 
|---|
| 14 | #include "log.hpp"
 | 
|---|
| 15 | #include "memoryallocator.hpp"
 | 
|---|
| 16 | #include "molecule.hpp"
 | 
|---|
| 17 | #include "World.hpp"
 | 
|---|
| 18 | 
 | 
|---|
| 19 | /************************************* Functions for class molecule *********************************/
 | 
|---|
| 20 | 
 | 
|---|
| 21 | 
 | 
|---|
| 22 | /** Centers the molecule in the box whose lengths are defined by vector \a *BoxLengths.
 | 
|---|
| 23 |  * \param *out output stream for debugging
 | 
|---|
| 24 |  */
 | 
|---|
| 25 | bool molecule::CenterInBox()
 | 
|---|
| 26 | {
 | 
|---|
| 27 |   bool status = true;
 | 
|---|
| 28 |   const Vector *Center = DetermineCenterOfAll();
 | 
|---|
| 29 |   double * const cell_size = World::getInstance().getDomain();
 | 
|---|
| 30 |   double *M = ReturnFullMatrixforSymmetric(cell_size);
 | 
|---|
| 31 |   double *Minv = InverseMatrix(M);
 | 
|---|
| 32 | 
 | 
|---|
| 33 |   // go through all atoms
 | 
|---|
| 34 |   ActOnAllVectors( &Vector::SubtractVector, *Center);
 | 
|---|
| 35 |   ActOnAllVectors( &Vector::WrapPeriodically, (const double *)M, (const double *)Minv);
 | 
|---|
| 36 | 
 | 
|---|
| 37 |   Free(&M);
 | 
|---|
| 38 |   Free(&Minv);
 | 
|---|
| 39 |   delete(Center);
 | 
|---|
| 40 |   return status;
 | 
|---|
| 41 | };
 | 
|---|
| 42 | 
 | 
|---|
| 43 | 
 | 
|---|
| 44 | /** Bounds the molecule in the box whose lengths are defined by vector \a *BoxLengths.
 | 
|---|
| 45 |  * \param *out output stream for debugging
 | 
|---|
| 46 |  */
 | 
|---|
| 47 | bool molecule::BoundInBox()
 | 
|---|
| 48 | {
 | 
|---|
| 49 |   bool status = true;
 | 
|---|
| 50 |   double * const cell_size = World::getInstance().getDomain();
 | 
|---|
| 51 |   double *M = ReturnFullMatrixforSymmetric(cell_size);
 | 
|---|
| 52 |   double *Minv = InverseMatrix(M);
 | 
|---|
| 53 | 
 | 
|---|
| 54 |   // go through all atoms
 | 
|---|
| 55 |   ActOnAllVectors( &Vector::WrapPeriodically, (const double *)M, (const double *)Minv);
 | 
|---|
| 56 | 
 | 
|---|
| 57 |   Free(&M);
 | 
|---|
| 58 |   Free(&Minv);
 | 
|---|
| 59 |   return status;
 | 
|---|
| 60 | };
 | 
|---|
| 61 | 
 | 
|---|
| 62 | /** Centers the edge of the atoms at (0,0,0).
 | 
|---|
| 63 |  * \param *out output stream for debugging
 | 
|---|
| 64 |  * \param *max coordinates of other edge, specifying box dimensions.
 | 
|---|
| 65 |  */
 | 
|---|
| 66 | void molecule::CenterEdge(Vector *max)
 | 
|---|
| 67 | {
 | 
|---|
| 68 |   Vector *min = new Vector;
 | 
|---|
| 69 | 
 | 
|---|
| 70 | //  Log() << Verbose(3) << "Begin of CenterEdge." << endl;
 | 
|---|
| 71 |   atom *ptr = start->next;  // start at first in list
 | 
|---|
| 72 |   if (ptr != end) {   //list not empty?
 | 
|---|
| 73 |     for (int i=NDIM;i--;) {
 | 
|---|
| 74 |       max->at(i) = ptr->x[i];
 | 
|---|
| 75 |       min->at(i) = ptr->x[i];
 | 
|---|
| 76 |     }
 | 
|---|
| 77 |     while (ptr->next != end) {  // continue with second if present
 | 
|---|
| 78 |       ptr = ptr->next;
 | 
|---|
| 79 |       //ptr->Output(1,1,out);
 | 
|---|
| 80 |       for (int i=NDIM;i--;) {
 | 
|---|
| 81 |         max->at(i) = (max->at(i) < ptr->x[i]) ? ptr->x[i] : max->at(i);
 | 
|---|
| 82 |         min->at(i) = (min->at(i) > ptr->x[i]) ? ptr->x[i] : min->at(i);
 | 
|---|
| 83 |       }
 | 
|---|
| 84 |     }
 | 
|---|
| 85 | //    Log() << Verbose(4) << "Maximum is ";
 | 
|---|
| 86 | //    max->Output(out);
 | 
|---|
| 87 | //    Log() << Verbose(0) << ", Minimum is ";
 | 
|---|
| 88 | //    min->Output(out);
 | 
|---|
| 89 | //    Log() << Verbose(0) << endl;
 | 
|---|
| 90 |     min->Scale(-1.);
 | 
|---|
| 91 |     (*max) += (*min);
 | 
|---|
| 92 |     Translate(min);
 | 
|---|
| 93 |     Center.Zero();
 | 
|---|
| 94 |   }
 | 
|---|
| 95 |   delete(min);
 | 
|---|
| 96 | //  Log() << Verbose(3) << "End of CenterEdge." << endl;
 | 
|---|
| 97 | };
 | 
|---|
| 98 | 
 | 
|---|
| 99 | /** Centers the center of the atoms at (0,0,0).
 | 
|---|
| 100 |  * \param *out output stream for debugging
 | 
|---|
| 101 |  * \param *center return vector for translation vector
 | 
|---|
| 102 |  */
 | 
|---|
| 103 | void molecule::CenterOrigin()
 | 
|---|
| 104 | {
 | 
|---|
| 105 |   int Num = 0;
 | 
|---|
| 106 |   atom *ptr = start;  // start at first in list
 | 
|---|
| 107 | 
 | 
|---|
| 108 |   Center.Zero();
 | 
|---|
| 109 | 
 | 
|---|
| 110 |   if (ptr->next != end) {   //list not empty?
 | 
|---|
| 111 |     while (ptr->next != end) {  // continue with second if present
 | 
|---|
| 112 |       ptr = ptr->next;
 | 
|---|
| 113 |       Num++;
 | 
|---|
| 114 |       Center += ptr->x;
 | 
|---|
| 115 |     }
 | 
|---|
| 116 |     Center.Scale(-1./Num); // divide through total number (and sign for direction)
 | 
|---|
| 117 |     Translate(&Center);
 | 
|---|
| 118 |     Center.Zero();
 | 
|---|
| 119 |   }
 | 
|---|
| 120 | };
 | 
|---|
| 121 | 
 | 
|---|
| 122 | /** Returns vector pointing to center of all atoms.
 | 
|---|
| 123 |  * \return pointer to center of all vector
 | 
|---|
| 124 |  */
 | 
|---|
| 125 | Vector * molecule::DetermineCenterOfAll() const
 | 
|---|
| 126 | {
 | 
|---|
| 127 |   atom *ptr = start->next;  // start at first in list
 | 
|---|
| 128 |   Vector *a = new Vector();
 | 
|---|
| 129 |   Vector tmp;
 | 
|---|
| 130 |   double Num = 0;
 | 
|---|
| 131 | 
 | 
|---|
| 132 |   a->Zero();
 | 
|---|
| 133 | 
 | 
|---|
| 134 |   if (ptr != end) {   //list not empty?
 | 
|---|
| 135 |     while (ptr->next != end) {  // continue with second if present
 | 
|---|
| 136 |       ptr = ptr->next;
 | 
|---|
| 137 |       Num += 1.;
 | 
|---|
| 138 |       tmp = ptr->x;
 | 
|---|
| 139 |       (*a) += tmp;
 | 
|---|
| 140 |     }
 | 
|---|
| 141 |     a->Scale(1./Num); // divide through total mass (and sign for direction)
 | 
|---|
| 142 |   }
 | 
|---|
| 143 |   return a;
 | 
|---|
| 144 | };
 | 
|---|
| 145 | 
 | 
|---|
| 146 | /** Returns vector pointing to center of gravity.
 | 
|---|
| 147 |  * \param *out output stream for debugging
 | 
|---|
| 148 |  * \return pointer to center of gravity vector
 | 
|---|
| 149 |  */
 | 
|---|
| 150 | Vector * molecule::DetermineCenterOfGravity()
 | 
|---|
| 151 | {
 | 
|---|
| 152 |   atom *ptr = start->next;  // start at first in list
 | 
|---|
| 153 |   Vector *a = new Vector();
 | 
|---|
| 154 |   Vector tmp;
 | 
|---|
| 155 |   double Num = 0;
 | 
|---|
| 156 | 
 | 
|---|
| 157 |   a->Zero();
 | 
|---|
| 158 | 
 | 
|---|
| 159 |   if (ptr != end) {   //list not empty?
 | 
|---|
| 160 |     while (ptr->next != end) {  // continue with second if present
 | 
|---|
| 161 |       ptr = ptr->next;
 | 
|---|
| 162 |       Num += ptr->type->mass;
 | 
|---|
| 163 |       tmp = ptr->type->mass * ptr->x;
 | 
|---|
| 164 |       (*a) += tmp;
 | 
|---|
| 165 |     }
 | 
|---|
| 166 |     a->Scale(-1./Num); // divide through total mass (and sign for direction)
 | 
|---|
| 167 |   }
 | 
|---|
| 168 | //  Log() << Verbose(1) << "Resulting center of gravity: ";
 | 
|---|
| 169 | //  a->Output(out);
 | 
|---|
| 170 | //  Log() << Verbose(0) << endl;
 | 
|---|
| 171 |   return a;
 | 
|---|
| 172 | };
 | 
|---|
| 173 | 
 | 
|---|
| 174 | /** Centers the center of gravity of the atoms at (0,0,0).
 | 
|---|
| 175 |  * \param *out output stream for debugging
 | 
|---|
| 176 |  * \param *center return vector for translation vector
 | 
|---|
| 177 |  */
 | 
|---|
| 178 | void molecule::CenterPeriodic()
 | 
|---|
| 179 | {
 | 
|---|
| 180 |   DeterminePeriodicCenter(Center);
 | 
|---|
| 181 | };
 | 
|---|
| 182 | 
 | 
|---|
| 183 | 
 | 
|---|
| 184 | /** Centers the center of gravity of the atoms at (0,0,0).
 | 
|---|
| 185 |  * \param *out output stream for debugging
 | 
|---|
| 186 |  * \param *center return vector for translation vector
 | 
|---|
| 187 |  */
 | 
|---|
| 188 | void molecule::CenterAtVector(Vector *newcenter)
 | 
|---|
| 189 | {
 | 
|---|
| 190 |   Center = *newcenter;
 | 
|---|
| 191 | };
 | 
|---|
| 192 | 
 | 
|---|
| 193 | 
 | 
|---|
| 194 | /** Scales all atoms by \a *factor.
 | 
|---|
| 195 |  * \param *factor pointer to scaling factor
 | 
|---|
| 196 |  *
 | 
|---|
| 197 |  * TODO: Is this realy what is meant, i.e.
 | 
|---|
| 198 |  * x=(x[0]*factor[0],x[1]*factor[1],x[2]*factor[2]) (current impl)
 | 
|---|
| 199 |  * or rather
 | 
|---|
| 200 |  * x=(**factor) * x (as suggested by comment)
 | 
|---|
| 201 |  */
 | 
|---|
| 202 | void molecule::Scale(const double ** const factor)
 | 
|---|
| 203 | {
 | 
|---|
| 204 |   atom *ptr = start;
 | 
|---|
| 205 | 
 | 
|---|
| 206 |   while (ptr->next != end) {
 | 
|---|
| 207 |     ptr = ptr->next;
 | 
|---|
| 208 |     for (int j=0;j<MDSteps;j++)
 | 
|---|
| 209 |       ptr->Trajectory.R.at(j).ScaleAll(*factor);
 | 
|---|
| 210 |     ptr->x.ScaleAll(*factor);
 | 
|---|
| 211 |   }
 | 
|---|
| 212 | };
 | 
|---|
| 213 | 
 | 
|---|
| 214 | /** Translate all atoms by given vector.
 | 
|---|
| 215 |  * \param trans[] translation vector.
 | 
|---|
| 216 |  */
 | 
|---|
| 217 | void molecule::Translate(const Vector *trans)
 | 
|---|
| 218 | {
 | 
|---|
| 219 |   atom *ptr = start;
 | 
|---|
| 220 | 
 | 
|---|
| 221 |   while (ptr->next != end) {
 | 
|---|
| 222 |     ptr = ptr->next;
 | 
|---|
| 223 |     for (int j=0;j<MDSteps;j++)
 | 
|---|
| 224 |       ptr->Trajectory.R.at(j) += (*trans);
 | 
|---|
| 225 |     ptr->x += (*trans);
 | 
|---|
| 226 |   }
 | 
|---|
| 227 | };
 | 
|---|
| 228 | 
 | 
|---|
| 229 | /** Translate the molecule periodically in the box.
 | 
|---|
| 230 |  * \param trans[] translation vector.
 | 
|---|
| 231 |  * TODO treatment of trajetories missing
 | 
|---|
| 232 |  */
 | 
|---|
| 233 | void molecule::TranslatePeriodically(const Vector *trans)
 | 
|---|
| 234 | {
 | 
|---|
| 235 |   double * const cell_size = World::getInstance().getDomain();
 | 
|---|
| 236 |   double *M = ReturnFullMatrixforSymmetric(cell_size);
 | 
|---|
| 237 |   double *Minv = InverseMatrix(M);
 | 
|---|
| 238 | 
 | 
|---|
| 239 |   // go through all atoms
 | 
|---|
| 240 |   ActOnAllVectors( &Vector::SubtractVector, *trans);
 | 
|---|
| 241 |   ActOnAllVectors( &Vector::WrapPeriodically, (const double *)M, (const double *)Minv);
 | 
|---|
| 242 | 
 | 
|---|
| 243 |   Free(&M);
 | 
|---|
| 244 |   Free(&Minv);
 | 
|---|
| 245 | };
 | 
|---|
| 246 | 
 | 
|---|
| 247 | 
 | 
|---|
| 248 | /** Mirrors all atoms against a given plane.
 | 
|---|
| 249 |  * \param n[] normal vector of mirror plane.
 | 
|---|
| 250 |  */
 | 
|---|
| 251 | void molecule::Mirror(const Vector *n)
 | 
|---|
| 252 | {
 | 
|---|
| 253 |   ActOnAllVectors( &Vector::Mirror, *n );
 | 
|---|
| 254 | };
 | 
|---|
| 255 | 
 | 
|---|
| 256 | /** Determines center of molecule (yet not considering atom masses).
 | 
|---|
| 257 |  * \param center reference to return vector
 | 
|---|
| 258 |  */
 | 
|---|
| 259 | void molecule::DeterminePeriodicCenter(Vector ¢er)
 | 
|---|
| 260 | {
 | 
|---|
| 261 |   atom *Walker = start;
 | 
|---|
| 262 |   double * const cell_size = World::getInstance().getDomain();
 | 
|---|
| 263 |   double *matrix = ReturnFullMatrixforSymmetric(cell_size);
 | 
|---|
| 264 |   double *inversematrix = InverseMatrix(cell_size);
 | 
|---|
| 265 |   double tmp;
 | 
|---|
| 266 |   bool flag;
 | 
|---|
| 267 |   Vector Testvector, Translationvector;
 | 
|---|
| 268 | 
 | 
|---|
| 269 |   do {
 | 
|---|
| 270 |     Center.Zero();
 | 
|---|
| 271 |     flag = true;
 | 
|---|
| 272 |     while (Walker->next != end) {
 | 
|---|
| 273 |       Walker = Walker->next;
 | 
|---|
| 274 | #ifdef ADDHYDROGEN
 | 
|---|
| 275 |       if (Walker->type->Z != 1) {
 | 
|---|
| 276 | #endif
 | 
|---|
| 277 |         Testvector = Walker->x;
 | 
|---|
| 278 |         Testvector.MatrixMultiplication(inversematrix);
 | 
|---|
| 279 |         Translationvector.Zero();
 | 
|---|
| 280 |         for (BondList::const_iterator Runner = Walker->ListOfBonds.begin(); Runner != Walker->ListOfBonds.end(); (++Runner)) {
 | 
|---|
| 281 |          if (Walker->nr < (*Runner)->GetOtherAtom(Walker)->nr) // otherwise we shift one to, the other fro and gain nothing
 | 
|---|
| 282 |             for (int j=0;j<NDIM;j++) {
 | 
|---|
| 283 |               tmp = Walker->x[j] - (*Runner)->GetOtherAtom(Walker)->x[j];
 | 
|---|
| 284 |               if ((fabs(tmp)) > BondDistance) {
 | 
|---|
| 285 |                 flag = false;
 | 
|---|
| 286 |                 DoLog(0) && (Log() << Verbose(0) << "Hit: atom " << Walker->Name << " in bond " << *(*Runner) << " has to be shifted due to " << tmp << "." << endl);
 | 
|---|
| 287 |                 if (tmp > 0)
 | 
|---|
| 288 |                   Translationvector[j] -= 1.;
 | 
|---|
| 289 |                 else
 | 
|---|
| 290 |                   Translationvector[j] += 1.;
 | 
|---|
| 291 |               }
 | 
|---|
| 292 |             }
 | 
|---|
| 293 |         }
 | 
|---|
| 294 |         Testvector += Translationvector;
 | 
|---|
| 295 |         Testvector.MatrixMultiplication(matrix);
 | 
|---|
| 296 |         Center += Testvector;
 | 
|---|
| 297 |         Log() << Verbose(1) << "vector is: " << Testvector << endl;
 | 
|---|
| 298 | #ifdef ADDHYDROGEN
 | 
|---|
| 299 |         // now also change all hydrogens
 | 
|---|
| 300 |         for (BondList::const_iterator Runner = Walker->ListOfBonds.begin(); Runner != Walker->ListOfBonds.end(); (++Runner)) {
 | 
|---|
| 301 |           if ((*Runner)->GetOtherAtom(Walker)->type->Z == 1) {
 | 
|---|
| 302 |             Testvector = (*Runner)->GetOtherAtom(Walker)->x;
 | 
|---|
| 303 |             Testvector.MatrixMultiplication(inversematrix);
 | 
|---|
| 304 |             Testvector += Translationvector;
 | 
|---|
| 305 |             Testvector.MatrixMultiplication(matrix);
 | 
|---|
| 306 |             Center += Testvector;
 | 
|---|
| 307 |             Log() << Verbose(1) << "Hydrogen vector is: " << Testvector << endl;
 | 
|---|
| 308 |           }
 | 
|---|
| 309 |         }
 | 
|---|
| 310 |       }
 | 
|---|
| 311 | #endif
 | 
|---|
| 312 |     }
 | 
|---|
| 313 |   } while (!flag);
 | 
|---|
| 314 |   Free(&matrix);
 | 
|---|
| 315 |   Free(&inversematrix);
 | 
|---|
| 316 | 
 | 
|---|
| 317 |   Center.Scale(1./(double)AtomCount);
 | 
|---|
| 318 | };
 | 
|---|
| 319 | 
 | 
|---|
| 320 | /** Transforms/Rotates the given molecule into its principal axis system.
 | 
|---|
| 321 |  * \param *out output stream for debugging
 | 
|---|
| 322 |  * \param DoRotate whether to rotate (true) or only to determine the PAS.
 | 
|---|
| 323 |  * TODO treatment of trajetories missing
 | 
|---|
| 324 |  */
 | 
|---|
| 325 | void molecule::PrincipalAxisSystem(bool DoRotate)
 | 
|---|
| 326 | {
 | 
|---|
| 327 |   atom *ptr = start;  // start at first in list
 | 
|---|
| 328 |   double InertiaTensor[NDIM*NDIM];
 | 
|---|
| 329 |   Vector *CenterOfGravity = DetermineCenterOfGravity();
 | 
|---|
| 330 | 
 | 
|---|
| 331 |   CenterPeriodic();
 | 
|---|
| 332 | 
 | 
|---|
| 333 |   // reset inertia tensor
 | 
|---|
| 334 |   for(int i=0;i<NDIM*NDIM;i++)
 | 
|---|
| 335 |     InertiaTensor[i] = 0.;
 | 
|---|
| 336 | 
 | 
|---|
| 337 |   // sum up inertia tensor
 | 
|---|
| 338 |   while (ptr->next != end) {
 | 
|---|
| 339 |     ptr = ptr->next;
 | 
|---|
| 340 |     Vector x = ptr->x;
 | 
|---|
| 341 |     //x.SubtractVector(CenterOfGravity);
 | 
|---|
| 342 |     InertiaTensor[0] += ptr->type->mass*(x[1]*x[1] + x[2]*x[2]);
 | 
|---|
| 343 |     InertiaTensor[1] += ptr->type->mass*(-x[0]*x[1]);
 | 
|---|
| 344 |     InertiaTensor[2] += ptr->type->mass*(-x[0]*x[2]);
 | 
|---|
| 345 |     InertiaTensor[3] += ptr->type->mass*(-x[1]*x[0]);
 | 
|---|
| 346 |     InertiaTensor[4] += ptr->type->mass*(x[0]*x[0] + x[2]*x[2]);
 | 
|---|
| 347 |     InertiaTensor[5] += ptr->type->mass*(-x[1]*x[2]);
 | 
|---|
| 348 |     InertiaTensor[6] += ptr->type->mass*(-x[2]*x[0]);
 | 
|---|
| 349 |     InertiaTensor[7] += ptr->type->mass*(-x[2]*x[1]);
 | 
|---|
| 350 |     InertiaTensor[8] += ptr->type->mass*(x[0]*x[0] + x[1]*x[1]);
 | 
|---|
| 351 |   }
 | 
|---|
| 352 |   // print InertiaTensor for debugging
 | 
|---|
| 353 |   DoLog(0) && (Log() << Verbose(0) << "The inertia tensor is:" << endl);
 | 
|---|
| 354 |   for(int i=0;i<NDIM;i++) {
 | 
|---|
| 355 |     for(int j=0;j<NDIM;j++)
 | 
|---|
| 356 |       DoLog(0) && (Log() << Verbose(0) << InertiaTensor[i*NDIM+j] << " ");
 | 
|---|
| 357 |     DoLog(0) && (Log() << Verbose(0) << endl);
 | 
|---|
| 358 |   }
 | 
|---|
| 359 |   DoLog(0) && (Log() << Verbose(0) << endl);
 | 
|---|
| 360 | 
 | 
|---|
| 361 |   // diagonalize to determine principal axis system
 | 
|---|
| 362 |   gsl_eigen_symmv_workspace *T = gsl_eigen_symmv_alloc(NDIM);
 | 
|---|
| 363 |   gsl_matrix_view m = gsl_matrix_view_array(InertiaTensor, NDIM, NDIM);
 | 
|---|
| 364 |   gsl_vector *eval = gsl_vector_alloc(NDIM);
 | 
|---|
| 365 |   gsl_matrix *evec = gsl_matrix_alloc(NDIM, NDIM);
 | 
|---|
| 366 |   gsl_eigen_symmv(&m.matrix, eval, evec, T);
 | 
|---|
| 367 |   gsl_eigen_symmv_free(T);
 | 
|---|
| 368 |   gsl_eigen_symmv_sort(eval, evec, GSL_EIGEN_SORT_ABS_DESC);
 | 
|---|
| 369 | 
 | 
|---|
| 370 |   for(int i=0;i<NDIM;i++) {
 | 
|---|
| 371 |     DoLog(1) && (Log() << Verbose(1) << "eigenvalue = " << gsl_vector_get(eval, i));
 | 
|---|
| 372 |     DoLog(0) && (Log() << Verbose(0) << ", eigenvector = (" << evec->data[i * evec->tda + 0] << "," << evec->data[i * evec->tda + 1] << "," << evec->data[i * evec->tda + 2] << ")" << endl);
 | 
|---|
| 373 |   }
 | 
|---|
| 374 | 
 | 
|---|
| 375 |   // check whether we rotate or not
 | 
|---|
| 376 |   if (DoRotate) {
 | 
|---|
| 377 |     DoLog(1) && (Log() << Verbose(1) << "Transforming molecule into PAS ... ");
 | 
|---|
| 378 |     // the eigenvectors specify the transformation matrix
 | 
|---|
| 379 |     ActOnAllVectors( &Vector::MatrixMultiplication, (const double *) evec->data );
 | 
|---|
| 380 |     DoLog(0) && (Log() << Verbose(0) << "done." << endl);
 | 
|---|
| 381 | 
 | 
|---|
| 382 |     // summing anew for debugging (resulting matrix has to be diagonal!)
 | 
|---|
| 383 |     // reset inertia tensor
 | 
|---|
| 384 |     for(int i=0;i<NDIM*NDIM;i++)
 | 
|---|
| 385 |       InertiaTensor[i] = 0.;
 | 
|---|
| 386 | 
 | 
|---|
| 387 |     // sum up inertia tensor
 | 
|---|
| 388 |     ptr = start;
 | 
|---|
| 389 |     while (ptr->next != end) {
 | 
|---|
| 390 |       ptr = ptr->next;
 | 
|---|
| 391 |       Vector x = ptr->x;
 | 
|---|
| 392 |       //x.SubtractVector(CenterOfGravity);
 | 
|---|
| 393 |       InertiaTensor[0] += ptr->type->mass*(x[1]*x[1] + x[2]*x[2]);
 | 
|---|
| 394 |       InertiaTensor[1] += ptr->type->mass*(-x[0]*x[1]);
 | 
|---|
| 395 |       InertiaTensor[2] += ptr->type->mass*(-x[0]*x[2]);
 | 
|---|
| 396 |       InertiaTensor[3] += ptr->type->mass*(-x[1]*x[0]);
 | 
|---|
| 397 |       InertiaTensor[4] += ptr->type->mass*(x[0]*x[0] + x[2]*x[2]);
 | 
|---|
| 398 |       InertiaTensor[5] += ptr->type->mass*(-x[1]*x[2]);
 | 
|---|
| 399 |       InertiaTensor[6] += ptr->type->mass*(-x[2]*x[0]);
 | 
|---|
| 400 |       InertiaTensor[7] += ptr->type->mass*(-x[2]*x[1]);
 | 
|---|
| 401 |       InertiaTensor[8] += ptr->type->mass*(x[0]*x[0] + x[1]*x[1]);
 | 
|---|
| 402 |     }
 | 
|---|
| 403 |     // print InertiaTensor for debugging
 | 
|---|
| 404 |     DoLog(0) && (Log() << Verbose(0) << "The inertia tensor is:" << endl);
 | 
|---|
| 405 |     for(int i=0;i<NDIM;i++) {
 | 
|---|
| 406 |       for(int j=0;j<NDIM;j++)
 | 
|---|
| 407 |         DoLog(0) && (Log() << Verbose(0) << InertiaTensor[i*NDIM+j] << " ");
 | 
|---|
| 408 |       DoLog(0) && (Log() << Verbose(0) << endl);
 | 
|---|
| 409 |     }
 | 
|---|
| 410 |     DoLog(0) && (Log() << Verbose(0) << endl);
 | 
|---|
| 411 |   }
 | 
|---|
| 412 | 
 | 
|---|
| 413 |   // free everything
 | 
|---|
| 414 |   delete(CenterOfGravity);
 | 
|---|
| 415 |   gsl_vector_free(eval);
 | 
|---|
| 416 |   gsl_matrix_free(evec);
 | 
|---|
| 417 | };
 | 
|---|
| 418 | 
 | 
|---|
| 419 | 
 | 
|---|
| 420 | /** Align all atoms in such a manner that given vector \a *n is along z axis.
 | 
|---|
| 421 |  * \param n[] alignment vector.
 | 
|---|
| 422 |  */
 | 
|---|
| 423 | void molecule::Align(Vector *n)
 | 
|---|
| 424 | {
 | 
|---|
| 425 |   atom *ptr = start;
 | 
|---|
| 426 |   double alpha, tmp;
 | 
|---|
| 427 |   Vector z_axis;
 | 
|---|
| 428 |   z_axis[0] = 0.;
 | 
|---|
| 429 |   z_axis[1] = 0.;
 | 
|---|
| 430 |   z_axis[2] = 1.;
 | 
|---|
| 431 | 
 | 
|---|
| 432 |   // rotate on z-x plane
 | 
|---|
| 433 |   DoLog(0) && (Log() << Verbose(0) << "Begin of Aligning all atoms." << endl);
 | 
|---|
| 434 |   alpha = atan(-n->at(0)/n->at(2));
 | 
|---|
| 435 |   DoLog(1) && (Log() << Verbose(1) << "Z-X-angle: " << alpha << " ... ");
 | 
|---|
| 436 |   while (ptr->next != end) {
 | 
|---|
| 437 |     ptr = ptr->next;
 | 
|---|
| 438 |     tmp = ptr->x[0];
 | 
|---|
| 439 |     ptr->x[0] =  cos(alpha) * tmp + sin(alpha) * ptr->x[2];
 | 
|---|
| 440 |     ptr->x[2] = -sin(alpha) * tmp + cos(alpha) * ptr->x[2];
 | 
|---|
| 441 |     for (int j=0;j<MDSteps;j++) {
 | 
|---|
| 442 |       tmp = ptr->Trajectory.R.at(j)[0];
 | 
|---|
| 443 |       ptr->Trajectory.R.at(j)[0] =  cos(alpha) * tmp + sin(alpha) * ptr->Trajectory.R.at(j)[2];
 | 
|---|
| 444 |       ptr->Trajectory.R.at(j)[2] = -sin(alpha) * tmp + cos(alpha) * ptr->Trajectory.R.at(j)[2];
 | 
|---|
| 445 |     }
 | 
|---|
| 446 |   }
 | 
|---|
| 447 |   // rotate n vector
 | 
|---|
| 448 |   tmp = n->at(0);
 | 
|---|
| 449 |   n->at(0) =  cos(alpha) * tmp +  sin(alpha) * n->at(2);
 | 
|---|
| 450 |   n->at(2) = -sin(alpha) * tmp +  cos(alpha) * n->at(2);
 | 
|---|
| 451 |   DoLog(1) && (Log() << Verbose(1) << "alignment vector after first rotation: " << n << endl);
 | 
|---|
| 452 | 
 | 
|---|
| 453 |   // rotate on z-y plane
 | 
|---|
| 454 |   ptr = start;
 | 
|---|
| 455 |   alpha = atan(-n->at(1)/n->at(2));
 | 
|---|
| 456 |   DoLog(1) && (Log() << Verbose(1) << "Z-Y-angle: " << alpha << " ... ");
 | 
|---|
| 457 |   while (ptr->next != end) {
 | 
|---|
| 458 |     ptr = ptr->next;
 | 
|---|
| 459 |     tmp = ptr->x[1];
 | 
|---|
| 460 |     ptr->x[1] =  cos(alpha) * tmp + sin(alpha) * ptr->x[2];
 | 
|---|
| 461 |     ptr->x[2] = -sin(alpha) * tmp + cos(alpha) * ptr->x[2];
 | 
|---|
| 462 |     for (int j=0;j<MDSteps;j++) {
 | 
|---|
| 463 |       tmp = ptr->Trajectory.R.at(j)[1];
 | 
|---|
| 464 |       ptr->Trajectory.R.at(j)[1] =  cos(alpha) * tmp + sin(alpha) * ptr->Trajectory.R.at(j)[2];
 | 
|---|
| 465 |       ptr->Trajectory.R.at(j)[2] = -sin(alpha) * tmp + cos(alpha) * ptr->Trajectory.R.at(j)[2];
 | 
|---|
| 466 |     }
 | 
|---|
| 467 |   }
 | 
|---|
| 468 |   // rotate n vector (for consistency check)
 | 
|---|
| 469 |   tmp = n->at(1);
 | 
|---|
| 470 |   n->at(1) =  cos(alpha) * tmp +  sin(alpha) * n->at(2);
 | 
|---|
| 471 |   n->at(2) = -sin(alpha) * tmp +  cos(alpha) * n->at(2);
 | 
|---|
| 472 | 
 | 
|---|
| 473 | 
 | 
|---|
| 474 |   DoLog(1) && (Log() << Verbose(1) << "alignment vector after second rotation: " << n << endl);
 | 
|---|
| 475 |   DoLog(0) && (Log() << Verbose(0) << "End of Aligning all atoms." << endl);
 | 
|---|
| 476 | };
 | 
|---|
| 477 | 
 | 
|---|
| 478 | 
 | 
|---|
| 479 | /** Calculates sum over least square distance to line hidden in \a *x.
 | 
|---|
| 480 |  * \param *x offset and direction vector
 | 
|---|
| 481 |  * \param *params pointer to lsq_params structure
 | 
|---|
| 482 |  * \return \f$ sum_i^N | y_i - (a + t_i b)|^2\f$
 | 
|---|
| 483 |  */
 | 
|---|
| 484 | double LeastSquareDistance (const gsl_vector * x, void * params)
 | 
|---|
| 485 | {
 | 
|---|
| 486 |   double res = 0, t;
 | 
|---|
| 487 |   Vector a,b,c,d;
 | 
|---|
| 488 |   struct lsq_params *par = (struct lsq_params *)params;
 | 
|---|
| 489 |   atom *ptr = par->mol->start;
 | 
|---|
| 490 | 
 | 
|---|
| 491 |   // initialize vectors
 | 
|---|
| 492 |   a[0] = gsl_vector_get(x,0);
 | 
|---|
| 493 |   a[1] = gsl_vector_get(x,1);
 | 
|---|
| 494 |   a[2] = gsl_vector_get(x,2);
 | 
|---|
| 495 |   b[0] = gsl_vector_get(x,3);
 | 
|---|
| 496 |   b[1] = gsl_vector_get(x,4);
 | 
|---|
| 497 |   b[2] = gsl_vector_get(x,5);
 | 
|---|
| 498 |   // go through all atoms
 | 
|---|
| 499 |   while (ptr != par->mol->end) {
 | 
|---|
| 500 |     ptr = ptr->next;
 | 
|---|
| 501 |     if (ptr->type == ((struct lsq_params *)params)->type) { // for specific type
 | 
|---|
| 502 |       c = ptr->x - a;
 | 
|---|
| 503 |       t = c.ScalarProduct(b);           // get direction parameter
 | 
|---|
| 504 |       d = t*b;       // and create vector
 | 
|---|
| 505 |       c -= d;   // ... yielding distance vector
 | 
|---|
| 506 |       res += d.ScalarProduct(d);        // add squared distance
 | 
|---|
| 507 |     }
 | 
|---|
| 508 |   }
 | 
|---|
| 509 |   return res;
 | 
|---|
| 510 | };
 | 
|---|
| 511 | 
 | 
|---|
| 512 | /** By minimizing the least square distance gains alignment vector.
 | 
|---|
| 513 |  * \bug this is not yet working properly it seems
 | 
|---|
| 514 |  */
 | 
|---|
| 515 | void molecule::GetAlignvector(struct lsq_params * par) const
 | 
|---|
| 516 | {
 | 
|---|
| 517 |     int np = 6;
 | 
|---|
| 518 | 
 | 
|---|
| 519 |    const gsl_multimin_fminimizer_type *T =
 | 
|---|
| 520 |      gsl_multimin_fminimizer_nmsimplex;
 | 
|---|
| 521 |    gsl_multimin_fminimizer *s = NULL;
 | 
|---|
| 522 |    gsl_vector *ss;
 | 
|---|
| 523 |    gsl_multimin_function minex_func;
 | 
|---|
| 524 | 
 | 
|---|
| 525 |    size_t iter = 0, i;
 | 
|---|
| 526 |    int status;
 | 
|---|
| 527 |    double size;
 | 
|---|
| 528 | 
 | 
|---|
| 529 |    /* Initial vertex size vector */
 | 
|---|
| 530 |    ss = gsl_vector_alloc (np);
 | 
|---|
| 531 | 
 | 
|---|
| 532 |    /* Set all step sizes to 1 */
 | 
|---|
| 533 |    gsl_vector_set_all (ss, 1.0);
 | 
|---|
| 534 | 
 | 
|---|
| 535 |    /* Starting point */
 | 
|---|
| 536 |    par->x = gsl_vector_alloc (np);
 | 
|---|
| 537 |    par->mol = this;
 | 
|---|
| 538 | 
 | 
|---|
| 539 |    gsl_vector_set (par->x, 0, 0.0);  // offset
 | 
|---|
| 540 |    gsl_vector_set (par->x, 1, 0.0);
 | 
|---|
| 541 |    gsl_vector_set (par->x, 2, 0.0);
 | 
|---|
| 542 |    gsl_vector_set (par->x, 3, 0.0);  // direction
 | 
|---|
| 543 |    gsl_vector_set (par->x, 4, 0.0);
 | 
|---|
| 544 |    gsl_vector_set (par->x, 5, 1.0);
 | 
|---|
| 545 | 
 | 
|---|
| 546 |    /* Initialize method and iterate */
 | 
|---|
| 547 |    minex_func.f = &LeastSquareDistance;
 | 
|---|
| 548 |    minex_func.n = np;
 | 
|---|
| 549 |    minex_func.params = (void *)par;
 | 
|---|
| 550 | 
 | 
|---|
| 551 |    s = gsl_multimin_fminimizer_alloc (T, np);
 | 
|---|
| 552 |    gsl_multimin_fminimizer_set (s, &minex_func, par->x, ss);
 | 
|---|
| 553 | 
 | 
|---|
| 554 |    do
 | 
|---|
| 555 |      {
 | 
|---|
| 556 |        iter++;
 | 
|---|
| 557 |        status = gsl_multimin_fminimizer_iterate(s);
 | 
|---|
| 558 | 
 | 
|---|
| 559 |        if (status)
 | 
|---|
| 560 |          break;
 | 
|---|
| 561 | 
 | 
|---|
| 562 |        size = gsl_multimin_fminimizer_size (s);
 | 
|---|
| 563 |        status = gsl_multimin_test_size (size, 1e-2);
 | 
|---|
| 564 | 
 | 
|---|
| 565 |        if (status == GSL_SUCCESS)
 | 
|---|
| 566 |          {
 | 
|---|
| 567 |            printf ("converged to minimum at\n");
 | 
|---|
| 568 |          }
 | 
|---|
| 569 | 
 | 
|---|
| 570 |        printf ("%5d ", (int)iter);
 | 
|---|
| 571 |        for (i = 0; i < (size_t)np; i++)
 | 
|---|
| 572 |          {
 | 
|---|
| 573 |            printf ("%10.3e ", gsl_vector_get (s->x, i));
 | 
|---|
| 574 |          }
 | 
|---|
| 575 |        printf ("f() = %7.3f size = %.3f\n", s->fval, size);
 | 
|---|
| 576 |      }
 | 
|---|
| 577 |    while (status == GSL_CONTINUE && iter < 100);
 | 
|---|
| 578 | 
 | 
|---|
| 579 |   for (i=0;i<(size_t)np;i++)
 | 
|---|
| 580 |     gsl_vector_set(par->x, i, gsl_vector_get(s->x, i));
 | 
|---|
| 581 |    //gsl_vector_free(par->x);
 | 
|---|
| 582 |    gsl_vector_free(ss);
 | 
|---|
| 583 |    gsl_multimin_fminimizer_free (s);
 | 
|---|
| 584 | };
 | 
|---|