ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1779
Committed: Mon Aug 20 17:51:39 2012 UTC (12 years, 8 months ago) by gezelter
File size: 32768 byte(s)
Log Message:
Prevent a memory access bug when dump files aren't configured correctly.

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 gezelter 1779 for (vector<Component*>::iterator i = components.begin();
92     i !=components.end(); ++i) {
93 gezelter 1528 molStamp = (*i)->getMoleculeStamp();
94     nMolWithSameStamp = (*i)->getNMol();
95 tim 770
96 gezelter 1528 addMoleculeStamp(molStamp, nMolWithSameStamp);
97    
98     //calculate atoms in molecules
99     nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp;
100    
101     //calculate atoms in cutoff groups
102     int nAtomsInGroups = 0;
103     int nCutoffGroupsInStamp = molStamp->getNCutoffGroups();
104    
105     for (int j=0; j < nCutoffGroupsInStamp; j++) {
106     cgStamp = molStamp->getCutoffGroupStamp(j);
107     nAtomsInGroups += cgStamp->getNMembers();
108 gezelter 507 }
109 gezelter 1528
110     nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
111    
112     nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;
113    
114     //calculate atoms in rigid bodies
115     int nAtomsInRigidBodies = 0;
116     int nRigidBodiesInStamp = molStamp->getNRigidBodies();
117    
118     for (int j=0; j < nRigidBodiesInStamp; j++) {
119     rbStamp = molStamp->getRigidBodyStamp(j);
120     nAtomsInRigidBodies += rbStamp->getNMembers();
121     }
122    
123     nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
124     nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;
125    
126     }
127    
128     //every free atom (atom does not belong to cutoff groups) is a cutoff
129     //group therefore the total number of cutoff groups in the system is
130     //equal to the total number of atoms minus number of atoms belong to
131     //cutoff group defined in meta-data file plus the number of cutoff
132     //groups defined in meta-data file
133 gezelter 1540
134 gezelter 1528 nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
135    
136     //every free atom (atom does not belong to rigid bodies) is an
137     //integrable object therefore the total number of integrable objects
138     //in the system is equal to the total number of atoms minus number of
139     //atoms belong to rigid body defined in meta-data file plus the number
140     //of rigid bodies defined in meta-data file
141     nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms
142     + nGlobalRigidBodies_;
143    
144     nGlobalMols_ = molStampIds_.size();
145     molToProcMap_.resize(nGlobalMols_);
146     }
147 chrisfen 645
148 gezelter 507 SimInfo::~SimInfo() {
149 gezelter 1528 map<int, Molecule*>::iterator i;
150 tim 398 for (i = molecules_.begin(); i != molecules_.end(); ++i) {
151 gezelter 507 delete i->second;
152 tim 398 }
153     molecules_.clear();
154 tim 490
155 gezelter 246 delete sman_;
156     delete simParams_;
157     delete forceField_;
158 gezelter 507 }
159 gezelter 2
160    
161 gezelter 507 bool SimInfo::addMolecule(Molecule* mol) {
162 gezelter 246 MoleculeIterator i;
163 gezelter 1528
164 gezelter 246 i = molecules_.find(mol->getGlobalIndex());
165     if (i == molecules_.end() ) {
166 gezelter 1528
167     molecules_.insert(make_pair(mol->getGlobalIndex(), mol));
168    
169 gezelter 507 nAtoms_ += mol->getNAtoms();
170     nBonds_ += mol->getNBonds();
171     nBends_ += mol->getNBends();
172     nTorsions_ += mol->getNTorsions();
173 gezelter 1277 nInversions_ += mol->getNInversions();
174 gezelter 507 nRigidBodies_ += mol->getNRigidBodies();
175     nIntegrableObjects_ += mol->getNIntegrableObjects();
176     nCutoffGroups_ += mol->getNCutoffGroups();
177     nConstraints_ += mol->getNConstraintPairs();
178 gezelter 1528
179 gezelter 1287 addInteractionPairs(mol);
180 gezelter 1528
181 gezelter 507 return true;
182 gezelter 246 } else {
183 gezelter 507 return false;
184 gezelter 246 }
185 gezelter 507 }
186 gezelter 1528
187 gezelter 507 bool SimInfo::removeMolecule(Molecule* mol) {
188 gezelter 246 MoleculeIterator i;
189     i = molecules_.find(mol->getGlobalIndex());
190 gezelter 2
191 gezelter 246 if (i != molecules_.end() ) {
192 gezelter 2
193 gezelter 507 assert(mol == i->second);
194 gezelter 246
195 gezelter 507 nAtoms_ -= mol->getNAtoms();
196     nBonds_ -= mol->getNBonds();
197     nBends_ -= mol->getNBends();
198     nTorsions_ -= mol->getNTorsions();
199 gezelter 1277 nInversions_ -= mol->getNInversions();
200 gezelter 507 nRigidBodies_ -= mol->getNRigidBodies();
201     nIntegrableObjects_ -= mol->getNIntegrableObjects();
202     nCutoffGroups_ -= mol->getNCutoffGroups();
203     nConstraints_ -= mol->getNConstraintPairs();
204 gezelter 2
205 gezelter 1287 removeInteractionPairs(mol);
206 gezelter 507 molecules_.erase(mol->getGlobalIndex());
207 gezelter 2
208 gezelter 507 delete mol;
209 gezelter 246
210 gezelter 507 return true;
211 gezelter 246 } else {
212 gezelter 507 return false;
213 gezelter 246 }
214 gezelter 507 }
215 gezelter 246
216    
217 gezelter 507 Molecule* SimInfo::beginMolecule(MoleculeIterator& i) {
218 gezelter 246 i = molecules_.begin();
219     return i == molecules_.end() ? NULL : i->second;
220 gezelter 507 }
221 gezelter 246
222 gezelter 507 Molecule* SimInfo::nextMolecule(MoleculeIterator& i) {
223 gezelter 246 ++i;
224     return i == molecules_.end() ? NULL : i->second;
225 gezelter 507 }
226 gezelter 2
227    
228 gezelter 507 void SimInfo::calcNdf() {
229 gezelter 1715 int ndf_local, nfq_local;
230 gezelter 246 MoleculeIterator i;
231 gezelter 1528 vector<StuntDouble*>::iterator j;
232 gezelter 1715 vector<Atom*>::iterator k;
233    
234 gezelter 246 Molecule* mol;
235 gezelter 1769 StuntDouble* sd;
236 gezelter 1715 Atom* atom;
237 gezelter 2
238 gezelter 246 ndf_local = 0;
239 gezelter 1715 nfq_local = 0;
240 gezelter 246
241     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
242 gezelter 2
243 gezelter 1769 for (sd = mol->beginIntegrableObject(j); sd != NULL;
244     sd = mol->nextIntegrableObject(j)) {
245    
246 gezelter 507 ndf_local += 3;
247 gezelter 2
248 gezelter 1769 if (sd->isDirectional()) {
249     if (sd->isLinear()) {
250 gezelter 507 ndf_local += 2;
251     } else {
252     ndf_local += 3;
253     }
254     }
255 tim 770 }
256 gezelter 1769
257 gezelter 1715 for (atom = mol->beginFluctuatingCharge(k); atom != NULL;
258     atom = mol->nextFluctuatingCharge(k)) {
259     if (atom->isFluctuatingCharge()) {
260     nfq_local++;
261     }
262     }
263 tim 770 }
264 gezelter 246
265 gezelter 1744 ndfLocal_ = ndf_local;
266    
267 gezelter 246 // n_constraints is local, so subtract them on each processor
268     ndf_local -= nConstraints_;
269    
270     #ifdef IS_MPI
271     MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
272 gezelter 1715 MPI_Allreduce(&nfq_local,&nGlobalFluctuatingCharges_,1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
273 gezelter 246 #else
274     ndf_ = ndf_local;
275 gezelter 1715 nGlobalFluctuatingCharges_ = nfq_local;
276 gezelter 246 #endif
277    
278     // nZconstraints_ is global, as are the 3 COM translations for the
279     // entire system:
280     ndf_ = ndf_ - 3 - nZconstraint_;
281    
282 gezelter 507 }
283 gezelter 2
284 gezelter 945 int SimInfo::getFdf() {
285     #ifdef IS_MPI
286     MPI_Allreduce(&fdf_local,&fdf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
287     #else
288     fdf_ = fdf_local;
289     #endif
290     return fdf_;
291     }
292 gezelter 1577
293     unsigned int SimInfo::getNLocalCutoffGroups(){
294     int nLocalCutoffAtoms = 0;
295     Molecule* mol;
296     MoleculeIterator mi;
297     CutoffGroup* cg;
298     Molecule::CutoffGroupIterator ci;
299 gezelter 945
300 gezelter 1577 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
301    
302     for (cg = mol->beginCutoffGroup(ci); cg != NULL;
303     cg = mol->nextCutoffGroup(ci)) {
304     nLocalCutoffAtoms += cg->getNumAtom();
305    
306     }
307     }
308    
309     return nAtoms_ - nLocalCutoffAtoms + nCutoffGroups_;
310     }
311    
312 gezelter 507 void SimInfo::calcNdfRaw() {
313 gezelter 246 int ndfRaw_local;
314 gezelter 2
315 gezelter 246 MoleculeIterator i;
316 gezelter 1528 vector<StuntDouble*>::iterator j;
317 gezelter 246 Molecule* mol;
318 gezelter 1769 StuntDouble* sd;
319 gezelter 246
320     // Raw degrees of freedom that we have to set
321     ndfRaw_local = 0;
322    
323     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
324    
325 gezelter 1769 for (sd = mol->beginIntegrableObject(j); sd != NULL;
326     sd = mol->nextIntegrableObject(j)) {
327    
328 gezelter 507 ndfRaw_local += 3;
329 gezelter 246
330 gezelter 1769 if (sd->isDirectional()) {
331     if (sd->isLinear()) {
332 gezelter 507 ndfRaw_local += 2;
333     } else {
334     ndfRaw_local += 3;
335     }
336     }
337 gezelter 246
338 gezelter 507 }
339 gezelter 246 }
340    
341     #ifdef IS_MPI
342     MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
343     #else
344     ndfRaw_ = ndfRaw_local;
345     #endif
346 gezelter 507 }
347 gezelter 2
348 gezelter 507 void SimInfo::calcNdfTrans() {
349 gezelter 246 int ndfTrans_local;
350 gezelter 2
351 gezelter 246 ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
352 gezelter 2
353    
354 gezelter 246 #ifdef IS_MPI
355     MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
356     #else
357     ndfTrans_ = ndfTrans_local;
358     #endif
359 gezelter 2
360 gezelter 246 ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
361    
362 gezelter 507 }
363 gezelter 2
364 gezelter 1287 void SimInfo::addInteractionPairs(Molecule* mol) {
365     ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
366 gezelter 1528 vector<Bond*>::iterator bondIter;
367     vector<Bend*>::iterator bendIter;
368     vector<Torsion*>::iterator torsionIter;
369     vector<Inversion*>::iterator inversionIter;
370 gezelter 246 Bond* bond;
371     Bend* bend;
372     Torsion* torsion;
373 gezelter 1277 Inversion* inversion;
374 gezelter 246 int a;
375     int b;
376     int c;
377     int d;
378 tim 749
379 gezelter 1287 // atomGroups can be used to add special interaction maps between
380     // groups of atoms that are in two separate rigid bodies.
381     // However, most site-site interactions between two rigid bodies
382     // are probably not special, just the ones between the physically
383     // bonded atoms. Interactions *within* a single rigid body should
384     // always be excluded. These are done at the bottom of this
385     // function.
386    
387 gezelter 1528 map<int, set<int> > atomGroups;
388 tim 749 Molecule::RigidBodyIterator rbIter;
389     RigidBody* rb;
390     Molecule::IntegrableObjectIterator ii;
391 gezelter 1769 StuntDouble* sd;
392 gezelter 246
393 gezelter 1769 for (sd = mol->beginIntegrableObject(ii); sd != NULL;
394     sd = mol->nextIntegrableObject(ii)) {
395 gezelter 1287
396 gezelter 1769 if (sd->isRigidBody()) {
397     rb = static_cast<RigidBody*>(sd);
398 gezelter 1528 vector<Atom*> atoms = rb->getAtoms();
399     set<int> rigidAtoms;
400 gezelter 1287 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
401     rigidAtoms.insert(atoms[i]->getGlobalIndex());
402     }
403     for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
404 gezelter 1528 atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
405 gezelter 1287 }
406 tim 749 } else {
407 gezelter 1528 set<int> oneAtomSet;
408 gezelter 1769 oneAtomSet.insert(sd->getGlobalIndex());
409     atomGroups.insert(map<int, set<int> >::value_type(sd->getGlobalIndex(), oneAtomSet));
410 tim 749 }
411     }
412 gezelter 1287
413     for (bond= mol->beginBond(bondIter); bond != NULL;
414     bond = mol->nextBond(bondIter)) {
415 tim 749
416 gezelter 1287 a = bond->getAtomA()->getGlobalIndex();
417     b = bond->getAtomB()->getGlobalIndex();
418 tim 749
419 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
420     oneTwoInteractions_.addPair(a, b);
421     } else {
422     excludedInteractions_.addPair(a, b);
423     }
424 gezelter 246 }
425 gezelter 2
426 gezelter 1287 for (bend= mol->beginBend(bendIter); bend != NULL;
427     bend = mol->nextBend(bendIter)) {
428    
429 gezelter 507 a = bend->getAtomA()->getGlobalIndex();
430     b = bend->getAtomB()->getGlobalIndex();
431     c = bend->getAtomC()->getGlobalIndex();
432 gezelter 1287
433     if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
434     oneTwoInteractions_.addPair(a, b);
435     oneTwoInteractions_.addPair(b, c);
436     } else {
437     excludedInteractions_.addPair(a, b);
438     excludedInteractions_.addPair(b, c);
439     }
440 gezelter 2
441 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
442     oneThreeInteractions_.addPair(a, c);
443     } else {
444     excludedInteractions_.addPair(a, c);
445     }
446 gezelter 246 }
447 gezelter 2
448 gezelter 1287 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
449     torsion = mol->nextTorsion(torsionIter)) {
450    
451 gezelter 507 a = torsion->getAtomA()->getGlobalIndex();
452     b = torsion->getAtomB()->getGlobalIndex();
453     c = torsion->getAtomC()->getGlobalIndex();
454 gezelter 1287 d = torsion->getAtomD()->getGlobalIndex();
455 cli2 1290
456 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
457     oneTwoInteractions_.addPair(a, b);
458     oneTwoInteractions_.addPair(b, c);
459     oneTwoInteractions_.addPair(c, d);
460     } else {
461     excludedInteractions_.addPair(a, b);
462     excludedInteractions_.addPair(b, c);
463     excludedInteractions_.addPair(c, d);
464     }
465 gezelter 2
466 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
467     oneThreeInteractions_.addPair(a, c);
468     oneThreeInteractions_.addPair(b, d);
469     } else {
470     excludedInteractions_.addPair(a, c);
471     excludedInteractions_.addPair(b, d);
472     }
473 tim 749
474 gezelter 1287 if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
475     oneFourInteractions_.addPair(a, d);
476     } else {
477     excludedInteractions_.addPair(a, d);
478     }
479 gezelter 2 }
480    
481 gezelter 1277 for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
482     inversion = mol->nextInversion(inversionIter)) {
483 gezelter 1287
484 gezelter 1277 a = inversion->getAtomA()->getGlobalIndex();
485     b = inversion->getAtomB()->getGlobalIndex();
486     c = inversion->getAtomC()->getGlobalIndex();
487     d = inversion->getAtomD()->getGlobalIndex();
488    
489 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
490     oneTwoInteractions_.addPair(a, b);
491     oneTwoInteractions_.addPair(a, c);
492     oneTwoInteractions_.addPair(a, d);
493     } else {
494     excludedInteractions_.addPair(a, b);
495     excludedInteractions_.addPair(a, c);
496     excludedInteractions_.addPair(a, d);
497     }
498 gezelter 1277
499 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
500     oneThreeInteractions_.addPair(b, c);
501     oneThreeInteractions_.addPair(b, d);
502     oneThreeInteractions_.addPair(c, d);
503     } else {
504     excludedInteractions_.addPair(b, c);
505     excludedInteractions_.addPair(b, d);
506     excludedInteractions_.addPair(c, d);
507     }
508 gezelter 1277 }
509    
510 gezelter 1287 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
511     rb = mol->nextRigidBody(rbIter)) {
512 gezelter 1528 vector<Atom*> atoms = rb->getAtoms();
513 gezelter 1287 for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
514     for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
515 gezelter 507 a = atoms[i]->getGlobalIndex();
516     b = atoms[j]->getGlobalIndex();
517 gezelter 1287 excludedInteractions_.addPair(a, b);
518 gezelter 507 }
519     }
520 tim 430 }
521    
522 gezelter 507 }
523 gezelter 246
524 gezelter 1287 void SimInfo::removeInteractionPairs(Molecule* mol) {
525     ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
526 gezelter 1528 vector<Bond*>::iterator bondIter;
527     vector<Bend*>::iterator bendIter;
528     vector<Torsion*>::iterator torsionIter;
529     vector<Inversion*>::iterator inversionIter;
530 gezelter 246 Bond* bond;
531     Bend* bend;
532     Torsion* torsion;
533 gezelter 1277 Inversion* inversion;
534 gezelter 246 int a;
535     int b;
536     int c;
537     int d;
538 tim 749
539 gezelter 1528 map<int, set<int> > atomGroups;
540 tim 749 Molecule::RigidBodyIterator rbIter;
541     RigidBody* rb;
542     Molecule::IntegrableObjectIterator ii;
543 gezelter 1769 StuntDouble* sd;
544 gezelter 246
545 gezelter 1769 for (sd = mol->beginIntegrableObject(ii); sd != NULL;
546     sd = mol->nextIntegrableObject(ii)) {
547 gezelter 1287
548 gezelter 1769 if (sd->isRigidBody()) {
549     rb = static_cast<RigidBody*>(sd);
550 gezelter 1528 vector<Atom*> atoms = rb->getAtoms();
551     set<int> rigidAtoms;
552 gezelter 1287 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
553     rigidAtoms.insert(atoms[i]->getGlobalIndex());
554     }
555     for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
556 gezelter 1528 atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
557 gezelter 1287 }
558 tim 749 } else {
559 gezelter 1528 set<int> oneAtomSet;
560 gezelter 1769 oneAtomSet.insert(sd->getGlobalIndex());
561     atomGroups.insert(map<int, set<int> >::value_type(sd->getGlobalIndex(), oneAtomSet));
562 tim 749 }
563     }
564    
565 gezelter 1287 for (bond= mol->beginBond(bondIter); bond != NULL;
566     bond = mol->nextBond(bondIter)) {
567    
568     a = bond->getAtomA()->getGlobalIndex();
569     b = bond->getAtomB()->getGlobalIndex();
570 tim 749
571 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
572     oneTwoInteractions_.removePair(a, b);
573     } else {
574     excludedInteractions_.removePair(a, b);
575     }
576 gezelter 2 }
577 gezelter 246
578 gezelter 1287 for (bend= mol->beginBend(bendIter); bend != NULL;
579     bend = mol->nextBend(bendIter)) {
580    
581 gezelter 507 a = bend->getAtomA()->getGlobalIndex();
582     b = bend->getAtomB()->getGlobalIndex();
583     c = bend->getAtomC()->getGlobalIndex();
584 gezelter 1287
585     if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
586     oneTwoInteractions_.removePair(a, b);
587     oneTwoInteractions_.removePair(b, c);
588     } else {
589     excludedInteractions_.removePair(a, b);
590     excludedInteractions_.removePair(b, c);
591     }
592 gezelter 246
593 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
594     oneThreeInteractions_.removePair(a, c);
595     } else {
596     excludedInteractions_.removePair(a, c);
597     }
598 gezelter 2 }
599 gezelter 246
600 gezelter 1287 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
601     torsion = mol->nextTorsion(torsionIter)) {
602    
603 gezelter 507 a = torsion->getAtomA()->getGlobalIndex();
604     b = torsion->getAtomB()->getGlobalIndex();
605     c = torsion->getAtomC()->getGlobalIndex();
606 gezelter 1287 d = torsion->getAtomD()->getGlobalIndex();
607    
608     if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
609     oneTwoInteractions_.removePair(a, b);
610     oneTwoInteractions_.removePair(b, c);
611     oneTwoInteractions_.removePair(c, d);
612     } else {
613     excludedInteractions_.removePair(a, b);
614     excludedInteractions_.removePair(b, c);
615     excludedInteractions_.removePair(c, d);
616     }
617 gezelter 246
618 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
619     oneThreeInteractions_.removePair(a, c);
620     oneThreeInteractions_.removePair(b, d);
621     } else {
622     excludedInteractions_.removePair(a, c);
623     excludedInteractions_.removePair(b, d);
624     }
625 tim 749
626 gezelter 1287 if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
627     oneFourInteractions_.removePair(a, d);
628     } else {
629     excludedInteractions_.removePair(a, d);
630     }
631     }
632 tim 749
633 gezelter 1287 for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
634     inversion = mol->nextInversion(inversionIter)) {
635 tim 749
636 gezelter 1277 a = inversion->getAtomA()->getGlobalIndex();
637     b = inversion->getAtomB()->getGlobalIndex();
638     c = inversion->getAtomC()->getGlobalIndex();
639     d = inversion->getAtomD()->getGlobalIndex();
640    
641 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
642     oneTwoInteractions_.removePair(a, b);
643     oneTwoInteractions_.removePair(a, c);
644     oneTwoInteractions_.removePair(a, d);
645     } else {
646     excludedInteractions_.removePair(a, b);
647     excludedInteractions_.removePair(a, c);
648     excludedInteractions_.removePair(a, d);
649     }
650 gezelter 1277
651 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
652     oneThreeInteractions_.removePair(b, c);
653     oneThreeInteractions_.removePair(b, d);
654     oneThreeInteractions_.removePair(c, d);
655     } else {
656     excludedInteractions_.removePair(b, c);
657     excludedInteractions_.removePair(b, d);
658     excludedInteractions_.removePair(c, d);
659     }
660 gezelter 1277 }
661    
662 gezelter 1287 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
663     rb = mol->nextRigidBody(rbIter)) {
664 gezelter 1528 vector<Atom*> atoms = rb->getAtoms();
665 gezelter 1287 for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
666     for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
667 gezelter 507 a = atoms[i]->getGlobalIndex();
668     b = atoms[j]->getGlobalIndex();
669 gezelter 1287 excludedInteractions_.removePair(a, b);
670 gezelter 507 }
671     }
672 tim 430 }
673 gezelter 1287
674 gezelter 507 }
675 gezelter 1287
676    
677 gezelter 507 void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) {
678 gezelter 246 int curStampId;
679 gezelter 1287
680 gezelter 246 //index from 0
681     curStampId = moleculeStamps_.size();
682 gezelter 2
683 gezelter 246 moleculeStamps_.push_back(molStamp);
684     molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
685 gezelter 507 }
686 gezelter 2
687 gezelter 1530
688     /**
689     * update
690     *
691 gezelter 1535 * Performs the global checks and variable settings after the
692     * objects have been created.
693 gezelter 1530 *
694     */
695 gezelter 1535 void SimInfo::update() {
696 gezelter 1530 setupSimVariables();
697 gezelter 246 calcNdf();
698     calcNdfRaw();
699     calcNdfTrans();
700 gezelter 507 }
701 gezelter 1528
702 gezelter 1535 /**
703     * getSimulatedAtomTypes
704     *
705     * Returns an STL set of AtomType* that are actually present in this
706     * simulation. Must query all processors to assemble this information.
707     *
708     */
709 gezelter 1528 set<AtomType*> SimInfo::getSimulatedAtomTypes() {
710 gezelter 246 SimInfo::MoleculeIterator mi;
711     Molecule* mol;
712     Molecule::AtomIterator ai;
713     Atom* atom;
714 gezelter 1528 set<AtomType*> atomTypes;
715    
716 gezelter 1588 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
717     for(atom = mol->beginAtom(ai); atom != NULL;
718     atom = mol->nextAtom(ai)) {
719 gezelter 507 atomTypes.insert(atom->getAtomType());
720 gezelter 1529 }
721     }
722 gezelter 1588
723 gezelter 1535 #ifdef IS_MPI
724    
725     // loop over the found atom types on this processor, and add their
726     // numerical idents to a vector:
727 gezelter 1588
728 gezelter 1535 vector<int> foundTypes;
729     set<AtomType*>::iterator i;
730     for (i = atomTypes.begin(); i != atomTypes.end(); ++i)
731     foundTypes.push_back( (*i)->getIdent() );
732    
733     // count_local holds the number of found types on this processor
734     int count_local = foundTypes.size();
735    
736 gezelter 1588 int nproc = MPI::COMM_WORLD.Get_size();
737 gezelter 1535
738 gezelter 1588 // we need arrays to hold the counts and displacement vectors for
739     // all processors
740     vector<int> counts(nproc, 0);
741     vector<int> disps(nproc, 0);
742 gezelter 1535
743 gezelter 1588 // fill the counts array
744     MPI::COMM_WORLD.Allgather(&count_local, 1, MPI::INT, &counts[0],
745     1, MPI::INT);
746    
747     // use the processor counts to compute the displacement array
748     disps[0] = 0;
749     int totalCount = counts[0];
750     for (int iproc = 1; iproc < nproc; iproc++) {
751     disps[iproc] = disps[iproc-1] + counts[iproc-1];
752     totalCount += counts[iproc];
753     }
754 gezelter 1535
755 gezelter 1588 // we need a (possibly redundant) set of all found types:
756     vector<int> ftGlobal(totalCount);
757 gezelter 2
758 gezelter 1588 // now spray out the foundTypes to all the other processors:
759 gezelter 1535 MPI::COMM_WORLD.Allgatherv(&foundTypes[0], count_local, MPI::INT,
760 gezelter 1588 &ftGlobal[0], &counts[0], &disps[0],
761     MPI::INT);
762 gezelter 1126
763 gezelter 1588 vector<int>::iterator j;
764    
765 gezelter 1535 // foundIdents is a stl set, so inserting an already found ident
766     // will have no effect.
767     set<int> foundIdents;
768 gezelter 1588
769 gezelter 1535 for (j = ftGlobal.begin(); j != ftGlobal.end(); ++j)
770     foundIdents.insert((*j));
771 gezelter 1528
772 gezelter 1535 // now iterate over the foundIdents and get the actual atom types
773     // that correspond to these:
774     set<int>::iterator it;
775 gezelter 1588 for (it = foundIdents.begin(); it != foundIdents.end(); ++it)
776 gezelter 1535 atomTypes.insert( forceField_->getAtomType((*it)) );
777    
778     #endif
779 gezelter 1588
780 gezelter 1535 return atomTypes;
781 gezelter 1528 }
782 chrisfen 611
783 gezelter 1534 void SimInfo::setupSimVariables() {
784     useAtomicVirial_ = simParams_->getUseAtomicVirial();
785 gezelter 1779 // we only call setAccumulateBoxDipole if the accumulateBoxDipole
786     // parameter is true
787 gezelter 1534 calcBoxDipole_ = false;
788     if ( simParams_->haveAccumulateBoxDipole() )
789     if ( simParams_->getAccumulateBoxDipole() ) {
790     calcBoxDipole_ = true;
791     }
792 gezelter 1583
793 gezelter 1528 set<AtomType*>::iterator i;
794     set<AtomType*> atomTypes;
795 gezelter 1534 atomTypes = getSimulatedAtomTypes();
796 gezelter 1767 bool usesElectrostatic = false;
797     bool usesMetallic = false;
798     bool usesDirectional = false;
799     bool usesFluctuatingCharges = false;
800 gezelter 246 //loop over all of the atom types
801     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
802 gezelter 1528 usesElectrostatic |= (*i)->isElectrostatic();
803     usesMetallic |= (*i)->isMetal();
804     usesDirectional |= (*i)->isDirectional();
805 gezelter 1715 usesFluctuatingCharges |= (*i)->isFluctuatingCharge();
806 gezelter 246 }
807 gezelter 1760
808 gezelter 1767 #ifdef IS_MPI
809     bool temp;
810 gezelter 1528 temp = usesDirectional;
811 gezelter 1767 MPI::COMM_WORLD.Allreduce(&temp, &usesDirectionalAtoms_, 1, MPI::BOOL,
812     MPI::LOR);
813    
814 gezelter 1528 temp = usesMetallic;
815 gezelter 1767 MPI::COMM_WORLD.Allreduce(&temp, &usesMetallicAtoms_, 1, MPI::BOOL,
816     MPI::LOR);
817 gezelter 1583
818 gezelter 1528 temp = usesElectrostatic;
819 gezelter 1767 MPI::COMM_WORLD.Allreduce(&temp, &usesElectrostaticAtoms_, 1, MPI::BOOL,
820     MPI::LOR);
821 gezelter 1715
822     temp = usesFluctuatingCharges;
823 gezelter 1767 MPI::COMM_WORLD.Allreduce(&temp, &usesFluctuatingCharges_, 1, MPI::BOOL,
824     MPI::LOR);
825 gezelter 1583 #else
826    
827     usesDirectionalAtoms_ = usesDirectional;
828     usesMetallicAtoms_ = usesMetallic;
829     usesElectrostaticAtoms_ = usesElectrostatic;
830 gezelter 1715 usesFluctuatingCharges_ = usesFluctuatingCharges;
831 gezelter 1583
832 gezelter 2 #endif
833 gezelter 1583
834     requiresPrepair_ = usesMetallicAtoms_ ? true : false;
835     requiresSkipCorrection_ = usesElectrostaticAtoms_ ? true : false;
836     requiresSelfCorrection_ = usesElectrostaticAtoms_ ? true : false;
837 gezelter 507 }
838 gezelter 2
839 gezelter 1549
840     vector<int> SimInfo::getGlobalAtomIndices() {
841     SimInfo::MoleculeIterator mi;
842     Molecule* mol;
843     Molecule::AtomIterator ai;
844     Atom* atom;
845    
846     vector<int> GlobalAtomIndices(getNAtoms(), 0);
847    
848     for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
849    
850     for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
851     GlobalAtomIndices[atom->getLocalIndex()] = atom->getGlobalIndex();
852     }
853     }
854     return GlobalAtomIndices;
855     }
856    
857    
858     vector<int> SimInfo::getGlobalGroupIndices() {
859     SimInfo::MoleculeIterator mi;
860     Molecule* mol;
861     Molecule::CutoffGroupIterator ci;
862     CutoffGroup* cg;
863    
864     vector<int> GlobalGroupIndices;
865    
866     for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
867    
868     //local index of cutoff group is trivial, it only depends on the
869     //order of travesing
870     for (cg = mol->beginCutoffGroup(ci); cg != NULL;
871     cg = mol->nextCutoffGroup(ci)) {
872     GlobalGroupIndices.push_back(cg->getGlobalIndex());
873     }
874     }
875     return GlobalGroupIndices;
876     }
877    
878    
879 gezelter 1569 void SimInfo::prepareTopology() {
880 gezelter 1287 int nExclude, nOneTwo, nOneThree, nOneFour;
881 gezelter 2
882 gezelter 246 //calculate mass ratio of cutoff group
883     SimInfo::MoleculeIterator mi;
884     Molecule* mol;
885     Molecule::CutoffGroupIterator ci;
886     CutoffGroup* cg;
887     Molecule::AtomIterator ai;
888     Atom* atom;
889 tim 963 RealType totalMass;
890 gezelter 246
891 gezelter 1581 /**
892     * The mass factor is the relative mass of an atom to the total
893     * mass of the cutoff group it belongs to. By default, all atoms
894     * are their own cutoff groups, and therefore have mass factors of
895     * 1. We need some special handling for massless atoms, which
896     * will be treated as carrying the entire mass of the cutoff
897     * group.
898     */
899 gezelter 1569 massFactors_.clear();
900 gezelter 1581 massFactors_.resize(getNAtoms(), 1.0);
901 gezelter 2
902 gezelter 246 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
903 gezelter 1569 for (cg = mol->beginCutoffGroup(ci); cg != NULL;
904     cg = mol->nextCutoffGroup(ci)) {
905 gezelter 2
906 gezelter 507 totalMass = cg->getMass();
907     for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
908 chrisfen 645 // Check for massless groups - set mfact to 1 if true
909 gezelter 1581 if (totalMass != 0)
910     massFactors_[atom->getLocalIndex()] = atom->getMass()/totalMass;
911 chrisfen 645 else
912 gezelter 1581 massFactors_[atom->getLocalIndex()] = 1.0;
913 gezelter 507 }
914     }
915 gezelter 246 }
916 gezelter 2
917 gezelter 1544 // Build the identArray_
918 gezelter 2
919 gezelter 1544 identArray_.clear();
920     identArray_.reserve(getNAtoms());
921 gezelter 246 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
922 gezelter 507 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
923 gezelter 1544 identArray_.push_back(atom->getIdent());
924 gezelter 507 }
925 gezelter 246 }
926    
927 gezelter 1569 //scan topology
928 gezelter 1287
929     nExclude = excludedInteractions_.getSize();
930     nOneTwo = oneTwoInteractions_.getSize();
931     nOneThree = oneThreeInteractions_.getSize();
932     nOneFour = oneFourInteractions_.getSize();
933    
934     int* excludeList = excludedInteractions_.getPairList();
935     int* oneTwoList = oneTwoInteractions_.getPairList();
936     int* oneThreeList = oneThreeInteractions_.getPairList();
937     int* oneFourList = oneFourInteractions_.getPairList();
938    
939 gezelter 1569 topologyDone_ = true;
940 gezelter 507 }
941 chrisfen 143
942 gezelter 507 void SimInfo::addProperty(GenericData* genData) {
943 gezelter 246 properties_.addProperty(genData);
944 gezelter 507 }
945 gezelter 2
946 gezelter 1528 void SimInfo::removeProperty(const string& propName) {
947 gezelter 246 properties_.removeProperty(propName);
948 gezelter 507 }
949 gezelter 2
950 gezelter 507 void SimInfo::clearProperties() {
951 gezelter 246 properties_.clearProperties();
952 gezelter 507 }
953 gezelter 2
954 gezelter 1528 vector<string> SimInfo::getPropertyNames() {
955 gezelter 246 return properties_.getPropertyNames();
956 gezelter 507 }
957 gezelter 246
958 gezelter 1528 vector<GenericData*> SimInfo::getProperties() {
959 gezelter 246 return properties_.getProperties();
960 gezelter 507 }
961 gezelter 2
962 gezelter 1528 GenericData* SimInfo::getPropertyByName(const string& propName) {
963 gezelter 246 return properties_.getPropertyByName(propName);
964 gezelter 507 }
965 gezelter 2
966 gezelter 507 void SimInfo::setSnapshotManager(SnapshotManager* sman) {
967 tim 432 if (sman_ == sman) {
968 gezelter 507 return;
969 tim 432 }
970     delete sman_;
971 gezelter 246 sman_ = sman;
972 gezelter 2
973 gezelter 246 Molecule* mol;
974     RigidBody* rb;
975     Atom* atom;
976 gezelter 1540 CutoffGroup* cg;
977 gezelter 246 SimInfo::MoleculeIterator mi;
978     Molecule::RigidBodyIterator rbIter;
979 gezelter 1540 Molecule::AtomIterator atomIter;
980     Molecule::CutoffGroupIterator cgIter;
981 gezelter 246
982     for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
983    
984 gezelter 1779 for (atom = mol->beginAtom(atomIter); atom != NULL;
985     atom = mol->nextAtom(atomIter)) {
986 gezelter 507 atom->setSnapshotManager(sman_);
987     }
988 gezelter 246
989 gezelter 1779 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
990     rb = mol->nextRigidBody(rbIter)) {
991 gezelter 507 rb->setSnapshotManager(sman_);
992     }
993 gezelter 1540
994 gezelter 1779 for (cg = mol->beginCutoffGroup(cgIter); cg != NULL;
995     cg = mol->nextCutoffGroup(cgIter)) {
996 gezelter 1540 cg->setSnapshotManager(sman_);
997     }
998 gezelter 246 }
999 gezelter 2
1000 gezelter 507 }
1001 gezelter 2
1002    
1003 gezelter 1528 ostream& operator <<(ostream& o, SimInfo& info) {
1004 gezelter 246
1005     return o;
1006 gezelter 507 }
1007 chuckv 555
1008 gezelter 1764
1009 tim 1024 StuntDouble* SimInfo::getIOIndexToIntegrableObject(int index) {
1010 gezelter 1779 if (index >= IOIndexToIntegrableObject.size()) {
1011     sprintf(painCave.errMsg,
1012     "SimInfo::getIOIndexToIntegrableObject Error: Integrable Object\n"
1013     "\tindex exceeds number of known objects!\n");
1014     painCave.isFatal = 1;
1015     simError();
1016     return NULL;
1017     } else
1018     return IOIndexToIntegrableObject.at(index);
1019 tim 1024 }
1020    
1021 gezelter 1528 void SimInfo::setIOIndexToIntegrableObject(const vector<StuntDouble*>& v) {
1022 tim 1024 IOIndexToIntegrableObject= v;
1023     }
1024    
1025 gezelter 1528 int SimInfo::getNGlobalConstraints() {
1026     int nGlobalConstraints;
1027     #ifdef IS_MPI
1028     MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM,
1029     MPI_COMM_WORLD);
1030     #else
1031     nGlobalConstraints = nConstraints_;
1032     #endif
1033     return nGlobalConstraints;
1034     }
1035    
1036 gezelter 1390 }//end namespace OpenMD
1037 gezelter 246

Properties

Name Value
svn:keywords Author Id Revision Date