ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1798
Committed: Thu Sep 13 14:10:11 2012 UTC (12 years, 7 months ago) by gezelter
File size: 32843 byte(s)
Log Message:
Merged trunk changes into the development branch

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

Properties

Name Value
svn:keywords Author Id Revision Date