ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1627
Committed: Tue Sep 13 22:05:04 2011 UTC (13 years, 7 months ago) by gezelter
File size: 38782 byte(s)
Log Message:
Splitting out ifstrstream into a header and an implementation.  This
means that much of the code that depends on it can be compiled only
once and the parallel I/O is concentrated into a few files.  To do
this, a number of files that relied on basic_ifstrstream to load the
mpi header had to be modified to manage their own headers.


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

Properties

Name Value
svn:keywords Author Id Revision Date