ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimInfo.cpp
Revision: 1930
Committed: Mon Aug 19 13:51:04 2013 UTC (11 years, 8 months ago) by gezelter
File size: 33121 byte(s)
Log Message:
region fixes, performance boosts

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

Properties

Name Value
svn:keywords Author Id Revision Date