ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimInfo.cpp
Revision: 1940
Committed: Fri Nov 1 19:31:41 2013 UTC (11 years, 6 months ago) by gezelter
File size: 33067 byte(s)
Log Message:
Angular velocity fix in RNEMD

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date