ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimInfo.cpp
Revision: 1953
Committed: Thu Dec 5 18:19:26 2013 UTC (11 years, 4 months ago) by gezelter
File size: 34331 byte(s)
Log Message:
Rewrote much of selection module, added a bond correlation function

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

Properties

Name Value
svn:keywords Author Id Revision Date