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