ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimInfo.cpp
Revision: 2071
Committed: Sat Mar 7 21:41:51 2015 UTC (10 years, 1 month ago) by gezelter
File size: 34682 byte(s)
Log Message:
Reducing the number of warnings when using g++ to compile.

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

Properties

Name Value
svn:keywords Author Id Revision Date