ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1767
Committed: Fri Jul 6 22:01:58 2012 UTC (12 years, 9 months ago) by gezelter
File size: 33156 byte(s)
Log Message:
Various fixes required to compile OpenMD with the MS Visual C++ compiler

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date