| 1 |
gezelter |
507 |
/* |
| 2 |
gezelter |
246 |
* 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. Acknowledgement of the program authors must be made in any |
| 10 |
|
|
* publication of scientific results based in part on use of the |
| 11 |
|
|
* program. An acceptable form of acknowledgement is citation of |
| 12 |
|
|
* the article in which the program was described (Matthew |
| 13 |
|
|
* A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher |
| 14 |
|
|
* J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented |
| 15 |
|
|
* Parallel Simulation Engine for Molecular Dynamics," |
| 16 |
|
|
* J. Comput. Chem. 26, pp. 252-271 (2005)) |
| 17 |
|
|
* |
| 18 |
|
|
* 2. Redistributions of source code must retain the above copyright |
| 19 |
|
|
* notice, this list of conditions and the following disclaimer. |
| 20 |
|
|
* |
| 21 |
|
|
* 3. Redistributions in binary form must reproduce the above copyright |
| 22 |
|
|
* notice, this list of conditions and the following disclaimer in the |
| 23 |
|
|
* documentation and/or other materials provided with the |
| 24 |
|
|
* distribution. |
| 25 |
|
|
* |
| 26 |
|
|
* This software is provided "AS IS," without a warranty of any |
| 27 |
|
|
* kind. All express or implied conditions, representations and |
| 28 |
|
|
* warranties, including any implied warranty of merchantability, |
| 29 |
|
|
* fitness for a particular purpose or non-infringement, are hereby |
| 30 |
|
|
* excluded. The University of Notre Dame and its licensors shall not |
| 31 |
|
|
* be liable for any damages suffered by licensee as a result of |
| 32 |
|
|
* using, modifying or distributing the software or its |
| 33 |
|
|
* derivatives. In no event will the University of Notre Dame or its |
| 34 |
|
|
* licensors be liable for any lost revenue, profit or data, or for |
| 35 |
|
|
* direct, indirect, special, consequential, incidental or punitive |
| 36 |
|
|
* damages, however caused and regardless of the theory of liability, |
| 37 |
|
|
* arising out of the use of or inability to use software, even if the |
| 38 |
|
|
* University of Notre Dame has been advised of the possibility of |
| 39 |
|
|
* such damages. |
| 40 |
|
|
*/ |
| 41 |
|
|
|
| 42 |
|
|
/** |
| 43 |
|
|
* @file SimInfo.cpp |
| 44 |
|
|
* @author tlin |
| 45 |
|
|
* @date 11/02/2004 |
| 46 |
|
|
* @version 1.0 |
| 47 |
|
|
*/ |
| 48 |
gezelter |
2 |
|
| 49 |
gezelter |
246 |
#include <algorithm> |
| 50 |
|
|
#include <set> |
| 51 |
gezelter |
2 |
|
| 52 |
tim |
3 |
#include "brains/SimInfo.hpp" |
| 53 |
gezelter |
246 |
#include "math/Vector3.hpp" |
| 54 |
|
|
#include "primitives/Molecule.hpp" |
| 55 |
gezelter |
586 |
#include "UseTheForce/fCutoffPolicy.h" |
| 56 |
chrisfen |
606 |
#include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h" |
| 57 |
gezelter |
246 |
#include "UseTheForce/doForces_interface.h" |
| 58 |
chrisfen |
610 |
#include "UseTheForce/DarkSide/electrostatic_interface.h" |
| 59 |
gezelter |
246 |
#include "UseTheForce/notifyCutoffs_interface.h" |
| 60 |
|
|
#include "utils/MemoryUtils.hpp" |
| 61 |
tim |
3 |
#include "utils/simError.h" |
| 62 |
tim |
316 |
#include "selection/SelectionManager.hpp" |
| 63 |
gezelter |
2 |
|
| 64 |
gezelter |
246 |
#ifdef IS_MPI |
| 65 |
|
|
#include "UseTheForce/mpiComponentPlan.h" |
| 66 |
|
|
#include "UseTheForce/DarkSide/simParallel_interface.h" |
| 67 |
|
|
#endif |
| 68 |
gezelter |
2 |
|
| 69 |
gezelter |
246 |
namespace oopse { |
| 70 |
gezelter |
2 |
|
| 71 |
gezelter |
507 |
SimInfo::SimInfo(MakeStamps* stamps, std::vector<std::pair<MoleculeStamp*, int> >& molStampPairs, |
| 72 |
|
|
ForceField* ff, Globals* simParams) : |
| 73 |
|
|
stamps_(stamps), forceField_(ff), simParams_(simParams), |
| 74 |
|
|
ndf_(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0), |
| 75 |
|
|
nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0), |
| 76 |
|
|
nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0), |
| 77 |
|
|
nAtoms_(0), nBonds_(0), nBends_(0), nTorsions_(0), nRigidBodies_(0), |
| 78 |
|
|
nIntegrableObjects_(0), nCutoffGroups_(0), nConstraints_(0), |
| 79 |
|
|
sman_(NULL), fortranInitialized_(false) { |
| 80 |
gezelter |
2 |
|
| 81 |
gezelter |
246 |
|
| 82 |
gezelter |
507 |
std::vector<std::pair<MoleculeStamp*, int> >::iterator i; |
| 83 |
|
|
MoleculeStamp* molStamp; |
| 84 |
|
|
int nMolWithSameStamp; |
| 85 |
|
|
int nCutoffAtoms = 0; // number of atoms belong to cutoff groups |
| 86 |
chrisfen |
645 |
int nGroups = 0; //total cutoff groups defined in meta-data file |
| 87 |
gezelter |
507 |
CutoffGroupStamp* cgStamp; |
| 88 |
|
|
RigidBodyStamp* rbStamp; |
| 89 |
|
|
int nRigidAtoms = 0; |
| 90 |
gezelter |
246 |
|
| 91 |
gezelter |
507 |
for (i = molStampPairs.begin(); i !=molStampPairs.end(); ++i) { |
| 92 |
gezelter |
246 |
molStamp = i->first; |
| 93 |
|
|
nMolWithSameStamp = i->second; |
| 94 |
|
|
|
| 95 |
|
|
addMoleculeStamp(molStamp, nMolWithSameStamp); |
| 96 |
gezelter |
2 |
|
| 97 |
gezelter |
246 |
//calculate atoms in molecules |
| 98 |
|
|
nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp; |
| 99 |
gezelter |
2 |
|
| 100 |
|
|
|
| 101 |
gezelter |
246 |
//calculate atoms in cutoff groups |
| 102 |
|
|
int nAtomsInGroups = 0; |
| 103 |
|
|
int nCutoffGroupsInStamp = molStamp->getNCutoffGroups(); |
| 104 |
|
|
|
| 105 |
|
|
for (int j=0; j < nCutoffGroupsInStamp; j++) { |
| 106 |
gezelter |
507 |
cgStamp = molStamp->getCutoffGroup(j); |
| 107 |
|
|
nAtomsInGroups += cgStamp->getNMembers(); |
| 108 |
gezelter |
246 |
} |
| 109 |
gezelter |
2 |
|
| 110 |
gezelter |
246 |
nGroups += nCutoffGroupsInStamp * nMolWithSameStamp; |
| 111 |
chrisfen |
645 |
|
| 112 |
gezelter |
246 |
nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp; |
| 113 |
gezelter |
2 |
|
| 114 |
gezelter |
246 |
//calculate atoms in rigid bodies |
| 115 |
|
|
int nAtomsInRigidBodies = 0; |
| 116 |
tim |
274 |
int nRigidBodiesInStamp = molStamp->getNRigidBodies(); |
| 117 |
gezelter |
246 |
|
| 118 |
|
|
for (int j=0; j < nRigidBodiesInStamp; j++) { |
| 119 |
gezelter |
507 |
rbStamp = molStamp->getRigidBody(j); |
| 120 |
|
|
nAtomsInRigidBodies += rbStamp->getNMembers(); |
| 121 |
gezelter |
246 |
} |
| 122 |
gezelter |
2 |
|
| 123 |
gezelter |
246 |
nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp; |
| 124 |
|
|
nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp; |
| 125 |
|
|
|
| 126 |
gezelter |
507 |
} |
| 127 |
chrisfen |
143 |
|
| 128 |
chrisfen |
645 |
//every free atom (atom does not belong to cutoff groups) is a cutoff |
| 129 |
|
|
//group therefore the total number of cutoff groups in the system is |
| 130 |
|
|
//equal to the total number of atoms minus number of atoms belong to |
| 131 |
|
|
//cutoff group defined in meta-data file plus the number of cutoff |
| 132 |
|
|
//groups defined in meta-data file |
| 133 |
gezelter |
507 |
nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups; |
| 134 |
gezelter |
2 |
|
| 135 |
chrisfen |
645 |
//every free atom (atom does not belong to rigid bodies) is an |
| 136 |
|
|
//integrable object therefore the total number of integrable objects |
| 137 |
|
|
//in the system is equal to the total number of atoms minus number of |
| 138 |
|
|
//atoms belong to rigid body defined in meta-data file plus the number |
| 139 |
|
|
//of rigid bodies defined in meta-data file |
| 140 |
|
|
nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms |
| 141 |
|
|
+ nGlobalRigidBodies_; |
| 142 |
|
|
|
| 143 |
gezelter |
507 |
nGlobalMols_ = molStampIds_.size(); |
| 144 |
gezelter |
2 |
|
| 145 |
gezelter |
246 |
#ifdef IS_MPI |
| 146 |
gezelter |
507 |
molToProcMap_.resize(nGlobalMols_); |
| 147 |
gezelter |
246 |
#endif |
| 148 |
tim |
292 |
|
| 149 |
gezelter |
507 |
} |
| 150 |
gezelter |
2 |
|
| 151 |
gezelter |
507 |
SimInfo::~SimInfo() { |
| 152 |
tim |
398 |
std::map<int, Molecule*>::iterator i; |
| 153 |
|
|
for (i = molecules_.begin(); i != molecules_.end(); ++i) { |
| 154 |
gezelter |
507 |
delete i->second; |
| 155 |
tim |
398 |
} |
| 156 |
|
|
molecules_.clear(); |
| 157 |
tim |
490 |
|
| 158 |
|
|
delete stamps_; |
| 159 |
gezelter |
246 |
delete sman_; |
| 160 |
|
|
delete simParams_; |
| 161 |
|
|
delete forceField_; |
| 162 |
gezelter |
507 |
} |
| 163 |
gezelter |
2 |
|
| 164 |
gezelter |
507 |
int SimInfo::getNGlobalConstraints() { |
| 165 |
gezelter |
246 |
int nGlobalConstraints; |
| 166 |
|
|
#ifdef IS_MPI |
| 167 |
|
|
MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM, |
| 168 |
|
|
MPI_COMM_WORLD); |
| 169 |
|
|
#else |
| 170 |
|
|
nGlobalConstraints = nConstraints_; |
| 171 |
|
|
#endif |
| 172 |
|
|
return nGlobalConstraints; |
| 173 |
gezelter |
507 |
} |
| 174 |
gezelter |
2 |
|
| 175 |
gezelter |
507 |
bool SimInfo::addMolecule(Molecule* mol) { |
| 176 |
gezelter |
246 |
MoleculeIterator i; |
| 177 |
gezelter |
2 |
|
| 178 |
gezelter |
246 |
i = molecules_.find(mol->getGlobalIndex()); |
| 179 |
|
|
if (i == molecules_.end() ) { |
| 180 |
gezelter |
2 |
|
| 181 |
gezelter |
507 |
molecules_.insert(std::make_pair(mol->getGlobalIndex(), mol)); |
| 182 |
gezelter |
246 |
|
| 183 |
gezelter |
507 |
nAtoms_ += mol->getNAtoms(); |
| 184 |
|
|
nBonds_ += mol->getNBonds(); |
| 185 |
|
|
nBends_ += mol->getNBends(); |
| 186 |
|
|
nTorsions_ += mol->getNTorsions(); |
| 187 |
|
|
nRigidBodies_ += mol->getNRigidBodies(); |
| 188 |
|
|
nIntegrableObjects_ += mol->getNIntegrableObjects(); |
| 189 |
|
|
nCutoffGroups_ += mol->getNCutoffGroups(); |
| 190 |
|
|
nConstraints_ += mol->getNConstraintPairs(); |
| 191 |
gezelter |
2 |
|
| 192 |
gezelter |
507 |
addExcludePairs(mol); |
| 193 |
gezelter |
246 |
|
| 194 |
gezelter |
507 |
return true; |
| 195 |
gezelter |
246 |
} else { |
| 196 |
gezelter |
507 |
return false; |
| 197 |
gezelter |
246 |
} |
| 198 |
gezelter |
507 |
} |
| 199 |
gezelter |
2 |
|
| 200 |
gezelter |
507 |
bool SimInfo::removeMolecule(Molecule* mol) { |
| 201 |
gezelter |
246 |
MoleculeIterator i; |
| 202 |
|
|
i = molecules_.find(mol->getGlobalIndex()); |
| 203 |
gezelter |
2 |
|
| 204 |
gezelter |
246 |
if (i != molecules_.end() ) { |
| 205 |
gezelter |
2 |
|
| 206 |
gezelter |
507 |
assert(mol == i->second); |
| 207 |
gezelter |
246 |
|
| 208 |
gezelter |
507 |
nAtoms_ -= mol->getNAtoms(); |
| 209 |
|
|
nBonds_ -= mol->getNBonds(); |
| 210 |
|
|
nBends_ -= mol->getNBends(); |
| 211 |
|
|
nTorsions_ -= mol->getNTorsions(); |
| 212 |
|
|
nRigidBodies_ -= mol->getNRigidBodies(); |
| 213 |
|
|
nIntegrableObjects_ -= mol->getNIntegrableObjects(); |
| 214 |
|
|
nCutoffGroups_ -= mol->getNCutoffGroups(); |
| 215 |
|
|
nConstraints_ -= mol->getNConstraintPairs(); |
| 216 |
gezelter |
2 |
|
| 217 |
gezelter |
507 |
removeExcludePairs(mol); |
| 218 |
|
|
molecules_.erase(mol->getGlobalIndex()); |
| 219 |
gezelter |
2 |
|
| 220 |
gezelter |
507 |
delete mol; |
| 221 |
gezelter |
246 |
|
| 222 |
gezelter |
507 |
return true; |
| 223 |
gezelter |
246 |
} else { |
| 224 |
gezelter |
507 |
return false; |
| 225 |
gezelter |
246 |
} |
| 226 |
|
|
|
| 227 |
|
|
|
| 228 |
gezelter |
507 |
} |
| 229 |
gezelter |
246 |
|
| 230 |
|
|
|
| 231 |
gezelter |
507 |
Molecule* SimInfo::beginMolecule(MoleculeIterator& i) { |
| 232 |
gezelter |
246 |
i = molecules_.begin(); |
| 233 |
|
|
return i == molecules_.end() ? NULL : i->second; |
| 234 |
gezelter |
507 |
} |
| 235 |
gezelter |
246 |
|
| 236 |
gezelter |
507 |
Molecule* SimInfo::nextMolecule(MoleculeIterator& i) { |
| 237 |
gezelter |
246 |
++i; |
| 238 |
|
|
return i == molecules_.end() ? NULL : i->second; |
| 239 |
gezelter |
507 |
} |
| 240 |
gezelter |
2 |
|
| 241 |
|
|
|
| 242 |
gezelter |
507 |
void SimInfo::calcNdf() { |
| 243 |
gezelter |
246 |
int ndf_local; |
| 244 |
|
|
MoleculeIterator i; |
| 245 |
|
|
std::vector<StuntDouble*>::iterator j; |
| 246 |
|
|
Molecule* mol; |
| 247 |
|
|
StuntDouble* integrableObject; |
| 248 |
gezelter |
2 |
|
| 249 |
gezelter |
246 |
ndf_local = 0; |
| 250 |
|
|
|
| 251 |
|
|
for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { |
| 252 |
gezelter |
507 |
for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL; |
| 253 |
|
|
integrableObject = mol->nextIntegrableObject(j)) { |
| 254 |
gezelter |
2 |
|
| 255 |
gezelter |
507 |
ndf_local += 3; |
| 256 |
gezelter |
2 |
|
| 257 |
gezelter |
507 |
if (integrableObject->isDirectional()) { |
| 258 |
|
|
if (integrableObject->isLinear()) { |
| 259 |
|
|
ndf_local += 2; |
| 260 |
|
|
} else { |
| 261 |
|
|
ndf_local += 3; |
| 262 |
|
|
} |
| 263 |
|
|
} |
| 264 |
gezelter |
246 |
|
| 265 |
gezelter |
507 |
}//end for (integrableObject) |
| 266 |
gezelter |
246 |
}// end for (mol) |
| 267 |
|
|
|
| 268 |
|
|
// n_constraints is local, so subtract them on each processor |
| 269 |
|
|
ndf_local -= nConstraints_; |
| 270 |
|
|
|
| 271 |
|
|
#ifdef IS_MPI |
| 272 |
|
|
MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD); |
| 273 |
|
|
#else |
| 274 |
|
|
ndf_ = ndf_local; |
| 275 |
|
|
#endif |
| 276 |
|
|
|
| 277 |
|
|
// nZconstraints_ is global, as are the 3 COM translations for the |
| 278 |
|
|
// entire system: |
| 279 |
|
|
ndf_ = ndf_ - 3 - nZconstraint_; |
| 280 |
|
|
|
| 281 |
gezelter |
507 |
} |
| 282 |
gezelter |
2 |
|
| 283 |
gezelter |
507 |
void SimInfo::calcNdfRaw() { |
| 284 |
gezelter |
246 |
int ndfRaw_local; |
| 285 |
gezelter |
2 |
|
| 286 |
gezelter |
246 |
MoleculeIterator i; |
| 287 |
|
|
std::vector<StuntDouble*>::iterator j; |
| 288 |
|
|
Molecule* mol; |
| 289 |
|
|
StuntDouble* integrableObject; |
| 290 |
|
|
|
| 291 |
|
|
// Raw degrees of freedom that we have to set |
| 292 |
|
|
ndfRaw_local = 0; |
| 293 |
|
|
|
| 294 |
|
|
for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { |
| 295 |
gezelter |
507 |
for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL; |
| 296 |
|
|
integrableObject = mol->nextIntegrableObject(j)) { |
| 297 |
gezelter |
246 |
|
| 298 |
gezelter |
507 |
ndfRaw_local += 3; |
| 299 |
gezelter |
246 |
|
| 300 |
gezelter |
507 |
if (integrableObject->isDirectional()) { |
| 301 |
|
|
if (integrableObject->isLinear()) { |
| 302 |
|
|
ndfRaw_local += 2; |
| 303 |
|
|
} else { |
| 304 |
|
|
ndfRaw_local += 3; |
| 305 |
|
|
} |
| 306 |
|
|
} |
| 307 |
gezelter |
246 |
|
| 308 |
gezelter |
507 |
} |
| 309 |
gezelter |
246 |
} |
| 310 |
|
|
|
| 311 |
|
|
#ifdef IS_MPI |
| 312 |
|
|
MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD); |
| 313 |
|
|
#else |
| 314 |
|
|
ndfRaw_ = ndfRaw_local; |
| 315 |
|
|
#endif |
| 316 |
gezelter |
507 |
} |
| 317 |
gezelter |
2 |
|
| 318 |
gezelter |
507 |
void SimInfo::calcNdfTrans() { |
| 319 |
gezelter |
246 |
int ndfTrans_local; |
| 320 |
gezelter |
2 |
|
| 321 |
gezelter |
246 |
ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_; |
| 322 |
gezelter |
2 |
|
| 323 |
|
|
|
| 324 |
gezelter |
246 |
#ifdef IS_MPI |
| 325 |
|
|
MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD); |
| 326 |
|
|
#else |
| 327 |
|
|
ndfTrans_ = ndfTrans_local; |
| 328 |
|
|
#endif |
| 329 |
gezelter |
2 |
|
| 330 |
gezelter |
246 |
ndfTrans_ = ndfTrans_ - 3 - nZconstraint_; |
| 331 |
|
|
|
| 332 |
gezelter |
507 |
} |
| 333 |
gezelter |
2 |
|
| 334 |
gezelter |
507 |
void SimInfo::addExcludePairs(Molecule* mol) { |
| 335 |
gezelter |
246 |
std::vector<Bond*>::iterator bondIter; |
| 336 |
|
|
std::vector<Bend*>::iterator bendIter; |
| 337 |
|
|
std::vector<Torsion*>::iterator torsionIter; |
| 338 |
|
|
Bond* bond; |
| 339 |
|
|
Bend* bend; |
| 340 |
|
|
Torsion* torsion; |
| 341 |
|
|
int a; |
| 342 |
|
|
int b; |
| 343 |
|
|
int c; |
| 344 |
|
|
int d; |
| 345 |
|
|
|
| 346 |
|
|
for (bond= mol->beginBond(bondIter); bond != NULL; bond = mol->nextBond(bondIter)) { |
| 347 |
gezelter |
507 |
a = bond->getAtomA()->getGlobalIndex(); |
| 348 |
|
|
b = bond->getAtomB()->getGlobalIndex(); |
| 349 |
|
|
exclude_.addPair(a, b); |
| 350 |
gezelter |
246 |
} |
| 351 |
gezelter |
2 |
|
| 352 |
gezelter |
246 |
for (bend= mol->beginBend(bendIter); bend != NULL; bend = mol->nextBend(bendIter)) { |
| 353 |
gezelter |
507 |
a = bend->getAtomA()->getGlobalIndex(); |
| 354 |
|
|
b = bend->getAtomB()->getGlobalIndex(); |
| 355 |
|
|
c = bend->getAtomC()->getGlobalIndex(); |
| 356 |
gezelter |
2 |
|
| 357 |
gezelter |
507 |
exclude_.addPair(a, b); |
| 358 |
|
|
exclude_.addPair(a, c); |
| 359 |
|
|
exclude_.addPair(b, c); |
| 360 |
gezelter |
246 |
} |
| 361 |
gezelter |
2 |
|
| 362 |
gezelter |
246 |
for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; torsion = mol->nextTorsion(torsionIter)) { |
| 363 |
gezelter |
507 |
a = torsion->getAtomA()->getGlobalIndex(); |
| 364 |
|
|
b = torsion->getAtomB()->getGlobalIndex(); |
| 365 |
|
|
c = torsion->getAtomC()->getGlobalIndex(); |
| 366 |
|
|
d = torsion->getAtomD()->getGlobalIndex(); |
| 367 |
gezelter |
2 |
|
| 368 |
gezelter |
507 |
exclude_.addPair(a, b); |
| 369 |
|
|
exclude_.addPair(a, c); |
| 370 |
|
|
exclude_.addPair(a, d); |
| 371 |
|
|
exclude_.addPair(b, c); |
| 372 |
|
|
exclude_.addPair(b, d); |
| 373 |
|
|
exclude_.addPair(c, d); |
| 374 |
gezelter |
2 |
} |
| 375 |
|
|
|
| 376 |
tim |
430 |
Molecule::RigidBodyIterator rbIter; |
| 377 |
|
|
RigidBody* rb; |
| 378 |
|
|
for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) { |
| 379 |
gezelter |
507 |
std::vector<Atom*> atoms = rb->getAtoms(); |
| 380 |
|
|
for (int i = 0; i < atoms.size() -1 ; ++i) { |
| 381 |
|
|
for (int j = i + 1; j < atoms.size(); ++j) { |
| 382 |
|
|
a = atoms[i]->getGlobalIndex(); |
| 383 |
|
|
b = atoms[j]->getGlobalIndex(); |
| 384 |
|
|
exclude_.addPair(a, b); |
| 385 |
|
|
} |
| 386 |
|
|
} |
| 387 |
tim |
430 |
} |
| 388 |
|
|
|
| 389 |
gezelter |
507 |
} |
| 390 |
gezelter |
246 |
|
| 391 |
gezelter |
507 |
void SimInfo::removeExcludePairs(Molecule* mol) { |
| 392 |
gezelter |
246 |
std::vector<Bond*>::iterator bondIter; |
| 393 |
|
|
std::vector<Bend*>::iterator bendIter; |
| 394 |
|
|
std::vector<Torsion*>::iterator torsionIter; |
| 395 |
|
|
Bond* bond; |
| 396 |
|
|
Bend* bend; |
| 397 |
|
|
Torsion* torsion; |
| 398 |
|
|
int a; |
| 399 |
|
|
int b; |
| 400 |
|
|
int c; |
| 401 |
|
|
int d; |
| 402 |
|
|
|
| 403 |
|
|
for (bond= mol->beginBond(bondIter); bond != NULL; bond = mol->nextBond(bondIter)) { |
| 404 |
gezelter |
507 |
a = bond->getAtomA()->getGlobalIndex(); |
| 405 |
|
|
b = bond->getAtomB()->getGlobalIndex(); |
| 406 |
|
|
exclude_.removePair(a, b); |
| 407 |
gezelter |
2 |
} |
| 408 |
gezelter |
246 |
|
| 409 |
|
|
for (bend= mol->beginBend(bendIter); bend != NULL; bend = mol->nextBend(bendIter)) { |
| 410 |
gezelter |
507 |
a = bend->getAtomA()->getGlobalIndex(); |
| 411 |
|
|
b = bend->getAtomB()->getGlobalIndex(); |
| 412 |
|
|
c = bend->getAtomC()->getGlobalIndex(); |
| 413 |
gezelter |
246 |
|
| 414 |
gezelter |
507 |
exclude_.removePair(a, b); |
| 415 |
|
|
exclude_.removePair(a, c); |
| 416 |
|
|
exclude_.removePair(b, c); |
| 417 |
gezelter |
2 |
} |
| 418 |
gezelter |
246 |
|
| 419 |
|
|
for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; torsion = mol->nextTorsion(torsionIter)) { |
| 420 |
gezelter |
507 |
a = torsion->getAtomA()->getGlobalIndex(); |
| 421 |
|
|
b = torsion->getAtomB()->getGlobalIndex(); |
| 422 |
|
|
c = torsion->getAtomC()->getGlobalIndex(); |
| 423 |
|
|
d = torsion->getAtomD()->getGlobalIndex(); |
| 424 |
gezelter |
246 |
|
| 425 |
gezelter |
507 |
exclude_.removePair(a, b); |
| 426 |
|
|
exclude_.removePair(a, c); |
| 427 |
|
|
exclude_.removePair(a, d); |
| 428 |
|
|
exclude_.removePair(b, c); |
| 429 |
|
|
exclude_.removePair(b, d); |
| 430 |
|
|
exclude_.removePair(c, d); |
| 431 |
gezelter |
246 |
} |
| 432 |
|
|
|
| 433 |
tim |
430 |
Molecule::RigidBodyIterator rbIter; |
| 434 |
|
|
RigidBody* rb; |
| 435 |
|
|
for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) { |
| 436 |
gezelter |
507 |
std::vector<Atom*> atoms = rb->getAtoms(); |
| 437 |
|
|
for (int i = 0; i < atoms.size() -1 ; ++i) { |
| 438 |
|
|
for (int j = i + 1; j < atoms.size(); ++j) { |
| 439 |
|
|
a = atoms[i]->getGlobalIndex(); |
| 440 |
|
|
b = atoms[j]->getGlobalIndex(); |
| 441 |
|
|
exclude_.removePair(a, b); |
| 442 |
|
|
} |
| 443 |
|
|
} |
| 444 |
tim |
430 |
} |
| 445 |
|
|
|
| 446 |
gezelter |
507 |
} |
| 447 |
gezelter |
2 |
|
| 448 |
|
|
|
| 449 |
gezelter |
507 |
void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) { |
| 450 |
gezelter |
246 |
int curStampId; |
| 451 |
gezelter |
2 |
|
| 452 |
gezelter |
246 |
//index from 0 |
| 453 |
|
|
curStampId = moleculeStamps_.size(); |
| 454 |
gezelter |
2 |
|
| 455 |
gezelter |
246 |
moleculeStamps_.push_back(molStamp); |
| 456 |
|
|
molStampIds_.insert(molStampIds_.end(), nmol, curStampId); |
| 457 |
gezelter |
507 |
} |
| 458 |
gezelter |
2 |
|
| 459 |
gezelter |
507 |
void SimInfo::update() { |
| 460 |
gezelter |
2 |
|
| 461 |
gezelter |
246 |
setupSimType(); |
| 462 |
gezelter |
2 |
|
| 463 |
gezelter |
246 |
#ifdef IS_MPI |
| 464 |
|
|
setupFortranParallel(); |
| 465 |
|
|
#endif |
| 466 |
gezelter |
2 |
|
| 467 |
gezelter |
246 |
setupFortranSim(); |
| 468 |
gezelter |
2 |
|
| 469 |
gezelter |
246 |
//setup fortran force field |
| 470 |
|
|
/** @deprecate */ |
| 471 |
|
|
int isError = 0; |
| 472 |
chrisfen |
598 |
|
| 473 |
chrisfen |
603 |
setupElectrostaticSummationMethod( isError ); |
| 474 |
chrisfen |
598 |
|
| 475 |
gezelter |
246 |
if(isError){ |
| 476 |
gezelter |
507 |
sprintf( painCave.errMsg, |
| 477 |
|
|
"ForceField error: There was an error initializing the forceField in fortran.\n" ); |
| 478 |
|
|
painCave.isFatal = 1; |
| 479 |
|
|
simError(); |
| 480 |
gezelter |
246 |
} |
| 481 |
gezelter |
2 |
|
| 482 |
gezelter |
246 |
|
| 483 |
|
|
setupCutoff(); |
| 484 |
gezelter |
2 |
|
| 485 |
gezelter |
246 |
calcNdf(); |
| 486 |
|
|
calcNdfRaw(); |
| 487 |
|
|
calcNdfTrans(); |
| 488 |
|
|
|
| 489 |
|
|
fortranInitialized_ = true; |
| 490 |
gezelter |
507 |
} |
| 491 |
gezelter |
2 |
|
| 492 |
gezelter |
507 |
std::set<AtomType*> SimInfo::getUniqueAtomTypes() { |
| 493 |
gezelter |
246 |
SimInfo::MoleculeIterator mi; |
| 494 |
|
|
Molecule* mol; |
| 495 |
|
|
Molecule::AtomIterator ai; |
| 496 |
|
|
Atom* atom; |
| 497 |
|
|
std::set<AtomType*> atomTypes; |
| 498 |
gezelter |
2 |
|
| 499 |
gezelter |
246 |
for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) { |
| 500 |
gezelter |
2 |
|
| 501 |
gezelter |
507 |
for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) { |
| 502 |
|
|
atomTypes.insert(atom->getAtomType()); |
| 503 |
|
|
} |
| 504 |
gezelter |
246 |
|
| 505 |
|
|
} |
| 506 |
gezelter |
2 |
|
| 507 |
gezelter |
246 |
return atomTypes; |
| 508 |
gezelter |
507 |
} |
| 509 |
gezelter |
2 |
|
| 510 |
gezelter |
507 |
void SimInfo::setupSimType() { |
| 511 |
gezelter |
246 |
std::set<AtomType*>::iterator i; |
| 512 |
|
|
std::set<AtomType*> atomTypes; |
| 513 |
|
|
atomTypes = getUniqueAtomTypes(); |
| 514 |
gezelter |
2 |
|
| 515 |
gezelter |
246 |
int useLennardJones = 0; |
| 516 |
|
|
int useElectrostatic = 0; |
| 517 |
|
|
int useEAM = 0; |
| 518 |
|
|
int useCharge = 0; |
| 519 |
|
|
int useDirectional = 0; |
| 520 |
|
|
int useDipole = 0; |
| 521 |
|
|
int useGayBerne = 0; |
| 522 |
|
|
int useSticky = 0; |
| 523 |
chrisfen |
523 |
int useStickyPower = 0; |
| 524 |
gezelter |
246 |
int useShape = 0; |
| 525 |
|
|
int useFLARB = 0; //it is not in AtomType yet |
| 526 |
|
|
int useDirectionalAtom = 0; |
| 527 |
|
|
int useElectrostatics = 0; |
| 528 |
|
|
//usePBC and useRF are from simParams |
| 529 |
tim |
665 |
int usePBC = simParams_->getUsePeriodicBoundaryConditions(); |
| 530 |
chrisfen |
611 |
int useRF; |
| 531 |
tim |
665 |
std::string myMethod; |
| 532 |
gezelter |
2 |
|
| 533 |
chrisfen |
611 |
// set the useRF logical |
| 534 |
tim |
665 |
useRF = 0; |
| 535 |
|
|
if (simParams_->haveElectrostaticSummationMethod()) { |
| 536 |
|
|
myMethod = simParams_->getElectrostaticSummationMethod(); |
| 537 |
|
|
if (myMethod == "REACTION_FIELD") |
| 538 |
|
|
useRF = 1; |
| 539 |
|
|
} |
| 540 |
chrisfen |
611 |
|
| 541 |
gezelter |
246 |
//loop over all of the atom types |
| 542 |
|
|
for (i = atomTypes.begin(); i != atomTypes.end(); ++i) { |
| 543 |
gezelter |
507 |
useLennardJones |= (*i)->isLennardJones(); |
| 544 |
|
|
useElectrostatic |= (*i)->isElectrostatic(); |
| 545 |
|
|
useEAM |= (*i)->isEAM(); |
| 546 |
|
|
useCharge |= (*i)->isCharge(); |
| 547 |
|
|
useDirectional |= (*i)->isDirectional(); |
| 548 |
|
|
useDipole |= (*i)->isDipole(); |
| 549 |
|
|
useGayBerne |= (*i)->isGayBerne(); |
| 550 |
|
|
useSticky |= (*i)->isSticky(); |
| 551 |
chrisfen |
523 |
useStickyPower |= (*i)->isStickyPower(); |
| 552 |
gezelter |
507 |
useShape |= (*i)->isShape(); |
| 553 |
gezelter |
246 |
} |
| 554 |
gezelter |
2 |
|
| 555 |
chrisfen |
523 |
if (useSticky || useStickyPower || useDipole || useGayBerne || useShape) { |
| 556 |
gezelter |
507 |
useDirectionalAtom = 1; |
| 557 |
gezelter |
246 |
} |
| 558 |
gezelter |
2 |
|
| 559 |
gezelter |
246 |
if (useCharge || useDipole) { |
| 560 |
gezelter |
507 |
useElectrostatics = 1; |
| 561 |
gezelter |
246 |
} |
| 562 |
gezelter |
2 |
|
| 563 |
gezelter |
246 |
#ifdef IS_MPI |
| 564 |
|
|
int temp; |
| 565 |
gezelter |
2 |
|
| 566 |
gezelter |
246 |
temp = usePBC; |
| 567 |
|
|
MPI_Allreduce(&temp, &usePBC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
| 568 |
gezelter |
2 |
|
| 569 |
gezelter |
246 |
temp = useDirectionalAtom; |
| 570 |
|
|
MPI_Allreduce(&temp, &useDirectionalAtom, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
| 571 |
gezelter |
2 |
|
| 572 |
gezelter |
246 |
temp = useLennardJones; |
| 573 |
|
|
MPI_Allreduce(&temp, &useLennardJones, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
| 574 |
gezelter |
2 |
|
| 575 |
gezelter |
246 |
temp = useElectrostatics; |
| 576 |
|
|
MPI_Allreduce(&temp, &useElectrostatics, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
| 577 |
gezelter |
2 |
|
| 578 |
gezelter |
246 |
temp = useCharge; |
| 579 |
|
|
MPI_Allreduce(&temp, &useCharge, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
| 580 |
gezelter |
2 |
|
| 581 |
gezelter |
246 |
temp = useDipole; |
| 582 |
|
|
MPI_Allreduce(&temp, &useDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
| 583 |
gezelter |
2 |
|
| 584 |
gezelter |
246 |
temp = useSticky; |
| 585 |
|
|
MPI_Allreduce(&temp, &useSticky, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
| 586 |
gezelter |
2 |
|
| 587 |
chrisfen |
523 |
temp = useStickyPower; |
| 588 |
|
|
MPI_Allreduce(&temp, &useStickyPower, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
| 589 |
|
|
|
| 590 |
gezelter |
246 |
temp = useGayBerne; |
| 591 |
|
|
MPI_Allreduce(&temp, &useGayBerne, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
| 592 |
gezelter |
2 |
|
| 593 |
gezelter |
246 |
temp = useEAM; |
| 594 |
|
|
MPI_Allreduce(&temp, &useEAM, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
| 595 |
gezelter |
2 |
|
| 596 |
gezelter |
246 |
temp = useShape; |
| 597 |
|
|
MPI_Allreduce(&temp, &useShape, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
| 598 |
|
|
|
| 599 |
|
|
temp = useFLARB; |
| 600 |
|
|
MPI_Allreduce(&temp, &useFLARB, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
| 601 |
|
|
|
| 602 |
chrisfen |
611 |
temp = useRF; |
| 603 |
|
|
MPI_Allreduce(&temp, &useRF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); |
| 604 |
|
|
|
| 605 |
gezelter |
2 |
#endif |
| 606 |
|
|
|
| 607 |
gezelter |
246 |
fInfo_.SIM_uses_PBC = usePBC; |
| 608 |
|
|
fInfo_.SIM_uses_DirectionalAtoms = useDirectionalAtom; |
| 609 |
|
|
fInfo_.SIM_uses_LennardJones = useLennardJones; |
| 610 |
|
|
fInfo_.SIM_uses_Electrostatics = useElectrostatics; |
| 611 |
|
|
fInfo_.SIM_uses_Charges = useCharge; |
| 612 |
|
|
fInfo_.SIM_uses_Dipoles = useDipole; |
| 613 |
|
|
fInfo_.SIM_uses_Sticky = useSticky; |
| 614 |
chrisfen |
523 |
fInfo_.SIM_uses_StickyPower = useStickyPower; |
| 615 |
gezelter |
246 |
fInfo_.SIM_uses_GayBerne = useGayBerne; |
| 616 |
|
|
fInfo_.SIM_uses_EAM = useEAM; |
| 617 |
|
|
fInfo_.SIM_uses_Shapes = useShape; |
| 618 |
|
|
fInfo_.SIM_uses_FLARB = useFLARB; |
| 619 |
chrisfen |
611 |
fInfo_.SIM_uses_RF = useRF; |
| 620 |
gezelter |
2 |
|
| 621 |
chrisfen |
611 |
if( fInfo_.SIM_uses_Dipoles && myMethod == "REACTION_FIELD") { |
| 622 |
gezelter |
2 |
|
| 623 |
gezelter |
507 |
if (simParams_->haveDielectric()) { |
| 624 |
|
|
fInfo_.dielect = simParams_->getDielectric(); |
| 625 |
|
|
} else { |
| 626 |
|
|
sprintf(painCave.errMsg, |
| 627 |
|
|
"SimSetup Error: No Dielectric constant was set.\n" |
| 628 |
|
|
"\tYou are trying to use Reaction Field without" |
| 629 |
|
|
"\tsetting a dielectric constant!\n"); |
| 630 |
|
|
painCave.isFatal = 1; |
| 631 |
|
|
simError(); |
| 632 |
|
|
} |
| 633 |
gezelter |
246 |
|
| 634 |
|
|
} else { |
| 635 |
gezelter |
507 |
fInfo_.dielect = 0.0; |
| 636 |
gezelter |
246 |
} |
| 637 |
|
|
|
| 638 |
gezelter |
507 |
} |
| 639 |
gezelter |
2 |
|
| 640 |
gezelter |
507 |
void SimInfo::setupFortranSim() { |
| 641 |
gezelter |
246 |
int isError; |
| 642 |
|
|
int nExclude; |
| 643 |
|
|
std::vector<int> fortranGlobalGroupMembership; |
| 644 |
|
|
|
| 645 |
|
|
nExclude = exclude_.getSize(); |
| 646 |
|
|
isError = 0; |
| 647 |
gezelter |
2 |
|
| 648 |
gezelter |
246 |
//globalGroupMembership_ is filled by SimCreator |
| 649 |
|
|
for (int i = 0; i < nGlobalAtoms_; i++) { |
| 650 |
gezelter |
507 |
fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1); |
| 651 |
gezelter |
246 |
} |
| 652 |
gezelter |
2 |
|
| 653 |
gezelter |
246 |
//calculate mass ratio of cutoff group |
| 654 |
|
|
std::vector<double> mfact; |
| 655 |
|
|
SimInfo::MoleculeIterator mi; |
| 656 |
|
|
Molecule* mol; |
| 657 |
|
|
Molecule::CutoffGroupIterator ci; |
| 658 |
|
|
CutoffGroup* cg; |
| 659 |
|
|
Molecule::AtomIterator ai; |
| 660 |
|
|
Atom* atom; |
| 661 |
|
|
double totalMass; |
| 662 |
|
|
|
| 663 |
|
|
//to avoid memory reallocation, reserve enough space for mfact |
| 664 |
|
|
mfact.reserve(getNCutoffGroups()); |
| 665 |
gezelter |
2 |
|
| 666 |
gezelter |
246 |
for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) { |
| 667 |
gezelter |
507 |
for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) { |
| 668 |
gezelter |
2 |
|
| 669 |
gezelter |
507 |
totalMass = cg->getMass(); |
| 670 |
|
|
for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) { |
| 671 |
chrisfen |
645 |
// Check for massless groups - set mfact to 1 if true |
| 672 |
|
|
if (totalMass != 0) |
| 673 |
|
|
mfact.push_back(atom->getMass()/totalMass); |
| 674 |
|
|
else |
| 675 |
|
|
mfact.push_back( 1.0 ); |
| 676 |
gezelter |
507 |
} |
| 677 |
gezelter |
2 |
|
| 678 |
gezelter |
507 |
} |
| 679 |
gezelter |
246 |
} |
| 680 |
gezelter |
2 |
|
| 681 |
gezelter |
246 |
//fill ident array of local atoms (it is actually ident of AtomType, it is so confusing !!!) |
| 682 |
|
|
std::vector<int> identArray; |
| 683 |
gezelter |
2 |
|
| 684 |
gezelter |
246 |
//to avoid memory reallocation, reserve enough space identArray |
| 685 |
|
|
identArray.reserve(getNAtoms()); |
| 686 |
|
|
|
| 687 |
|
|
for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) { |
| 688 |
gezelter |
507 |
for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) { |
| 689 |
|
|
identArray.push_back(atom->getIdent()); |
| 690 |
|
|
} |
| 691 |
gezelter |
246 |
} |
| 692 |
gezelter |
2 |
|
| 693 |
gezelter |
246 |
//fill molMembershipArray |
| 694 |
|
|
//molMembershipArray is filled by SimCreator |
| 695 |
|
|
std::vector<int> molMembershipArray(nGlobalAtoms_); |
| 696 |
|
|
for (int i = 0; i < nGlobalAtoms_; i++) { |
| 697 |
gezelter |
507 |
molMembershipArray[i] = globalMolMembership_[i] + 1; |
| 698 |
gezelter |
246 |
} |
| 699 |
|
|
|
| 700 |
|
|
//setup fortran simulation |
| 701 |
|
|
int nGlobalExcludes = 0; |
| 702 |
|
|
int* globalExcludes = NULL; |
| 703 |
|
|
int* excludeList = exclude_.getExcludeList(); |
| 704 |
|
|
setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray[0], &nExclude, excludeList , |
| 705 |
gezelter |
507 |
&nGlobalExcludes, globalExcludes, &molMembershipArray[0], |
| 706 |
|
|
&mfact[0], &nCutoffGroups_, &fortranGlobalGroupMembership[0], &isError); |
| 707 |
gezelter |
2 |
|
| 708 |
gezelter |
246 |
if( isError ){ |
| 709 |
gezelter |
2 |
|
| 710 |
gezelter |
507 |
sprintf( painCave.errMsg, |
| 711 |
|
|
"There was an error setting the simulation information in fortran.\n" ); |
| 712 |
|
|
painCave.isFatal = 1; |
| 713 |
|
|
painCave.severity = OOPSE_ERROR; |
| 714 |
|
|
simError(); |
| 715 |
gezelter |
246 |
} |
| 716 |
|
|
|
| 717 |
|
|
#ifdef IS_MPI |
| 718 |
|
|
sprintf( checkPointMsg, |
| 719 |
gezelter |
507 |
"succesfully sent the simulation information to fortran.\n"); |
| 720 |
gezelter |
246 |
MPIcheckPoint(); |
| 721 |
|
|
#endif // is_mpi |
| 722 |
gezelter |
507 |
} |
| 723 |
gezelter |
2 |
|
| 724 |
|
|
|
| 725 |
gezelter |
246 |
#ifdef IS_MPI |
| 726 |
gezelter |
507 |
void SimInfo::setupFortranParallel() { |
| 727 |
gezelter |
246 |
|
| 728 |
|
|
//SimInfo is responsible for creating localToGlobalAtomIndex and localToGlobalGroupIndex |
| 729 |
|
|
std::vector<int> localToGlobalAtomIndex(getNAtoms(), 0); |
| 730 |
|
|
std::vector<int> localToGlobalCutoffGroupIndex; |
| 731 |
|
|
SimInfo::MoleculeIterator mi; |
| 732 |
|
|
Molecule::AtomIterator ai; |
| 733 |
|
|
Molecule::CutoffGroupIterator ci; |
| 734 |
|
|
Molecule* mol; |
| 735 |
|
|
Atom* atom; |
| 736 |
|
|
CutoffGroup* cg; |
| 737 |
|
|
mpiSimData parallelData; |
| 738 |
|
|
int isError; |
| 739 |
gezelter |
2 |
|
| 740 |
gezelter |
246 |
for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) { |
| 741 |
gezelter |
2 |
|
| 742 |
gezelter |
507 |
//local index(index in DataStorge) of atom is important |
| 743 |
|
|
for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) { |
| 744 |
|
|
localToGlobalAtomIndex[atom->getLocalIndex()] = atom->getGlobalIndex() + 1; |
| 745 |
|
|
} |
| 746 |
gezelter |
2 |
|
| 747 |
gezelter |
507 |
//local index of cutoff group is trivial, it only depends on the order of travesing |
| 748 |
|
|
for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) { |
| 749 |
|
|
localToGlobalCutoffGroupIndex.push_back(cg->getGlobalIndex() + 1); |
| 750 |
|
|
} |
| 751 |
gezelter |
246 |
|
| 752 |
|
|
} |
| 753 |
gezelter |
2 |
|
| 754 |
gezelter |
246 |
//fill up mpiSimData struct |
| 755 |
|
|
parallelData.nMolGlobal = getNGlobalMolecules(); |
| 756 |
|
|
parallelData.nMolLocal = getNMolecules(); |
| 757 |
|
|
parallelData.nAtomsGlobal = getNGlobalAtoms(); |
| 758 |
|
|
parallelData.nAtomsLocal = getNAtoms(); |
| 759 |
|
|
parallelData.nGroupsGlobal = getNGlobalCutoffGroups(); |
| 760 |
|
|
parallelData.nGroupsLocal = getNCutoffGroups(); |
| 761 |
|
|
parallelData.myNode = worldRank; |
| 762 |
|
|
MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors)); |
| 763 |
gezelter |
2 |
|
| 764 |
gezelter |
246 |
//pass mpiSimData struct and index arrays to fortran |
| 765 |
|
|
setFsimParallel(¶llelData, &(parallelData.nAtomsLocal), |
| 766 |
|
|
&localToGlobalAtomIndex[0], &(parallelData.nGroupsLocal), |
| 767 |
|
|
&localToGlobalCutoffGroupIndex[0], &isError); |
| 768 |
gezelter |
2 |
|
| 769 |
gezelter |
246 |
if (isError) { |
| 770 |
gezelter |
507 |
sprintf(painCave.errMsg, |
| 771 |
|
|
"mpiRefresh errror: fortran didn't like something we gave it.\n"); |
| 772 |
|
|
painCave.isFatal = 1; |
| 773 |
|
|
simError(); |
| 774 |
gezelter |
246 |
} |
| 775 |
gezelter |
2 |
|
| 776 |
gezelter |
246 |
sprintf(checkPointMsg, " mpiRefresh successful.\n"); |
| 777 |
|
|
MPIcheckPoint(); |
| 778 |
gezelter |
2 |
|
| 779 |
|
|
|
| 780 |
gezelter |
507 |
} |
| 781 |
chrisfen |
143 |
|
| 782 |
gezelter |
246 |
#endif |
| 783 |
chrisfen |
143 |
|
| 784 |
gezelter |
507 |
double SimInfo::calcMaxCutoffRadius() { |
| 785 |
chrisfen |
143 |
|
| 786 |
|
|
|
| 787 |
gezelter |
246 |
std::set<AtomType*> atomTypes; |
| 788 |
|
|
std::set<AtomType*>::iterator i; |
| 789 |
|
|
std::vector<double> cutoffRadius; |
| 790 |
gezelter |
2 |
|
| 791 |
gezelter |
246 |
//get the unique atom types |
| 792 |
|
|
atomTypes = getUniqueAtomTypes(); |
| 793 |
|
|
|
| 794 |
|
|
//query the max cutoff radius among these atom types |
| 795 |
|
|
for (i = atomTypes.begin(); i != atomTypes.end(); ++i) { |
| 796 |
gezelter |
507 |
cutoffRadius.push_back(forceField_->getRcutFromAtomType(*i)); |
| 797 |
gezelter |
246 |
} |
| 798 |
|
|
|
| 799 |
|
|
double maxCutoffRadius = *(std::max_element(cutoffRadius.begin(), cutoffRadius.end())); |
| 800 |
gezelter |
2 |
#ifdef IS_MPI |
| 801 |
gezelter |
246 |
//pick the max cutoff radius among the processors |
| 802 |
gezelter |
2 |
#endif |
| 803 |
|
|
|
| 804 |
gezelter |
246 |
return maxCutoffRadius; |
| 805 |
gezelter |
507 |
} |
| 806 |
gezelter |
246 |
|
| 807 |
gezelter |
507 |
void SimInfo::getCutoff(double& rcut, double& rsw) { |
| 808 |
gezelter |
2 |
|
| 809 |
gezelter |
246 |
if (fInfo_.SIM_uses_Charges | fInfo_.SIM_uses_Dipoles | fInfo_.SIM_uses_RF) { |
| 810 |
|
|
|
| 811 |
tim |
665 |
if (!simParams_->haveCutoffRadius()){ |
| 812 |
gezelter |
507 |
sprintf(painCave.errMsg, |
| 813 |
gezelter |
246 |
"SimCreator Warning: No value was set for the cutoffRadius.\n" |
| 814 |
|
|
"\tOOPSE will use a default value of 15.0 angstroms" |
| 815 |
|
|
"\tfor the cutoffRadius.\n"); |
| 816 |
gezelter |
507 |
painCave.isFatal = 0; |
| 817 |
|
|
simError(); |
| 818 |
|
|
rcut = 15.0; |
| 819 |
|
|
} else{ |
| 820 |
tim |
665 |
rcut = simParams_->getCutoffRadius(); |
| 821 |
gezelter |
507 |
} |
| 822 |
gezelter |
246 |
|
| 823 |
tim |
665 |
if (!simParams_->haveSwitchingRadius()){ |
| 824 |
gezelter |
507 |
sprintf(painCave.errMsg, |
| 825 |
gezelter |
246 |
"SimCreator Warning: No value was set for switchingRadius.\n" |
| 826 |
|
|
"\tOOPSE will use a default value of\n" |
| 827 |
|
|
"\t0.95 * cutoffRadius for the switchingRadius\n"); |
| 828 |
gezelter |
507 |
painCave.isFatal = 0; |
| 829 |
|
|
simError(); |
| 830 |
|
|
rsw = 0.95 * rcut; |
| 831 |
|
|
} else{ |
| 832 |
tim |
665 |
rsw = simParams_->getSwitchingRadius(); |
| 833 |
gezelter |
507 |
} |
| 834 |
gezelter |
246 |
|
| 835 |
|
|
} else { |
| 836 |
gezelter |
507 |
// if charge, dipole or reaction field is not used and the cutofff radius is not specified in |
| 837 |
|
|
//meta-data file, the maximum cutoff radius calculated from forcefiled will be used |
| 838 |
gezelter |
246 |
|
| 839 |
tim |
665 |
if (simParams_->haveCutoffRadius()) { |
| 840 |
|
|
rcut = simParams_->getCutoffRadius(); |
| 841 |
gezelter |
507 |
} else { |
| 842 |
|
|
//set cutoff radius to the maximum cutoff radius based on atom types in the whole system |
| 843 |
|
|
rcut = calcMaxCutoffRadius(); |
| 844 |
|
|
} |
| 845 |
gezelter |
246 |
|
| 846 |
tim |
665 |
if (simParams_->haveSwitchingRadius()) { |
| 847 |
|
|
rsw = simParams_->getSwitchingRadius(); |
| 848 |
gezelter |
507 |
} else { |
| 849 |
|
|
rsw = rcut; |
| 850 |
|
|
} |
| 851 |
gezelter |
246 |
|
| 852 |
|
|
} |
| 853 |
gezelter |
507 |
} |
| 854 |
tim |
326 |
|
| 855 |
gezelter |
586 |
void SimInfo::setupCutoff() { |
| 856 |
tim |
326 |
getCutoff(rcut_, rsw_); |
| 857 |
gezelter |
246 |
double rnblist = rcut_ + 1; // skin of neighbor list |
| 858 |
|
|
|
| 859 |
|
|
//Pass these cutoff radius etc. to fortran. This function should be called once and only once |
| 860 |
gezelter |
586 |
|
| 861 |
|
|
int cp = TRADITIONAL_CUTOFF_POLICY; |
| 862 |
|
|
if (simParams_->haveCutoffPolicy()) { |
| 863 |
|
|
std::string myPolicy = simParams_->getCutoffPolicy(); |
| 864 |
tim |
665 |
toUpper(myPolicy); |
| 865 |
gezelter |
586 |
if (myPolicy == "MIX") { |
| 866 |
|
|
cp = MIX_CUTOFF_POLICY; |
| 867 |
|
|
} else { |
| 868 |
|
|
if (myPolicy == "MAX") { |
| 869 |
|
|
cp = MAX_CUTOFF_POLICY; |
| 870 |
|
|
} else { |
| 871 |
|
|
if (myPolicy == "TRADITIONAL") { |
| 872 |
|
|
cp = TRADITIONAL_CUTOFF_POLICY; |
| 873 |
|
|
} else { |
| 874 |
|
|
// throw error |
| 875 |
|
|
sprintf( painCave.errMsg, |
| 876 |
|
|
"SimInfo error: Unknown cutoffPolicy. (Input file specified %s .)\n\tcutoffPolicy must be one of: \"Mix\", \"Max\", or \"Traditional\".", myPolicy.c_str() ); |
| 877 |
|
|
painCave.isFatal = 1; |
| 878 |
|
|
simError(); |
| 879 |
|
|
} |
| 880 |
|
|
} |
| 881 |
|
|
} |
| 882 |
|
|
} |
| 883 |
chuckv |
629 |
|
| 884 |
|
|
|
| 885 |
|
|
if (simParams_->haveSkinThickness()) { |
| 886 |
|
|
double skinThickness = simParams_->getSkinThickness(); |
| 887 |
|
|
} |
| 888 |
|
|
|
| 889 |
gezelter |
586 |
notifyFortranCutoffs(&rcut_, &rsw_, &rnblist, &cp); |
| 890 |
chrisfen |
610 |
// also send cutoff notification to electrostatics |
| 891 |
|
|
setElectrostaticCutoffRadius(&rcut_); |
| 892 |
gezelter |
507 |
} |
| 893 |
gezelter |
2 |
|
| 894 |
chrisfen |
603 |
void SimInfo::setupElectrostaticSummationMethod( int isError ) { |
| 895 |
chrisfen |
598 |
|
| 896 |
|
|
int errorOut; |
| 897 |
chrisfen |
603 |
int esm = NONE; |
| 898 |
chrisfen |
598 |
double alphaVal; |
| 899 |
chrisfen |
610 |
double dielectric; |
| 900 |
chrisfen |
598 |
|
| 901 |
|
|
errorOut = isError; |
| 902 |
chrisfen |
610 |
alphaVal = simParams_->getDampingAlpha(); |
| 903 |
|
|
dielectric = simParams_->getDielectric(); |
| 904 |
chrisfen |
598 |
|
| 905 |
chrisfen |
603 |
if (simParams_->haveElectrostaticSummationMethod()) { |
| 906 |
chrisfen |
604 |
std::string myMethod = simParams_->getElectrostaticSummationMethod(); |
| 907 |
tim |
665 |
toUpper(myMethod); |
| 908 |
chrisfen |
603 |
if (myMethod == "NONE") { |
| 909 |
|
|
esm = NONE; |
| 910 |
chrisfen |
598 |
} else { |
| 911 |
chrisfen |
603 |
if (myMethod == "UNDAMPED_WOLF") { |
| 912 |
|
|
esm = UNDAMPED_WOLF; |
| 913 |
chrisfen |
598 |
} else { |
| 914 |
chrisfen |
603 |
if (myMethod == "DAMPED_WOLF") { |
| 915 |
chrisfen |
604 |
esm = DAMPED_WOLF; |
| 916 |
chrisfen |
598 |
if (!simParams_->haveDampingAlpha()) { |
| 917 |
|
|
//throw error |
| 918 |
|
|
sprintf( painCave.errMsg, |
| 919 |
chrisfen |
610 |
"SimInfo warning: dampingAlpha was not specified in the input file. A default value of %f (1/ang) will be used for the Damped Wolf Method.", alphaVal); |
| 920 |
chrisfen |
598 |
painCave.isFatal = 0; |
| 921 |
|
|
simError(); |
| 922 |
|
|
} |
| 923 |
|
|
} else { |
| 924 |
chrisfen |
603 |
if (myMethod == "REACTION_FIELD") { |
| 925 |
|
|
esm = REACTION_FIELD; |
| 926 |
chrisfen |
598 |
} else { |
| 927 |
|
|
// throw error |
| 928 |
|
|
sprintf( painCave.errMsg, |
| 929 |
chrisfen |
603 |
"SimInfo error: Unknown electrostaticSummationMethod. (Input file specified %s .)\n\telectrostaticSummationMethod must be one of: \"none\", \"undamped_wolf\", \"damped_wolf\", or \"reaction_field\".", myMethod.c_str() ); |
| 930 |
chrisfen |
598 |
painCave.isFatal = 1; |
| 931 |
|
|
simError(); |
| 932 |
|
|
} |
| 933 |
|
|
} |
| 934 |
|
|
} |
| 935 |
|
|
} |
| 936 |
|
|
} |
| 937 |
chrisfen |
610 |
// let's pass some summation method variables to fortran |
| 938 |
|
|
setElectrostaticSummationMethod( &esm ); |
| 939 |
|
|
setDampedWolfAlpha( &alphaVal ); |
| 940 |
|
|
setReactionFieldDielectric( &dielectric ); |
| 941 |
|
|
initFortranFF( &esm, &errorOut ); |
| 942 |
chrisfen |
598 |
} |
| 943 |
|
|
|
| 944 |
gezelter |
507 |
void SimInfo::addProperty(GenericData* genData) { |
| 945 |
gezelter |
246 |
properties_.addProperty(genData); |
| 946 |
gezelter |
507 |
} |
| 947 |
gezelter |
2 |
|
| 948 |
gezelter |
507 |
void SimInfo::removeProperty(const std::string& propName) { |
| 949 |
gezelter |
246 |
properties_.removeProperty(propName); |
| 950 |
gezelter |
507 |
} |
| 951 |
gezelter |
2 |
|
| 952 |
gezelter |
507 |
void SimInfo::clearProperties() { |
| 953 |
gezelter |
246 |
properties_.clearProperties(); |
| 954 |
gezelter |
507 |
} |
| 955 |
gezelter |
2 |
|
| 956 |
gezelter |
507 |
std::vector<std::string> SimInfo::getPropertyNames() { |
| 957 |
gezelter |
246 |
return properties_.getPropertyNames(); |
| 958 |
gezelter |
507 |
} |
| 959 |
gezelter |
246 |
|
| 960 |
gezelter |
507 |
std::vector<GenericData*> SimInfo::getProperties() { |
| 961 |
gezelter |
246 |
return properties_.getProperties(); |
| 962 |
gezelter |
507 |
} |
| 963 |
gezelter |
2 |
|
| 964 |
gezelter |
507 |
GenericData* SimInfo::getPropertyByName(const std::string& propName) { |
| 965 |
gezelter |
246 |
return properties_.getPropertyByName(propName); |
| 966 |
gezelter |
507 |
} |
| 967 |
gezelter |
2 |
|
| 968 |
gezelter |
507 |
void SimInfo::setSnapshotManager(SnapshotManager* sman) { |
| 969 |
tim |
432 |
if (sman_ == sman) { |
| 970 |
gezelter |
507 |
return; |
| 971 |
tim |
432 |
} |
| 972 |
|
|
delete sman_; |
| 973 |
gezelter |
246 |
sman_ = sman; |
| 974 |
gezelter |
2 |
|
| 975 |
gezelter |
246 |
Molecule* mol; |
| 976 |
|
|
RigidBody* rb; |
| 977 |
|
|
Atom* atom; |
| 978 |
|
|
SimInfo::MoleculeIterator mi; |
| 979 |
|
|
Molecule::RigidBodyIterator rbIter; |
| 980 |
|
|
Molecule::AtomIterator atomIter;; |
| 981 |
|
|
|
| 982 |
|
|
for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) { |
| 983 |
|
|
|
| 984 |
gezelter |
507 |
for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) { |
| 985 |
|
|
atom->setSnapshotManager(sman_); |
| 986 |
|
|
} |
| 987 |
gezelter |
246 |
|
| 988 |
gezelter |
507 |
for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) { |
| 989 |
|
|
rb->setSnapshotManager(sman_); |
| 990 |
|
|
} |
| 991 |
gezelter |
246 |
} |
| 992 |
gezelter |
2 |
|
| 993 |
gezelter |
507 |
} |
| 994 |
gezelter |
2 |
|
| 995 |
gezelter |
507 |
Vector3d SimInfo::getComVel(){ |
| 996 |
gezelter |
246 |
SimInfo::MoleculeIterator i; |
| 997 |
|
|
Molecule* mol; |
| 998 |
gezelter |
2 |
|
| 999 |
gezelter |
246 |
Vector3d comVel(0.0); |
| 1000 |
|
|
double totalMass = 0.0; |
| 1001 |
gezelter |
2 |
|
| 1002 |
gezelter |
246 |
|
| 1003 |
|
|
for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { |
| 1004 |
gezelter |
507 |
double mass = mol->getMass(); |
| 1005 |
|
|
totalMass += mass; |
| 1006 |
|
|
comVel += mass * mol->getComVel(); |
| 1007 |
gezelter |
246 |
} |
| 1008 |
gezelter |
2 |
|
| 1009 |
gezelter |
246 |
#ifdef IS_MPI |
| 1010 |
|
|
double tmpMass = totalMass; |
| 1011 |
|
|
Vector3d tmpComVel(comVel); |
| 1012 |
|
|
MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
| 1013 |
|
|
MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
| 1014 |
|
|
#endif |
| 1015 |
|
|
|
| 1016 |
|
|
comVel /= totalMass; |
| 1017 |
|
|
|
| 1018 |
|
|
return comVel; |
| 1019 |
gezelter |
507 |
} |
| 1020 |
gezelter |
2 |
|
| 1021 |
gezelter |
507 |
Vector3d SimInfo::getCom(){ |
| 1022 |
gezelter |
246 |
SimInfo::MoleculeIterator i; |
| 1023 |
|
|
Molecule* mol; |
| 1024 |
gezelter |
2 |
|
| 1025 |
gezelter |
246 |
Vector3d com(0.0); |
| 1026 |
|
|
double totalMass = 0.0; |
| 1027 |
|
|
|
| 1028 |
|
|
for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { |
| 1029 |
gezelter |
507 |
double mass = mol->getMass(); |
| 1030 |
|
|
totalMass += mass; |
| 1031 |
|
|
com += mass * mol->getCom(); |
| 1032 |
gezelter |
246 |
} |
| 1033 |
gezelter |
2 |
|
| 1034 |
|
|
#ifdef IS_MPI |
| 1035 |
gezelter |
246 |
double tmpMass = totalMass; |
| 1036 |
|
|
Vector3d tmpCom(com); |
| 1037 |
|
|
MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
| 1038 |
|
|
MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
| 1039 |
gezelter |
2 |
#endif |
| 1040 |
|
|
|
| 1041 |
gezelter |
246 |
com /= totalMass; |
| 1042 |
gezelter |
2 |
|
| 1043 |
gezelter |
246 |
return com; |
| 1044 |
gezelter |
2 |
|
| 1045 |
gezelter |
507 |
} |
| 1046 |
gezelter |
246 |
|
| 1047 |
gezelter |
507 |
std::ostream& operator <<(std::ostream& o, SimInfo& info) { |
| 1048 |
gezelter |
246 |
|
| 1049 |
|
|
return o; |
| 1050 |
gezelter |
507 |
} |
| 1051 |
chuckv |
555 |
|
| 1052 |
|
|
|
| 1053 |
|
|
/* |
| 1054 |
|
|
Returns center of mass and center of mass velocity in one function call. |
| 1055 |
|
|
*/ |
| 1056 |
|
|
|
| 1057 |
|
|
void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){ |
| 1058 |
|
|
SimInfo::MoleculeIterator i; |
| 1059 |
|
|
Molecule* mol; |
| 1060 |
|
|
|
| 1061 |
|
|
|
| 1062 |
|
|
double totalMass = 0.0; |
| 1063 |
|
|
|
| 1064 |
gezelter |
246 |
|
| 1065 |
chuckv |
555 |
for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { |
| 1066 |
|
|
double mass = mol->getMass(); |
| 1067 |
|
|
totalMass += mass; |
| 1068 |
|
|
com += mass * mol->getCom(); |
| 1069 |
|
|
comVel += mass * mol->getComVel(); |
| 1070 |
|
|
} |
| 1071 |
|
|
|
| 1072 |
|
|
#ifdef IS_MPI |
| 1073 |
|
|
double tmpMass = totalMass; |
| 1074 |
|
|
Vector3d tmpCom(com); |
| 1075 |
|
|
Vector3d tmpComVel(comVel); |
| 1076 |
|
|
MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
| 1077 |
|
|
MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
| 1078 |
|
|
MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
| 1079 |
|
|
#endif |
| 1080 |
|
|
|
| 1081 |
|
|
com /= totalMass; |
| 1082 |
|
|
comVel /= totalMass; |
| 1083 |
|
|
} |
| 1084 |
|
|
|
| 1085 |
|
|
/* |
| 1086 |
|
|
Return intertia tensor for entire system and angular momentum Vector. |
| 1087 |
chuckv |
557 |
|
| 1088 |
|
|
|
| 1089 |
|
|
[ Ixx -Ixy -Ixz ] |
| 1090 |
|
|
J =| -Iyx Iyy -Iyz | |
| 1091 |
|
|
[ -Izx -Iyz Izz ] |
| 1092 |
chuckv |
555 |
*/ |
| 1093 |
|
|
|
| 1094 |
|
|
void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){ |
| 1095 |
|
|
|
| 1096 |
|
|
|
| 1097 |
|
|
double xx = 0.0; |
| 1098 |
|
|
double yy = 0.0; |
| 1099 |
|
|
double zz = 0.0; |
| 1100 |
|
|
double xy = 0.0; |
| 1101 |
|
|
double xz = 0.0; |
| 1102 |
|
|
double yz = 0.0; |
| 1103 |
|
|
Vector3d com(0.0); |
| 1104 |
|
|
Vector3d comVel(0.0); |
| 1105 |
|
|
|
| 1106 |
|
|
getComAll(com, comVel); |
| 1107 |
|
|
|
| 1108 |
|
|
SimInfo::MoleculeIterator i; |
| 1109 |
|
|
Molecule* mol; |
| 1110 |
|
|
|
| 1111 |
|
|
Vector3d thisq(0.0); |
| 1112 |
|
|
Vector3d thisv(0.0); |
| 1113 |
|
|
|
| 1114 |
|
|
double thisMass = 0.0; |
| 1115 |
|
|
|
| 1116 |
|
|
|
| 1117 |
|
|
|
| 1118 |
|
|
|
| 1119 |
|
|
for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { |
| 1120 |
|
|
|
| 1121 |
|
|
thisq = mol->getCom()-com; |
| 1122 |
|
|
thisv = mol->getComVel()-comVel; |
| 1123 |
|
|
thisMass = mol->getMass(); |
| 1124 |
|
|
// Compute moment of intertia coefficients. |
| 1125 |
|
|
xx += thisq[0]*thisq[0]*thisMass; |
| 1126 |
|
|
yy += thisq[1]*thisq[1]*thisMass; |
| 1127 |
|
|
zz += thisq[2]*thisq[2]*thisMass; |
| 1128 |
|
|
|
| 1129 |
|
|
// compute products of intertia |
| 1130 |
|
|
xy += thisq[0]*thisq[1]*thisMass; |
| 1131 |
|
|
xz += thisq[0]*thisq[2]*thisMass; |
| 1132 |
|
|
yz += thisq[1]*thisq[2]*thisMass; |
| 1133 |
|
|
|
| 1134 |
|
|
angularMomentum += cross( thisq, thisv ) * thisMass; |
| 1135 |
|
|
|
| 1136 |
|
|
} |
| 1137 |
|
|
|
| 1138 |
|
|
|
| 1139 |
|
|
inertiaTensor(0,0) = yy + zz; |
| 1140 |
|
|
inertiaTensor(0,1) = -xy; |
| 1141 |
|
|
inertiaTensor(0,2) = -xz; |
| 1142 |
|
|
inertiaTensor(1,0) = -xy; |
| 1143 |
chuckv |
557 |
inertiaTensor(1,1) = xx + zz; |
| 1144 |
chuckv |
555 |
inertiaTensor(1,2) = -yz; |
| 1145 |
|
|
inertiaTensor(2,0) = -xz; |
| 1146 |
|
|
inertiaTensor(2,1) = -yz; |
| 1147 |
|
|
inertiaTensor(2,2) = xx + yy; |
| 1148 |
|
|
|
| 1149 |
|
|
#ifdef IS_MPI |
| 1150 |
|
|
Mat3x3d tmpI(inertiaTensor); |
| 1151 |
|
|
Vector3d tmpAngMom; |
| 1152 |
|
|
MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
| 1153 |
|
|
MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
| 1154 |
|
|
#endif |
| 1155 |
|
|
|
| 1156 |
|
|
return; |
| 1157 |
|
|
} |
| 1158 |
|
|
|
| 1159 |
|
|
//Returns the angular momentum of the system |
| 1160 |
|
|
Vector3d SimInfo::getAngularMomentum(){ |
| 1161 |
|
|
|
| 1162 |
|
|
Vector3d com(0.0); |
| 1163 |
|
|
Vector3d comVel(0.0); |
| 1164 |
|
|
Vector3d angularMomentum(0.0); |
| 1165 |
|
|
|
| 1166 |
|
|
getComAll(com,comVel); |
| 1167 |
|
|
|
| 1168 |
|
|
SimInfo::MoleculeIterator i; |
| 1169 |
|
|
Molecule* mol; |
| 1170 |
|
|
|
| 1171 |
chuckv |
557 |
Vector3d thisr(0.0); |
| 1172 |
|
|
Vector3d thisp(0.0); |
| 1173 |
chuckv |
555 |
|
| 1174 |
chuckv |
557 |
double thisMass; |
| 1175 |
chuckv |
555 |
|
| 1176 |
|
|
for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) { |
| 1177 |
chuckv |
557 |
thisMass = mol->getMass(); |
| 1178 |
|
|
thisr = mol->getCom()-com; |
| 1179 |
|
|
thisp = (mol->getComVel()-comVel)*thisMass; |
| 1180 |
chuckv |
555 |
|
| 1181 |
chuckv |
557 |
angularMomentum += cross( thisr, thisp ); |
| 1182 |
|
|
|
| 1183 |
chuckv |
555 |
} |
| 1184 |
|
|
|
| 1185 |
|
|
#ifdef IS_MPI |
| 1186 |
|
|
Vector3d tmpAngMom; |
| 1187 |
|
|
MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
| 1188 |
|
|
#endif |
| 1189 |
|
|
|
| 1190 |
|
|
return angularMomentum; |
| 1191 |
|
|
} |
| 1192 |
|
|
|
| 1193 |
|
|
|
| 1194 |
gezelter |
246 |
}//end namespace oopse |
| 1195 |
|
|
|