ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/SHAPES/RigidBody.cpp
Revision: 1281
Committed: Mon Jun 21 13:38:55 2004 UTC (20 years ago) by chrisfen
File size: 12995 byte(s)
Log Message:
Now calculates energies and builds the grid files.  There are still bugs, but
it compiles and runs...

File Contents

# User Rev Content
1 gezelter 1271 #include <math.h>
2 chrisfen 1281 #include <iostream>
3 gezelter 1271 #include "RigidBody.hpp"
4     #include "VDWAtom.hpp"
5     #include "MatVec3.h"
6    
7     RigidBody::RigidBody() {
8     is_linear = false;
9     linear_axis = -1;
10     momIntTol = 1e-6;
11     }
12    
13     RigidBody::~RigidBody() {
14     }
15    
16     void RigidBody::addAtom(VDWAtom* at) {
17    
18     vec3 coords;
19    
20     myAtoms.push_back(at);
21    
22     at->getPos(coords.vec);
23     refCoords.push_back(coords);
24     }
25    
26     void RigidBody::getPos(double theP[3]){
27     for (int i = 0; i < 3 ; i++)
28     theP[i] = pos[i];
29     }
30    
31     void RigidBody::setPos(double theP[3]){
32     for (int i = 0; i < 3 ; i++)
33     pos[i] = theP[i];
34     }
35    
36    
37     void RigidBody::setEuler( double phi, double theta, double psi ){
38    
39     A[0][0] = (cos(phi) * cos(psi)) - (sin(phi) * cos(theta) * sin(psi));
40     A[0][1] = (sin(phi) * cos(psi)) + (cos(phi) * cos(theta) * sin(psi));
41     A[0][2] = sin(theta) * sin(psi);
42    
43     A[1][0] = -(cos(phi) * sin(psi)) - (sin(phi) * cos(theta) * cos(psi));
44     A[1][1] = -(sin(phi) * sin(psi)) + (cos(phi) * cos(theta) * cos(psi));
45     A[1][2] = sin(theta) * cos(psi);
46    
47     A[2][0] = sin(phi) * sin(theta);
48     A[2][1] = -cos(phi) * sin(theta);
49     A[2][2] = cos(theta);
50    
51     }
52    
53     void RigidBody::getQ( double q[4] ){
54    
55     double t, s;
56     double ad1, ad2, ad3;
57    
58     t = A[0][0] + A[1][1] + A[2][2] + 1.0;
59     if( t > 0.0 ){
60    
61     s = 0.5 / sqrt( t );
62     q[0] = 0.25 / s;
63     q[1] = (A[1][2] - A[2][1]) * s;
64     q[2] = (A[2][0] - A[0][2]) * s;
65     q[3] = (A[0][1] - A[1][0]) * s;
66     }
67     else{
68    
69     ad1 = fabs( A[0][0] );
70     ad2 = fabs( A[1][1] );
71     ad3 = fabs( A[2][2] );
72    
73     if( ad1 >= ad2 && ad1 >= ad3 ){
74    
75     s = 2.0 * sqrt( 1.0 + A[0][0] - A[1][1] - A[2][2] );
76     q[0] = (A[1][2] + A[2][1]) / s;
77     q[1] = 0.5 / s;
78     q[2] = (A[0][1] + A[1][0]) / s;
79     q[3] = (A[0][2] + A[2][0]) / s;
80     }
81     else if( ad2 >= ad1 && ad2 >= ad3 ){
82    
83     s = sqrt( 1.0 + A[1][1] - A[0][0] - A[2][2] ) * 2.0;
84     q[0] = (A[0][2] + A[2][0]) / s;
85     q[1] = (A[0][1] + A[1][0]) / s;
86     q[2] = 0.5 / s;
87     q[3] = (A[1][2] + A[2][1]) / s;
88     }
89     else{
90    
91     s = sqrt( 1.0 + A[2][2] - A[0][0] - A[1][1] ) * 2.0;
92     q[0] = (A[0][1] + A[1][0]) / s;
93     q[1] = (A[0][2] + A[2][0]) / s;
94     q[2] = (A[1][2] + A[2][1]) / s;
95     q[3] = 0.5 / s;
96     }
97     }
98     }
99    
100     void RigidBody::setQ( double the_q[4] ){
101    
102     double q0Sqr, q1Sqr, q2Sqr, q3Sqr;
103    
104     q0Sqr = the_q[0] * the_q[0];
105     q1Sqr = the_q[1] * the_q[1];
106     q2Sqr = the_q[2] * the_q[2];
107     q3Sqr = the_q[3] * the_q[3];
108    
109     A[0][0] = q0Sqr + q1Sqr - q2Sqr - q3Sqr;
110     A[0][1] = 2.0 * ( the_q[1] * the_q[2] + the_q[0] * the_q[3] );
111     A[0][2] = 2.0 * ( the_q[1] * the_q[3] - the_q[0] * the_q[2] );
112    
113     A[1][0] = 2.0 * ( the_q[1] * the_q[2] - the_q[0] * the_q[3] );
114     A[1][1] = q0Sqr - q1Sqr + q2Sqr - q3Sqr;
115     A[1][2] = 2.0 * ( the_q[2] * the_q[3] + the_q[0] * the_q[1] );
116    
117     A[2][0] = 2.0 * ( the_q[1] * the_q[3] + the_q[0] * the_q[2] );
118     A[2][1] = 2.0 * ( the_q[2] * the_q[3] - the_q[0] * the_q[1] );
119     A[2][2] = q0Sqr - q1Sqr -q2Sqr +q3Sqr;
120    
121     }
122    
123     void RigidBody::getA( double the_A[3][3] ){
124    
125     for (int i = 0; i < 3; i++)
126     for (int j = 0; j < 3; j++)
127     the_A[i][j] = A[i][j];
128    
129     }
130    
131     void RigidBody::setA( double the_A[3][3] ){
132    
133     for (int i = 0; i < 3; i++)
134     for (int j = 0; j < 3; j++)
135     A[i][j] = the_A[i][j];
136    
137     }
138    
139     void RigidBody::getI( double the_I[3][3] ){
140    
141     for (int i = 0; i < 3; i++)
142     for (int j = 0; j < 3; j++)
143     the_I[i][j] = I[i][j];
144    
145     }
146    
147     void RigidBody::lab2Body( double r[3] ){
148    
149     double rl[3]; // the lab frame vector
150    
151     rl[0] = r[0];
152     rl[1] = r[1];
153     rl[2] = r[2];
154    
155     r[0] = (A[0][0] * rl[0]) + (A[0][1] * rl[1]) + (A[0][2] * rl[2]);
156     r[1] = (A[1][0] * rl[0]) + (A[1][1] * rl[1]) + (A[1][2] * rl[2]);
157     r[2] = (A[2][0] * rl[0]) + (A[2][1] * rl[1]) + (A[2][2] * rl[2]);
158    
159     }
160    
161     void RigidBody::body2Lab( double r[3] ){
162    
163     double rb[3]; // the body frame vector
164    
165     rb[0] = r[0];
166     rb[1] = r[1];
167     rb[2] = r[2];
168    
169     r[0] = (A[0][0] * rb[0]) + (A[1][0] * rb[1]) + (A[2][0] * rb[2]);
170     r[1] = (A[0][1] * rb[0]) + (A[1][1] * rb[1]) + (A[2][1] * rb[2]);
171     r[2] = (A[0][2] * rb[0]) + (A[1][2] * rb[1]) + (A[2][2] * rb[2]);
172    
173     }
174    
175     void RigidBody::calcRefCoords( ) {
176    
177 chrisfen 1279 int i, j, it, n_linear_coords, pAxis, maxAxis, midAxis;
178 gezelter 1271 double mtmp;
179     vec3 apos;
180     double refCOM[3];
181     vec3 ptmp;
182     double Itmp[3][3];
183 chrisfen 1279 double pAxisMat[3][3], pAxisRotMat[3][3];
184 gezelter 1271 double evals[3];
185     double r, r2, len;
186 chrisfen 1279 double iMat[3][3];
187 chrisfen 1281 double test[3];
188    
189 gezelter 1271 // First, find the center of mass:
190    
191     mass = 0.0;
192     for (j=0; j<3; j++)
193     refCOM[j] = 0.0;
194    
195     for (i = 0; i < myAtoms.size(); i++) {
196     mtmp = myAtoms[i]->getMass();
197     mass += mtmp;
198    
199     apos = refCoords[i];
200     for(j = 0; j < 3; j++) {
201     refCOM[j] += apos[j]*mtmp;
202     }
203     }
204    
205     for(j = 0; j < 3; j++)
206     refCOM[j] /= mass;
207    
208     // Next, move the origin of the reference coordinate system to the COM:
209    
210     for (i = 0; i < myAtoms.size(); i++) {
211     apos = refCoords[i];
212     for (j=0; j < 3; j++) {
213     apos[j] = apos[j] - refCOM[j];
214     }
215     refCoords[i] = apos;
216     }
217    
218     // Moment of Inertia calculation
219    
220     for (i = 0; i < 3; i++)
221     for (j = 0; j < 3; j++)
222     Itmp[i][j] = 0.0;
223    
224     for (it = 0; it < myAtoms.size(); it++) {
225    
226     mtmp = myAtoms[it]->getMass();
227     ptmp = refCoords[it];
228     r= norm3(ptmp.vec);
229     r2 = r*r;
230    
231     for (i = 0; i < 3; i++) {
232     for (j = 0; j < 3; j++) {
233    
234     if (i==j) Itmp[i][j] += mtmp * r2;
235    
236     Itmp[i][j] -= mtmp * ptmp.vec[i]*ptmp.vec[j];
237     }
238     }
239     }
240    
241     diagonalize3x3(Itmp, evals, sU);
242    
243     // zero out I and then fill the diagonals with the moments of inertia:
244    
245     n_linear_coords = 0;
246    
247     for (i = 0; i < 3; i++) {
248     for (j = 0; j < 3; j++) {
249     I[i][j] = 0.0;
250     }
251     I[i][i] = evals[i];
252    
253     if (fabs(evals[i]) < momIntTol) {
254     is_linear = true;
255     n_linear_coords++;
256     linear_axis = i;
257     }
258     }
259    
260     if (n_linear_coords > 1) {
261     printf(
262     "RigidBody error.\n"
263 chrisfen 1281 "\tSHAPES found more than one axis in this rigid body with a vanishing \n"
264 gezelter 1271 "\tmoment of inertia. This can happen in one of three ways:\n"
265     "\t 1) Only one atom was specified, or \n"
266     "\t 2) All atoms were specified at the same location, or\n"
267     "\t 3) The programmers did something stupid.\n"
268     "\tIt is silly to use a rigid body to describe this situation. Be smarter.\n"
269     );
270     exit(-1);
271     }
272 chrisfen 1281
273 gezelter 1271 // renormalize column vectors:
274    
275     for (i=0; i < 3; i++) {
276     len = 0.0;
277     for (j = 0; j < 3; j++) {
278     len += sU[i][j]*sU[i][j];
279     }
280     len = sqrt(len);
281     for (j = 0; j < 3; j++) {
282     sU[i][j] /= len;
283     }
284     }
285 chrisfen 1281
286     //sort and reorder the moment axes
287    
288     // The only problem below is for molecules like C60 with 3 nearly identical
289     // non-zero moments of inertia. In this case it doesn't really matter which is
290     // the principal axis, so they get assigned nearly randomly depending on the
291     // floating point comparison between eigenvalues
292     if (! is_linear) {
293     pAxis = 0;
294     maxAxis = 0;
295    
296     for (i = 0; i < 3; i++) {
297     if (evals[i] < evals[pAxis]) pAxis = i;
298     if (evals[i] > evals[maxAxis]) maxAxis = i;
299     }
300    
301     midAxis = 0;
302     for (i=0; i < 3; i++) {
303     if (pAxis != i && maxAxis != i) midAxis = i;
304     }
305     } else {
306     pAxis = linear_axis;
307     // linear molecules have one zero moment of inertia and two identical
308     // moments of inertia. In this case, it doesn't matter which is chosen
309     // as mid and which is max, so just permute from the pAxis:
310     midAxis = (pAxis + 1)%3;
311     maxAxis = (pAxis + 2)%3;
312     }
313    
314     //let z be the smallest and x be the largest eigenvalue axes
315     for (i=0; i<3; i++){
316     pAxisMat[i][2] = sU[i][pAxis];
317     pAxisMat[i][1] = sU[i][midAxis];
318     pAxisMat[i][0] = sU[i][maxAxis];
319     }
320    
321     //calculate the proper rotation matrix
322     transposeMat3(pAxisMat, pAxisRotMat);
323    
324     for (i=0; i<myAtoms.size(); i++){
325     getAtomPos(test, i);
326     printf("%d\t%d\t%d\n",test[0],test[1],test[2]);
327     }
328    
329     //rotate the rigid body to the principle axis frame
330     for (i = 0; i < myAtoms.size(); i++) {
331     matVecMul3(pAxisRotMat, refCoords[i].vec, refCoords[i].vec);
332     myAtoms[i]->setPos(refCoords[i].vec);
333     }
334    
335     for (i=0; i<myAtoms.size(); i++){
336     getAtomPos(test,i);
337     printf("%d\t%d\t%d\n",test[0],test[1],test[2]);
338     }
339    
340     identityMat3(iMat);
341     setA(iMat);
342 gezelter 1271 }
343    
344 chrisfen 1276 void RigidBody::doEulerToRotMat(double euler[3], double myA[3][3] ){
345 gezelter 1271
346     double phi, theta, psi;
347    
348     phi = euler[0];
349     theta = euler[1];
350     psi = euler[2];
351    
352     myA[0][0] = (cos(phi) * cos(psi)) - (sin(phi) * cos(theta) * sin(psi));
353     myA[0][1] = (sin(phi) * cos(psi)) + (cos(phi) * cos(theta) * sin(psi));
354     myA[0][2] = sin(theta) * sin(psi);
355    
356     myA[1][0] = -(cos(phi) * sin(psi)) - (sin(phi) * cos(theta) * cos(psi));
357     myA[1][1] = -(sin(phi) * sin(psi)) + (cos(phi) * cos(theta) * cos(psi));
358     myA[1][2] = sin(theta) * cos(psi);
359    
360     myA[2][0] = sin(phi) * sin(theta);
361     myA[2][1] = -cos(phi) * sin(theta);
362     myA[2][2] = cos(theta);
363    
364     }
365    
366     void RigidBody::updateAtoms() {
367     int i, j;
368     vec3 ref;
369     double apos[3];
370    
371     for (i = 0; i < myAtoms.size(); i++) {
372    
373     ref = refCoords[i];
374    
375     body2Lab(ref.vec);
376    
377     for (j = 0; j<3; j++)
378     apos[j] = pos[j] + ref.vec[j];
379    
380     myAtoms[i]->setPos(apos);
381    
382     }
383     }
384    
385     /**
386     * getEulerAngles computes a set of Euler angle values consistent
387     * with an input rotation matrix. They are returned in the following
388     * order:
389     * myEuler[0] = phi;
390     * myEuler[1] = theta;
391     * myEuler[2] = psi;
392     */
393     void RigidBody::getEulerAngles(double myEuler[3]) {
394    
395     // We use so-called "x-convention", which is the most common
396     // definition. In this convention, the rotation given by Euler
397     // angles (phi, theta, psi), where the first rotation is by an angle
398     // phi about the z-axis, the second is by an angle theta (0 <= theta
399     // <= 180) about the x-axis, and the third is by an angle psi about
400     // the z-axis (again).
401    
402    
403     double phi,theta,psi,eps;
404 chrisfen 1276 double ctheta;
405     double stheta;
406    
407 gezelter 1271 // set the tolerance for Euler angles and rotation elements
408    
409     eps = 1.0e-8;
410    
411     theta = acos(min(1.0,max(-1.0,A[2][2])));
412     ctheta = A[2][2];
413     stheta = sqrt(1.0 - ctheta * ctheta);
414    
415     // when sin(theta) is close to 0, we need to consider the
416     // possibility of a singularity. In this case, we can assign an
417     // arbitary value to phi (or psi), and then determine the psi (or
418     // phi) or vice-versa. We'll assume that phi always gets the
419     // rotation, and psi is 0 in cases of singularity. we use atan2
420     // instead of atan, since atan2 will give us -Pi to Pi. Since 0 <=
421     // theta <= 180, sin(theta) will be always non-negative. Therefore,
422     // it never changes the sign of both of the parameters passed to
423     // atan2.
424    
425     if (fabs(stheta) <= eps){
426     psi = 0.0;
427     phi = atan2(-A[1][0], A[0][0]);
428     }
429     // we only have one unique solution
430     else{
431     phi = atan2(A[2][0], -A[2][1]);
432     psi = atan2(A[0][2], A[1][2]);
433     }
434    
435     //wrap phi and psi, make sure they are in the range from 0 to 2*Pi
436     //if (phi < 0)
437     // phi += M_PI;
438    
439     //if (psi < 0)
440     // psi += M_PI;
441    
442     myEuler[0] = phi;
443     myEuler[1] = theta;
444     myEuler[2] = psi;
445    
446     return;
447     }
448    
449     double RigidBody::max(double x, double y) {
450     return (x > y) ? x : y;
451     }
452    
453     double RigidBody::min(double x, double y) {
454     return (x > y) ? y : x;
455     }
456    
457 chrisfen 1276 double RigidBody::findMaxExtent(){
458     int i;
459     double refAtomPos[3];
460     double maxExtent;
461     double tempExtent;
462    
463     //zero the extent variables
464     maxExtent = 0.0;
465     tempExtent = 0.0;
466     for (i=0; i<3; i++)
467     refAtomPos[i] = 0.0;
468    
469     //loop over all atoms
470     for (i=0; i<myAtoms.size(); i++){
471     getAtomRefCoor(refAtomPos, i);
472     tempExtent = sqrt(refAtomPos[0]*refAtomPos[0] + refAtomPos[1]*refAtomPos[1]
473     + refAtomPos[2]*refAtomPos[2]);
474     if (tempExtent > maxExtent)
475     maxExtent = tempExtent;
476     }
477     return maxExtent;
478     }
479    
480 gezelter 1271 void RigidBody::findCOM() {
481    
482     size_t i;
483     int j;
484     double mtmp;
485     double ptmp[3];
486 chrisfen 1276
487 gezelter 1271 for(j = 0; j < 3; j++) {
488     pos[j] = 0.0;
489     }
490     mass = 0.0;
491    
492     for (i = 0; i < myAtoms.size(); i++) {
493    
494     mtmp = myAtoms[i]->getMass();
495     myAtoms[i]->getPos(ptmp);
496    
497     mass += mtmp;
498    
499     for(j = 0; j < 3; j++) {
500     pos[j] += ptmp[j]*mtmp;
501     }
502    
503     }
504    
505     for(j = 0; j < 3; j++) {
506     pos[j] /= mass;
507     }
508    
509     }
510    
511     void RigidBody::getAtomPos(double theP[3], int index){
512     vec3 ref;
513    
514     if (index >= myAtoms.size())
515     printf( "%d is an invalid index, current rigid body contains "
516     "%d atoms\n", index, myAtoms.size());
517    
518     ref = refCoords[index];
519     body2Lab(ref.vec);
520    
521     theP[0] = pos[0] + ref[0];
522     theP[1] = pos[1] + ref[1];
523     theP[2] = pos[2] + ref[2];
524     }
525    
526    
527     void RigidBody::getAtomRefCoor(double pos[3], int index){
528     vec3 ref;
529    
530     ref = refCoords[index];
531     pos[0] = ref[0];
532     pos[1] = ref[1];
533     pos[2] = ref[2];
534    
535     }
536 chrisfen 1281
537     double RigidBody::getAtomRpar(int index){
538    
539     return myAtoms[index]->getRpar();
540    
541     }
542    
543     double RigidBody::getAtomEps(int index){
544    
545     return myAtoms[index]->getEps();
546    
547     }