| 1 | /* | 
| 2 | * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. | 
| 3 | * | 
| 4 | * The University of Notre Dame grants you ("Licensee") a | 
| 5 | * non-exclusive, royalty free, license to use, modify and | 
| 6 | * redistribute this software in source and binary code form, provided | 
| 7 | * that the following conditions are met: | 
| 8 | * | 
| 9 | * 1. Redistributions of source code must retain the above copyright | 
| 10 | *    notice, this list of conditions and the following disclaimer. | 
| 11 | * | 
| 12 | * 2. Redistributions in binary form must reproduce the above copyright | 
| 13 | *    notice, this list of conditions and the following disclaimer in the | 
| 14 | *    documentation and/or other materials provided with the | 
| 15 | *    distribution. | 
| 16 | * | 
| 17 | * This software is provided "AS IS," without a warranty of any | 
| 18 | * kind. All express or implied conditions, representations and | 
| 19 | * warranties, including any implied warranty of merchantability, | 
| 20 | * fitness for a particular purpose or non-infringement, are hereby | 
| 21 | * excluded.  The University of Notre Dame and its licensors shall not | 
| 22 | * be liable for any damages suffered by licensee as a result of | 
| 23 | * using, modifying or distributing the software or its | 
| 24 | * derivatives. In no event will the University of Notre Dame or its | 
| 25 | * licensors be liable for any lost revenue, profit or data, or for | 
| 26 | * direct, indirect, special, consequential, incidental or punitive | 
| 27 | * damages, however caused and regardless of the theory of liability, | 
| 28 | * arising out of the use of or inability to use software, even if the | 
| 29 | * University of Notre Dame has been advised of the possibility of | 
| 30 | * such damages. | 
| 31 | * | 
| 32 | * SUPPORT OPEN SCIENCE!  If you use OpenMD or its source code in your | 
| 33 | * research, please cite the appropriate papers when you publish your | 
| 34 | * work.  Good starting points are: | 
| 35 | * | 
| 36 | * [1]  Meineke, et al., J. Comp. Chem. 26, 252-271 (2005). | 
| 37 | * [2]  Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006). | 
| 38 | * [3]  Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008). | 
| 39 | * [4]  Vardeman & Gezelter, in progress (2009). | 
| 40 | */ | 
| 41 |  | 
| 42 | /** | 
| 43 | * @file SimInfo.cpp | 
| 44 | * @author    tlin | 
| 45 | * @date  11/02/2004 | 
| 46 | * @version 1.0 | 
| 47 | */ | 
| 48 |  | 
| 49 | #include <algorithm> | 
| 50 | #include <set> | 
| 51 | #include <map> | 
| 52 |  | 
| 53 | #include "brains/SimInfo.hpp" | 
| 54 | #include "math/Vector3.hpp" | 
| 55 | #include "primitives/Molecule.hpp" | 
| 56 | #include "primitives/StuntDouble.hpp" | 
| 57 | #include "UseTheForce/fCutoffPolicy.h" | 
| 58 | #include "UseTheForce/DarkSide/fSwitchingFunctionType.h" | 
| 59 | #include "UseTheForce/doForces_interface.h" | 
| 60 | #include "UseTheForce/DarkSide/neighborLists_interface.h" | 
| 61 | #include "UseTheForce/DarkSide/switcheroo_interface.h" | 
| 62 | #include "utils/MemoryUtils.hpp" | 
| 63 | #include "utils/simError.h" | 
| 64 | #include "selection/SelectionManager.hpp" | 
| 65 | #include "io/ForceFieldOptions.hpp" | 
| 66 | #include "UseTheForce/ForceField.hpp" | 
| 67 |  | 
| 68 |  | 
| 69 | #ifdef IS_MPI | 
| 70 | #include "UseTheForce/mpiComponentPlan.h" | 
| 71 | #include "UseTheForce/DarkSide/simParallel_interface.h" | 
| 72 | #endif | 
| 73 |  | 
| 74 | namespace OpenMD { | 
| 75 | std::set<int> getRigidSet(int index, std::map<int, std::set<int> >& container) { | 
| 76 | std::map<int, std::set<int> >::iterator i = container.find(index); | 
| 77 | std::set<int> result; | 
| 78 | if (i != container.end()) { | 
| 79 | result = i->second; | 
| 80 | } | 
| 81 |  | 
| 82 | return result; | 
| 83 | } | 
| 84 |  | 
| 85 | SimInfo::SimInfo(ForceField* ff, Globals* simParams) : | 
| 86 | forceField_(ff), simParams_(simParams), | 
| 87 | ndf_(0), fdf_local(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0), | 
| 88 | nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0), | 
| 89 | nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0), | 
| 90 | nAtoms_(0), nBonds_(0),  nBends_(0), nTorsions_(0), nInversions_(0), | 
| 91 | nRigidBodies_(0), nIntegrableObjects_(0), nCutoffGroups_(0), | 
| 92 | nConstraints_(0), sman_(NULL), fortranInitialized_(false), | 
| 93 | calcBoxDipole_(false), useAtomicVirial_(true) { | 
| 94 |  | 
| 95 |  | 
| 96 | MoleculeStamp* molStamp; | 
| 97 | int nMolWithSameStamp; | 
| 98 | int nCutoffAtoms = 0; // number of atoms belong to cutoff groups | 
| 99 | int nGroups = 0;      //total cutoff groups defined in meta-data file | 
| 100 | CutoffGroupStamp* cgStamp; | 
| 101 | RigidBodyStamp* rbStamp; | 
| 102 | int nRigidAtoms = 0; | 
| 103 |  | 
| 104 | std::vector<Component*> components = simParams->getComponents(); | 
| 105 |  | 
| 106 | for (std::vector<Component*>::iterator i = components.begin(); i !=components.end(); ++i) { | 
| 107 | molStamp = (*i)->getMoleculeStamp(); | 
| 108 | nMolWithSameStamp = (*i)->getNMol(); | 
| 109 |  | 
| 110 | addMoleculeStamp(molStamp, nMolWithSameStamp); | 
| 111 |  | 
| 112 | //calculate atoms in molecules | 
| 113 | nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp; | 
| 114 |  | 
| 115 | //calculate atoms in cutoff groups | 
| 116 | int nAtomsInGroups = 0; | 
| 117 | int nCutoffGroupsInStamp = molStamp->getNCutoffGroups(); | 
| 118 |  | 
| 119 | for (int j=0; j < nCutoffGroupsInStamp; j++) { | 
| 120 | cgStamp = molStamp->getCutoffGroupStamp(j); | 
| 121 | nAtomsInGroups += cgStamp->getNMembers(); | 
| 122 | } | 
| 123 |  | 
| 124 | nGroups += nCutoffGroupsInStamp * nMolWithSameStamp; | 
| 125 |  | 
| 126 | nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp; | 
| 127 |  | 
| 128 | //calculate atoms in rigid bodies | 
| 129 | int nAtomsInRigidBodies = 0; | 
| 130 | int nRigidBodiesInStamp = molStamp->getNRigidBodies(); | 
| 131 |  | 
| 132 | for (int j=0; j < nRigidBodiesInStamp; j++) { | 
| 133 | rbStamp = molStamp->getRigidBodyStamp(j); | 
| 134 | nAtomsInRigidBodies += rbStamp->getNMembers(); | 
| 135 | } | 
| 136 |  | 
| 137 | nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp; | 
| 138 | nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp; | 
| 139 |  | 
| 140 | } | 
| 141 |  | 
| 142 | //every free atom (atom does not belong to cutoff groups) is a cutoff | 
| 143 | //group therefore the total number of cutoff groups in the system is | 
| 144 | //equal to the total number of atoms minus number of atoms belong to | 
| 145 | //cutoff group defined in meta-data file plus the number of cutoff | 
| 146 | //groups defined in meta-data file | 
| 147 | nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups; | 
| 148 |  | 
| 149 | //every free atom (atom does not belong to rigid bodies) is an | 
| 150 | //integrable object therefore the total number of integrable objects | 
| 151 | //in the system is equal to the total number of atoms minus number of | 
| 152 | //atoms belong to rigid body defined in meta-data file plus the number | 
| 153 | //of rigid bodies defined in meta-data file | 
| 154 | nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms | 
| 155 | + nGlobalRigidBodies_; | 
| 156 |  | 
| 157 | nGlobalMols_ = molStampIds_.size(); | 
| 158 | molToProcMap_.resize(nGlobalMols_); | 
| 159 | } | 
| 160 |  | 
| 161 | SimInfo::~SimInfo() { | 
| 162 | std::map<int, Molecule*>::iterator i; | 
| 163 | for (i = molecules_.begin(); i != molecules_.end(); ++i) { | 
| 164 | delete i->second; | 
| 165 | } | 
| 166 | molecules_.clear(); | 
| 167 |  | 
| 168 | delete sman_; | 
| 169 | delete simParams_; | 
| 170 | delete forceField_; | 
| 171 | } | 
| 172 |  | 
| 173 | int SimInfo::getNGlobalConstraints() { | 
| 174 | int nGlobalConstraints; | 
| 175 | #ifdef IS_MPI | 
| 176 | MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM, | 
| 177 | MPI_COMM_WORLD); | 
| 178 | #else | 
| 179 | nGlobalConstraints =  nConstraints_; | 
| 180 | #endif | 
| 181 | return nGlobalConstraints; | 
| 182 | } | 
| 183 |  | 
| 184 | bool SimInfo::addMolecule(Molecule* mol) { | 
| 185 | MoleculeIterator i; | 
| 186 |  | 
| 187 | i = molecules_.find(mol->getGlobalIndex()); | 
| 188 | if (i == molecules_.end() ) { | 
| 189 |  | 
| 190 | molecules_.insert(std::make_pair(mol->getGlobalIndex(), mol)); | 
| 191 |  | 
| 192 | nAtoms_ += mol->getNAtoms(); | 
| 193 | nBonds_ += mol->getNBonds(); | 
| 194 | nBends_ += mol->getNBends(); | 
| 195 | nTorsions_ += mol->getNTorsions(); | 
| 196 | nInversions_ += mol->getNInversions(); | 
| 197 | nRigidBodies_ += mol->getNRigidBodies(); | 
| 198 | nIntegrableObjects_ += mol->getNIntegrableObjects(); | 
| 199 | nCutoffGroups_ += mol->getNCutoffGroups(); | 
| 200 | nConstraints_ += mol->getNConstraintPairs(); | 
| 201 |  | 
| 202 | addInteractionPairs(mol); | 
| 203 |  | 
| 204 | return true; | 
| 205 | } else { | 
| 206 | return false; | 
| 207 | } | 
| 208 | } | 
| 209 |  | 
| 210 | bool SimInfo::removeMolecule(Molecule* mol) { | 
| 211 | MoleculeIterator i; | 
| 212 | i = molecules_.find(mol->getGlobalIndex()); | 
| 213 |  | 
| 214 | if (i != molecules_.end() ) { | 
| 215 |  | 
| 216 | assert(mol == i->second); | 
| 217 |  | 
| 218 | nAtoms_ -= mol->getNAtoms(); | 
| 219 | nBonds_ -= mol->getNBonds(); | 
| 220 | nBends_ -= mol->getNBends(); | 
| 221 | nTorsions_ -= mol->getNTorsions(); | 
| 222 | nInversions_ -= mol->getNInversions(); | 
| 223 | nRigidBodies_ -= mol->getNRigidBodies(); | 
| 224 | nIntegrableObjects_ -= mol->getNIntegrableObjects(); | 
| 225 | nCutoffGroups_ -= mol->getNCutoffGroups(); | 
| 226 | nConstraints_ -= mol->getNConstraintPairs(); | 
| 227 |  | 
| 228 | removeInteractionPairs(mol); | 
| 229 | molecules_.erase(mol->getGlobalIndex()); | 
| 230 |  | 
| 231 | delete mol; | 
| 232 |  | 
| 233 | return true; | 
| 234 | } else { | 
| 235 | return false; | 
| 236 | } | 
| 237 |  | 
| 238 |  | 
| 239 | } | 
| 240 |  | 
| 241 |  | 
| 242 | Molecule* SimInfo::beginMolecule(MoleculeIterator& i) { | 
| 243 | i = molecules_.begin(); | 
| 244 | return i == molecules_.end() ? NULL : i->second; | 
| 245 | } | 
| 246 |  | 
| 247 | Molecule* SimInfo::nextMolecule(MoleculeIterator& i) { | 
| 248 | ++i; | 
| 249 | return i == molecules_.end() ? NULL : i->second; | 
| 250 | } | 
| 251 |  | 
| 252 |  | 
| 253 | void SimInfo::calcNdf() { | 
| 254 | int ndf_local; | 
| 255 | MoleculeIterator i; | 
| 256 | std::vector<StuntDouble*>::iterator j; | 
| 257 | Molecule* mol; | 
| 258 | StuntDouble* integrableObject; | 
| 259 |  | 
| 260 | ndf_local = 0; | 
| 261 |  | 
| 262 | for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { | 
| 263 | for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL; | 
| 264 | integrableObject = mol->nextIntegrableObject(j)) { | 
| 265 |  | 
| 266 | ndf_local += 3; | 
| 267 |  | 
| 268 | if (integrableObject->isDirectional()) { | 
| 269 | if (integrableObject->isLinear()) { | 
| 270 | ndf_local += 2; | 
| 271 | } else { | 
| 272 | ndf_local += 3; | 
| 273 | } | 
| 274 | } | 
| 275 |  | 
| 276 | } | 
| 277 | } | 
| 278 |  | 
| 279 | // n_constraints is local, so subtract them on each processor | 
| 280 | ndf_local -= nConstraints_; | 
| 281 |  | 
| 282 | #ifdef IS_MPI | 
| 283 | MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD); | 
| 284 | #else | 
| 285 | ndf_ = ndf_local; | 
| 286 | #endif | 
| 287 |  | 
| 288 | // nZconstraints_ is global, as are the 3 COM translations for the | 
| 289 | // entire system: | 
| 290 | ndf_ = ndf_ - 3 - nZconstraint_; | 
| 291 |  | 
| 292 | } | 
| 293 |  | 
| 294 | int SimInfo::getFdf() { | 
| 295 | #ifdef IS_MPI | 
| 296 | MPI_Allreduce(&fdf_local,&fdf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD); | 
| 297 | #else | 
| 298 | fdf_ = fdf_local; | 
| 299 | #endif | 
| 300 | return fdf_; | 
| 301 | } | 
| 302 |  | 
| 303 | void SimInfo::calcNdfRaw() { | 
| 304 | int ndfRaw_local; | 
| 305 |  | 
| 306 | MoleculeIterator i; | 
| 307 | std::vector<StuntDouble*>::iterator j; | 
| 308 | Molecule* mol; | 
| 309 | StuntDouble* integrableObject; | 
| 310 |  | 
| 311 | // Raw degrees of freedom that we have to set | 
| 312 | ndfRaw_local = 0; | 
| 313 |  | 
| 314 | for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { | 
| 315 | for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL; | 
| 316 | integrableObject = mol->nextIntegrableObject(j)) { | 
| 317 |  | 
| 318 | ndfRaw_local += 3; | 
| 319 |  | 
| 320 | if (integrableObject->isDirectional()) { | 
| 321 | if (integrableObject->isLinear()) { | 
| 322 | ndfRaw_local += 2; | 
| 323 | } else { | 
| 324 | ndfRaw_local += 3; | 
| 325 | } | 
| 326 | } | 
| 327 |  | 
| 328 | } | 
| 329 | } | 
| 330 |  | 
| 331 | #ifdef IS_MPI | 
| 332 | MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD); | 
| 333 | #else | 
| 334 | ndfRaw_ = ndfRaw_local; | 
| 335 | #endif | 
| 336 | } | 
| 337 |  | 
| 338 | void SimInfo::calcNdfTrans() { | 
| 339 | int ndfTrans_local; | 
| 340 |  | 
| 341 | ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_; | 
| 342 |  | 
| 343 |  | 
| 344 | #ifdef IS_MPI | 
| 345 | MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD); | 
| 346 | #else | 
| 347 | ndfTrans_ = ndfTrans_local; | 
| 348 | #endif | 
| 349 |  | 
| 350 | ndfTrans_ = ndfTrans_ - 3 - nZconstraint_; | 
| 351 |  | 
| 352 | } | 
| 353 |  | 
| 354 | void SimInfo::addInteractionPairs(Molecule* mol) { | 
| 355 | ForceFieldOptions& options_ = forceField_->getForceFieldOptions(); | 
| 356 | std::vector<Bond*>::iterator bondIter; | 
| 357 | std::vector<Bend*>::iterator bendIter; | 
| 358 | std::vector<Torsion*>::iterator torsionIter; | 
| 359 | std::vector<Inversion*>::iterator inversionIter; | 
| 360 | Bond* bond; | 
| 361 | Bend* bend; | 
| 362 | Torsion* torsion; | 
| 363 | Inversion* inversion; | 
| 364 | int a; | 
| 365 | int b; | 
| 366 | int c; | 
| 367 | int d; | 
| 368 |  | 
| 369 | // atomGroups can be used to add special interaction maps between | 
| 370 | // groups of atoms that are in two separate rigid bodies. | 
| 371 | // However, most site-site interactions between two rigid bodies | 
| 372 | // are probably not special, just the ones between the physically | 
| 373 | // bonded atoms.  Interactions *within* a single rigid body should | 
| 374 | // always be excluded.  These are done at the bottom of this | 
| 375 | // function. | 
| 376 |  | 
| 377 | std::map<int, std::set<int> > atomGroups; | 
| 378 | Molecule::RigidBodyIterator rbIter; | 
| 379 | RigidBody* rb; | 
| 380 | Molecule::IntegrableObjectIterator ii; | 
| 381 | StuntDouble* integrableObject; | 
| 382 |  | 
| 383 | for (integrableObject = mol->beginIntegrableObject(ii); | 
| 384 | integrableObject != NULL; | 
| 385 | integrableObject = mol->nextIntegrableObject(ii)) { | 
| 386 |  | 
| 387 | if (integrableObject->isRigidBody()) { | 
| 388 | rb = static_cast<RigidBody*>(integrableObject); | 
| 389 | std::vector<Atom*> atoms = rb->getAtoms(); | 
| 390 | std::set<int> rigidAtoms; | 
| 391 | for (int i = 0; i < static_cast<int>(atoms.size()); ++i) { | 
| 392 | rigidAtoms.insert(atoms[i]->getGlobalIndex()); | 
| 393 | } | 
| 394 | for (int i = 0; i < static_cast<int>(atoms.size()); ++i) { | 
| 395 | atomGroups.insert(std::map<int, std::set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms)); | 
| 396 | } | 
| 397 | } else { | 
| 398 | std::set<int> oneAtomSet; | 
| 399 | oneAtomSet.insert(integrableObject->getGlobalIndex()); | 
| 400 | atomGroups.insert(std::map<int, std::set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet)); | 
| 401 | } | 
| 402 | } | 
| 403 |  | 
| 404 | for (bond= mol->beginBond(bondIter); bond != NULL; | 
| 405 | bond = mol->nextBond(bondIter)) { | 
| 406 |  | 
| 407 | a = bond->getAtomA()->getGlobalIndex(); | 
| 408 | b = bond->getAtomB()->getGlobalIndex(); | 
| 409 |  | 
| 410 | if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) { | 
| 411 | oneTwoInteractions_.addPair(a, b); | 
| 412 | } else { | 
| 413 | excludedInteractions_.addPair(a, b); | 
| 414 | } | 
| 415 | } | 
| 416 |  | 
| 417 | for (bend= mol->beginBend(bendIter); bend != NULL; | 
| 418 | bend = mol->nextBend(bendIter)) { | 
| 419 |  | 
| 420 | a = bend->getAtomA()->getGlobalIndex(); | 
| 421 | b = bend->getAtomB()->getGlobalIndex(); | 
| 422 | c = bend->getAtomC()->getGlobalIndex(); | 
| 423 |  | 
| 424 | if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) { | 
| 425 | oneTwoInteractions_.addPair(a, b); | 
| 426 | oneTwoInteractions_.addPair(b, c); | 
| 427 | } else { | 
| 428 | excludedInteractions_.addPair(a, b); | 
| 429 | excludedInteractions_.addPair(b, c); | 
| 430 | } | 
| 431 |  | 
| 432 | if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) { | 
| 433 | oneThreeInteractions_.addPair(a, c); | 
| 434 | } else { | 
| 435 | excludedInteractions_.addPair(a, c); | 
| 436 | } | 
| 437 | } | 
| 438 |  | 
| 439 | for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; | 
| 440 | torsion = mol->nextTorsion(torsionIter)) { | 
| 441 |  | 
| 442 | a = torsion->getAtomA()->getGlobalIndex(); | 
| 443 | b = torsion->getAtomB()->getGlobalIndex(); | 
| 444 | c = torsion->getAtomC()->getGlobalIndex(); | 
| 445 | d = torsion->getAtomD()->getGlobalIndex(); | 
| 446 |  | 
| 447 | if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) { | 
| 448 | oneTwoInteractions_.addPair(a, b); | 
| 449 | oneTwoInteractions_.addPair(b, c); | 
| 450 | oneTwoInteractions_.addPair(c, d); | 
| 451 | } else { | 
| 452 | excludedInteractions_.addPair(a, b); | 
| 453 | excludedInteractions_.addPair(b, c); | 
| 454 | excludedInteractions_.addPair(c, d); | 
| 455 | } | 
| 456 |  | 
| 457 | if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) { | 
| 458 | oneThreeInteractions_.addPair(a, c); | 
| 459 | oneThreeInteractions_.addPair(b, d); | 
| 460 | } else { | 
| 461 | excludedInteractions_.addPair(a, c); | 
| 462 | excludedInteractions_.addPair(b, d); | 
| 463 | } | 
| 464 |  | 
| 465 | if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) { | 
| 466 | oneFourInteractions_.addPair(a, d); | 
| 467 | } else { | 
| 468 | excludedInteractions_.addPair(a, d); | 
| 469 | } | 
| 470 | } | 
| 471 |  | 
| 472 | for (inversion= mol->beginInversion(inversionIter); inversion != NULL; | 
| 473 | inversion = mol->nextInversion(inversionIter)) { | 
| 474 |  | 
| 475 | a = inversion->getAtomA()->getGlobalIndex(); | 
| 476 | b = inversion->getAtomB()->getGlobalIndex(); | 
| 477 | c = inversion->getAtomC()->getGlobalIndex(); | 
| 478 | d = inversion->getAtomD()->getGlobalIndex(); | 
| 479 |  | 
| 480 | if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) { | 
| 481 | oneTwoInteractions_.addPair(a, b); | 
| 482 | oneTwoInteractions_.addPair(a, c); | 
| 483 | oneTwoInteractions_.addPair(a, d); | 
| 484 | } else { | 
| 485 | excludedInteractions_.addPair(a, b); | 
| 486 | excludedInteractions_.addPair(a, c); | 
| 487 | excludedInteractions_.addPair(a, d); | 
| 488 | } | 
| 489 |  | 
| 490 | if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) { | 
| 491 | oneThreeInteractions_.addPair(b, c); | 
| 492 | oneThreeInteractions_.addPair(b, d); | 
| 493 | oneThreeInteractions_.addPair(c, d); | 
| 494 | } else { | 
| 495 | excludedInteractions_.addPair(b, c); | 
| 496 | excludedInteractions_.addPair(b, d); | 
| 497 | excludedInteractions_.addPair(c, d); | 
| 498 | } | 
| 499 | } | 
| 500 |  | 
| 501 | for (rb = mol->beginRigidBody(rbIter); rb != NULL; | 
| 502 | rb = mol->nextRigidBody(rbIter)) { | 
| 503 | std::vector<Atom*> atoms = rb->getAtoms(); | 
| 504 | for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) { | 
| 505 | for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) { | 
| 506 | a = atoms[i]->getGlobalIndex(); | 
| 507 | b = atoms[j]->getGlobalIndex(); | 
| 508 | excludedInteractions_.addPair(a, b); | 
| 509 | } | 
| 510 | } | 
| 511 | } | 
| 512 |  | 
| 513 | } | 
| 514 |  | 
| 515 | void SimInfo::removeInteractionPairs(Molecule* mol) { | 
| 516 | ForceFieldOptions& options_ = forceField_->getForceFieldOptions(); | 
| 517 | std::vector<Bond*>::iterator bondIter; | 
| 518 | std::vector<Bend*>::iterator bendIter; | 
| 519 | std::vector<Torsion*>::iterator torsionIter; | 
| 520 | std::vector<Inversion*>::iterator inversionIter; | 
| 521 | Bond* bond; | 
| 522 | Bend* bend; | 
| 523 | Torsion* torsion; | 
| 524 | Inversion* inversion; | 
| 525 | int a; | 
| 526 | int b; | 
| 527 | int c; | 
| 528 | int d; | 
| 529 |  | 
| 530 | std::map<int, std::set<int> > atomGroups; | 
| 531 | Molecule::RigidBodyIterator rbIter; | 
| 532 | RigidBody* rb; | 
| 533 | Molecule::IntegrableObjectIterator ii; | 
| 534 | StuntDouble* integrableObject; | 
| 535 |  | 
| 536 | for (integrableObject = mol->beginIntegrableObject(ii); | 
| 537 | integrableObject != NULL; | 
| 538 | integrableObject = mol->nextIntegrableObject(ii)) { | 
| 539 |  | 
| 540 | if (integrableObject->isRigidBody()) { | 
| 541 | rb = static_cast<RigidBody*>(integrableObject); | 
| 542 | std::vector<Atom*> atoms = rb->getAtoms(); | 
| 543 | std::set<int> rigidAtoms; | 
| 544 | for (int i = 0; i < static_cast<int>(atoms.size()); ++i) { | 
| 545 | rigidAtoms.insert(atoms[i]->getGlobalIndex()); | 
| 546 | } | 
| 547 | for (int i = 0; i < static_cast<int>(atoms.size()); ++i) { | 
| 548 | atomGroups.insert(std::map<int, std::set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms)); | 
| 549 | } | 
| 550 | } else { | 
| 551 | std::set<int> oneAtomSet; | 
| 552 | oneAtomSet.insert(integrableObject->getGlobalIndex()); | 
| 553 | atomGroups.insert(std::map<int, std::set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet)); | 
| 554 | } | 
| 555 | } | 
| 556 |  | 
| 557 | for (bond= mol->beginBond(bondIter); bond != NULL; | 
| 558 | bond = mol->nextBond(bondIter)) { | 
| 559 |  | 
| 560 | a = bond->getAtomA()->getGlobalIndex(); | 
| 561 | b = bond->getAtomB()->getGlobalIndex(); | 
| 562 |  | 
| 563 | if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) { | 
| 564 | oneTwoInteractions_.removePair(a, b); | 
| 565 | } else { | 
| 566 | excludedInteractions_.removePair(a, b); | 
| 567 | } | 
| 568 | } | 
| 569 |  | 
| 570 | for (bend= mol->beginBend(bendIter); bend != NULL; | 
| 571 | bend = mol->nextBend(bendIter)) { | 
| 572 |  | 
| 573 | a = bend->getAtomA()->getGlobalIndex(); | 
| 574 | b = bend->getAtomB()->getGlobalIndex(); | 
| 575 | c = bend->getAtomC()->getGlobalIndex(); | 
| 576 |  | 
| 577 | if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) { | 
| 578 | oneTwoInteractions_.removePair(a, b); | 
| 579 | oneTwoInteractions_.removePair(b, c); | 
| 580 | } else { | 
| 581 | excludedInteractions_.removePair(a, b); | 
| 582 | excludedInteractions_.removePair(b, c); | 
| 583 | } | 
| 584 |  | 
| 585 | if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) { | 
| 586 | oneThreeInteractions_.removePair(a, c); | 
| 587 | } else { | 
| 588 | excludedInteractions_.removePair(a, c); | 
| 589 | } | 
| 590 | } | 
| 591 |  | 
| 592 | for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; | 
| 593 | torsion = mol->nextTorsion(torsionIter)) { | 
| 594 |  | 
| 595 | a = torsion->getAtomA()->getGlobalIndex(); | 
| 596 | b = torsion->getAtomB()->getGlobalIndex(); | 
| 597 | c = torsion->getAtomC()->getGlobalIndex(); | 
| 598 | d = torsion->getAtomD()->getGlobalIndex(); | 
| 599 |  | 
| 600 | if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) { | 
| 601 | oneTwoInteractions_.removePair(a, b); | 
| 602 | oneTwoInteractions_.removePair(b, c); | 
| 603 | oneTwoInteractions_.removePair(c, d); | 
| 604 | } else { | 
| 605 | excludedInteractions_.removePair(a, b); | 
| 606 | excludedInteractions_.removePair(b, c); | 
| 607 | excludedInteractions_.removePair(c, d); | 
| 608 | } | 
| 609 |  | 
| 610 | if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) { | 
| 611 | oneThreeInteractions_.removePair(a, c); | 
| 612 | oneThreeInteractions_.removePair(b, d); | 
| 613 | } else { | 
| 614 | excludedInteractions_.removePair(a, c); | 
| 615 | excludedInteractions_.removePair(b, d); | 
| 616 | } | 
| 617 |  | 
| 618 | if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) { | 
| 619 | oneFourInteractions_.removePair(a, d); | 
| 620 | } else { | 
| 621 | excludedInteractions_.removePair(a, d); | 
| 622 | } | 
| 623 | } | 
| 624 |  | 
| 625 | for (inversion= mol->beginInversion(inversionIter); inversion != NULL; | 
| 626 | inversion = mol->nextInversion(inversionIter)) { | 
| 627 |  | 
| 628 | a = inversion->getAtomA()->getGlobalIndex(); | 
| 629 | b = inversion->getAtomB()->getGlobalIndex(); | 
| 630 | c = inversion->getAtomC()->getGlobalIndex(); | 
| 631 | d = inversion->getAtomD()->getGlobalIndex(); | 
| 632 |  | 
| 633 | if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) { | 
| 634 | oneTwoInteractions_.removePair(a, b); | 
| 635 | oneTwoInteractions_.removePair(a, c); | 
| 636 | oneTwoInteractions_.removePair(a, d); | 
| 637 | } else { | 
| 638 | excludedInteractions_.removePair(a, b); | 
| 639 | excludedInteractions_.removePair(a, c); | 
| 640 | excludedInteractions_.removePair(a, d); | 
| 641 | } | 
| 642 |  | 
| 643 | if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) { | 
| 644 | oneThreeInteractions_.removePair(b, c); | 
| 645 | oneThreeInteractions_.removePair(b, d); | 
| 646 | oneThreeInteractions_.removePair(c, d); | 
| 647 | } else { | 
| 648 | excludedInteractions_.removePair(b, c); | 
| 649 | excludedInteractions_.removePair(b, d); | 
| 650 | excludedInteractions_.removePair(c, d); | 
| 651 | } | 
| 652 | } | 
| 653 |  | 
| 654 | for (rb = mol->beginRigidBody(rbIter); rb != NULL; | 
| 655 | rb = mol->nextRigidBody(rbIter)) { | 
| 656 | std::vector<Atom*> atoms = rb->getAtoms(); | 
| 657 | for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) { | 
| 658 | for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) { | 
| 659 | a = atoms[i]->getGlobalIndex(); | 
| 660 | b = atoms[j]->getGlobalIndex(); | 
| 661 | excludedInteractions_.removePair(a, b); | 
| 662 | } | 
| 663 | } | 
| 664 | } | 
| 665 |  | 
| 666 | } | 
| 667 |  | 
| 668 |  | 
| 669 | void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) { | 
| 670 | int curStampId; | 
| 671 |  | 
| 672 | //index from 0 | 
| 673 | curStampId = moleculeStamps_.size(); | 
| 674 |  | 
| 675 | moleculeStamps_.push_back(molStamp); | 
| 676 | molStampIds_.insert(molStampIds_.end(), nmol, curStampId); | 
| 677 | } | 
| 678 |  | 
| 679 | void SimInfo::update() { | 
| 680 |  | 
| 681 | setupSimType(); | 
| 682 |  | 
| 683 | #ifdef IS_MPI | 
| 684 | setupFortranParallel(); | 
| 685 | #endif | 
| 686 |  | 
| 687 | setupFortranSim(); | 
| 688 |  | 
| 689 | //setup fortran force field | 
| 690 | /** @deprecate */ | 
| 691 | int isError = 0; | 
| 692 |  | 
| 693 | setupCutoff(); | 
| 694 |  | 
| 695 | setupElectrostaticSummationMethod( isError ); | 
| 696 | setupSwitchingFunction(); | 
| 697 | setupAccumulateBoxDipole(); | 
| 698 |  | 
| 699 | if(isError){ | 
| 700 | sprintf( painCave.errMsg, | 
| 701 | "ForceField error: There was an error initializing the forceField in fortran.\n" ); | 
| 702 | painCave.isFatal = 1; | 
| 703 | simError(); | 
| 704 | } | 
| 705 |  | 
| 706 | calcNdf(); | 
| 707 | calcNdfRaw(); | 
| 708 | calcNdfTrans(); | 
| 709 |  | 
| 710 | fortranInitialized_ = true; | 
| 711 | } | 
| 712 |  | 
| 713 | std::set<AtomType*> SimInfo::getUniqueAtomTypes() { | 
| 714 | SimInfo::MoleculeIterator mi; | 
| 715 | Molecule* mol; | 
| 716 | Molecule::AtomIterator ai; | 
| 717 | Atom* atom; | 
| 718 | std::set<AtomType*> atomTypes; | 
| 719 |  | 
| 720 | for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) { | 
| 721 |  | 
| 722 | for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) { | 
| 723 | atomTypes.insert(atom->getAtomType()); | 
| 724 | } | 
| 725 |  | 
| 726 | } | 
| 727 |  | 
| 728 | return atomTypes; | 
| 729 | } | 
| 730 |  | 
| 731 | void SimInfo::setupSimType() { | 
| 732 | std::set<AtomType*>::iterator i; | 
| 733 | std::set<AtomType*> atomTypes; | 
| 734 | atomTypes = getUniqueAtomTypes(); | 
| 735 |  | 
| 736 | int useLennardJones = 0; | 
| 737 | int useElectrostatic = 0; | 
| 738 | int useEAM = 0; | 
| 739 | int useSC = 0; | 
| 740 | int useCharge = 0; | 
| 741 | int useDirectional = 0; | 
| 742 | int useDipole = 0; | 
| 743 | int useGayBerne = 0; | 
| 744 | int useSticky = 0; | 
| 745 | int useStickyPower = 0; | 
| 746 | int useShape = 0; | 
| 747 | int useFLARB = 0; //it is not in AtomType yet | 
| 748 | int useDirectionalAtom = 0; | 
| 749 | int useElectrostatics = 0; | 
| 750 | //usePBC and useRF are from simParams | 
| 751 | int usePBC = simParams_->getUsePeriodicBoundaryConditions(); | 
| 752 | int useRF; | 
| 753 | int useSF; | 
| 754 | int useSP; | 
| 755 | int useBoxDipole; | 
| 756 |  | 
| 757 | std::string myMethod; | 
| 758 |  | 
| 759 | // set the useRF logical | 
| 760 | useRF = 0; | 
| 761 | useSF = 0; | 
| 762 | useSP = 0; | 
| 763 | useBoxDipole = 0; | 
| 764 |  | 
| 765 | if (simParams_->haveElectrostaticSummationMethod()) { | 
| 766 | std::string myMethod = simParams_->getElectrostaticSummationMethod(); | 
| 767 | toUpper(myMethod); | 
| 768 | if (myMethod == "REACTION_FIELD"){ | 
| 769 | useRF = 1; | 
| 770 | } else if (myMethod == "SHIFTED_FORCE"){ | 
| 771 | useSF = 1; | 
| 772 | } else if (myMethod == "SHIFTED_POTENTIAL"){ | 
| 773 | useSP = 1; | 
| 774 | } | 
| 775 | } | 
| 776 |  | 
| 777 | if (simParams_->haveAccumulateBoxDipole()) | 
| 778 | if (simParams_->getAccumulateBoxDipole()) | 
| 779 | useBoxDipole = 1; | 
| 780 |  | 
| 781 | useAtomicVirial_ = simParams_->getUseAtomicVirial(); | 
| 782 |  | 
| 783 | //loop over all of the atom types | 
| 784 | for (i = atomTypes.begin(); i != atomTypes.end(); ++i) { | 
| 785 | useLennardJones |= (*i)->isLennardJones(); | 
| 786 | useElectrostatic |= (*i)->isElectrostatic(); | 
| 787 | useEAM |= (*i)->isEAM(); | 
| 788 | useSC |= (*i)->isSC(); | 
| 789 | useCharge |= (*i)->isCharge(); | 
| 790 | useDirectional |= (*i)->isDirectional(); | 
| 791 | useDipole |= (*i)->isDipole(); | 
| 792 | useGayBerne |= (*i)->isGayBerne(); | 
| 793 | useSticky |= (*i)->isSticky(); | 
| 794 | useStickyPower |= (*i)->isStickyPower(); | 
| 795 | useShape |= (*i)->isShape(); | 
| 796 | } | 
| 797 |  | 
| 798 | if (useSticky || useStickyPower || useDipole || useGayBerne || useShape) { | 
| 799 | useDirectionalAtom = 1; | 
| 800 | } | 
| 801 |  | 
| 802 | if (useCharge || useDipole) { | 
| 803 | useElectrostatics = 1; | 
| 804 | } | 
| 805 |  | 
| 806 | #ifdef IS_MPI | 
| 807 | int temp; | 
| 808 |  | 
| 809 | temp = usePBC; | 
| 810 | MPI_Allreduce(&temp, &usePBC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 811 |  | 
| 812 | temp = useDirectionalAtom; | 
| 813 | MPI_Allreduce(&temp, &useDirectionalAtom, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 814 |  | 
| 815 | temp = useLennardJones; | 
| 816 | MPI_Allreduce(&temp, &useLennardJones, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 817 |  | 
| 818 | temp = useElectrostatics; | 
| 819 | MPI_Allreduce(&temp, &useElectrostatics, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 820 |  | 
| 821 | temp = useCharge; | 
| 822 | MPI_Allreduce(&temp, &useCharge, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 823 |  | 
| 824 | temp = useDipole; | 
| 825 | MPI_Allreduce(&temp, &useDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 826 |  | 
| 827 | temp = useSticky; | 
| 828 | MPI_Allreduce(&temp, &useSticky, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 829 |  | 
| 830 | temp = useStickyPower; | 
| 831 | MPI_Allreduce(&temp, &useStickyPower, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 832 |  | 
| 833 | temp = useGayBerne; | 
| 834 | MPI_Allreduce(&temp, &useGayBerne, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 835 |  | 
| 836 | temp = useEAM; | 
| 837 | MPI_Allreduce(&temp, &useEAM, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 838 |  | 
| 839 | temp = useSC; | 
| 840 | MPI_Allreduce(&temp, &useSC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 841 |  | 
| 842 | temp = useShape; | 
| 843 | MPI_Allreduce(&temp, &useShape, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 844 |  | 
| 845 | temp = useFLARB; | 
| 846 | MPI_Allreduce(&temp, &useFLARB, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 847 |  | 
| 848 | temp = useRF; | 
| 849 | MPI_Allreduce(&temp, &useRF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 850 |  | 
| 851 | temp = useSF; | 
| 852 | MPI_Allreduce(&temp, &useSF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 853 |  | 
| 854 | temp = useSP; | 
| 855 | MPI_Allreduce(&temp, &useSP, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 856 |  | 
| 857 | temp = useBoxDipole; | 
| 858 | MPI_Allreduce(&temp, &useBoxDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 859 |  | 
| 860 | temp = useAtomicVirial_; | 
| 861 | MPI_Allreduce(&temp, &useAtomicVirial_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); | 
| 862 |  | 
| 863 | #endif | 
| 864 | fInfo_.SIM_uses_PBC = usePBC; | 
| 865 | fInfo_.SIM_uses_DirectionalAtoms = useDirectionalAtom; | 
| 866 | fInfo_.SIM_uses_LennardJones = useLennardJones; | 
| 867 | fInfo_.SIM_uses_Electrostatics = useElectrostatics; | 
| 868 | fInfo_.SIM_uses_Charges = useCharge; | 
| 869 | fInfo_.SIM_uses_Dipoles = useDipole; | 
| 870 | fInfo_.SIM_uses_Sticky = useSticky; | 
| 871 | fInfo_.SIM_uses_StickyPower = useStickyPower; | 
| 872 | fInfo_.SIM_uses_GayBerne = useGayBerne; | 
| 873 | fInfo_.SIM_uses_EAM = useEAM; | 
| 874 | fInfo_.SIM_uses_SC = useSC; | 
| 875 | fInfo_.SIM_uses_Shapes = useShape; | 
| 876 | fInfo_.SIM_uses_FLARB = useFLARB; | 
| 877 | fInfo_.SIM_uses_RF = useRF; | 
| 878 | fInfo_.SIM_uses_SF = useSF; | 
| 879 | fInfo_.SIM_uses_SP = useSP; | 
| 880 | fInfo_.SIM_uses_BoxDipole = useBoxDipole; | 
| 881 | fInfo_.SIM_uses_AtomicVirial = useAtomicVirial_; | 
| 882 | } | 
| 883 |  | 
| 884 | void SimInfo::setupFortranSim() { | 
| 885 | int isError; | 
| 886 | int nExclude, nOneTwo, nOneThree, nOneFour; | 
| 887 | std::vector<int> fortranGlobalGroupMembership; | 
| 888 |  | 
| 889 | isError = 0; | 
| 890 |  | 
| 891 | //globalGroupMembership_ is filled by SimCreator | 
| 892 | for (int i = 0; i < nGlobalAtoms_; i++) { | 
| 893 | fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1); | 
| 894 | } | 
| 895 |  | 
| 896 | //calculate mass ratio of cutoff group | 
| 897 | std::vector<RealType> mfact; | 
| 898 | SimInfo::MoleculeIterator mi; | 
| 899 | Molecule* mol; | 
| 900 | Molecule::CutoffGroupIterator ci; | 
| 901 | CutoffGroup* cg; | 
| 902 | Molecule::AtomIterator ai; | 
| 903 | Atom* atom; | 
| 904 | RealType totalMass; | 
| 905 |  | 
| 906 | //to avoid memory reallocation, reserve enough space for mfact | 
| 907 | mfact.reserve(getNCutoffGroups()); | 
| 908 |  | 
| 909 | for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) { | 
| 910 | for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) { | 
| 911 |  | 
| 912 | totalMass = cg->getMass(); | 
| 913 | for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) { | 
| 914 | // Check for massless groups - set mfact to 1 if true | 
| 915 | if (totalMass != 0) | 
| 916 | mfact.push_back(atom->getMass()/totalMass); | 
| 917 | else | 
| 918 | mfact.push_back( 1.0 ); | 
| 919 | } | 
| 920 | } | 
| 921 | } | 
| 922 |  | 
| 923 | //fill ident array of local atoms (it is actually ident of AtomType, it is so confusing !!!) | 
| 924 | std::vector<int> identArray; | 
| 925 |  | 
| 926 | //to avoid memory reallocation, reserve enough space identArray | 
| 927 | identArray.reserve(getNAtoms()); | 
| 928 |  | 
| 929 | for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) { | 
| 930 | for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) { | 
| 931 | identArray.push_back(atom->getIdent()); | 
| 932 | } | 
| 933 | } | 
| 934 |  | 
| 935 | //fill molMembershipArray | 
| 936 | //molMembershipArray is filled by SimCreator | 
| 937 | std::vector<int> molMembershipArray(nGlobalAtoms_); | 
| 938 | for (int i = 0; i < nGlobalAtoms_; i++) { | 
| 939 | molMembershipArray[i] = globalMolMembership_[i] + 1; | 
| 940 | } | 
| 941 |  | 
| 942 | //setup fortran simulation | 
| 943 |  | 
| 944 | nExclude = excludedInteractions_.getSize(); | 
| 945 | nOneTwo = oneTwoInteractions_.getSize(); | 
| 946 | nOneThree = oneThreeInteractions_.getSize(); | 
| 947 | nOneFour = oneFourInteractions_.getSize(); | 
| 948 |  | 
| 949 | int* excludeList = excludedInteractions_.getPairList(); | 
| 950 | int* oneTwoList = oneTwoInteractions_.getPairList(); | 
| 951 | int* oneThreeList = oneThreeInteractions_.getPairList(); | 
| 952 | int* oneFourList = oneFourInteractions_.getPairList(); | 
| 953 |  | 
| 954 | setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray[0], | 
| 955 | &nExclude, excludeList, | 
| 956 | &nOneTwo, oneTwoList, | 
| 957 | &nOneThree, oneThreeList, | 
| 958 | &nOneFour, oneFourList, | 
| 959 | &molMembershipArray[0], &mfact[0], &nCutoffGroups_, | 
| 960 | &fortranGlobalGroupMembership[0], &isError); | 
| 961 |  | 
| 962 | if( isError ){ | 
| 963 |  | 
| 964 | sprintf( painCave.errMsg, | 
| 965 | "There was an error setting the simulation information in fortran.\n" ); | 
| 966 | painCave.isFatal = 1; | 
| 967 | painCave.severity = OPENMD_ERROR; | 
| 968 | simError(); | 
| 969 | } | 
| 970 |  | 
| 971 |  | 
| 972 | sprintf( checkPointMsg, | 
| 973 | "succesfully sent the simulation information to fortran.\n"); | 
| 974 |  | 
| 975 | errorCheckPoint(); | 
| 976 |  | 
| 977 | // Setup number of neighbors in neighbor list if present | 
| 978 | if (simParams_->haveNeighborListNeighbors()) { | 
| 979 | int nlistNeighbors = simParams_->getNeighborListNeighbors(); | 
| 980 | setNeighbors(&nlistNeighbors); | 
| 981 | } | 
| 982 |  | 
| 983 |  | 
| 984 | } | 
| 985 |  | 
| 986 |  | 
| 987 | void SimInfo::setupFortranParallel() { | 
| 988 | #ifdef IS_MPI | 
| 989 | //SimInfo is responsible for creating localToGlobalAtomIndex and localToGlobalGroupIndex | 
| 990 | std::vector<int> localToGlobalAtomIndex(getNAtoms(), 0); | 
| 991 | std::vector<int> localToGlobalCutoffGroupIndex; | 
| 992 | SimInfo::MoleculeIterator mi; | 
| 993 | Molecule::AtomIterator ai; | 
| 994 | Molecule::CutoffGroupIterator ci; | 
| 995 | Molecule* mol; | 
| 996 | Atom* atom; | 
| 997 | CutoffGroup* cg; | 
| 998 | mpiSimData parallelData; | 
| 999 | int isError; | 
| 1000 |  | 
| 1001 | for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) { | 
| 1002 |  | 
| 1003 | //local index(index in DataStorge) of atom is important | 
| 1004 | for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) { | 
| 1005 | localToGlobalAtomIndex[atom->getLocalIndex()] = atom->getGlobalIndex() + 1; | 
| 1006 | } | 
| 1007 |  | 
| 1008 | //local index of cutoff group is trivial, it only depends on the order of travesing | 
| 1009 | for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) { | 
| 1010 | localToGlobalCutoffGroupIndex.push_back(cg->getGlobalIndex() + 1); | 
| 1011 | } | 
| 1012 |  | 
| 1013 | } | 
| 1014 |  | 
| 1015 | //fill up mpiSimData struct | 
| 1016 | parallelData.nMolGlobal = getNGlobalMolecules(); | 
| 1017 | parallelData.nMolLocal = getNMolecules(); | 
| 1018 | parallelData.nAtomsGlobal = getNGlobalAtoms(); | 
| 1019 | parallelData.nAtomsLocal = getNAtoms(); | 
| 1020 | parallelData.nGroupsGlobal = getNGlobalCutoffGroups(); | 
| 1021 | parallelData.nGroupsLocal = getNCutoffGroups(); | 
| 1022 | parallelData.myNode = worldRank; | 
| 1023 | MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors)); | 
| 1024 |  | 
| 1025 | //pass mpiSimData struct and index arrays to fortran | 
| 1026 | setFsimParallel(¶llelData, &(parallelData.nAtomsLocal), | 
| 1027 | &localToGlobalAtomIndex[0],  &(parallelData.nGroupsLocal), | 
| 1028 | &localToGlobalCutoffGroupIndex[0], &isError); | 
| 1029 |  | 
| 1030 | if (isError) { | 
| 1031 | sprintf(painCave.errMsg, | 
| 1032 | "mpiRefresh errror: fortran didn't like something we gave it.\n"); | 
| 1033 | painCave.isFatal = 1; | 
| 1034 | simError(); | 
| 1035 | } | 
| 1036 |  | 
| 1037 | sprintf(checkPointMsg, " mpiRefresh successful.\n"); | 
| 1038 | errorCheckPoint(); | 
| 1039 |  | 
| 1040 | #endif | 
| 1041 | } | 
| 1042 |  | 
| 1043 | void SimInfo::setupCutoff() { | 
| 1044 |  | 
| 1045 | ForceFieldOptions& forceFieldOptions_ = forceField_->getForceFieldOptions(); | 
| 1046 |  | 
| 1047 | // Check the cutoff policy | 
| 1048 | int cp =  TRADITIONAL_CUTOFF_POLICY; // Set to traditional by default | 
| 1049 |  | 
| 1050 | // Set LJ shifting bools to false | 
| 1051 | ljsp_ = 0; | 
| 1052 | ljsf_ = 0; | 
| 1053 |  | 
| 1054 | std::string myPolicy; | 
| 1055 | if (forceFieldOptions_.haveCutoffPolicy()){ | 
| 1056 | myPolicy = forceFieldOptions_.getCutoffPolicy(); | 
| 1057 | }else if (simParams_->haveCutoffPolicy()) { | 
| 1058 | myPolicy = simParams_->getCutoffPolicy(); | 
| 1059 | } | 
| 1060 |  | 
| 1061 | if (!myPolicy.empty()){ | 
| 1062 | toUpper(myPolicy); | 
| 1063 | if (myPolicy == "MIX") { | 
| 1064 | cp = MIX_CUTOFF_POLICY; | 
| 1065 | } else { | 
| 1066 | if (myPolicy == "MAX") { | 
| 1067 | cp = MAX_CUTOFF_POLICY; | 
| 1068 | } else { | 
| 1069 | if (myPolicy == "TRADITIONAL") { | 
| 1070 | cp = TRADITIONAL_CUTOFF_POLICY; | 
| 1071 | } else { | 
| 1072 | // throw error | 
| 1073 | sprintf( painCave.errMsg, | 
| 1074 | "SimInfo error: Unknown cutoffPolicy. (Input file specified %s .)\n\tcutoffPolicy must be one of: \"Mix\", \"Max\", or \"Traditional\".", myPolicy.c_str() ); | 
| 1075 | painCave.isFatal = 1; | 
| 1076 | simError(); | 
| 1077 | } | 
| 1078 | } | 
| 1079 | } | 
| 1080 | } | 
| 1081 | notifyFortranCutoffPolicy(&cp); | 
| 1082 |  | 
| 1083 | // Check the Skin Thickness for neighborlists | 
| 1084 | RealType skin; | 
| 1085 | if (simParams_->haveSkinThickness()) { | 
| 1086 | skin = simParams_->getSkinThickness(); | 
| 1087 | notifyFortranSkinThickness(&skin); | 
| 1088 | } | 
| 1089 |  | 
| 1090 | // Check if the cutoff was set explicitly: | 
| 1091 | if (simParams_->haveCutoffRadius()) { | 
| 1092 | rcut_ = simParams_->getCutoffRadius(); | 
| 1093 | if (simParams_->haveSwitchingRadius()) { | 
| 1094 | rsw_  = simParams_->getSwitchingRadius(); | 
| 1095 | } else { | 
| 1096 | if (fInfo_.SIM_uses_Charges | | 
| 1097 | fInfo_.SIM_uses_Dipoles | | 
| 1098 | fInfo_.SIM_uses_RF) { | 
| 1099 |  | 
| 1100 | rsw_ = 0.85 * rcut_; | 
| 1101 | sprintf(painCave.errMsg, | 
| 1102 | "SimCreator Warning: No value was set for the switchingRadius.\n" | 
| 1103 | "\tOpenMD will use a default value of 85 percent of the cutoffRadius.\n" | 
| 1104 | "\tswitchingRadius = %f. for this simulation\n", rsw_); | 
| 1105 | painCave.isFatal = 0; | 
| 1106 | simError(); | 
| 1107 | } else { | 
| 1108 | rsw_ = rcut_; | 
| 1109 | sprintf(painCave.errMsg, | 
| 1110 | "SimCreator Warning: No value was set for the switchingRadius.\n" | 
| 1111 | "\tOpenMD will use the same value as the cutoffRadius.\n" | 
| 1112 | "\tswitchingRadius = %f. for this simulation\n", rsw_); | 
| 1113 | painCave.isFatal = 0; | 
| 1114 | simError(); | 
| 1115 | } | 
| 1116 | } | 
| 1117 |  | 
| 1118 | if (simParams_->haveElectrostaticSummationMethod()) { | 
| 1119 | std::string myMethod = simParams_->getElectrostaticSummationMethod(); | 
| 1120 | toUpper(myMethod); | 
| 1121 |  | 
| 1122 | if (myMethod == "SHIFTED_POTENTIAL") { | 
| 1123 | ljsp_ = 1; | 
| 1124 | } else if (myMethod == "SHIFTED_FORCE") { | 
| 1125 | ljsf_ = 1; | 
| 1126 | } | 
| 1127 | } | 
| 1128 |  | 
| 1129 | notifyFortranCutoffs(&rcut_, &rsw_, &ljsp_, &ljsf_); | 
| 1130 |  | 
| 1131 | } else { | 
| 1132 |  | 
| 1133 | // For electrostatic atoms, we'll assume a large safe value: | 
| 1134 | if (fInfo_.SIM_uses_Charges | fInfo_.SIM_uses_Dipoles | fInfo_.SIM_uses_RF) { | 
| 1135 | sprintf(painCave.errMsg, | 
| 1136 | "SimCreator Warning: No value was set for the cutoffRadius.\n" | 
| 1137 | "\tOpenMD will use a default value of 15.0 angstroms" | 
| 1138 | "\tfor the cutoffRadius.\n"); | 
| 1139 | painCave.isFatal = 0; | 
| 1140 | simError(); | 
| 1141 | rcut_ = 15.0; | 
| 1142 |  | 
| 1143 | if (simParams_->haveElectrostaticSummationMethod()) { | 
| 1144 | std::string myMethod = simParams_->getElectrostaticSummationMethod(); | 
| 1145 | toUpper(myMethod); | 
| 1146 |  | 
| 1147 | // For the time being, we're tethering the LJ shifted behavior to the | 
| 1148 | // electrostaticSummationMethod keyword options | 
| 1149 | if (myMethod == "SHIFTED_POTENTIAL") { | 
| 1150 | ljsp_ = 1; | 
| 1151 | } else if (myMethod == "SHIFTED_FORCE") { | 
| 1152 | ljsf_ = 1; | 
| 1153 | } | 
| 1154 | if (myMethod == "SHIFTED_POTENTIAL" || myMethod == "SHIFTED_FORCE") { | 
| 1155 | if (simParams_->haveSwitchingRadius()){ | 
| 1156 | sprintf(painCave.errMsg, | 
| 1157 | "SimInfo Warning: A value was set for the switchingRadius\n" | 
| 1158 | "\teven though the electrostaticSummationMethod was\n" | 
| 1159 | "\tset to %s\n", myMethod.c_str()); | 
| 1160 | painCave.isFatal = 1; | 
| 1161 | simError(); | 
| 1162 | } | 
| 1163 | } | 
| 1164 | } | 
| 1165 |  | 
| 1166 | if (simParams_->haveSwitchingRadius()){ | 
| 1167 | rsw_ = simParams_->getSwitchingRadius(); | 
| 1168 | } else { | 
| 1169 | sprintf(painCave.errMsg, | 
| 1170 | "SimCreator Warning: No value was set for switchingRadius.\n" | 
| 1171 | "\tOpenMD will use a default value of\n" | 
| 1172 | "\t0.85 * cutoffRadius for the switchingRadius\n"); | 
| 1173 | painCave.isFatal = 0; | 
| 1174 | simError(); | 
| 1175 | rsw_ = 0.85 * rcut_; | 
| 1176 | } | 
| 1177 |  | 
| 1178 | Electrostatic::setElectrostaticCutoffRadius(rcut_, rsw_); | 
| 1179 | notifyFortranCutoffs(&rcut_, &rsw_, &ljsp_, &ljsf_); | 
| 1180 |  | 
| 1181 | } else { | 
| 1182 | // We didn't set rcut explicitly, and we don't have electrostatic atoms, so | 
| 1183 | // We'll punt and let fortran figure out the cutoffs later. | 
| 1184 |  | 
| 1185 | notifyFortranYouAreOnYourOwn(); | 
| 1186 |  | 
| 1187 | } | 
| 1188 | } | 
| 1189 | } | 
| 1190 |  | 
| 1191 | void SimInfo::setupElectrostaticSummationMethod( int isError ) { | 
| 1192 |  | 
| 1193 | int errorOut; | 
| 1194 | ElectrostaticSummationMethod esm = NONE; | 
| 1195 | ElectrostaticScreeningMethod sm = UNDAMPED; | 
| 1196 | RealType alphaVal; | 
| 1197 | RealType dielectric; | 
| 1198 |  | 
| 1199 | errorOut = isError; | 
| 1200 |  | 
| 1201 | if (simParams_->haveElectrostaticSummationMethod()) { | 
| 1202 | std::string myMethod = simParams_->getElectrostaticSummationMethod(); | 
| 1203 | toUpper(myMethod); | 
| 1204 | if (myMethod == "NONE") { | 
| 1205 | esm = NONE; | 
| 1206 | } else { | 
| 1207 | if (myMethod == "SWITCHING_FUNCTION") { | 
| 1208 | esm = SWITCHING_FUNCTION; | 
| 1209 | } else { | 
| 1210 | if (myMethod == "SHIFTED_POTENTIAL") { | 
| 1211 | esm = SHIFTED_POTENTIAL; | 
| 1212 | } else { | 
| 1213 | if (myMethod == "SHIFTED_FORCE") { | 
| 1214 | esm = SHIFTED_FORCE; | 
| 1215 | } else { | 
| 1216 | if (myMethod == "REACTION_FIELD") { | 
| 1217 | esm = REACTION_FIELD; | 
| 1218 | dielectric = simParams_->getDielectric(); | 
| 1219 | if (!simParams_->haveDielectric()) { | 
| 1220 | // throw warning | 
| 1221 | sprintf( painCave.errMsg, | 
| 1222 | "SimInfo warning: dielectric was not specified in the input file\n\tfor the reaction field correction method.\n" | 
| 1223 | "\tA default value of %f will be used for the dielectric.\n", dielectric); | 
| 1224 | painCave.isFatal = 0; | 
| 1225 | simError(); | 
| 1226 | } | 
| 1227 | } else { | 
| 1228 | // throw error | 
| 1229 | sprintf( painCave.errMsg, | 
| 1230 | "SimInfo error: Unknown electrostaticSummationMethod.\n" | 
| 1231 | "\t(Input file specified %s .)\n" | 
| 1232 | "\telectrostaticSummationMethod must be one of: \"none\",\n" | 
| 1233 | "\t\"shifted_potential\", \"shifted_force\", or \n" | 
| 1234 | "\t\"reaction_field\".\n", myMethod.c_str() ); | 
| 1235 | painCave.isFatal = 1; | 
| 1236 | simError(); | 
| 1237 | } | 
| 1238 | } | 
| 1239 | } | 
| 1240 | } | 
| 1241 | } | 
| 1242 | } | 
| 1243 |  | 
| 1244 | if (simParams_->haveElectrostaticScreeningMethod()) { | 
| 1245 | std::string myScreen = simParams_->getElectrostaticScreeningMethod(); | 
| 1246 | toUpper(myScreen); | 
| 1247 | if (myScreen == "UNDAMPED") { | 
| 1248 | sm = UNDAMPED; | 
| 1249 | } else { | 
| 1250 | if (myScreen == "DAMPED") { | 
| 1251 | sm = DAMPED; | 
| 1252 | if (!simParams_->haveDampingAlpha()) { | 
| 1253 | // first set a cutoff dependent alpha value | 
| 1254 | // we assume alpha depends linearly with rcut from 0 to 20.5 ang | 
| 1255 | alphaVal = 0.5125 - rcut_* 0.025; | 
| 1256 | // for values rcut > 20.5, alpha is zero | 
| 1257 | if (alphaVal < 0) alphaVal = 0; | 
| 1258 |  | 
| 1259 | // throw warning | 
| 1260 | sprintf( painCave.errMsg, | 
| 1261 | "SimInfo warning: dampingAlpha was not specified in the input file.\n" | 
| 1262 | "\tA default value of %f (1/ang) will be used for the cutoff of\n\t%f (ang).\n", alphaVal, rcut_); | 
| 1263 | painCave.isFatal = 0; | 
| 1264 | simError(); | 
| 1265 | } else { | 
| 1266 | alphaVal = simParams_->getDampingAlpha(); | 
| 1267 | } | 
| 1268 |  | 
| 1269 | } else { | 
| 1270 | // throw error | 
| 1271 | sprintf( painCave.errMsg, | 
| 1272 | "SimInfo error: Unknown electrostaticScreeningMethod.\n" | 
| 1273 | "\t(Input file specified %s .)\n" | 
| 1274 | "\telectrostaticScreeningMethod must be one of: \"undamped\"\n" | 
| 1275 | "or \"damped\".\n", myScreen.c_str() ); | 
| 1276 | painCave.isFatal = 1; | 
| 1277 | simError(); | 
| 1278 | } | 
| 1279 | } | 
| 1280 | } | 
| 1281 |  | 
| 1282 |  | 
| 1283 | Electrostatic::setElectrostaticSummationMethod( esm ); | 
| 1284 | Electrostatic::setElectrostaticScreeningMethod( sm ); | 
| 1285 | Electrostatic::setDampingAlpha( alphaVal ); | 
| 1286 | Electrostatic::setReactionFieldDielectric( dielectric ); | 
| 1287 | initFortranFF( &errorOut ); | 
| 1288 | } | 
| 1289 |  | 
| 1290 | void SimInfo::setupSwitchingFunction() { | 
| 1291 | int ft = CUBIC; | 
| 1292 |  | 
| 1293 | if (simParams_->haveSwitchingFunctionType()) { | 
| 1294 | std::string funcType = simParams_->getSwitchingFunctionType(); | 
| 1295 | toUpper(funcType); | 
| 1296 | if (funcType == "CUBIC") { | 
| 1297 | ft = CUBIC; | 
| 1298 | } else { | 
| 1299 | if (funcType == "FIFTH_ORDER_POLYNOMIAL") { | 
| 1300 | ft = FIFTH_ORDER_POLY; | 
| 1301 | } else { | 
| 1302 | // throw error | 
| 1303 | sprintf( painCave.errMsg, | 
| 1304 | "SimInfo error: Unknown switchingFunctionType. (Input file specified %s .)\n\tswitchingFunctionType must be one of: \"cubic\" or \"fifth_order_polynomial\".", funcType.c_str() ); | 
| 1305 | painCave.isFatal = 1; | 
| 1306 | simError(); | 
| 1307 | } | 
| 1308 | } | 
| 1309 | } | 
| 1310 |  | 
| 1311 | // send switching function notification to switcheroo | 
| 1312 | setFunctionType(&ft); | 
| 1313 |  | 
| 1314 | } | 
| 1315 |  | 
| 1316 | void SimInfo::setupAccumulateBoxDipole() { | 
| 1317 |  | 
| 1318 | // we only call setAccumulateBoxDipole if the accumulateBoxDipole parameter is true | 
| 1319 | if ( simParams_->haveAccumulateBoxDipole() ) | 
| 1320 | if ( simParams_->getAccumulateBoxDipole() ) { | 
| 1321 | setAccumulateBoxDipole(); | 
| 1322 | calcBoxDipole_ = true; | 
| 1323 | } | 
| 1324 |  | 
| 1325 | } | 
| 1326 |  | 
| 1327 | void SimInfo::addProperty(GenericData* genData) { | 
| 1328 | properties_.addProperty(genData); | 
| 1329 | } | 
| 1330 |  | 
| 1331 | void SimInfo::removeProperty(const std::string& propName) { | 
| 1332 | properties_.removeProperty(propName); | 
| 1333 | } | 
| 1334 |  | 
| 1335 | void SimInfo::clearProperties() { | 
| 1336 | properties_.clearProperties(); | 
| 1337 | } | 
| 1338 |  | 
| 1339 | std::vector<std::string> SimInfo::getPropertyNames() { | 
| 1340 | return properties_.getPropertyNames(); | 
| 1341 | } | 
| 1342 |  | 
| 1343 | std::vector<GenericData*> SimInfo::getProperties() { | 
| 1344 | return properties_.getProperties(); | 
| 1345 | } | 
| 1346 |  | 
| 1347 | GenericData* SimInfo::getPropertyByName(const std::string& propName) { | 
| 1348 | return properties_.getPropertyByName(propName); | 
| 1349 | } | 
| 1350 |  | 
| 1351 | void SimInfo::setSnapshotManager(SnapshotManager* sman) { | 
| 1352 | if (sman_ == sman) { | 
| 1353 | return; | 
| 1354 | } | 
| 1355 | delete sman_; | 
| 1356 | sman_ = sman; | 
| 1357 |  | 
| 1358 | Molecule* mol; | 
| 1359 | RigidBody* rb; | 
| 1360 | Atom* atom; | 
| 1361 | SimInfo::MoleculeIterator mi; | 
| 1362 | Molecule::RigidBodyIterator rbIter; | 
| 1363 | Molecule::AtomIterator atomIter;; | 
| 1364 |  | 
| 1365 | for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) { | 
| 1366 |  | 
| 1367 | for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) { | 
| 1368 | atom->setSnapshotManager(sman_); | 
| 1369 | } | 
| 1370 |  | 
| 1371 | for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) { | 
| 1372 | rb->setSnapshotManager(sman_); | 
| 1373 | } | 
| 1374 | } | 
| 1375 |  | 
| 1376 | } | 
| 1377 |  | 
| 1378 | Vector3d SimInfo::getComVel(){ | 
| 1379 | SimInfo::MoleculeIterator i; | 
| 1380 | Molecule* mol; | 
| 1381 |  | 
| 1382 | Vector3d comVel(0.0); | 
| 1383 | RealType totalMass = 0.0; | 
| 1384 |  | 
| 1385 |  | 
| 1386 | for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { | 
| 1387 | RealType mass = mol->getMass(); | 
| 1388 | totalMass += mass; | 
| 1389 | comVel += mass * mol->getComVel(); | 
| 1390 | } | 
| 1391 |  | 
| 1392 | #ifdef IS_MPI | 
| 1393 | RealType tmpMass = totalMass; | 
| 1394 | Vector3d tmpComVel(comVel); | 
| 1395 | MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); | 
| 1396 | MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); | 
| 1397 | #endif | 
| 1398 |  | 
| 1399 | comVel /= totalMass; | 
| 1400 |  | 
| 1401 | return comVel; | 
| 1402 | } | 
| 1403 |  | 
| 1404 | Vector3d SimInfo::getCom(){ | 
| 1405 | SimInfo::MoleculeIterator i; | 
| 1406 | Molecule* mol; | 
| 1407 |  | 
| 1408 | Vector3d com(0.0); | 
| 1409 | RealType totalMass = 0.0; | 
| 1410 |  | 
| 1411 | for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { | 
| 1412 | RealType mass = mol->getMass(); | 
| 1413 | totalMass += mass; | 
| 1414 | com += mass * mol->getCom(); | 
| 1415 | } | 
| 1416 |  | 
| 1417 | #ifdef IS_MPI | 
| 1418 | RealType tmpMass = totalMass; | 
| 1419 | Vector3d tmpCom(com); | 
| 1420 | MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); | 
| 1421 | MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); | 
| 1422 | #endif | 
| 1423 |  | 
| 1424 | com /= totalMass; | 
| 1425 |  | 
| 1426 | return com; | 
| 1427 |  | 
| 1428 | } | 
| 1429 |  | 
| 1430 | std::ostream& operator <<(std::ostream& o, SimInfo& info) { | 
| 1431 |  | 
| 1432 | return o; | 
| 1433 | } | 
| 1434 |  | 
| 1435 |  | 
| 1436 | /* | 
| 1437 | Returns center of mass and center of mass velocity in one function call. | 
| 1438 | */ | 
| 1439 |  | 
| 1440 | void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){ | 
| 1441 | SimInfo::MoleculeIterator i; | 
| 1442 | Molecule* mol; | 
| 1443 |  | 
| 1444 |  | 
| 1445 | RealType totalMass = 0.0; | 
| 1446 |  | 
| 1447 |  | 
| 1448 | for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { | 
| 1449 | RealType mass = mol->getMass(); | 
| 1450 | totalMass += mass; | 
| 1451 | com += mass * mol->getCom(); | 
| 1452 | comVel += mass * mol->getComVel(); | 
| 1453 | } | 
| 1454 |  | 
| 1455 | #ifdef IS_MPI | 
| 1456 | RealType tmpMass = totalMass; | 
| 1457 | Vector3d tmpCom(com); | 
| 1458 | Vector3d tmpComVel(comVel); | 
| 1459 | MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); | 
| 1460 | MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); | 
| 1461 | MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); | 
| 1462 | #endif | 
| 1463 |  | 
| 1464 | com /= totalMass; | 
| 1465 | comVel /= totalMass; | 
| 1466 | } | 
| 1467 |  | 
| 1468 | /* | 
| 1469 | Return intertia tensor for entire system and angular momentum Vector. | 
| 1470 |  | 
| 1471 |  | 
| 1472 | [  Ixx -Ixy  -Ixz ] | 
| 1473 | J =| -Iyx  Iyy  -Iyz | | 
| 1474 | [ -Izx -Iyz   Izz ] | 
| 1475 | */ | 
| 1476 |  | 
| 1477 | void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){ | 
| 1478 |  | 
| 1479 |  | 
| 1480 | RealType xx = 0.0; | 
| 1481 | RealType yy = 0.0; | 
| 1482 | RealType zz = 0.0; | 
| 1483 | RealType xy = 0.0; | 
| 1484 | RealType xz = 0.0; | 
| 1485 | RealType yz = 0.0; | 
| 1486 | Vector3d com(0.0); | 
| 1487 | Vector3d comVel(0.0); | 
| 1488 |  | 
| 1489 | getComAll(com, comVel); | 
| 1490 |  | 
| 1491 | SimInfo::MoleculeIterator i; | 
| 1492 | Molecule* mol; | 
| 1493 |  | 
| 1494 | Vector3d thisq(0.0); | 
| 1495 | Vector3d thisv(0.0); | 
| 1496 |  | 
| 1497 | RealType thisMass = 0.0; | 
| 1498 |  | 
| 1499 |  | 
| 1500 |  | 
| 1501 |  | 
| 1502 | for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { | 
| 1503 |  | 
| 1504 | thisq = mol->getCom()-com; | 
| 1505 | thisv = mol->getComVel()-comVel; | 
| 1506 | thisMass = mol->getMass(); | 
| 1507 | // Compute moment of intertia coefficients. | 
| 1508 | xx += thisq[0]*thisq[0]*thisMass; | 
| 1509 | yy += thisq[1]*thisq[1]*thisMass; | 
| 1510 | zz += thisq[2]*thisq[2]*thisMass; | 
| 1511 |  | 
| 1512 | // compute products of intertia | 
| 1513 | xy += thisq[0]*thisq[1]*thisMass; | 
| 1514 | xz += thisq[0]*thisq[2]*thisMass; | 
| 1515 | yz += thisq[1]*thisq[2]*thisMass; | 
| 1516 |  | 
| 1517 | angularMomentum += cross( thisq, thisv ) * thisMass; | 
| 1518 |  | 
| 1519 | } | 
| 1520 |  | 
| 1521 |  | 
| 1522 | inertiaTensor(0,0) = yy + zz; | 
| 1523 | inertiaTensor(0,1) = -xy; | 
| 1524 | inertiaTensor(0,2) = -xz; | 
| 1525 | inertiaTensor(1,0) = -xy; | 
| 1526 | inertiaTensor(1,1) = xx + zz; | 
| 1527 | inertiaTensor(1,2) = -yz; | 
| 1528 | inertiaTensor(2,0) = -xz; | 
| 1529 | inertiaTensor(2,1) = -yz; | 
| 1530 | inertiaTensor(2,2) = xx + yy; | 
| 1531 |  | 
| 1532 | #ifdef IS_MPI | 
| 1533 | Mat3x3d tmpI(inertiaTensor); | 
| 1534 | Vector3d tmpAngMom; | 
| 1535 | MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); | 
| 1536 | MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); | 
| 1537 | #endif | 
| 1538 |  | 
| 1539 | return; | 
| 1540 | } | 
| 1541 |  | 
| 1542 | //Returns the angular momentum of the system | 
| 1543 | Vector3d SimInfo::getAngularMomentum(){ | 
| 1544 |  | 
| 1545 | Vector3d com(0.0); | 
| 1546 | Vector3d comVel(0.0); | 
| 1547 | Vector3d angularMomentum(0.0); | 
| 1548 |  | 
| 1549 | getComAll(com,comVel); | 
| 1550 |  | 
| 1551 | SimInfo::MoleculeIterator i; | 
| 1552 | Molecule* mol; | 
| 1553 |  | 
| 1554 | Vector3d thisr(0.0); | 
| 1555 | Vector3d thisp(0.0); | 
| 1556 |  | 
| 1557 | RealType thisMass; | 
| 1558 |  | 
| 1559 | for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { | 
| 1560 | thisMass = mol->getMass(); | 
| 1561 | thisr = mol->getCom()-com; | 
| 1562 | thisp = (mol->getComVel()-comVel)*thisMass; | 
| 1563 |  | 
| 1564 | angularMomentum += cross( thisr, thisp ); | 
| 1565 |  | 
| 1566 | } | 
| 1567 |  | 
| 1568 | #ifdef IS_MPI | 
| 1569 | Vector3d tmpAngMom; | 
| 1570 | MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD); | 
| 1571 | #endif | 
| 1572 |  | 
| 1573 | return angularMomentum; | 
| 1574 | } | 
| 1575 |  | 
| 1576 | StuntDouble* SimInfo::getIOIndexToIntegrableObject(int index) { | 
| 1577 | return IOIndexToIntegrableObject.at(index); | 
| 1578 | } | 
| 1579 |  | 
| 1580 | void SimInfo::setIOIndexToIntegrableObject(const std::vector<StuntDouble*>& v) { | 
| 1581 | IOIndexToIntegrableObject= v; | 
| 1582 | } | 
| 1583 |  | 
| 1584 | /* Returns the Volume of the simulation based on a ellipsoid with semi-axes | 
| 1585 | based on the radius of gyration V=4/3*Pi*R_1*R_2*R_3 | 
| 1586 | where R_i are related to the principle inertia moments R_i = sqrt(C*I_i/N), this reduces to | 
| 1587 | V = 4/3*Pi*(C/N)^3/2*sqrt(det(I)). See S.E. Baltazar et. al. Comp. Mat. Sci. 37 (2006) 526-536. | 
| 1588 | */ | 
| 1589 | void SimInfo::getGyrationalVolume(RealType &volume){ | 
| 1590 | Mat3x3d intTensor; | 
| 1591 | RealType det; | 
| 1592 | Vector3d dummyAngMom; | 
| 1593 | RealType sysconstants; | 
| 1594 | RealType geomCnst; | 
| 1595 |  | 
| 1596 | geomCnst = 3.0/2.0; | 
| 1597 | /* Get the inertial tensor and angular momentum for free*/ | 
| 1598 | getInertiaTensor(intTensor,dummyAngMom); | 
| 1599 |  | 
| 1600 | det = intTensor.determinant(); | 
| 1601 | sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_; | 
| 1602 | volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(det); | 
| 1603 | return; | 
| 1604 | } | 
| 1605 |  | 
| 1606 | void SimInfo::getGyrationalVolume(RealType &volume, RealType &detI){ | 
| 1607 | Mat3x3d intTensor; | 
| 1608 | Vector3d dummyAngMom; | 
| 1609 | RealType sysconstants; | 
| 1610 | RealType geomCnst; | 
| 1611 |  | 
| 1612 | geomCnst = 3.0/2.0; | 
| 1613 | /* Get the inertial tensor and angular momentum for free*/ | 
| 1614 | getInertiaTensor(intTensor,dummyAngMom); | 
| 1615 |  | 
| 1616 | detI = intTensor.determinant(); | 
| 1617 | sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_; | 
| 1618 | volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(detI); | 
| 1619 | return; | 
| 1620 | } | 
| 1621 | /* | 
| 1622 | void SimInfo::setStuntDoubleFromGlobalIndex(std::vector<StuntDouble*> v) { | 
| 1623 | assert( v.size() == nAtoms_ + nRigidBodies_); | 
| 1624 | sdByGlobalIndex_ = v; | 
| 1625 | } | 
| 1626 |  | 
| 1627 | StuntDouble* SimInfo::getStuntDoubleFromGlobalIndex(int index) { | 
| 1628 | //assert(index < nAtoms_ + nRigidBodies_); | 
| 1629 | return sdByGlobalIndex_.at(index); | 
| 1630 | } | 
| 1631 | */ | 
| 1632 | }//end namespace OpenMD | 
| 1633 |  |