49#include "utils/Constants.hpp"
50#include "utils/simError.h"
54 RigidBody::RigidBody() :
55 StuntDouble(otRigidBody, &Snapshot::rigidbodyData), inertiaTensor_(0.0) {}
58 ((snapshotMan_->getPrevSnapshot())->*storage_).aMat[localIndex_] = a;
60 for (
unsigned int i = 0; i < atoms_.size(); ++i) {
61 if (atoms_[i]->isDirectional()) {
62 atoms_[i]->setPrevA(refOrients_[i].transpose() * a);
68 ((snapshotMan_->getCurrentSnapshot())->*storage_).aMat[localIndex_] = a;
70 for (
unsigned int i = 0; i < atoms_.size(); ++i) {
71 if (atoms_[i]->isDirectional()) {
72 atoms_[i]->setA(refOrients_[i].transpose() * a);
77 void RigidBody::setA(
const RotMat3x3d& a,
int snapshotNo) {
78 ((snapshotMan_->getSnapshot(snapshotNo))->*storage_).aMat[localIndex_] = a;
80 for (
unsigned int i = 0; i < atoms_.size(); ++i) {
81 if (atoms_[i]->isDirectional()) {
82 atoms_[i]->setA(refOrients_[i].transpose() * a, snapshotNo);
87 Mat3x3d RigidBody::getI() {
return inertiaTensor_; }
89 std::vector<RealType> RigidBody::getGrad() {
90 std::vector<RealType> grad(6, 0.0);
95 RealType cphi, sphi, ctheta, stheta;
103 myEuler = getA().toEulerAngles();
113 if (fabs(stheta) < 1.0E-9) { stheta = 1.0E-9; }
115 ephi[0] = -sphi * ctheta / stheta;
116 ephi[1] = cphi * ctheta / stheta;
123 epsi[0] = sphi / stheta;
124 epsi[1] = -cphi / stheta;
128 for (
int j = 0; j < 3; j++)
131 for (
int j = 0; j < 3; j++) {
132 grad[3] -= torque[j] * ephi[j];
133 grad[4] -= torque[j] * etheta[j];
134 grad[5] -= torque[j] * epsi[j];
143 void RigidBody::calcRefCoords() {
147 for (std::size_t i = 0; i < atoms_.size(); ++i) {
148 mtmp = atoms_[i]->getMass();
150 refCOM += refCoords_[i] * mtmp;
152 refCOM_ = refCOM / mass_;
155 for (std::size_t i = 0; i < atoms_.size(); ++i) {
156 refCoords_[i] -= refCOM_;
161 for (std::size_t i = 0; i < atoms_.size(); i++) {
163 mtmp = atoms_[i]->getMass();
164 IAtom -= outProduct(refCoords_[i], refCoords_[i]) * mtmp;
165 RealType r2 = refCoords_[i].lengthSquare();
166 IAtom(0, 0) += mtmp * r2;
167 IAtom(1, 1) += mtmp * r2;
168 IAtom(2, 2) += mtmp * r2;
172 if (atoms_[i]->isDirectional()) {
173 Itmp += refOrients_[i].transpose() * atoms_[i]->getI() * refOrients_[i];
179 Mat3x3d::diagonalize(Itmp, evals, sU_);
182 inertiaTensor_(0, 0) = evals[0];
183 inertiaTensor_(1, 1) = evals[1];
184 inertiaTensor_(2, 2) = evals[2];
187 for (
int i = 0; i < 3; i++) {
188 if (fabs(evals[i]) < OpenMD::epsilon) {
195 if (nLinearAxis > 1) {
197 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
199 "\tOpenMD found more than one axis in this rigid body with a "
201 "\tmoment of inertia. This can happen in one of three ways:\n"
202 "\t 1) Only one atom was specified, or \n"
203 "\t 2) All atoms were specified at the same location, or\n"
204 "\t 3) The programmers did something stupid.\n"
205 "\tIt is silly to use a rigid body to describe this situation. Be "
207 painCave.isFatal = 1;
212 void RigidBody::calcForcesAndTorques() {
225 ((snapshotMan_->getCurrentSnapshot())->*storage_).getStorageLayout();
227 for (
unsigned int i = 0; i < atoms_.size(); i++) {
228 atype = atoms_[i]->getAtomType();
230 afrc = atoms_[i]->getFrc();
231 apos = atoms_[i]->getPos();
236 trq[0] += rpos[1] * afrc[2] - rpos[2] * afrc[1];
237 trq[1] += rpos[2] * afrc[0] - rpos[0] * afrc[2];
238 trq[2] += rpos[0] * afrc[1] - rpos[1] * afrc[0];
243 if (atoms_[i]->isDirectional()) {
244 atrq = atoms_[i]->getTrq();
248 if ((sl & DataStorage::dslElectricField) && (atype->isElectrostatic())) {
249 ef += atoms_[i]->getElectricField();
256 if (sl & DataStorage::dslElectricField) {
258 setElectricField(ef);
262 Mat3x3d RigidBody::calcForcesAndTorquesAndVirial() {
278 ((snapshotMan_->getCurrentSnapshot())->*storage_).getStorageLayout();
280 for (
unsigned int i = 0; i < atoms_.size(); i++) {
281 atype = atoms_[i]->getAtomType();
283 afrc = atoms_[i]->getFrc();
284 apos = atoms_[i]->getPos();
289 trq[0] += rpos[1] * afrc[2] - rpos[2] * afrc[1];
290 trq[1] += rpos[2] * afrc[0] - rpos[0] * afrc[2];
291 trq[2] += rpos[0] * afrc[1] - rpos[1] * afrc[0];
296 if (atoms_[i]->isDirectional()) {
297 atrq = atoms_[i]->getTrq();
301 if ((sl & DataStorage::dslElectricField) && (atype->isElectrostatic())) {
302 ef += atoms_[i]->getElectricField();
306 tau_(0, 0) -= rpos[0] * afrc[0];
307 tau_(0, 1) -= rpos[0] * afrc[1];
308 tau_(0, 2) -= rpos[0] * afrc[2];
309 tau_(1, 0) -= rpos[1] * afrc[0];
310 tau_(1, 1) -= rpos[1] * afrc[1];
311 tau_(1, 2) -= rpos[1] * afrc[2];
312 tau_(2, 0) -= rpos[2] * afrc[0];
313 tau_(2, 1) -= rpos[2] * afrc[1];
314 tau_(2, 2) -= rpos[2] * afrc[2];
319 if (sl & DataStorage::dslElectricField) {
321 setElectricField(ef);
327 void RigidBody::updateAtoms() {
335 for (i = 0; i < atoms_.size(); i++) {
336 ref = body2Lab(refCoords_[i]);
340 atoms_[i]->setPos(apos);
342 if (atoms_[i]->isDirectional()) {
344 dAtom->
setA(refOrients_[i].transpose() * a);
349 void RigidBody::updateAtoms(
int frame) {
357 for (i = 0; i < atoms_.size(); i++) {
358 ref = body2Lab(refCoords_[i], frame);
362 atoms_[i]->setPos(apos, frame);
364 if (atoms_[i]->isDirectional()) {
366 dAtom->
setA(refOrients_[i].transpose() * a, frame);
371 void RigidBody::updateAtomVel() {
374 Vector3d ji = getJ();
378 skewMat(0, 1) = ji[2] / I(2, 2);
379 skewMat(0, 2) = -ji[1] / I(1, 1);
381 skewMat(1, 0) = -ji[2] / I(2, 2);
383 skewMat(1, 2) = ji[0] / I(0, 0);
385 skewMat(2, 0) = ji[1] / I(1, 1);
386 skewMat(2, 1) = -ji[0] / I(0, 0);
389 Mat3x3d mat = (getA() * skewMat).transpose();
390 Vector3d rbVel = getVel();
393 for (
unsigned int i = 0; i < refCoords_.size(); ++i) {
394 atoms_[i]->setVel(rbVel + mat * refCoords_[i]);
398 void RigidBody::updateAtomVel(
int frame) {
402 Vector3d ji = getJ(frame);
406 skewMat(0, 1) = ji[2] / I(2, 2);
407 skewMat(0, 2) = -ji[1] / I(1, 1);
409 skewMat(1, 0) = -ji[2] / I(2, 2);
411 skewMat(1, 2) = ji[0] / I(0, 0);
413 skewMat(2, 0) = ji[1] / I(1, 1);
414 skewMat(2, 1) = -ji[0] / I(0, 0);
417 Mat3x3d mat = (getA(frame) * skewMat).transpose();
418 Vector3d rbVel = getVel(frame);
421 for (
unsigned int i = 0; i < refCoords_.size(); ++i) {
422 atoms_[i]->setVel(rbVel + mat * refCoords_[i], frame);
426 bool RigidBody::getAtomPos(
Vector3d& pos,
unsigned int index) {
427 if (index < atoms_.size()) {
428 Vector3d ref = body2Lab(refCoords_[index]);
429 pos = getPos() + ref;
433 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
434 "%d is an invalid index. The current rigid body contans %lu atoms.\n",
435 index, atoms_.size());
436 painCave.isFatal = 0;
443 std::vector<Atom*>::iterator i;
444 i = std::find(atoms_.begin(), atoms_.end(), atom);
445 if (i != atoms_.end()) {
447 Vector3d ref = body2Lab(refCoords_[i - atoms_.begin()]);
448 pos = getPos() + ref;
451 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
452 "Atom %d does not belong to rigid body %d.\n",
454 painCave.isFatal = 0;
459 bool RigidBody::getAtomVel(
Vector3d& vel,
unsigned int index) {
462 if (index < atoms_.size()) {
471 skewMat(0, 1) = ji[2] / I(2, 2);
472 skewMat(0, 2) = -ji[1] / I(1, 1);
474 skewMat(1, 0) = -ji[2] / I(2, 2);
476 skewMat(1, 2) = ji[0] / I(0, 0);
478 skewMat(2, 0) = ji[1] / I(1, 1);
479 skewMat(2, 1) = -ji[0] / I(0, 0);
482 velRot = (getA() * skewMat).transpose() * ref;
484 vel = getVel() + velRot;
489 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
490 "Index %d is an invalid index, the current rigid body contains %lu "
492 index, atoms_.size());
493 painCave.isFatal = 0;
500 std::vector<Atom*>::iterator i;
501 i = std::find(atoms_.begin(), atoms_.end(), atom);
502 if (i != atoms_.end()) {
503 return getAtomVel(vel, i - atoms_.begin());
505 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
506 "Atom %d does not belong to rigid body %d.\n",
508 painCave.isFatal = 0;
514 bool RigidBody::getAtomRefCoor(
Vector3d& coor,
unsigned int index) {
515 if (index < atoms_.size()) {
516 coor = refCoords_[index];
520 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
521 "Index %d is an invalid index, the current rigid body contains %lu "
523 index, atoms_.size());
524 painCave.isFatal = 0;
531 std::vector<Atom*>::iterator i;
532 i = std::find(atoms_.begin(), atoms_.end(), atom);
533 if (i != atoms_.end()) {
535 coor = refCoords_[i - atoms_.begin()];
538 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
539 "Atom %d does not belong to rigid body %d.\n",
541 painCave.isFatal = 0;
551 atoms_.push_back(at);
553 if (!ats->havePosition()) {
554 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
556 "\tAtom %s does not have a position specified.\n"
557 "\tThis means RigidBody cannot set up reference coordinates.\n",
558 ats->getType().c_str());
559 painCave.isFatal = 1;
563 coords[0] = ats->getPosX();
564 coords[1] = ats->getPosY();
565 coords[2] = ats->getPosZ();
567 refCoords_.push_back(coords);
569 RotMat3x3d identMat = RotMat3x3d::identity();
571 if (at->isDirectional()) {
572 if (!ats->haveOrientation()) {
574 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
576 "\tAtom %s does not have an orientation specified.\n"
577 "\tThis means RigidBody cannot set up reference orientations.\n",
578 ats->getType().c_str());
579 painCave.isFatal = 1;
583 euler[0] = ats->getEulerPhi() * Constants::PI / 180.0;
584 euler[1] = ats->getEulerTheta() * Constants::PI / 180.0;
585 euler[2] = ats->getEulerPsi() * Constants::PI / 180.0;
587 RotMat3x3d Atmp(euler);
588 refOrients_.push_back(Atmp);
591 refOrients_.push_back(identMat);
AtomType is what OpenMD looks to for unchanging data about an atom.
virtual void setA(const RotMat3x3d &a)
Sets the current rotation matrix of this stuntdouble.
int getGlobalIndex()
Returns the global index of this stuntDouble.
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.