| 1 | #include <math.h> | 
| 2 | #include "RigidBody.hpp" | 
| 3 | #include "DirectionalAtom.hpp" | 
| 4 | #include "simError.h" | 
| 5 | #include "MatVec3.h" | 
| 6 |  | 
| 7 | RigidBody::RigidBody() : StuntDouble() { | 
| 8 | objType = OT_RIGIDBODY; | 
| 9 | is_linear = false; | 
| 10 | linear_axis =  -1; | 
| 11 | momIntTol = 1e-6; | 
| 12 | } | 
| 13 |  | 
| 14 | RigidBody::~RigidBody() { | 
| 15 | } | 
| 16 |  | 
| 17 | void RigidBody::addAtom(Atom* at, AtomStamp* ats) { | 
| 18 |  | 
| 19 | vec3 coords; | 
| 20 | vec3 euler; | 
| 21 | mat3x3 Atmp; | 
| 22 |  | 
| 23 | myAtoms.push_back(at); | 
| 24 |  | 
| 25 | if( !ats->havePosition() ){ | 
| 26 | sprintf( painCave.errMsg, | 
| 27 | "RigidBody error.\n" | 
| 28 | "\tAtom %s does not have a position specified.\n" | 
| 29 | "\tThis means RigidBody cannot set up reference coordinates.\n", | 
| 30 | ats->getType() ); | 
| 31 | painCave.isFatal = 1; | 
| 32 | simError(); | 
| 33 | } | 
| 34 |  | 
| 35 | coords[0] = ats->getPosX(); | 
| 36 | coords[1] = ats->getPosY(); | 
| 37 | coords[2] = ats->getPosZ(); | 
| 38 |  | 
| 39 | refCoords.push_back(coords); | 
| 40 |  | 
| 41 | if (at->isDirectional()) { | 
| 42 |  | 
| 43 | if( !ats->haveOrientation() ){ | 
| 44 | sprintf( painCave.errMsg, | 
| 45 | "RigidBody error.\n" | 
| 46 | "\tAtom %s does not have an orientation specified.\n" | 
| 47 | "\tThis means RigidBody cannot set up reference orientations.\n", | 
| 48 | ats->getType() ); | 
| 49 | painCave.isFatal = 1; | 
| 50 | simError(); | 
| 51 | } | 
| 52 |  | 
| 53 | euler[0] = ats->getEulerPhi(); | 
| 54 | euler[1] = ats->getEulerTheta(); | 
| 55 | euler[2] = ats->getEulerPsi(); | 
| 56 |  | 
| 57 | doEulerToRotMat(euler, Atmp); | 
| 58 |  | 
| 59 | refOrients.push_back(Atmp); | 
| 60 |  | 
| 61 | } | 
| 62 | } | 
| 63 |  | 
| 64 | void RigidBody::getPos(double theP[3]){ | 
| 65 | for (int i = 0; i < 3 ; i++) | 
| 66 | theP[i] = pos[i]; | 
| 67 | } | 
| 68 |  | 
| 69 | void RigidBody::setPos(double theP[3]){ | 
| 70 | for (int i = 0; i < 3 ; i++) | 
| 71 | pos[i] = theP[i]; | 
| 72 | } | 
| 73 |  | 
| 74 | void RigidBody::getVel(double theV[3]){ | 
| 75 | for (int i = 0; i < 3 ; i++) | 
| 76 | theV[i] = vel[i]; | 
| 77 | } | 
| 78 |  | 
| 79 | void RigidBody::setVel(double theV[3]){ | 
| 80 | for (int i = 0; i < 3 ; i++) | 
| 81 | vel[i] = theV[i]; | 
| 82 | } | 
| 83 |  | 
| 84 | void RigidBody::getFrc(double theF[3]){ | 
| 85 | for (int i = 0; i < 3 ; i++) | 
| 86 | theF[i] = frc[i]; | 
| 87 | } | 
| 88 |  | 
| 89 | void RigidBody::addFrc(double theF[3]){ | 
| 90 | for (int i = 0; i < 3 ; i++) | 
| 91 | frc[i] += theF[i]; | 
| 92 | } | 
| 93 |  | 
| 94 | void RigidBody::zeroForces() { | 
| 95 |  | 
| 96 | for (int i = 0; i < 3; i++) { | 
| 97 | frc[i] = 0.0; | 
| 98 | trq[i] = 0.0; | 
| 99 | } | 
| 100 |  | 
| 101 | } | 
| 102 |  | 
| 103 | void RigidBody::setEuler( double phi, double theta, double psi ){ | 
| 104 |  | 
| 105 | A[0][0] = (cos(phi) * cos(psi)) - (sin(phi) * cos(theta) * sin(psi)); | 
| 106 | A[0][1] = (sin(phi) * cos(psi)) + (cos(phi) * cos(theta) * sin(psi)); | 
| 107 | A[0][2] = sin(theta) * sin(psi); | 
| 108 |  | 
| 109 | A[1][0] = -(cos(phi) * sin(psi)) - (sin(phi) * cos(theta) * cos(psi)); | 
| 110 | A[1][1] = -(sin(phi) * sin(psi)) + (cos(phi) * cos(theta) * cos(psi)); | 
| 111 | A[1][2] = sin(theta) * cos(psi); | 
| 112 |  | 
| 113 | A[2][0] = sin(phi) * sin(theta); | 
| 114 | A[2][1] = -cos(phi) * sin(theta); | 
| 115 | A[2][2] = cos(theta); | 
| 116 |  | 
| 117 | } | 
| 118 |  | 
| 119 | void RigidBody::getQ( double q[4] ){ | 
| 120 |  | 
| 121 | double t, s; | 
| 122 | double ad1, ad2, ad3; | 
| 123 |  | 
| 124 | t = A[0][0] + A[1][1] + A[2][2] + 1.0; | 
| 125 | if( t > 0.0 ){ | 
| 126 |  | 
| 127 | s = 0.5 / sqrt( t ); | 
| 128 | q[0] = 0.25 / s; | 
| 129 | q[1] = (A[1][2] - A[2][1]) * s; | 
| 130 | q[2] = (A[2][0] - A[0][2]) * s; | 
| 131 | q[3] = (A[0][1] - A[1][0]) * s; | 
| 132 | } | 
| 133 | else{ | 
| 134 |  | 
| 135 | ad1 = fabs( A[0][0] ); | 
| 136 | ad2 = fabs( A[1][1] ); | 
| 137 | ad3 = fabs( A[2][2] ); | 
| 138 |  | 
| 139 | if( ad1 >= ad2 && ad1 >= ad3 ){ | 
| 140 |  | 
| 141 | s = 2.0 * sqrt( 1.0 + A[0][0] - A[1][1] - A[2][2] ); | 
| 142 | q[0] = (A[1][2] + A[2][1]) / s; | 
| 143 | q[1] = 0.5 / s; | 
| 144 | q[2] = (A[0][1] + A[1][0]) / s; | 
| 145 | q[3] = (A[0][2] + A[2][0]) / s; | 
| 146 | } | 
| 147 | else if( ad2 >= ad1 && ad2 >= ad3 ){ | 
| 148 |  | 
| 149 | s = sqrt( 1.0 + A[1][1] - A[0][0] - A[2][2] ) * 2.0; | 
| 150 | q[0] = (A[0][2] + A[2][0]) / s; | 
| 151 | q[1] = (A[0][1] + A[1][0]) / s; | 
| 152 | q[2] = 0.5 / s; | 
| 153 | q[3] = (A[1][2] + A[2][1]) / s; | 
| 154 | } | 
| 155 | else{ | 
| 156 |  | 
| 157 | s = sqrt( 1.0 + A[2][2] - A[0][0] - A[1][1] ) * 2.0; | 
| 158 | q[0] = (A[0][1] + A[1][0]) / s; | 
| 159 | q[1] = (A[0][2] + A[2][0]) / s; | 
| 160 | q[2] = (A[1][2] + A[2][1]) / s; | 
| 161 | q[3] = 0.5 / s; | 
| 162 | } | 
| 163 | } | 
| 164 | } | 
| 165 |  | 
| 166 | void RigidBody::setQ( double the_q[4] ){ | 
| 167 |  | 
| 168 | double q0Sqr, q1Sqr, q2Sqr, q3Sqr; | 
| 169 |  | 
| 170 | q0Sqr = the_q[0] * the_q[0]; | 
| 171 | q1Sqr = the_q[1] * the_q[1]; | 
| 172 | q2Sqr = the_q[2] * the_q[2]; | 
| 173 | q3Sqr = the_q[3] * the_q[3]; | 
| 174 |  | 
| 175 | A[0][0] = q0Sqr + q1Sqr - q2Sqr - q3Sqr; | 
| 176 | A[0][1] = 2.0 * ( the_q[1] * the_q[2] + the_q[0] * the_q[3] ); | 
| 177 | A[0][2] = 2.0 * ( the_q[1] * the_q[3] - the_q[0] * the_q[2] ); | 
| 178 |  | 
| 179 | A[1][0] = 2.0 * ( the_q[1] * the_q[2] - the_q[0] * the_q[3] ); | 
| 180 | A[1][1] = q0Sqr - q1Sqr + q2Sqr - q3Sqr; | 
| 181 | A[1][2] = 2.0 * ( the_q[2] * the_q[3] + the_q[0] * the_q[1] ); | 
| 182 |  | 
| 183 | A[2][0] = 2.0 * ( the_q[1] * the_q[3] + the_q[0] * the_q[2] ); | 
| 184 | A[2][1] = 2.0 * ( the_q[2] * the_q[3] - the_q[0] * the_q[1] ); | 
| 185 | A[2][2] = q0Sqr - q1Sqr -q2Sqr +q3Sqr; | 
| 186 |  | 
| 187 | } | 
| 188 |  | 
| 189 | void RigidBody::getA( double the_A[3][3] ){ | 
| 190 |  | 
| 191 | for (int i = 0; i < 3; i++) | 
| 192 | for (int j = 0; j < 3; j++) | 
| 193 | the_A[i][j] = A[i][j]; | 
| 194 |  | 
| 195 | } | 
| 196 |  | 
| 197 | void RigidBody::setA( double the_A[3][3] ){ | 
| 198 |  | 
| 199 | for (int i = 0; i < 3; i++) | 
| 200 | for (int j = 0; j < 3; j++) | 
| 201 | A[i][j] = the_A[i][j]; | 
| 202 |  | 
| 203 | } | 
| 204 |  | 
| 205 | void RigidBody::getJ( double theJ[3] ){ | 
| 206 |  | 
| 207 | for (int i = 0; i < 3; i++) | 
| 208 | theJ[i] = ji[i]; | 
| 209 |  | 
| 210 | } | 
| 211 |  | 
| 212 | void RigidBody::setJ( double theJ[3] ){ | 
| 213 |  | 
| 214 | for (int i = 0; i < 3; i++) | 
| 215 | ji[i] = theJ[i]; | 
| 216 |  | 
| 217 | } | 
| 218 |  | 
| 219 | void RigidBody::getTrq(double theT[3]){ | 
| 220 | for (int i = 0; i < 3 ; i++) | 
| 221 | theT[i] = trq[i]; | 
| 222 | } | 
| 223 |  | 
| 224 | void RigidBody::addTrq(double theT[3]){ | 
| 225 | for (int i = 0; i < 3 ; i++) | 
| 226 | trq[i] += theT[i]; | 
| 227 | } | 
| 228 |  | 
| 229 | void RigidBody::getI( double the_I[3][3] ){ | 
| 230 |  | 
| 231 | for (int i = 0; i < 3; i++) | 
| 232 | for (int j = 0; j < 3; j++) | 
| 233 | the_I[i][j] = I[i][j]; | 
| 234 |  | 
| 235 | } | 
| 236 |  | 
| 237 | void RigidBody::lab2Body( double r[3] ){ | 
| 238 |  | 
| 239 | double rl[3]; // the lab frame vector | 
| 240 |  | 
| 241 | rl[0] = r[0]; | 
| 242 | rl[1] = r[1]; | 
| 243 | rl[2] = r[2]; | 
| 244 |  | 
| 245 | r[0] = (A[0][0] * rl[0]) + (A[0][1] * rl[1]) + (A[0][2] * rl[2]); | 
| 246 | r[1] = (A[1][0] * rl[0]) + (A[1][1] * rl[1]) + (A[1][2] * rl[2]); | 
| 247 | r[2] = (A[2][0] * rl[0]) + (A[2][1] * rl[1]) + (A[2][2] * rl[2]); | 
| 248 |  | 
| 249 | } | 
| 250 |  | 
| 251 | void RigidBody::body2Lab( double r[3] ){ | 
| 252 |  | 
| 253 | double rb[3]; // the body frame vector | 
| 254 |  | 
| 255 | rb[0] = r[0]; | 
| 256 | rb[1] = r[1]; | 
| 257 | rb[2] = r[2]; | 
| 258 |  | 
| 259 | r[0] = (A[0][0] * rb[0]) + (A[1][0] * rb[1]) + (A[2][0] * rb[2]); | 
| 260 | r[1] = (A[0][1] * rb[0]) + (A[1][1] * rb[1]) + (A[2][1] * rb[2]); | 
| 261 | r[2] = (A[0][2] * rb[0]) + (A[1][2] * rb[1]) + (A[2][2] * rb[2]); | 
| 262 |  | 
| 263 | } | 
| 264 |  | 
| 265 | void RigidBody::calcRefCoords( ) { | 
| 266 |  | 
| 267 | int i,j,k, it, n_linear_coords; | 
| 268 | double mtmp; | 
| 269 | vec3 apos; | 
| 270 | double refCOM[3]; | 
| 271 | vec3 ptmp; | 
| 272 | double Itmp[3][3]; | 
| 273 | double evals[3]; | 
| 274 | double evects[3][3]; | 
| 275 | double r, r2, len; | 
| 276 |  | 
| 277 | // First, find the center of mass: | 
| 278 |  | 
| 279 | mass = 0.0; | 
| 280 | for (j=0; j<3; j++) | 
| 281 | refCOM[j] = 0.0; | 
| 282 |  | 
| 283 | for (i = 0; i < myAtoms.size(); i++) { | 
| 284 | mtmp = myAtoms[i]->getMass(); | 
| 285 | mass += mtmp; | 
| 286 |  | 
| 287 | apos = refCoords[i]; | 
| 288 |  | 
| 289 | for(j = 0; j < 3; j++) { | 
| 290 | refCOM[j] += apos[j]*mtmp; | 
| 291 | } | 
| 292 | } | 
| 293 |  | 
| 294 | for(j = 0; j < 3; j++) | 
| 295 | refCOM[j] /= mass; | 
| 296 |  | 
| 297 | // Next, move the origin of the reference coordinate system to the COM: | 
| 298 |  | 
| 299 | for (i = 0; i < myAtoms.size(); i++) { | 
| 300 | apos = refCoords[i]; | 
| 301 | for (j=0; j < 3; j++) { | 
| 302 | apos[j] = apos[j] - refCOM[j]; | 
| 303 | } | 
| 304 | refCoords[i] = apos; | 
| 305 | } | 
| 306 |  | 
| 307 | // Moment of Inertia calculation | 
| 308 |  | 
| 309 | for (i = 0; i < 3; i++) | 
| 310 | for (j = 0; j < 3; j++) | 
| 311 | Itmp[i][j] = 0.0; | 
| 312 |  | 
| 313 | for (it = 0; it < myAtoms.size(); it++) { | 
| 314 |  | 
| 315 | mtmp = myAtoms[it]->getMass(); | 
| 316 | ptmp = refCoords[it]; | 
| 317 | r= norm3(ptmp.vec); | 
| 318 | r2 = r*r; | 
| 319 |  | 
| 320 | for (i = 0; i < 3; i++) { | 
| 321 | for (j = 0; j < 3; j++) { | 
| 322 |  | 
| 323 | if (i==j) Itmp[i][j] += mtmp * r2; | 
| 324 |  | 
| 325 | Itmp[i][j] -= mtmp * ptmp.vec[i]*ptmp.vec[j]; | 
| 326 | } | 
| 327 | } | 
| 328 | } | 
| 329 |  | 
| 330 | diagonalize3x3(Itmp, evals, sU); | 
| 331 |  | 
| 332 | // zero out I and then fill the diagonals with the moments of inertia: | 
| 333 |  | 
| 334 | n_linear_coords = 0; | 
| 335 |  | 
| 336 | for (i = 0; i < 3; i++) { | 
| 337 | for (j = 0; j < 3; j++) { | 
| 338 | I[i][j] = 0.0; | 
| 339 | } | 
| 340 | I[i][i] = evals[i]; | 
| 341 |  | 
| 342 | if (fabs(evals[i]) < momIntTol) { | 
| 343 | is_linear = true; | 
| 344 | n_linear_coords++; | 
| 345 | linear_axis = i; | 
| 346 | } | 
| 347 | } | 
| 348 |  | 
| 349 | if (n_linear_coords > 1) { | 
| 350 | sprintf( painCave.errMsg, | 
| 351 | "RigidBody error.\n" | 
| 352 | "\tOOPSE found more than one axis in this rigid body with a vanishing \n" | 
| 353 | "\tmoment of inertia.  This can happen in one of three ways:\n" | 
| 354 | "\t 1) Only one atom was specified, or \n" | 
| 355 | "\t 2) All atoms were specified at the same location, or\n" | 
| 356 | "\t 3) The programmers did something stupid.\n" | 
| 357 | "\tIt is silly to use a rigid body to describe this situation.  Be smarter.\n" | 
| 358 | ); | 
| 359 | painCave.isFatal = 1; | 
| 360 | simError(); | 
| 361 | } | 
| 362 |  | 
| 363 | // renormalize column vectors: | 
| 364 |  | 
| 365 | for (i=0; i < 3; i++) { | 
| 366 | len = 0.0; | 
| 367 | for (j = 0; j < 3; j++) { | 
| 368 | len += sU[i][j]*sU[i][j]; | 
| 369 | } | 
| 370 | len = sqrt(len); | 
| 371 | for (j = 0; j < 3; j++) { | 
| 372 | sU[i][j] /= len; | 
| 373 | } | 
| 374 | } | 
| 375 | } | 
| 376 |  | 
| 377 | void RigidBody::doEulerToRotMat(vec3 &euler, mat3x3 &myA ){ | 
| 378 |  | 
| 379 | double phi, theta, psi; | 
| 380 |  | 
| 381 | phi = euler[0]; | 
| 382 | theta = euler[1]; | 
| 383 | psi = euler[2]; | 
| 384 |  | 
| 385 | myA[0][0] = (cos(phi) * cos(psi)) - (sin(phi) * cos(theta) * sin(psi)); | 
| 386 | myA[0][1] = (sin(phi) * cos(psi)) + (cos(phi) * cos(theta) * sin(psi)); | 
| 387 | myA[0][2] = sin(theta) * sin(psi); | 
| 388 |  | 
| 389 | myA[1][0] = -(cos(phi) * sin(psi)) - (sin(phi) * cos(theta) * cos(psi)); | 
| 390 | myA[1][1] = -(sin(phi) * sin(psi)) + (cos(phi) * cos(theta) * cos(psi)); | 
| 391 | myA[1][2] = sin(theta) * cos(psi); | 
| 392 |  | 
| 393 | myA[2][0] = sin(phi) * sin(theta); | 
| 394 | myA[2][1] = -cos(phi) * sin(theta); | 
| 395 | myA[2][2] = cos(theta); | 
| 396 |  | 
| 397 | } | 
| 398 |  | 
| 399 | void RigidBody::calcForcesAndTorques() { | 
| 400 |  | 
| 401 | // Convert Atomic forces and torques to total forces and torques: | 
| 402 | int i, j; | 
| 403 | double apos[3]; | 
| 404 | double afrc[3]; | 
| 405 | double atrq[3]; | 
| 406 | double rpos[3]; | 
| 407 |  | 
| 408 | zeroForces(); | 
| 409 |  | 
| 410 | for (i = 0; i < myAtoms.size(); i++) { | 
| 411 |  | 
| 412 | myAtoms[i]->getPos(apos); | 
| 413 | myAtoms[i]->getFrc(afrc); | 
| 414 |  | 
| 415 | for (j=0; j<3; j++) { | 
| 416 | rpos[j] = apos[j] - pos[j]; | 
| 417 | frc[j] += afrc[j]; | 
| 418 | } | 
| 419 |  | 
| 420 | trq[0] += rpos[1]*afrc[2] - rpos[2]*afrc[1]; | 
| 421 | trq[1] += rpos[2]*afrc[0] - rpos[0]*afrc[2]; | 
| 422 | trq[2] += rpos[0]*afrc[1] - rpos[1]*afrc[0]; | 
| 423 |  | 
| 424 | // If the atom has a torque associated with it, then we also need to | 
| 425 | // migrate the torques onto the center of mass: | 
| 426 |  | 
| 427 | if (myAtoms[i]->isDirectional()) { | 
| 428 |  | 
| 429 | myAtoms[i]->getTrq(atrq); | 
| 430 |  | 
| 431 | for (j=0; j<3; j++) | 
| 432 | trq[j] += atrq[j]; | 
| 433 | } | 
| 434 | } | 
| 435 |  | 
| 436 | // Convert Torque to Body-fixed coordinates: | 
| 437 | // (Actually, on second thought, don't.  Integrator does this now.) | 
| 438 | // lab2Body(trq); | 
| 439 |  | 
| 440 | } | 
| 441 |  | 
| 442 | void RigidBody::updateAtoms() { | 
| 443 | int i, j; | 
| 444 | vec3 ref; | 
| 445 | double apos[3]; | 
| 446 | DirectionalAtom* dAtom; | 
| 447 |  | 
| 448 | for (i = 0; i < myAtoms.size(); i++) { | 
| 449 |  | 
| 450 | ref = refCoords[i]; | 
| 451 |  | 
| 452 | body2Lab(ref.vec); | 
| 453 |  | 
| 454 | for (j = 0; j<3; j++) | 
| 455 | apos[j] = pos[j] + ref.vec[j]; | 
| 456 |  | 
| 457 | myAtoms[i]->setPos(apos); | 
| 458 |  | 
| 459 | if (myAtoms[i]->isDirectional()) { | 
| 460 |  | 
| 461 | dAtom = (DirectionalAtom *) myAtoms[i]; | 
| 462 | dAtom->rotateBy( A ); | 
| 463 |  | 
| 464 | } | 
| 465 | } | 
| 466 | } | 
| 467 |  | 
| 468 | void RigidBody::getGrad( double grad[6] ) { | 
| 469 |  | 
| 470 | double myEuler[3]; | 
| 471 | double phi, theta, psi; | 
| 472 | double cphi, sphi, ctheta, stheta; | 
| 473 | double ephi[3]; | 
| 474 | double etheta[3]; | 
| 475 | double epsi[3]; | 
| 476 |  | 
| 477 | this->getEulerAngles(myEuler); | 
| 478 |  | 
| 479 | phi = myEuler[0]; | 
| 480 | theta = myEuler[1]; | 
| 481 | psi = myEuler[2]; | 
| 482 |  | 
| 483 | cphi = cos(phi); | 
| 484 | sphi = sin(phi); | 
| 485 | ctheta = cos(theta); | 
| 486 | stheta = sin(theta); | 
| 487 |  | 
| 488 | // get unit vectors along the phi, theta and psi rotation axes | 
| 489 |  | 
| 490 | ephi[0] = 0.0; | 
| 491 | ephi[1] = 0.0; | 
| 492 | ephi[2] = 1.0; | 
| 493 |  | 
| 494 | etheta[0] = cphi; | 
| 495 | etheta[1] = sphi; | 
| 496 | etheta[2] = 0.0; | 
| 497 |  | 
| 498 | epsi[0] = stheta * cphi; | 
| 499 | epsi[1] = stheta * sphi; | 
| 500 | epsi[2] = ctheta; | 
| 501 |  | 
| 502 | for (int j = 0 ; j<3; j++) | 
| 503 | grad[j] = frc[j]; | 
| 504 |  | 
| 505 | grad[3] = 0.0; | 
| 506 | grad[4] = 0.0; | 
| 507 | grad[5] = 0.0; | 
| 508 |  | 
| 509 | for (int j = 0; j < 3; j++ ) { | 
| 510 |  | 
| 511 | grad[3] += trq[j]*ephi[j]; | 
| 512 | grad[4] += trq[j]*etheta[j]; | 
| 513 | grad[5] += trq[j]*epsi[j]; | 
| 514 |  | 
| 515 | } | 
| 516 |  | 
| 517 | } | 
| 518 |  | 
| 519 | /** | 
| 520 | * getEulerAngles computes a set of Euler angle values consistent | 
| 521 | * with an input rotation matrix.  They are returned in the following | 
| 522 | * order: | 
| 523 | *  myEuler[0] = phi; | 
| 524 | *  myEuler[1] = theta; | 
| 525 | *  myEuler[2] = psi; | 
| 526 | */ | 
| 527 | void RigidBody::getEulerAngles(double myEuler[3]) { | 
| 528 |  | 
| 529 | // We use so-called "x-convention", which is the most common | 
| 530 | // definition.  In this convention, the rotation given by Euler | 
| 531 | // angles (phi, theta, psi), where the first rotation is by an angle | 
| 532 | // phi about the z-axis, the second is by an angle theta (0 <= theta | 
| 533 | // <= 180) about the x-axis, and the third is by an angle psi about | 
| 534 | // the z-axis (again). | 
| 535 |  | 
| 536 |  | 
| 537 | double phi,theta,psi,eps; | 
| 538 | double pi; | 
| 539 | double cphi,ctheta,cpsi; | 
| 540 | double sphi,stheta,spsi; | 
| 541 | double b[3]; | 
| 542 | int flip[3]; | 
| 543 |  | 
| 544 | // set the tolerance for Euler angles and rotation elements | 
| 545 |  | 
| 546 | eps = 1.0e-8; | 
| 547 |  | 
| 548 | theta = acos(min(1.0,max(-1.0,A[2][2]))); | 
| 549 | ctheta = A[2][2]; | 
| 550 | stheta = sqrt(1.0 - ctheta * ctheta); | 
| 551 |  | 
| 552 | // when sin(theta) is close to 0, we need to consider the | 
| 553 | // possibility of a singularity. In this case, we can assign an | 
| 554 | // arbitary value to phi (or psi), and then determine the psi (or | 
| 555 | // phi) or vice-versa.  We'll assume that phi always gets the | 
| 556 | // rotation, and psi is 0 in cases of singularity.  we use atan2 | 
| 557 | // instead of atan, since atan2 will give us -Pi to Pi.  Since 0 <= | 
| 558 | // theta <= 180, sin(theta) will be always non-negative. Therefore, | 
| 559 | // it never changes the sign of both of the parameters passed to | 
| 560 | // atan2. | 
| 561 |  | 
| 562 | if (fabs(stheta) <= eps){ | 
| 563 | psi = 0.0; | 
| 564 | phi = atan2(-A[1][0], A[0][0]); | 
| 565 | } | 
| 566 | // we only have one unique solution | 
| 567 | else{ | 
| 568 | phi = atan2(A[2][0], -A[2][1]); | 
| 569 | psi = atan2(A[0][2], A[1][2]); | 
| 570 | } | 
| 571 |  | 
| 572 | //wrap phi and psi, make sure they are in the range from 0 to 2*Pi | 
| 573 | //if (phi < 0) | 
| 574 | //  phi += M_PI; | 
| 575 |  | 
| 576 | //if (psi < 0) | 
| 577 | //  psi += M_PI; | 
| 578 |  | 
| 579 | myEuler[0] = phi; | 
| 580 | myEuler[1] = theta; | 
| 581 | myEuler[2] = psi; | 
| 582 |  | 
| 583 | return; | 
| 584 | } | 
| 585 |  | 
| 586 | double RigidBody::max(double x, double  y) { | 
| 587 | return (x > y) ? x : y; | 
| 588 | } | 
| 589 |  | 
| 590 | double RigidBody::min(double x, double  y) { | 
| 591 | return (x > y) ? y : x; | 
| 592 | } | 
| 593 |  | 
| 594 | void RigidBody::findCOM() { | 
| 595 |  | 
| 596 | size_t i; | 
| 597 | int j; | 
| 598 | double mtmp; | 
| 599 | double ptmp[3]; | 
| 600 | double vtmp[3]; | 
| 601 |  | 
| 602 | for(j = 0; j < 3; j++) { | 
| 603 | pos[j] = 0.0; | 
| 604 | vel[j] = 0.0; | 
| 605 | } | 
| 606 | mass = 0.0; | 
| 607 |  | 
| 608 | for (i = 0; i < myAtoms.size(); i++) { | 
| 609 |  | 
| 610 | mtmp = myAtoms[i]->getMass(); | 
| 611 | myAtoms[i]->getPos(ptmp); | 
| 612 | myAtoms[i]->getVel(vtmp); | 
| 613 |  | 
| 614 | mass += mtmp; | 
| 615 |  | 
| 616 | for(j = 0; j < 3; j++) { | 
| 617 | pos[j] += ptmp[j]*mtmp; | 
| 618 | vel[j] += vtmp[j]*mtmp; | 
| 619 | } | 
| 620 |  | 
| 621 | } | 
| 622 |  | 
| 623 | for(j = 0; j < 3; j++) { | 
| 624 | pos[j] /= mass; | 
| 625 | vel[j] /= mass; | 
| 626 | } | 
| 627 |  | 
| 628 | } | 
| 629 |  | 
| 630 | void RigidBody::accept(BaseVisitor* v){ | 
| 631 | vector<Atom*>::iterator atomIter; | 
| 632 | v->visit(this); | 
| 633 |  | 
| 634 | //for(atomIter = myAtoms.begin(); atomIter != myAtoms.end(); ++atomIter) | 
| 635 | //  (*atomIter)->accept(v); | 
| 636 | } |