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