ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1665
Committed: Tue Nov 22 20:38:56 2011 UTC (13 years, 5 months ago) by gezelter
File size: 38848 byte(s)
Log Message:
updated copyright notices

File Contents

# User Rev Content
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 gezelter 1390 * 1. Redistributions of source code must retain the above copyright
10 gezelter 246 * notice, this list of conditions and the following disclaimer.
11     *
12 gezelter 1390 * 2. Redistributions in binary form must reproduce the above copyright
13 gezelter 246 * 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 gezelter 1390 *
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 gezelter 1665 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40     * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 gezelter 246 */
42    
43     /**
44     * @file SimInfo.cpp
45     * @author tlin
46     * @date 11/02/2004
47     * @version 1.0
48     */
49 gezelter 2
50 gezelter 246 #include <algorithm>
51     #include <set>
52 tim 749 #include <map>
53 gezelter 2
54 tim 3 #include "brains/SimInfo.hpp"
55 gezelter 246 #include "math/Vector3.hpp"
56     #include "primitives/Molecule.hpp"
57 tim 1024 #include "primitives/StuntDouble.hpp"
58 gezelter 246 #include "utils/MemoryUtils.hpp"
59 tim 3 #include "utils/simError.h"
60 tim 316 #include "selection/SelectionManager.hpp"
61 chuckv 834 #include "io/ForceFieldOptions.hpp"
62     #include "UseTheForce/ForceField.hpp"
63 gezelter 1530 #include "nonbonded/SwitchingFunction.hpp"
64 gezelter 1627 #ifdef IS_MPI
65     #include <mpi.h>
66     #endif
67 gezelter 2
68 gezelter 1528 using namespace std;
69 gezelter 1390 namespace OpenMD {
70 tim 749
71 tim 770 SimInfo::SimInfo(ForceField* ff, Globals* simParams) :
72     forceField_(ff), simParams_(simParams),
73 gezelter 945 ndf_(0), fdf_local(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0),
74 gezelter 507 nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0),
75     nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
76 gezelter 1277 nAtoms_(0), nBonds_(0), nBends_(0), nTorsions_(0), nInversions_(0),
77     nRigidBodies_(0), nIntegrableObjects_(0), nCutoffGroups_(0),
78 gezelter 1569 nConstraints_(0), sman_(NULL), topologyDone_(false),
79 gezelter 1528 calcBoxDipole_(false), useAtomicVirial_(true) {
80    
81     MoleculeStamp* molStamp;
82     int nMolWithSameStamp;
83     int nCutoffAtoms = 0; // number of atoms belong to cutoff groups
84     int nGroups = 0; //total cutoff groups defined in meta-data file
85     CutoffGroupStamp* cgStamp;
86     RigidBodyStamp* rbStamp;
87     int nRigidAtoms = 0;
88    
89     vector<Component*> components = simParams->getComponents();
90    
91     for (vector<Component*>::iterator i = components.begin(); i !=components.end(); ++i) {
92     molStamp = (*i)->getMoleculeStamp();
93     nMolWithSameStamp = (*i)->getNMol();
94 tim 770
95 gezelter 1528 addMoleculeStamp(molStamp, nMolWithSameStamp);
96    
97     //calculate atoms in molecules
98     nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp;
99    
100     //calculate atoms in cutoff groups
101     int nAtomsInGroups = 0;
102     int nCutoffGroupsInStamp = molStamp->getNCutoffGroups();
103    
104     for (int j=0; j < nCutoffGroupsInStamp; j++) {
105     cgStamp = molStamp->getCutoffGroupStamp(j);
106     nAtomsInGroups += cgStamp->getNMembers();
107 gezelter 507 }
108 gezelter 1528
109     nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
110    
111     nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;
112    
113     //calculate atoms in rigid bodies
114     int nAtomsInRigidBodies = 0;
115     int nRigidBodiesInStamp = molStamp->getNRigidBodies();
116    
117     for (int j=0; j < nRigidBodiesInStamp; j++) {
118     rbStamp = molStamp->getRigidBodyStamp(j);
119     nAtomsInRigidBodies += rbStamp->getNMembers();
120     }
121    
122     nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
123     nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;
124    
125     }
126    
127     //every free atom (atom does not belong to cutoff groups) is a cutoff
128     //group therefore the total number of cutoff groups in the system is
129     //equal to the total number of atoms minus number of atoms belong to
130     //cutoff group defined in meta-data file plus the number of cutoff
131     //groups defined in meta-data file
132 gezelter 1540
133 gezelter 1528 nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
134    
135     //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     nGlobalMols_ = molStampIds_.size();
144     molToProcMap_.resize(nGlobalMols_);
145     }
146 chrisfen 645
147 gezelter 507 SimInfo::~SimInfo() {
148 gezelter 1528 map<int, Molecule*>::iterator i;
149 tim 398 for (i = molecules_.begin(); i != molecules_.end(); ++i) {
150 gezelter 507 delete i->second;
151 tim 398 }
152     molecules_.clear();
153 tim 490
154 gezelter 246 delete sman_;
155     delete simParams_;
156     delete forceField_;
157 gezelter 507 }
158 gezelter 2
159    
160 gezelter 507 bool SimInfo::addMolecule(Molecule* mol) {
161 gezelter 246 MoleculeIterator i;
162 gezelter 1528
163 gezelter 246 i = molecules_.find(mol->getGlobalIndex());
164     if (i == molecules_.end() ) {
165 gezelter 1528
166     molecules_.insert(make_pair(mol->getGlobalIndex(), mol));
167    
168 gezelter 507 nAtoms_ += mol->getNAtoms();
169     nBonds_ += mol->getNBonds();
170     nBends_ += mol->getNBends();
171     nTorsions_ += mol->getNTorsions();
172 gezelter 1277 nInversions_ += mol->getNInversions();
173 gezelter 507 nRigidBodies_ += mol->getNRigidBodies();
174     nIntegrableObjects_ += mol->getNIntegrableObjects();
175     nCutoffGroups_ += mol->getNCutoffGroups();
176     nConstraints_ += mol->getNConstraintPairs();
177 gezelter 1528
178 gezelter 1287 addInteractionPairs(mol);
179 gezelter 1528
180 gezelter 507 return true;
181 gezelter 246 } else {
182 gezelter 507 return false;
183 gezelter 246 }
184 gezelter 507 }
185 gezelter 1528
186 gezelter 507 bool SimInfo::removeMolecule(Molecule* mol) {
187 gezelter 246 MoleculeIterator i;
188     i = molecules_.find(mol->getGlobalIndex());
189 gezelter 2
190 gezelter 246 if (i != molecules_.end() ) {
191 gezelter 2
192 gezelter 507 assert(mol == i->second);
193 gezelter 246
194 gezelter 507 nAtoms_ -= mol->getNAtoms();
195     nBonds_ -= mol->getNBonds();
196     nBends_ -= mol->getNBends();
197     nTorsions_ -= mol->getNTorsions();
198 gezelter 1277 nInversions_ -= mol->getNInversions();
199 gezelter 507 nRigidBodies_ -= mol->getNRigidBodies();
200     nIntegrableObjects_ -= mol->getNIntegrableObjects();
201     nCutoffGroups_ -= mol->getNCutoffGroups();
202     nConstraints_ -= mol->getNConstraintPairs();
203 gezelter 2
204 gezelter 1287 removeInteractionPairs(mol);
205 gezelter 507 molecules_.erase(mol->getGlobalIndex());
206 gezelter 2
207 gezelter 507 delete mol;
208 gezelter 246
209 gezelter 507 return true;
210 gezelter 246 } else {
211 gezelter 507 return false;
212 gezelter 246 }
213 gezelter 507 }
214 gezelter 246
215    
216 gezelter 507 Molecule* SimInfo::beginMolecule(MoleculeIterator& i) {
217 gezelter 246 i = molecules_.begin();
218     return i == molecules_.end() ? NULL : i->second;
219 gezelter 507 }
220 gezelter 246
221 gezelter 507 Molecule* SimInfo::nextMolecule(MoleculeIterator& i) {
222 gezelter 246 ++i;
223     return i == molecules_.end() ? NULL : i->second;
224 gezelter 507 }
225 gezelter 2
226    
227 gezelter 507 void SimInfo::calcNdf() {
228 gezelter 246 int ndf_local;
229     MoleculeIterator i;
230 gezelter 1528 vector<StuntDouble*>::iterator j;
231 gezelter 246 Molecule* mol;
232     StuntDouble* integrableObject;
233 gezelter 2
234 gezelter 246 ndf_local = 0;
235    
236     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
237 gezelter 507 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
238     integrableObject = mol->nextIntegrableObject(j)) {
239 gezelter 2
240 gezelter 507 ndf_local += 3;
241 gezelter 2
242 gezelter 507 if (integrableObject->isDirectional()) {
243     if (integrableObject->isLinear()) {
244     ndf_local += 2;
245     } else {
246     ndf_local += 3;
247     }
248     }
249 gezelter 246
250 tim 770 }
251     }
252 gezelter 246
253     // n_constraints is local, so subtract them on each processor
254     ndf_local -= nConstraints_;
255    
256     #ifdef IS_MPI
257     MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
258     #else
259     ndf_ = ndf_local;
260     #endif
261    
262     // nZconstraints_ is global, as are the 3 COM translations for the
263     // entire system:
264     ndf_ = ndf_ - 3 - nZconstraint_;
265    
266 gezelter 507 }
267 gezelter 2
268 gezelter 945 int SimInfo::getFdf() {
269     #ifdef IS_MPI
270     MPI_Allreduce(&fdf_local,&fdf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
271     #else
272     fdf_ = fdf_local;
273     #endif
274     return fdf_;
275     }
276 gezelter 1577
277     unsigned int SimInfo::getNLocalCutoffGroups(){
278     int nLocalCutoffAtoms = 0;
279     Molecule* mol;
280     MoleculeIterator mi;
281     CutoffGroup* cg;
282     Molecule::CutoffGroupIterator ci;
283 gezelter 945
284 gezelter 1577 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
285    
286     for (cg = mol->beginCutoffGroup(ci); cg != NULL;
287     cg = mol->nextCutoffGroup(ci)) {
288     nLocalCutoffAtoms += cg->getNumAtom();
289    
290     }
291     }
292    
293     return nAtoms_ - nLocalCutoffAtoms + nCutoffGroups_;
294     }
295    
296 gezelter 507 void SimInfo::calcNdfRaw() {
297 gezelter 246 int ndfRaw_local;
298 gezelter 2
299 gezelter 246 MoleculeIterator i;
300 gezelter 1528 vector<StuntDouble*>::iterator j;
301 gezelter 246 Molecule* mol;
302     StuntDouble* integrableObject;
303    
304     // Raw degrees of freedom that we have to set
305     ndfRaw_local = 0;
306    
307     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
308 gezelter 507 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
309     integrableObject = mol->nextIntegrableObject(j)) {
310 gezelter 246
311 gezelter 507 ndfRaw_local += 3;
312 gezelter 246
313 gezelter 507 if (integrableObject->isDirectional()) {
314     if (integrableObject->isLinear()) {
315     ndfRaw_local += 2;
316     } else {
317     ndfRaw_local += 3;
318     }
319     }
320 gezelter 246
321 gezelter 507 }
322 gezelter 246 }
323    
324     #ifdef IS_MPI
325     MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
326     #else
327     ndfRaw_ = ndfRaw_local;
328     #endif
329 gezelter 507 }
330 gezelter 2
331 gezelter 507 void SimInfo::calcNdfTrans() {
332 gezelter 246 int ndfTrans_local;
333 gezelter 2
334 gezelter 246 ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
335 gezelter 2
336    
337 gezelter 246 #ifdef IS_MPI
338     MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
339     #else
340     ndfTrans_ = ndfTrans_local;
341     #endif
342 gezelter 2
343 gezelter 246 ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
344    
345 gezelter 507 }
346 gezelter 2
347 gezelter 1287 void SimInfo::addInteractionPairs(Molecule* mol) {
348     ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
349 gezelter 1528 vector<Bond*>::iterator bondIter;
350     vector<Bend*>::iterator bendIter;
351     vector<Torsion*>::iterator torsionIter;
352     vector<Inversion*>::iterator inversionIter;
353 gezelter 246 Bond* bond;
354     Bend* bend;
355     Torsion* torsion;
356 gezelter 1277 Inversion* inversion;
357 gezelter 246 int a;
358     int b;
359     int c;
360     int d;
361 tim 749
362 gezelter 1287 // atomGroups can be used to add special interaction maps between
363     // groups of atoms that are in two separate rigid bodies.
364     // However, most site-site interactions between two rigid bodies
365     // are probably not special, just the ones between the physically
366     // bonded atoms. Interactions *within* a single rigid body should
367     // always be excluded. These are done at the bottom of this
368     // function.
369    
370 gezelter 1528 map<int, set<int> > atomGroups;
371 tim 749 Molecule::RigidBodyIterator rbIter;
372     RigidBody* rb;
373     Molecule::IntegrableObjectIterator ii;
374     StuntDouble* integrableObject;
375 gezelter 246
376 gezelter 1287 for (integrableObject = mol->beginIntegrableObject(ii);
377     integrableObject != NULL;
378     integrableObject = mol->nextIntegrableObject(ii)) {
379    
380 tim 749 if (integrableObject->isRigidBody()) {
381 gezelter 1287 rb = static_cast<RigidBody*>(integrableObject);
382 gezelter 1528 vector<Atom*> atoms = rb->getAtoms();
383     set<int> rigidAtoms;
384 gezelter 1287 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
385     rigidAtoms.insert(atoms[i]->getGlobalIndex());
386     }
387     for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
388 gezelter 1528 atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
389 gezelter 1287 }
390 tim 749 } else {
391 gezelter 1528 set<int> oneAtomSet;
392 tim 749 oneAtomSet.insert(integrableObject->getGlobalIndex());
393 gezelter 1528 atomGroups.insert(map<int, set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));
394 tim 749 }
395     }
396 gezelter 1287
397     for (bond= mol->beginBond(bondIter); bond != NULL;
398     bond = mol->nextBond(bondIter)) {
399 tim 749
400 gezelter 1287 a = bond->getAtomA()->getGlobalIndex();
401     b = bond->getAtomB()->getGlobalIndex();
402 tim 749
403 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
404     oneTwoInteractions_.addPair(a, b);
405     } else {
406     excludedInteractions_.addPair(a, b);
407     }
408 gezelter 246 }
409 gezelter 2
410 gezelter 1287 for (bend= mol->beginBend(bendIter); bend != NULL;
411     bend = mol->nextBend(bendIter)) {
412    
413 gezelter 507 a = bend->getAtomA()->getGlobalIndex();
414     b = bend->getAtomB()->getGlobalIndex();
415     c = bend->getAtomC()->getGlobalIndex();
416 gezelter 1287
417     if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
418     oneTwoInteractions_.addPair(a, b);
419     oneTwoInteractions_.addPair(b, c);
420     } else {
421     excludedInteractions_.addPair(a, b);
422     excludedInteractions_.addPair(b, c);
423     }
424 gezelter 2
425 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
426     oneThreeInteractions_.addPair(a, c);
427     } else {
428     excludedInteractions_.addPair(a, c);
429     }
430 gezelter 246 }
431 gezelter 2
432 gezelter 1287 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
433     torsion = mol->nextTorsion(torsionIter)) {
434    
435 gezelter 507 a = torsion->getAtomA()->getGlobalIndex();
436     b = torsion->getAtomB()->getGlobalIndex();
437     c = torsion->getAtomC()->getGlobalIndex();
438 gezelter 1287 d = torsion->getAtomD()->getGlobalIndex();
439 cli2 1290
440 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
441     oneTwoInteractions_.addPair(a, b);
442     oneTwoInteractions_.addPair(b, c);
443     oneTwoInteractions_.addPair(c, d);
444     } else {
445     excludedInteractions_.addPair(a, b);
446     excludedInteractions_.addPair(b, c);
447     excludedInteractions_.addPair(c, d);
448     }
449 gezelter 2
450 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
451     oneThreeInteractions_.addPair(a, c);
452     oneThreeInteractions_.addPair(b, d);
453     } else {
454     excludedInteractions_.addPair(a, c);
455     excludedInteractions_.addPair(b, d);
456     }
457 tim 749
458 gezelter 1287 if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
459     oneFourInteractions_.addPair(a, d);
460     } else {
461     excludedInteractions_.addPair(a, d);
462     }
463 gezelter 2 }
464    
465 gezelter 1277 for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
466     inversion = mol->nextInversion(inversionIter)) {
467 gezelter 1287
468 gezelter 1277 a = inversion->getAtomA()->getGlobalIndex();
469     b = inversion->getAtomB()->getGlobalIndex();
470     c = inversion->getAtomC()->getGlobalIndex();
471     d = inversion->getAtomD()->getGlobalIndex();
472    
473 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
474     oneTwoInteractions_.addPair(a, b);
475     oneTwoInteractions_.addPair(a, c);
476     oneTwoInteractions_.addPair(a, d);
477     } else {
478     excludedInteractions_.addPair(a, b);
479     excludedInteractions_.addPair(a, c);
480     excludedInteractions_.addPair(a, d);
481     }
482 gezelter 1277
483 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
484     oneThreeInteractions_.addPair(b, c);
485     oneThreeInteractions_.addPair(b, d);
486     oneThreeInteractions_.addPair(c, d);
487     } else {
488     excludedInteractions_.addPair(b, c);
489     excludedInteractions_.addPair(b, d);
490     excludedInteractions_.addPair(c, d);
491     }
492 gezelter 1277 }
493    
494 gezelter 1287 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
495     rb = mol->nextRigidBody(rbIter)) {
496 gezelter 1528 vector<Atom*> atoms = rb->getAtoms();
497 gezelter 1287 for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
498     for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
499 gezelter 507 a = atoms[i]->getGlobalIndex();
500     b = atoms[j]->getGlobalIndex();
501 gezelter 1287 excludedInteractions_.addPair(a, b);
502 gezelter 507 }
503     }
504 tim 430 }
505    
506 gezelter 507 }
507 gezelter 246
508 gezelter 1287 void SimInfo::removeInteractionPairs(Molecule* mol) {
509     ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
510 gezelter 1528 vector<Bond*>::iterator bondIter;
511     vector<Bend*>::iterator bendIter;
512     vector<Torsion*>::iterator torsionIter;
513     vector<Inversion*>::iterator inversionIter;
514 gezelter 246 Bond* bond;
515     Bend* bend;
516     Torsion* torsion;
517 gezelter 1277 Inversion* inversion;
518 gezelter 246 int a;
519     int b;
520     int c;
521     int d;
522 tim 749
523 gezelter 1528 map<int, set<int> > atomGroups;
524 tim 749 Molecule::RigidBodyIterator rbIter;
525     RigidBody* rb;
526     Molecule::IntegrableObjectIterator ii;
527     StuntDouble* integrableObject;
528 gezelter 246
529 gezelter 1287 for (integrableObject = mol->beginIntegrableObject(ii);
530     integrableObject != NULL;
531     integrableObject = mol->nextIntegrableObject(ii)) {
532    
533 tim 749 if (integrableObject->isRigidBody()) {
534 gezelter 1287 rb = static_cast<RigidBody*>(integrableObject);
535 gezelter 1528 vector<Atom*> atoms = rb->getAtoms();
536     set<int> rigidAtoms;
537 gezelter 1287 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
538     rigidAtoms.insert(atoms[i]->getGlobalIndex());
539     }
540     for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
541 gezelter 1528 atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
542 gezelter 1287 }
543 tim 749 } else {
544 gezelter 1528 set<int> oneAtomSet;
545 tim 749 oneAtomSet.insert(integrableObject->getGlobalIndex());
546 gezelter 1528 atomGroups.insert(map<int, set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));
547 tim 749 }
548     }
549    
550 gezelter 1287 for (bond= mol->beginBond(bondIter); bond != NULL;
551     bond = mol->nextBond(bondIter)) {
552    
553     a = bond->getAtomA()->getGlobalIndex();
554     b = bond->getAtomB()->getGlobalIndex();
555 tim 749
556 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
557     oneTwoInteractions_.removePair(a, b);
558     } else {
559     excludedInteractions_.removePair(a, b);
560     }
561 gezelter 2 }
562 gezelter 246
563 gezelter 1287 for (bend= mol->beginBend(bendIter); bend != NULL;
564     bend = mol->nextBend(bendIter)) {
565    
566 gezelter 507 a = bend->getAtomA()->getGlobalIndex();
567     b = bend->getAtomB()->getGlobalIndex();
568     c = bend->getAtomC()->getGlobalIndex();
569 gezelter 1287
570     if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
571     oneTwoInteractions_.removePair(a, b);
572     oneTwoInteractions_.removePair(b, c);
573     } else {
574     excludedInteractions_.removePair(a, b);
575     excludedInteractions_.removePair(b, c);
576     }
577 gezelter 246
578 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
579     oneThreeInteractions_.removePair(a, c);
580     } else {
581     excludedInteractions_.removePair(a, c);
582     }
583 gezelter 2 }
584 gezelter 246
585 gezelter 1287 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
586     torsion = mol->nextTorsion(torsionIter)) {
587    
588 gezelter 507 a = torsion->getAtomA()->getGlobalIndex();
589     b = torsion->getAtomB()->getGlobalIndex();
590     c = torsion->getAtomC()->getGlobalIndex();
591 gezelter 1287 d = torsion->getAtomD()->getGlobalIndex();
592    
593     if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
594     oneTwoInteractions_.removePair(a, b);
595     oneTwoInteractions_.removePair(b, c);
596     oneTwoInteractions_.removePair(c, d);
597     } else {
598     excludedInteractions_.removePair(a, b);
599     excludedInteractions_.removePair(b, c);
600     excludedInteractions_.removePair(c, d);
601     }
602 gezelter 246
603 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
604     oneThreeInteractions_.removePair(a, c);
605     oneThreeInteractions_.removePair(b, d);
606     } else {
607     excludedInteractions_.removePair(a, c);
608     excludedInteractions_.removePair(b, d);
609     }
610 tim 749
611 gezelter 1287 if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
612     oneFourInteractions_.removePair(a, d);
613     } else {
614     excludedInteractions_.removePair(a, d);
615     }
616     }
617 tim 749
618 gezelter 1287 for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
619     inversion = mol->nextInversion(inversionIter)) {
620 tim 749
621 gezelter 1277 a = inversion->getAtomA()->getGlobalIndex();
622     b = inversion->getAtomB()->getGlobalIndex();
623     c = inversion->getAtomC()->getGlobalIndex();
624     d = inversion->getAtomD()->getGlobalIndex();
625    
626 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
627     oneTwoInteractions_.removePair(a, b);
628     oneTwoInteractions_.removePair(a, c);
629     oneTwoInteractions_.removePair(a, d);
630     } else {
631     excludedInteractions_.removePair(a, b);
632     excludedInteractions_.removePair(a, c);
633     excludedInteractions_.removePair(a, d);
634     }
635 gezelter 1277
636 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
637     oneThreeInteractions_.removePair(b, c);
638     oneThreeInteractions_.removePair(b, d);
639     oneThreeInteractions_.removePair(c, d);
640     } else {
641     excludedInteractions_.removePair(b, c);
642     excludedInteractions_.removePair(b, d);
643     excludedInteractions_.removePair(c, d);
644     }
645 gezelter 1277 }
646    
647 gezelter 1287 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
648     rb = mol->nextRigidBody(rbIter)) {
649 gezelter 1528 vector<Atom*> atoms = rb->getAtoms();
650 gezelter 1287 for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
651     for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
652 gezelter 507 a = atoms[i]->getGlobalIndex();
653     b = atoms[j]->getGlobalIndex();
654 gezelter 1287 excludedInteractions_.removePair(a, b);
655 gezelter 507 }
656     }
657 tim 430 }
658 gezelter 1287
659 gezelter 507 }
660 gezelter 1287
661    
662 gezelter 507 void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) {
663 gezelter 246 int curStampId;
664 gezelter 1287
665 gezelter 246 //index from 0
666     curStampId = moleculeStamps_.size();
667 gezelter 2
668 gezelter 246 moleculeStamps_.push_back(molStamp);
669     molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
670 gezelter 507 }
671 gezelter 2
672 gezelter 1530
673     /**
674     * update
675     *
676 gezelter 1535 * Performs the global checks and variable settings after the
677     * objects have been created.
678 gezelter 1530 *
679     */
680 gezelter 1535 void SimInfo::update() {
681 gezelter 1530 setupSimVariables();
682 gezelter 246 calcNdf();
683     calcNdfRaw();
684     calcNdfTrans();
685 gezelter 507 }
686 gezelter 1528
687 gezelter 1535 /**
688     * getSimulatedAtomTypes
689     *
690     * Returns an STL set of AtomType* that are actually present in this
691     * simulation. Must query all processors to assemble this information.
692     *
693     */
694 gezelter 1528 set<AtomType*> SimInfo::getSimulatedAtomTypes() {
695 gezelter 246 SimInfo::MoleculeIterator mi;
696     Molecule* mol;
697     Molecule::AtomIterator ai;
698     Atom* atom;
699 gezelter 1528 set<AtomType*> atomTypes;
700    
701 gezelter 1588 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
702     for(atom = mol->beginAtom(ai); atom != NULL;
703     atom = mol->nextAtom(ai)) {
704 gezelter 507 atomTypes.insert(atom->getAtomType());
705 gezelter 1529 }
706     }
707 gezelter 1588
708 gezelter 1535 #ifdef IS_MPI
709    
710     // loop over the found atom types on this processor, and add their
711     // numerical idents to a vector:
712 gezelter 1588
713 gezelter 1535 vector<int> foundTypes;
714     set<AtomType*>::iterator i;
715     for (i = atomTypes.begin(); i != atomTypes.end(); ++i)
716     foundTypes.push_back( (*i)->getIdent() );
717    
718     // count_local holds the number of found types on this processor
719     int count_local = foundTypes.size();
720    
721 gezelter 1588 int nproc = MPI::COMM_WORLD.Get_size();
722 gezelter 1535
723 gezelter 1588 // we need arrays to hold the counts and displacement vectors for
724     // all processors
725     vector<int> counts(nproc, 0);
726     vector<int> disps(nproc, 0);
727 gezelter 1535
728 gezelter 1588 // fill the counts array
729     MPI::COMM_WORLD.Allgather(&count_local, 1, MPI::INT, &counts[0],
730     1, MPI::INT);
731    
732     // use the processor counts to compute the displacement array
733     disps[0] = 0;
734     int totalCount = counts[0];
735     for (int iproc = 1; iproc < nproc; iproc++) {
736     disps[iproc] = disps[iproc-1] + counts[iproc-1];
737     totalCount += counts[iproc];
738     }
739 gezelter 1535
740 gezelter 1588 // we need a (possibly redundant) set of all found types:
741     vector<int> ftGlobal(totalCount);
742 gezelter 2
743 gezelter 1588 // now spray out the foundTypes to all the other processors:
744 gezelter 1535 MPI::COMM_WORLD.Allgatherv(&foundTypes[0], count_local, MPI::INT,
745 gezelter 1588 &ftGlobal[0], &counts[0], &disps[0],
746     MPI::INT);
747 gezelter 1126
748 gezelter 1588 vector<int>::iterator j;
749    
750 gezelter 1535 // foundIdents is a stl set, so inserting an already found ident
751     // will have no effect.
752     set<int> foundIdents;
753 gezelter 1588
754 gezelter 1535 for (j = ftGlobal.begin(); j != ftGlobal.end(); ++j)
755     foundIdents.insert((*j));
756 gezelter 1528
757 gezelter 1535 // now iterate over the foundIdents and get the actual atom types
758     // that correspond to these:
759     set<int>::iterator it;
760 gezelter 1588 for (it = foundIdents.begin(); it != foundIdents.end(); ++it)
761 gezelter 1535 atomTypes.insert( forceField_->getAtomType((*it)) );
762    
763     #endif
764 gezelter 1588
765 gezelter 1535 return atomTypes;
766 gezelter 1528 }
767 chrisfen 611
768 gezelter 1534 void SimInfo::setupSimVariables() {
769     useAtomicVirial_ = simParams_->getUseAtomicVirial();
770     // we only call setAccumulateBoxDipole if the accumulateBoxDipole parameter is true
771     calcBoxDipole_ = false;
772     if ( simParams_->haveAccumulateBoxDipole() )
773     if ( simParams_->getAccumulateBoxDipole() ) {
774     calcBoxDipole_ = true;
775     }
776 gezelter 1583
777 gezelter 1528 set<AtomType*>::iterator i;
778     set<AtomType*> atomTypes;
779 gezelter 1534 atomTypes = getSimulatedAtomTypes();
780 gezelter 1528 int usesElectrostatic = 0;
781     int usesMetallic = 0;
782     int usesDirectional = 0;
783 gezelter 246 //loop over all of the atom types
784     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
785 gezelter 1528 usesElectrostatic |= (*i)->isElectrostatic();
786     usesMetallic |= (*i)->isMetal();
787     usesDirectional |= (*i)->isDirectional();
788 gezelter 246 }
789 gezelter 1583
790 gezelter 246 #ifdef IS_MPI
791     int temp;
792 gezelter 1528 temp = usesDirectional;
793     MPI_Allreduce(&temp, &usesDirectionalAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
794 gezelter 1583
795 gezelter 1528 temp = usesMetallic;
796     MPI_Allreduce(&temp, &usesMetallicAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
797 gezelter 1583
798 gezelter 1528 temp = usesElectrostatic;
799     MPI_Allreduce(&temp, &usesElectrostaticAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
800 gezelter 1583 #else
801    
802     usesDirectionalAtoms_ = usesDirectional;
803     usesMetallicAtoms_ = usesMetallic;
804     usesElectrostaticAtoms_ = usesElectrostatic;
805    
806 gezelter 2 #endif
807 gezelter 1583
808     requiresPrepair_ = usesMetallicAtoms_ ? true : false;
809     requiresSkipCorrection_ = usesElectrostaticAtoms_ ? true : false;
810     requiresSelfCorrection_ = usesElectrostaticAtoms_ ? true : false;
811 gezelter 507 }
812 gezelter 2
813 gezelter 1549
814     vector<int> SimInfo::getGlobalAtomIndices() {
815     SimInfo::MoleculeIterator mi;
816     Molecule* mol;
817     Molecule::AtomIterator ai;
818     Atom* atom;
819    
820     vector<int> GlobalAtomIndices(getNAtoms(), 0);
821    
822     for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
823    
824     for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
825     GlobalAtomIndices[atom->getLocalIndex()] = atom->getGlobalIndex();
826     }
827     }
828     return GlobalAtomIndices;
829     }
830    
831    
832     vector<int> SimInfo::getGlobalGroupIndices() {
833     SimInfo::MoleculeIterator mi;
834     Molecule* mol;
835     Molecule::CutoffGroupIterator ci;
836     CutoffGroup* cg;
837    
838     vector<int> GlobalGroupIndices;
839    
840     for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
841    
842     //local index of cutoff group is trivial, it only depends on the
843     //order of travesing
844     for (cg = mol->beginCutoffGroup(ci); cg != NULL;
845     cg = mol->nextCutoffGroup(ci)) {
846     GlobalGroupIndices.push_back(cg->getGlobalIndex());
847     }
848     }
849     return GlobalGroupIndices;
850     }
851    
852    
853 gezelter 1569 void SimInfo::prepareTopology() {
854 gezelter 1287 int nExclude, nOneTwo, nOneThree, nOneFour;
855 gezelter 2
856 gezelter 246 //calculate mass ratio of cutoff group
857     SimInfo::MoleculeIterator mi;
858     Molecule* mol;
859     Molecule::CutoffGroupIterator ci;
860     CutoffGroup* cg;
861     Molecule::AtomIterator ai;
862     Atom* atom;
863 tim 963 RealType totalMass;
864 gezelter 246
865 gezelter 1581 /**
866     * The mass factor is the relative mass of an atom to the total
867     * mass of the cutoff group it belongs to. By default, all atoms
868     * are their own cutoff groups, and therefore have mass factors of
869     * 1. We need some special handling for massless atoms, which
870     * will be treated as carrying the entire mass of the cutoff
871     * group.
872     */
873 gezelter 1569 massFactors_.clear();
874 gezelter 1581 massFactors_.resize(getNAtoms(), 1.0);
875 gezelter 2
876 gezelter 246 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
877 gezelter 1569 for (cg = mol->beginCutoffGroup(ci); cg != NULL;
878     cg = mol->nextCutoffGroup(ci)) {
879 gezelter 2
880 gezelter 507 totalMass = cg->getMass();
881     for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
882 chrisfen 645 // Check for massless groups - set mfact to 1 if true
883 gezelter 1581 if (totalMass != 0)
884     massFactors_[atom->getLocalIndex()] = atom->getMass()/totalMass;
885 chrisfen 645 else
886 gezelter 1581 massFactors_[atom->getLocalIndex()] = 1.0;
887 gezelter 507 }
888     }
889 gezelter 246 }
890 gezelter 2
891 gezelter 1544 // Build the identArray_
892 gezelter 2
893 gezelter 1544 identArray_.clear();
894     identArray_.reserve(getNAtoms());
895 gezelter 246 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
896 gezelter 507 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
897 gezelter 1544 identArray_.push_back(atom->getIdent());
898 gezelter 507 }
899 gezelter 246 }
900    
901 gezelter 1569 //scan topology
902 gezelter 1287
903     nExclude = excludedInteractions_.getSize();
904     nOneTwo = oneTwoInteractions_.getSize();
905     nOneThree = oneThreeInteractions_.getSize();
906     nOneFour = oneFourInteractions_.getSize();
907    
908     int* excludeList = excludedInteractions_.getPairList();
909     int* oneTwoList = oneTwoInteractions_.getPairList();
910     int* oneThreeList = oneThreeInteractions_.getPairList();
911     int* oneFourList = oneFourInteractions_.getPairList();
912    
913 gezelter 1569 topologyDone_ = true;
914 gezelter 507 }
915 chrisfen 143
916 gezelter 507 void SimInfo::addProperty(GenericData* genData) {
917 gezelter 246 properties_.addProperty(genData);
918 gezelter 507 }
919 gezelter 2
920 gezelter 1528 void SimInfo::removeProperty(const string& propName) {
921 gezelter 246 properties_.removeProperty(propName);
922 gezelter 507 }
923 gezelter 2
924 gezelter 507 void SimInfo::clearProperties() {
925 gezelter 246 properties_.clearProperties();
926 gezelter 507 }
927 gezelter 2
928 gezelter 1528 vector<string> SimInfo::getPropertyNames() {
929 gezelter 246 return properties_.getPropertyNames();
930 gezelter 507 }
931 gezelter 246
932 gezelter 1528 vector<GenericData*> SimInfo::getProperties() {
933 gezelter 246 return properties_.getProperties();
934 gezelter 507 }
935 gezelter 2
936 gezelter 1528 GenericData* SimInfo::getPropertyByName(const string& propName) {
937 gezelter 246 return properties_.getPropertyByName(propName);
938 gezelter 507 }
939 gezelter 2
940 gezelter 507 void SimInfo::setSnapshotManager(SnapshotManager* sman) {
941 tim 432 if (sman_ == sman) {
942 gezelter 507 return;
943 tim 432 }
944     delete sman_;
945 gezelter 246 sman_ = sman;
946 gezelter 2
947 gezelter 246 Molecule* mol;
948     RigidBody* rb;
949     Atom* atom;
950 gezelter 1540 CutoffGroup* cg;
951 gezelter 246 SimInfo::MoleculeIterator mi;
952     Molecule::RigidBodyIterator rbIter;
953 gezelter 1540 Molecule::AtomIterator atomIter;
954     Molecule::CutoffGroupIterator cgIter;
955 gezelter 246
956     for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
957    
958 gezelter 507 for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) {
959     atom->setSnapshotManager(sman_);
960     }
961 gezelter 246
962 gezelter 507 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
963     rb->setSnapshotManager(sman_);
964     }
965 gezelter 1540
966     for (cg = mol->beginCutoffGroup(cgIter); cg != NULL; cg = mol->nextCutoffGroup(cgIter)) {
967     cg->setSnapshotManager(sman_);
968     }
969 gezelter 246 }
970 gezelter 2
971 gezelter 507 }
972 gezelter 2
973 gezelter 507 Vector3d SimInfo::getComVel(){
974 gezelter 246 SimInfo::MoleculeIterator i;
975     Molecule* mol;
976 gezelter 2
977 gezelter 246 Vector3d comVel(0.0);
978 tim 963 RealType totalMass = 0.0;
979 gezelter 2
980 gezelter 246
981     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
982 tim 963 RealType mass = mol->getMass();
983 gezelter 507 totalMass += mass;
984     comVel += mass * mol->getComVel();
985 gezelter 246 }
986 gezelter 2
987 gezelter 246 #ifdef IS_MPI
988 tim 963 RealType tmpMass = totalMass;
989 gezelter 246 Vector3d tmpComVel(comVel);
990 tim 963 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
991     MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
992 gezelter 246 #endif
993    
994     comVel /= totalMass;
995    
996     return comVel;
997 gezelter 507 }
998 gezelter 2
999 gezelter 507 Vector3d SimInfo::getCom(){
1000 gezelter 246 SimInfo::MoleculeIterator i;
1001     Molecule* mol;
1002 gezelter 2
1003 gezelter 246 Vector3d com(0.0);
1004 tim 963 RealType totalMass = 0.0;
1005 gezelter 246
1006     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1007 tim 963 RealType mass = mol->getMass();
1008 gezelter 507 totalMass += mass;
1009     com += mass * mol->getCom();
1010 gezelter 246 }
1011 gezelter 2
1012     #ifdef IS_MPI
1013 tim 963 RealType tmpMass = totalMass;
1014 gezelter 246 Vector3d tmpCom(com);
1015 tim 963 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1016     MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1017 gezelter 2 #endif
1018    
1019 gezelter 246 com /= totalMass;
1020 gezelter 2
1021 gezelter 246 return com;
1022 gezelter 2
1023 gezelter 507 }
1024 gezelter 246
1025 gezelter 1528 ostream& operator <<(ostream& o, SimInfo& info) {
1026 gezelter 246
1027     return o;
1028 gezelter 507 }
1029 chuckv 555
1030    
1031     /*
1032     Returns center of mass and center of mass velocity in one function call.
1033     */
1034    
1035     void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){
1036     SimInfo::MoleculeIterator i;
1037     Molecule* mol;
1038    
1039    
1040 tim 963 RealType totalMass = 0.0;
1041 chuckv 555
1042 gezelter 246
1043 chuckv 555 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1044 tim 963 RealType mass = mol->getMass();
1045 chuckv 555 totalMass += mass;
1046     com += mass * mol->getCom();
1047     comVel += mass * mol->getComVel();
1048     }
1049    
1050     #ifdef IS_MPI
1051 tim 963 RealType tmpMass = totalMass;
1052 chuckv 555 Vector3d tmpCom(com);
1053     Vector3d tmpComVel(comVel);
1054 tim 963 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1055     MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1056     MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1057 chuckv 555 #endif
1058    
1059     com /= totalMass;
1060     comVel /= totalMass;
1061     }
1062    
1063     /*
1064     Return intertia tensor for entire system and angular momentum Vector.
1065 chuckv 557
1066    
1067     [ Ixx -Ixy -Ixz ]
1068 gezelter 1505 J =| -Iyx Iyy -Iyz |
1069 chuckv 557 [ -Izx -Iyz Izz ]
1070 chuckv 555 */
1071    
1072     void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){
1073    
1074    
1075 tim 963 RealType xx = 0.0;
1076     RealType yy = 0.0;
1077     RealType zz = 0.0;
1078     RealType xy = 0.0;
1079     RealType xz = 0.0;
1080     RealType yz = 0.0;
1081 chuckv 555 Vector3d com(0.0);
1082     Vector3d comVel(0.0);
1083    
1084     getComAll(com, comVel);
1085    
1086     SimInfo::MoleculeIterator i;
1087     Molecule* mol;
1088    
1089     Vector3d thisq(0.0);
1090     Vector3d thisv(0.0);
1091    
1092 tim 963 RealType thisMass = 0.0;
1093 chuckv 555
1094    
1095    
1096    
1097     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1098    
1099     thisq = mol->getCom()-com;
1100     thisv = mol->getComVel()-comVel;
1101     thisMass = mol->getMass();
1102     // Compute moment of intertia coefficients.
1103     xx += thisq[0]*thisq[0]*thisMass;
1104     yy += thisq[1]*thisq[1]*thisMass;
1105     zz += thisq[2]*thisq[2]*thisMass;
1106    
1107     // compute products of intertia
1108     xy += thisq[0]*thisq[1]*thisMass;
1109     xz += thisq[0]*thisq[2]*thisMass;
1110     yz += thisq[1]*thisq[2]*thisMass;
1111    
1112     angularMomentum += cross( thisq, thisv ) * thisMass;
1113    
1114     }
1115    
1116    
1117     inertiaTensor(0,0) = yy + zz;
1118     inertiaTensor(0,1) = -xy;
1119     inertiaTensor(0,2) = -xz;
1120     inertiaTensor(1,0) = -xy;
1121 chuckv 557 inertiaTensor(1,1) = xx + zz;
1122 chuckv 555 inertiaTensor(1,2) = -yz;
1123     inertiaTensor(2,0) = -xz;
1124     inertiaTensor(2,1) = -yz;
1125     inertiaTensor(2,2) = xx + yy;
1126    
1127     #ifdef IS_MPI
1128     Mat3x3d tmpI(inertiaTensor);
1129     Vector3d tmpAngMom;
1130 tim 963 MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1131     MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1132 chuckv 555 #endif
1133    
1134     return;
1135     }
1136    
1137     //Returns the angular momentum of the system
1138     Vector3d SimInfo::getAngularMomentum(){
1139    
1140     Vector3d com(0.0);
1141     Vector3d comVel(0.0);
1142     Vector3d angularMomentum(0.0);
1143    
1144     getComAll(com,comVel);
1145    
1146     SimInfo::MoleculeIterator i;
1147     Molecule* mol;
1148    
1149 chuckv 557 Vector3d thisr(0.0);
1150     Vector3d thisp(0.0);
1151 chuckv 555
1152 tim 963 RealType thisMass;
1153 chuckv 555
1154     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1155 chuckv 557 thisMass = mol->getMass();
1156     thisr = mol->getCom()-com;
1157     thisp = (mol->getComVel()-comVel)*thisMass;
1158 chuckv 555
1159 chuckv 557 angularMomentum += cross( thisr, thisp );
1160    
1161 chuckv 555 }
1162    
1163     #ifdef IS_MPI
1164     Vector3d tmpAngMom;
1165 tim 963 MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1166 chuckv 555 #endif
1167    
1168     return angularMomentum;
1169     }
1170    
1171 tim 1024 StuntDouble* SimInfo::getIOIndexToIntegrableObject(int index) {
1172     return IOIndexToIntegrableObject.at(index);
1173     }
1174    
1175 gezelter 1528 void SimInfo::setIOIndexToIntegrableObject(const vector<StuntDouble*>& v) {
1176 tim 1024 IOIndexToIntegrableObject= v;
1177     }
1178    
1179 chuckv 1103 /* Returns the Volume of the simulation based on a ellipsoid with semi-axes
1180     based on the radius of gyration V=4/3*Pi*R_1*R_2*R_3
1181     where R_i are related to the principle inertia moments R_i = sqrt(C*I_i/N), this reduces to
1182     V = 4/3*Pi*(C/N)^3/2*sqrt(det(I)). See S.E. Baltazar et. al. Comp. Mat. Sci. 37 (2006) 526-536.
1183     */
1184     void SimInfo::getGyrationalVolume(RealType &volume){
1185     Mat3x3d intTensor;
1186     RealType det;
1187     Vector3d dummyAngMom;
1188     RealType sysconstants;
1189     RealType geomCnst;
1190    
1191     geomCnst = 3.0/2.0;
1192     /* Get the inertial tensor and angular momentum for free*/
1193     getInertiaTensor(intTensor,dummyAngMom);
1194    
1195     det = intTensor.determinant();
1196     sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1197     volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(det);
1198     return;
1199     }
1200    
1201     void SimInfo::getGyrationalVolume(RealType &volume, RealType &detI){
1202     Mat3x3d intTensor;
1203     Vector3d dummyAngMom;
1204     RealType sysconstants;
1205     RealType geomCnst;
1206    
1207     geomCnst = 3.0/2.0;
1208     /* Get the inertial tensor and angular momentum for free*/
1209     getInertiaTensor(intTensor,dummyAngMom);
1210    
1211     detI = intTensor.determinant();
1212     sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1213     volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(detI);
1214     return;
1215     }
1216 tim 1024 /*
1217 gezelter 1528 void SimInfo::setStuntDoubleFromGlobalIndex(vector<StuntDouble*> v) {
1218 tim 1024 assert( v.size() == nAtoms_ + nRigidBodies_);
1219     sdByGlobalIndex_ = v;
1220     }
1221    
1222     StuntDouble* SimInfo::getStuntDoubleFromGlobalIndex(int index) {
1223     //assert(index < nAtoms_ + nRigidBodies_);
1224     return sdByGlobalIndex_.at(index);
1225     }
1226     */
1227 gezelter 1528 int SimInfo::getNGlobalConstraints() {
1228     int nGlobalConstraints;
1229     #ifdef IS_MPI
1230     MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM,
1231     MPI_COMM_WORLD);
1232     #else
1233     nGlobalConstraints = nConstraints_;
1234     #endif
1235     return nGlobalConstraints;
1236     }
1237    
1238 gezelter 1390 }//end namespace OpenMD
1239 gezelter 246

Properties

Name Value
svn:keywords Author Id Revision Date