ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimInfo.cpp
Revision: 1987
Committed: Thu Apr 17 19:07:31 2014 UTC (11 years ago) by gezelter
File size: 34358 byte(s)
Log Message:
Removing vestiges of deprecated MPI:: namespace C++ MPI bindings

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date