ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/SHAPES/RigidBody.cpp
Revision: 1279
Committed: Fri Jun 18 19:01:40 2004 UTC (20 years ago) by chrisfen
File size: 12661 byte(s)
Log Message:
Built principle axis determination and molecule rotation into calcRefCoords in RigidBody.cpp

File Contents

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