ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1744
Committed: Tue Jun 5 18:07:08 2012 UTC (12 years, 10 months ago) by gezelter
File size: 39659 byte(s)
Log Message:
Fixes for minimization

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

Properties

Name Value
svn:keywords Author Id Revision Date