ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1769
Committed: Mon Jul 9 14:15:52 2012 UTC (12 years, 9 months ago) by gezelter
File size: 32747 byte(s)
Log Message:
Code readability updates.

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

Properties

Name Value
svn:keywords Author Id Revision Date