ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimInfo.cpp
(Generate patch)

Comparing trunk/src/brains/SimInfo.cpp (file contents):
Revision 1386 by gezelter, Fri Oct 23 18:41:09 2009 UTC vs.
Revision 1987 by gezelter, Thu Apr 17 19:07:31 2014 UTC

# Line 6 | Line 6
6   * redistribute this software in source and binary code form, provided
7   * that the following conditions are met:
8   *
9 < * 1. Acknowledgement of the program authors must be made in any
10 < *    publication of scientific results based in part on use of the
11 < *    program.  An acceptable form of acknowledgement is citation of
12 < *    the article in which the program was described (Matthew
13 < *    A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14 < *    J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15 < *    Parallel Simulation Engine for Molecular Dynamics,"
16 < *    J. Comput. Chem. 26, pp. 252-271 (2005))
17 < *
18 < * 2. Redistributions of source code must retain the above copyright
9 > * 1. Redistributions of source code must retain the above copyright
10   *    notice, this list of conditions and the following disclaimer.
11   *
12 < * 3. Redistributions in binary form must reproduce the above copyright
12 > * 2. Redistributions in binary form must reproduce the above copyright
13   *    notice, this list of conditions and the following disclaimer in the
14   *    documentation and/or other materials provided with the
15   *    distribution.
# Line 37 | Line 28
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 + *
32 + * SUPPORT OPEN SCIENCE!  If you use OpenMD or its source code in your
33 + * research, please cite the appropriate papers when you publish your
34 + * work.  Good starting points are:
35 + *                                                                      
36 + * [1]  Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).            
37 + * [2]  Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).          
38 + * [3]  Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).          
39 + * [4]  Kuang & Gezelter,  J. Chem. Phys. 133, 164101 (2010).
40 + * [5]  Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41   */
42  
43   /**
# Line 46 | Line 47
47   * @version 1.0
48   */
49  
50 + #ifdef IS_MPI
51 + #include <mpi.h>
52 + #endif
53   #include <algorithm>
54   #include <set>
55   #include <map>
# Line 54 | Line 58
58   #include "math/Vector3.hpp"
59   #include "primitives/Molecule.hpp"
60   #include "primitives/StuntDouble.hpp"
57 #include "UseTheForce/fCutoffPolicy.h"
58 #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
59 #include "UseTheForce/DarkSide/fElectrostaticScreeningMethod.h"
60 #include "UseTheForce/DarkSide/fSwitchingFunctionType.h"
61 #include "UseTheForce/doForces_interface.h"
62 #include "UseTheForce/DarkSide/neighborLists_interface.h"
63 #include "UseTheForce/DarkSide/electrostatic_interface.h"
64 #include "UseTheForce/DarkSide/switcheroo_interface.h"
61   #include "utils/MemoryUtils.hpp"
62   #include "utils/simError.h"
63   #include "selection/SelectionManager.hpp"
64   #include "io/ForceFieldOptions.hpp"
65 < #include "UseTheForce/ForceField.hpp"
65 > #include "brains/ForceField.hpp"
66 > #include "nonbonded/SwitchingFunction.hpp"
67  
68 <
69 < #ifdef IS_MPI
73 < #include "UseTheForce/mpiComponentPlan.h"
74 < #include "UseTheForce/DarkSide/simParallel_interface.h"
75 < #endif
76 <
77 < namespace oopse {
78 <  std::set<int> getRigidSet(int index, std::map<int, std::set<int> >& container) {
79 <    std::map<int, std::set<int> >::iterator i = container.find(index);
80 <    std::set<int> result;
81 <    if (i != container.end()) {
82 <        result = i->second;
83 <    }
84 <
85 <    return result;
86 <  }
68 > using namespace std;
69 > namespace OpenMD {
70    
71    SimInfo::SimInfo(ForceField* ff, Globals* simParams) :
72      forceField_(ff), simParams_(simParams),
73      ndf_(0), fdf_local(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0),
74      nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0),
75 <    nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
76 <    nAtoms_(0), nBonds_(0),  nBends_(0), nTorsions_(0), nInversions_(0),
75 >    nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
76 >    nGlobalFluctuatingCharges_(0), nGlobalBonds_(0), nGlobalBends_(0),
77 >    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), sman_(NULL), fortranInitialized_(false),
81 <    calcBoxDipole_(false), useAtomicVirial_(true) {
80 >    nConstraints_(0), nFluctuatingCharges_(0), sman_(NULL),
81 >    topologyDone_(false), calcBoxDipole_(false), useAtomicVirial_(true),
82 >    hasNGlobalConstraints_(false) {    
83 >    
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 >      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 <
99 <      MoleculeStamp* molStamp;
100 <      int nMolWithSameStamp;
101 <      int nCutoffAtoms = 0; // number of atoms belong to cutoff groups
102 <      int nGroups = 0;      //total cutoff groups defined in meta-data file
103 <      CutoffGroupStamp* cgStamp;    
104 <      RigidBodyStamp* rbStamp;
105 <      int nRigidAtoms = 0;
106 <
107 <      std::vector<Component*> components = simParams->getComponents();
104 >      nMolWithSameStamp = (*i)->getNMol();
105        
106 <      for (std::vector<Component*>::iterator i = components.begin(); i !=components.end(); ++i) {
107 <        molStamp = (*i)->getMoleculeStamp();
108 <        nMolWithSameStamp = (*i)->getNMol();
109 <        
110 <        addMoleculeStamp(molStamp, nMolWithSameStamp);
111 <
112 <        //calculate atoms in molecules
113 <        nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp;  
114 <
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();
125 <        }
126 <
127 <        nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
128 <
129 <        nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;            
130 <
131 <        //calculate atoms in rigid bodies
132 <        int nAtomsInRigidBodies = 0;
133 <        int nRigidBodiesInStamp = molStamp->getNRigidBodies();
134 <        
135 <        for (int j=0; j < nRigidBodiesInStamp; j++) {
136 <          rbStamp = molStamp->getRigidBodyStamp(j);
137 <          nAtomsInRigidBodies += rbStamp->getNMembers();
138 <        }
139 <
140 <        nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
141 <        nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;            
142 <        
106 >      addMoleculeStamp(molStamp, nMolWithSameStamp);
107 >      
108 >      //calculate atoms in molecules
109 >      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 >      
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        }
123 <
124 <      //every free atom (atom does not belong to cutoff groups) is a cutoff
125 <      //group therefore the total number of cutoff groups in the system is
126 <      //equal to the total number of atoms minus number of atoms belong to
127 <      //cutoff group defined in meta-data file plus the number of cutoff
128 <      //groups defined in meta-data file
129 <      nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
130 <
131 <      //every free atom (atom does not belong to rigid bodies) is an
132 <      //integrable object therefore the total number of integrable objects
133 <      //in the system is equal to the total number of atoms minus number of
134 <      //atoms belong to rigid body defined in meta-data file plus the number
135 <      //of rigid bodies defined in meta-data file
136 <      nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms
137 <                                                + nGlobalRigidBodies_;
138 <  
139 <      nGlobalMols_ = molStampIds_.size();
161 <      molToProcMap_.resize(nGlobalMols_);
123 >      
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  
148 +    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 +  
162    SimInfo::~SimInfo() {
163 <    std::map<int, Molecule*>::iterator i;
163 >    map<int, Molecule*>::iterator i;
164      for (i = molecules_.begin(); i != molecules_.end(); ++i) {
165        delete i->second;
166      }
# Line 173 | Line 171 | namespace oopse {
171      delete forceField_;
172    }
173  
176  int SimInfo::getNGlobalConstraints() {
177    int nGlobalConstraints;
178 #ifdef IS_MPI
179    MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM,
180                  MPI_COMM_WORLD);    
181 #else
182    nGlobalConstraints =  nConstraints_;
183 #endif
184    return nGlobalConstraints;
185  }
174  
175    bool SimInfo::addMolecule(Molecule* mol) {
176      MoleculeIterator i;
177 <
177 >    
178      i = molecules_.find(mol->getGlobalIndex());
179      if (i == molecules_.end() ) {
180 <
181 <      molecules_.insert(std::make_pair(mol->getGlobalIndex(), mol));
182 <        
180 >      
181 >      molecules_.insert(make_pair(mol->getGlobalIndex(), mol));
182 >      
183        nAtoms_ += mol->getNAtoms();
184        nBonds_ += mol->getNBonds();
185        nBends_ += mol->getNBends();
# Line 201 | Line 189 | namespace oopse {
189        nIntegrableObjects_ += mol->getNIntegrableObjects();
190        nCutoffGroups_ += mol->getNCutoffGroups();
191        nConstraints_ += mol->getNConstraintPairs();
192 <
192 >      
193        addInteractionPairs(mol);
194 <  
194 >      
195        return true;
196      } else {
197        return false;
198      }
199    }
200 <
200 >  
201    bool SimInfo::removeMolecule(Molecule* mol) {
202      MoleculeIterator i;
203      i = molecules_.find(mol->getGlobalIndex());
# Line 237 | Line 225 | namespace oopse {
225      } else {
226        return false;
227      }
240
241
228    }    
229  
230          
# Line 254 | Line 240 | namespace oopse {
240  
241  
242    void SimInfo::calcNdf() {
243 <    int ndf_local;
243 >    int ndf_local, nfq_local;
244      MoleculeIterator i;
245 <    std::vector<StuntDouble*>::iterator j;
245 >    vector<StuntDouble*>::iterator j;
246 >    vector<Atom*>::iterator k;
247 >
248      Molecule* mol;
249 <    StuntDouble* integrableObject;
249 >    StuntDouble* sd;
250 >    Atom* atom;
251  
252      ndf_local = 0;
253 +    nfq_local = 0;
254      
255      for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
266      for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
267           integrableObject = mol->nextIntegrableObject(j)) {
256  
257 +      for (sd = mol->beginIntegrableObject(j); sd != NULL;
258 +           sd = mol->nextIntegrableObject(j)) {
259 +
260          ndf_local += 3;
261  
262 <        if (integrableObject->isDirectional()) {
263 <          if (integrableObject->isLinear()) {
262 >        if (sd->isDirectional()) {
263 >          if (sd->isLinear()) {
264              ndf_local += 2;
265            } else {
266              ndf_local += 3;
267            }
268          }
278            
269        }
270 +
271 +      for (atom = mol->beginFluctuatingCharge(k); atom != NULL;
272 +           atom = mol->nextFluctuatingCharge(k)) {
273 +        if (atom->isFluctuatingCharge()) {
274 +          nfq_local++;
275 +        }
276 +      }
277      }
278      
279 +    ndfLocal_ = ndf_local;
280 +
281      // n_constraints is local, so subtract them on each processor
282      ndf_local -= nConstraints_;
283  
284   #ifdef IS_MPI
285 <    MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
285 >    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   #else
289      ndf_ = ndf_local;
290 +    nGlobalFluctuatingCharges_ = nfq_local;
291   #endif
292  
293      // nZconstraints_ is global, as are the 3 COM translations for the
# Line 296 | Line 298 | namespace oopse {
298  
299    int SimInfo::getFdf() {
300   #ifdef IS_MPI
301 <    MPI_Allreduce(&fdf_local,&fdf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
301 >    MPI_Allreduce(&fdf_local, &fdf_, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
302   #else
303      fdf_ = fdf_local;
304   #endif
305      return fdf_;
306    }
307 +  
308 +  unsigned int SimInfo::getNLocalCutoffGroups(){
309 +    int nLocalCutoffAtoms = 0;
310 +    Molecule* mol;
311 +    MoleculeIterator mi;
312 +    CutoffGroup* cg;
313 +    Molecule::CutoffGroupIterator ci;
314      
315 +    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    void SimInfo::calcNdfRaw() {
328      int ndfRaw_local;
329  
330      MoleculeIterator i;
331 <    std::vector<StuntDouble*>::iterator j;
331 >    vector<StuntDouble*>::iterator j;
332      Molecule* mol;
333 <    StuntDouble* integrableObject;
333 >    StuntDouble* sd;
334  
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)) {
318      for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
319           integrableObject = mol->nextIntegrableObject(j)) {
339  
340 +      for (sd = mol->beginIntegrableObject(j); sd != NULL;
341 +           sd = mol->nextIntegrableObject(j)) {
342 +
343          ndfRaw_local += 3;
344  
345 <        if (integrableObject->isDirectional()) {
346 <          if (integrableObject->isLinear()) {
345 >        if (sd->isDirectional()) {
346 >          if (sd->isLinear()) {
347              ndfRaw_local += 2;
348            } else {
349              ndfRaw_local += 3;
# Line 332 | Line 354 | namespace oopse {
354      }
355      
356   #ifdef IS_MPI
357 <    MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
357 >    MPI_Allreduce(&ndfRaw_local, &ndfRaw_, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
358   #else
359      ndfRaw_ = ndfRaw_local;
360   #endif
# Line 342 | Line 364 | namespace oopse {
364      int ndfTrans_local;
365  
366      ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
345
367  
368   #ifdef IS_MPI
369 <    MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
369 >    MPI_Allreduce(&ndfTrans_local, &ndfTrans_, 1, MPI_INT, MPI_SUM,
370 >                  MPI_COMM_WORLD);
371   #else
372      ndfTrans_ = ndfTrans_local;
373   #endif
374  
375      ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
354
376    }
377  
378    void SimInfo::addInteractionPairs(Molecule* mol) {
379      ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
380 <    std::vector<Bond*>::iterator bondIter;
381 <    std::vector<Bend*>::iterator bendIter;
382 <    std::vector<Torsion*>::iterator torsionIter;
383 <    std::vector<Inversion*>::iterator inversionIter;
380 >    vector<Bond*>::iterator bondIter;
381 >    vector<Bend*>::iterator bendIter;
382 >    vector<Torsion*>::iterator torsionIter;
383 >    vector<Inversion*>::iterator inversionIter;
384      Bond* bond;
385      Bend* bend;
386      Torsion* torsion;
# Line 377 | Line 398 | namespace oopse {
398      // always be excluded.  These are done at the bottom of this
399      // function.
400  
401 <    std::map<int, std::set<int> > atomGroups;
401 >    map<int, set<int> > atomGroups;
402      Molecule::RigidBodyIterator rbIter;
403      RigidBody* rb;
404      Molecule::IntegrableObjectIterator ii;
405 <    StuntDouble* integrableObject;
405 >    StuntDouble* sd;
406      
407 <    for (integrableObject = mol->beginIntegrableObject(ii);
408 <         integrableObject != NULL;
388 <         integrableObject = mol->nextIntegrableObject(ii)) {
407 >    for (sd = mol->beginIntegrableObject(ii); sd != NULL;
408 >         sd = mol->nextIntegrableObject(ii)) {
409        
410 <      if (integrableObject->isRigidBody()) {
411 <        rb = static_cast<RigidBody*>(integrableObject);
412 <        std::vector<Atom*> atoms = rb->getAtoms();
413 <        std::set<int> rigidAtoms;
410 >      if (sd->isRigidBody()) {
411 >        rb = static_cast<RigidBody*>(sd);
412 >        vector<Atom*> atoms = rb->getAtoms();
413 >        set<int> rigidAtoms;
414          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 <          atomGroups.insert(std::map<int, std::set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
418 >          atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
419          }      
420        } else {
421 <        std::set<int> oneAtomSet;
422 <        oneAtomSet.insert(integrableObject->getGlobalIndex());
423 <        atomGroups.insert(std::map<int, std::set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));        
421 >        set<int> oneAtomSet;
422 >        oneAtomSet.insert(sd->getGlobalIndex());
423 >        atomGroups.insert(map<int, set<int> >::value_type(sd->getGlobalIndex(), oneAtomSet));        
424        }
425      }  
426 +
427            
428      for (bond= mol->beginBond(bondIter); bond != NULL;
429           bond = mol->nextBond(bondIter)) {
430  
431        a = bond->getAtomA()->getGlobalIndex();
432        b = bond->getAtomB()->getGlobalIndex();  
433 <    
433 >
434        if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
435          oneTwoInteractions_.addPair(a, b);
436        } else {
# Line 503 | Line 524 | namespace oopse {
524  
525      for (rb = mol->beginRigidBody(rbIter); rb != NULL;
526           rb = mol->nextRigidBody(rbIter)) {
527 <      std::vector<Atom*> atoms = rb->getAtoms();
527 >      vector<Atom*> atoms = rb->getAtoms();
528        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            a = atoms[i]->getGlobalIndex();
# Line 517 | Line 538 | namespace oopse {
538  
539    void SimInfo::removeInteractionPairs(Molecule* mol) {
540      ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
541 <    std::vector<Bond*>::iterator bondIter;
542 <    std::vector<Bend*>::iterator bendIter;
543 <    std::vector<Torsion*>::iterator torsionIter;
544 <    std::vector<Inversion*>::iterator inversionIter;
541 >    vector<Bond*>::iterator bondIter;
542 >    vector<Bend*>::iterator bendIter;
543 >    vector<Torsion*>::iterator torsionIter;
544 >    vector<Inversion*>::iterator inversionIter;
545      Bond* bond;
546      Bend* bend;
547      Torsion* torsion;
# Line 530 | Line 551 | namespace oopse {
551      int c;
552      int d;
553  
554 <    std::map<int, std::set<int> > atomGroups;
554 >    map<int, set<int> > atomGroups;
555      Molecule::RigidBodyIterator rbIter;
556      RigidBody* rb;
557      Molecule::IntegrableObjectIterator ii;
558 <    StuntDouble* integrableObject;
558 >    StuntDouble* sd;
559      
560 <    for (integrableObject = mol->beginIntegrableObject(ii);
561 <         integrableObject != NULL;
541 <         integrableObject = mol->nextIntegrableObject(ii)) {
560 >    for (sd = mol->beginIntegrableObject(ii); sd != NULL;
561 >         sd = mol->nextIntegrableObject(ii)) {
562        
563 <      if (integrableObject->isRigidBody()) {
564 <        rb = static_cast<RigidBody*>(integrableObject);
565 <        std::vector<Atom*> atoms = rb->getAtoms();
566 <        std::set<int> rigidAtoms;
563 >      if (sd->isRigidBody()) {
564 >        rb = static_cast<RigidBody*>(sd);
565 >        vector<Atom*> atoms = rb->getAtoms();
566 >        set<int> rigidAtoms;
567          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 <          atomGroups.insert(std::map<int, std::set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
571 >          atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
572          }      
573        } else {
574 <        std::set<int> oneAtomSet;
575 <        oneAtomSet.insert(integrableObject->getGlobalIndex());
576 <        atomGroups.insert(std::map<int, std::set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));        
574 >        set<int> oneAtomSet;
575 >        oneAtomSet.insert(sd->getGlobalIndex());
576 >        atomGroups.insert(map<int, set<int> >::value_type(sd->getGlobalIndex(), oneAtomSet));        
577        }
578      }  
579  
# Line 656 | Line 676 | namespace oopse {
676  
677      for (rb = mol->beginRigidBody(rbIter); rb != NULL;
678           rb = mol->nextRigidBody(rbIter)) {
679 <      std::vector<Atom*> atoms = rb->getAtoms();
679 >      vector<Atom*> atoms = rb->getAtoms();
680        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            a = atoms[i]->getGlobalIndex();
# Line 679 | Line 699 | namespace oopse {
699      molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
700    }
701  
682  void SimInfo::update() {
702  
703 <    setupSimType();
704 <
705 < #ifdef IS_MPI
706 <    setupFortranParallel();
707 < #endif
708 <
709 <    setupFortranSim();
710 <
711 <    //setup fortran force field
712 <    /** @deprecate */    
694 <    int isError = 0;
695 <    
696 <    setupCutoff();
697 <    
698 <    setupElectrostaticSummationMethod( isError );
699 <    setupSwitchingFunction();
700 <    setupAccumulateBoxDipole();
701 <
702 <    if(isError){
703 <      sprintf( painCave.errMsg,
704 <               "ForceField error: There was an error initializing the forceField in fortran.\n" );
705 <      painCave.isFatal = 1;
706 <      simError();
707 <    }
708 <
703 >  /**
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 >    calcNConstraints();
713      calcNdf();
714      calcNdfRaw();
715      calcNdfTrans();
712
713    fortranInitialized_ = true;
716    }
717 <
718 <  std::set<AtomType*> SimInfo::getUniqueAtomTypes() {
717 >  
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      SimInfo::MoleculeIterator mi;
727      Molecule* mol;
728      Molecule::AtomIterator ai;
729      Atom* atom;
730 <    std::set<AtomType*> atomTypes;
731 <
730 >    set<AtomType*> atomTypes;
731 >    
732      for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
733 <
734 <      for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
733 >      for(atom = mol->beginAtom(ai); atom != NULL;
734 >          atom = mol->nextAtom(ai)) {
735          atomTypes.insert(atom->getAtomType());
736 <      }
737 <        
729 <    }
730 <
731 <    return atomTypes;        
732 <  }
733 <
734 <  void SimInfo::setupSimType() {
735 <    std::set<AtomType*>::iterator i;
736 <    std::set<AtomType*> atomTypes;
737 <    atomTypes = getUniqueAtomTypes();
736 >      }      
737 >    }    
738      
739 <    int useLennardJones = 0;
740 <    int useElectrostatic = 0;
741 <    int useEAM = 0;
742 <    int useSC = 0;
743 <    int useCharge = 0;
744 <    int useDirectional = 0;
745 <    int useDipole = 0;
746 <    int useGayBerne = 0;
747 <    int useSticky = 0;
748 <    int useStickyPower = 0;
749 <    int useShape = 0;
750 <    int useFLARB = 0; //it is not in AtomType yet
751 <    int useDirectionalAtom = 0;    
752 <    int useElectrostatics = 0;
753 <    //usePBC and useRF are from simParams
754 <    int usePBC = simParams_->getUsePeriodicBoundaryConditions();
755 <    int useRF;
756 <    int useSF;
757 <    int useSP;
758 <    int useBoxDipole;
739 > #ifdef IS_MPI
740  
741 <    std::string myMethod;
742 <
762 <    // set the useRF logical
763 <    useRF = 0;
764 <    useSF = 0;
765 <    useSP = 0;
766 <    useBoxDipole = 0;
767 <
768 <
769 <    if (simParams_->haveElectrostaticSummationMethod()) {
770 <      std::string myMethod = simParams_->getElectrostaticSummationMethod();
771 <      toUpper(myMethod);
772 <      if (myMethod == "REACTION_FIELD"){
773 <        useRF = 1;
774 <      } else if (myMethod == "SHIFTED_FORCE"){
775 <        useSF = 1;
776 <      } else if (myMethod == "SHIFTED_POTENTIAL"){
777 <        useSP = 1;
778 <      }
779 <    }
741 >    // loop over the found atom types on this processor, and add their
742 >    // numerical idents to a vector:
743      
744 <    if (simParams_->haveAccumulateBoxDipole())
745 <      if (simParams_->getAccumulateBoxDipole())
746 <        useBoxDipole = 1;
744 >    vector<int> foundTypes;
745 >    set<AtomType*>::iterator i;
746 >    for (i = atomTypes.begin(); i != atomTypes.end(); ++i)
747 >      foundTypes.push_back( (*i)->getIdent() );
748  
749 <    useAtomicVirial_ = simParams_->getUseAtomicVirial();
749 >    // count_local holds the number of found types on this processor
750 >    int count_local = foundTypes.size();
751  
752 <    //loop over all of the atom types
753 <    for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
789 <      useLennardJones |= (*i)->isLennardJones();
790 <      useElectrostatic |= (*i)->isElectrostatic();
791 <      useEAM |= (*i)->isEAM();
792 <      useSC |= (*i)->isSC();
793 <      useCharge |= (*i)->isCharge();
794 <      useDirectional |= (*i)->isDirectional();
795 <      useDipole |= (*i)->isDipole();
796 <      useGayBerne |= (*i)->isGayBerne();
797 <      useSticky |= (*i)->isSticky();
798 <      useStickyPower |= (*i)->isStickyPower();
799 <      useShape |= (*i)->isShape();
800 <    }
752 >    int nproc;
753 >    MPI_Comm_size( MPI_COMM_WORLD, &nproc);
754  
755 <    if (useSticky || useStickyPower || useDipole || useGayBerne || useShape) {
756 <      useDirectionalAtom = 1;
757 <    }
755 >    // 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  
760 <    if (useCharge || useDipole) {
761 <      useElectrostatics = 1;
760 >    // fill the counts array
761 >    MPI_Allgather(&count_local, 1, MPI_INT, &counts[0],
762 >                  1, MPI_INT, MPI_COMM_WORLD);
763 >  
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      }
771  
772 < #ifdef IS_MPI    
773 <    int temp;
772 >    // 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 >    MPI_Allgatherv(&foundTypes[0], count_local, MPI_INT,
777 >                   &ftGlobal[0], &counts[0], &disps[0],
778 >                   MPI_INT, MPI_COMM_WORLD);
779  
780 <    temp = usePBC;
814 <    MPI_Allreduce(&temp, &usePBC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
780 >    vector<int>::iterator j;
781  
782 <    temp = useDirectionalAtom;
783 <    MPI_Allreduce(&temp, &useDirectionalAtom, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
782 >    // foundIdents is a stl set, so inserting an already found ident
783 >    // will have no effect.
784 >    set<int> foundIdents;
785  
786 <    temp = useLennardJones;
787 <    MPI_Allreduce(&temp, &useLennardJones, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
786 >    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  
797 <    temp = useElectrostatics;
798 <    MPI_Allreduce(&temp, &useElectrostatics, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
797 >    return atomTypes;        
798 >  }
799  
825    temp = useCharge;
826    MPI_Allreduce(&temp, &useCharge, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
800  
801 <    temp = useDipole;
802 <    MPI_Allreduce(&temp, &useDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
801 >  int getGlobalCountOfType(AtomType* atype) {
802 >    /*
803 >    set<AtomType*> atypes = getSimulatedAtomTypes();
804 >    map<AtomType*, int> counts_;
805  
806 <    temp = useSticky;
807 <    MPI_Allreduce(&temp, &useSticky, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
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 <    temp = useStickyPower;
817 <    MPI_Allreduce(&temp, &useStickyPower, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
816 >  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 <    temp = useGayBerne;
827 <    MPI_Allreduce(&temp, &useGayBerne, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
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  
841 <    temp = useEAM;
842 <    MPI_Allreduce(&temp, &useEAM, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
841 > #ifdef IS_MPI
842 >    int temp;
843  
844 <    temp = useSC;
845 <    MPI_Allreduce(&temp, &useSC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
844 >    temp = usesDirectional;
845 >    MPI_Allreduce(MPI_IN_PLACE, &temp, 1, MPI_INT,  MPI_LOR, MPI_COMM_WORLD);
846 >    usesDirectionalAtoms_ = (temp == 0) ? false : true;
847      
848 <    temp = useShape;
849 <    MPI_Allreduce(&temp, &useShape, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);  
848 >    temp = usesMetallic;
849 >    MPI_Allreduce(MPI_IN_PLACE, &temp, 1, MPI_INT,  MPI_LOR, MPI_COMM_WORLD);
850 >    usesMetallicAtoms_ = (temp == 0) ? false : true;
851  
852 <    temp = useFLARB;
853 <    MPI_Allreduce(&temp, &useFLARB, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
852 >    temp = usesElectrostatic;
853 >    MPI_Allreduce(MPI_IN_PLACE, &temp, 1, MPI_INT,  MPI_LOR, MPI_COMM_WORLD);
854 >    usesElectrostaticAtoms_ = (temp == 0) ? false : true;
855  
856 <    temp = useRF;
857 <    MPI_Allreduce(&temp, &useRF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
856 >    temp = usesFluctuatingCharges;
857 >    MPI_Allreduce(MPI_IN_PLACE, &temp, 1, MPI_INT,  MPI_LOR, MPI_COMM_WORLD);
858 >    usesFluctuatingCharges_ = (temp == 0) ? false : true;
859 > #else
860  
861 <    temp = useSF;
862 <    MPI_Allreduce(&temp, &useSF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);  
861 >    usesDirectionalAtoms_ = usesDirectional;
862 >    usesMetallicAtoms_ = usesMetallic;
863 >    usesElectrostaticAtoms_ = usesElectrostatic;
864 >    usesFluctuatingCharges_ = usesFluctuatingCharges;
865  
866 <    temp = useSP;
867 <    MPI_Allreduce(&temp, &useSP, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
866 > #endif
867 >    
868 >    requiresPrepair_ = usesMetallicAtoms_ ? true : false;
869 >    requiresSkipCorrection_ = usesElectrostaticAtoms_ ? true : false;
870 >    requiresSelfCorrection_ = usesElectrostaticAtoms_ ? true : false;    
871 >  }
872  
861    temp = useBoxDipole;
862    MPI_Allreduce(&temp, &useBoxDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
873  
874 <    temp = useAtomicVirial_;
875 <    MPI_Allreduce(&temp, &useAtomicVirial_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
874 >  vector<int> SimInfo::getGlobalAtomIndices() {
875 >    SimInfo::MoleculeIterator mi;
876 >    Molecule* mol;
877 >    Molecule::AtomIterator ai;
878 >    Atom* atom;
879  
880 < #endif
881 <
882 <    fInfo_.SIM_uses_PBC = usePBC;    
883 <    fInfo_.SIM_uses_DirectionalAtoms = useDirectionalAtom;
884 <    fInfo_.SIM_uses_LennardJones = useLennardJones;
885 <    fInfo_.SIM_uses_Electrostatics = useElectrostatics;    
886 <    fInfo_.SIM_uses_Charges = useCharge;
887 <    fInfo_.SIM_uses_Dipoles = useDipole;
888 <    fInfo_.SIM_uses_Sticky = useSticky;
876 <    fInfo_.SIM_uses_StickyPower = useStickyPower;
877 <    fInfo_.SIM_uses_GayBerne = useGayBerne;
878 <    fInfo_.SIM_uses_EAM = useEAM;
879 <    fInfo_.SIM_uses_SC = useSC;
880 <    fInfo_.SIM_uses_Shapes = useShape;
881 <    fInfo_.SIM_uses_FLARB = useFLARB;
882 <    fInfo_.SIM_uses_RF = useRF;
883 <    fInfo_.SIM_uses_SF = useSF;
884 <    fInfo_.SIM_uses_SP = useSP;
885 <    fInfo_.SIM_uses_BoxDipole = useBoxDipole;
886 <    fInfo_.SIM_uses_AtomicVirial = useAtomicVirial_;
880 >    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  
889  void SimInfo::setupFortranSim() {
890    int isError;
891    int nExclude, nOneTwo, nOneThree, nOneFour;
892    std::vector<int> fortranGlobalGroupMembership;
893    
894    isError = 0;
891  
892 <    //globalGroupMembership_ is filled by SimCreator    
893 <    for (int i = 0; i < nGlobalAtoms_; i++) {
894 <      fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
892 >  vector<int> SimInfo::getGlobalGroupIndices() {
893 >    SimInfo::MoleculeIterator mi;
894 >    Molecule* mol;
895 >    Molecule::CutoffGroupIterator ci;
896 >    CutoffGroup* cg;
897 >
898 >    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  
912 +
913 +  void SimInfo::prepareTopology() {
914 +
915      //calculate mass ratio of cutoff group
902    std::vector<RealType> mfact;
916      SimInfo::MoleculeIterator mi;
917      Molecule* mol;
918      Molecule::CutoffGroupIterator ci;
# Line 908 | Line 921 | namespace oopse {
921      Atom* atom;
922      RealType totalMass;
923  
924 <    //to avoid memory reallocation, reserve enough space for mfact
925 <    mfact.reserve(getNCutoffGroups());
924 >    /**
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      
935      for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {        
936 <      for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
936 >      for (cg = mol->beginCutoffGroup(ci); cg != NULL;
937 >           cg = mol->nextCutoffGroup(ci)) {
938  
939          totalMass = cg->getMass();
940          for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
941            // Check for massless groups - set mfact to 1 if true
942 <          if (totalMass != 0)
943 <            mfact.push_back(atom->getMass()/totalMass);
942 >          if (totalMass != 0)
943 >            massFactors_[atom->getLocalIndex()] = atom->getMass()/totalMass;
944            else
945 <            mfact.push_back( 1.0 );
945 >            massFactors_[atom->getLocalIndex()] = 1.0;
946          }
947        }      
948      }
949  
950 <    //fill ident array of local atoms (it is actually ident of AtomType, it is so confusing !!!)
929 <    std::vector<int> identArray;
950 >    // Build the identArray_ and regions_
951  
952 <    //to avoid memory reallocation, reserve enough space identArray
953 <    identArray.reserve(getNAtoms());
954 <    
955 <    for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {        
952 >    identArray_.clear();
953 >    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        for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
960 <        identArray.push_back(atom->getIdent());
960 >        identArray_.push_back(atom->getIdent());
961 >        regions_.push_back(reg);
962        }
963      }    
964 +      
965 +    topologyDone_ = true;
966 +  }
967  
968 <    //fill molMembershipArray
969 <    //molMembershipArray is filled by SimCreator    
970 <    std::vector<int> molMembershipArray(nGlobalAtoms_);
943 <    for (int i = 0; i < nGlobalAtoms_; i++) {
944 <      molMembershipArray[i] = globalMolMembership_[i] + 1;
945 <    }
946 <    
947 <    //setup fortran simulation
968 >  void SimInfo::addProperty(GenericData* genData) {
969 >    properties_.addProperty(genData);  
970 >  }
971  
972 <    nExclude = excludedInteractions_.getSize();
973 <    nOneTwo = oneTwoInteractions_.getSize();
974 <    nOneThree = oneThreeInteractions_.getSize();
952 <    nOneFour = oneFourInteractions_.getSize();
972 >  void SimInfo::removeProperty(const string& propName) {
973 >    properties_.removeProperty(propName);  
974 >  }
975  
976 <    int* excludeList = excludedInteractions_.getPairList();
977 <    int* oneTwoList = oneTwoInteractions_.getPairList();
978 <    int* oneThreeList = oneThreeInteractions_.getPairList();
957 <    int* oneFourList = oneFourInteractions_.getPairList();
976 >  void SimInfo::clearProperties() {
977 >    properties_.clearProperties();
978 >  }
979  
980 <    setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray[0],
981 <                   &nExclude, excludeList,
982 <                   &nOneTwo, oneTwoList,
962 <                   &nOneThree, oneThreeList,
963 <                   &nOneFour, oneFourList,
964 <                   &molMembershipArray[0], &mfact[0], &nCutoffGroups_,
965 <                   &fortranGlobalGroupMembership[0], &isError);
966 <    
967 <    if( isError ){
980 >  vector<string> SimInfo::getPropertyNames() {
981 >    return properties_.getPropertyNames();  
982 >  }
983        
984 <      sprintf( painCave.errMsg,
985 <               "There was an error setting the simulation information in fortran.\n" );
986 <      painCave.isFatal = 1;
972 <      painCave.severity = OOPSE_ERROR;
973 <      simError();
974 <    }
975 <    
976 <    
977 <    sprintf( checkPointMsg,
978 <             "succesfully sent the simulation information to fortran.\n");
979 <    
980 <    errorCheckPoint();
981 <    
982 <    // Setup number of neighbors in neighbor list if present
983 <    if (simParams_->haveNeighborListNeighbors()) {
984 <      int nlistNeighbors = simParams_->getNeighborListNeighbors();
985 <      setNeighbors(&nlistNeighbors);
986 <    }
987 <  
984 >  vector<GenericData*> SimInfo::getProperties() {
985 >    return properties_.getProperties();
986 >  }
987  
988 +  GenericData* SimInfo::getPropertyByName(const string& propName) {
989 +    return properties_.getPropertyByName(propName);
990    }
991  
992 +  void SimInfo::setSnapshotManager(SnapshotManager* sman) {
993 +    if (sman_ == sman) {
994 +      return;
995 +    }    
996 +    delete sman_;
997 +    sman_ = sman;
998  
992  void SimInfo::setupFortranParallel() {
993 #ifdef IS_MPI    
994    //SimInfo is responsible for creating localToGlobalAtomIndex and localToGlobalGroupIndex
995    std::vector<int> localToGlobalAtomIndex(getNAtoms(), 0);
996    std::vector<int> localToGlobalCutoffGroupIndex;
999      SimInfo::MoleculeIterator mi;
1000      Molecule::AtomIterator ai;
1001 <    Molecule::CutoffGroupIterator ci;
1001 >    Molecule::RigidBodyIterator rbIter;
1002 >    Molecule::CutoffGroupIterator cgIter;
1003 >    Molecule::BondIterator bondIter;
1004 >    Molecule::BendIterator bendIter;
1005 >    Molecule::TorsionIterator torsionIter;
1006 >    Molecule::InversionIterator inversionIter;
1007 >
1008      Molecule* mol;
1009      Atom* atom;
1010 +    RigidBody* rb;
1011      CutoffGroup* cg;
1012 <    mpiSimData parallelData;
1013 <    int isError;
1012 >    Bond* bond;
1013 >    Bend* bend;
1014 >    Torsion* torsion;
1015 >    Inversion* inversion;    
1016  
1017 <    for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) {
1007 <
1008 <      //local index(index in DataStorge) of atom is important
1009 <      for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
1010 <        localToGlobalAtomIndex[atom->getLocalIndex()] = atom->getGlobalIndex() + 1;
1011 <      }
1012 <
1013 <      //local index of cutoff group is trivial, it only depends on the order of travesing
1014 <      for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
1015 <        localToGlobalCutoffGroupIndex.push_back(cg->getGlobalIndex() + 1);
1016 <      }        
1017 >    for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
1018          
1019 <    }
1020 <
1021 <    //fill up mpiSimData struct
1022 <    parallelData.nMolGlobal = getNGlobalMolecules();
1023 <    parallelData.nMolLocal = getNMolecules();
1024 <    parallelData.nAtomsGlobal = getNGlobalAtoms();
1025 <    parallelData.nAtomsLocal = getNAtoms();
1025 <    parallelData.nGroupsGlobal = getNGlobalCutoffGroups();
1026 <    parallelData.nGroupsLocal = getNCutoffGroups();
1027 <    parallelData.myNode = worldRank;
1028 <    MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors));
1029 <
1030 <    //pass mpiSimData struct and index arrays to fortran
1031 <    setFsimParallel(&parallelData, &(parallelData.nAtomsLocal),
1032 <                    &localToGlobalAtomIndex[0],  &(parallelData.nGroupsLocal),
1033 <                    &localToGlobalCutoffGroupIndex[0], &isError);
1034 <
1035 <    if (isError) {
1036 <      sprintf(painCave.errMsg,
1037 <              "mpiRefresh errror: fortran didn't like something we gave it.\n");
1038 <      painCave.isFatal = 1;
1039 <      simError();
1040 <    }
1041 <
1042 <    sprintf(checkPointMsg, " mpiRefresh successful.\n");
1043 <    errorCheckPoint();
1044 <
1045 < #endif
1046 <  }
1047 <
1048 <  void SimInfo::setupCutoff() {          
1049 <    
1050 <    ForceFieldOptions& forceFieldOptions_ = forceField_->getForceFieldOptions();
1051 <
1052 <    // Check the cutoff policy
1053 <    int cp =  TRADITIONAL_CUTOFF_POLICY; // Set to traditional by default
1054 <
1055 <    // Set LJ shifting bools to false
1056 <    ljsp_ = 0;
1057 <    ljsf_ = 0;
1058 <
1059 <    std::string myPolicy;
1060 <    if (forceFieldOptions_.haveCutoffPolicy()){
1061 <      myPolicy = forceFieldOptions_.getCutoffPolicy();
1062 <    }else if (simParams_->haveCutoffPolicy()) {
1063 <      myPolicy = simParams_->getCutoffPolicy();
1064 <    }
1065 <
1066 <    if (!myPolicy.empty()){
1067 <      toUpper(myPolicy);
1068 <      if (myPolicy == "MIX") {
1069 <        cp = MIX_CUTOFF_POLICY;
1070 <      } else {
1071 <        if (myPolicy == "MAX") {
1072 <          cp = MAX_CUTOFF_POLICY;
1073 <        } else {
1074 <          if (myPolicy == "TRADITIONAL") {            
1075 <            cp = TRADITIONAL_CUTOFF_POLICY;
1076 <          } else {
1077 <            // throw error        
1078 <            sprintf( painCave.errMsg,
1079 <                     "SimInfo error: Unknown cutoffPolicy. (Input file specified %s .)\n\tcutoffPolicy must be one of: \"Mix\", \"Max\", or \"Traditional\".", myPolicy.c_str() );
1080 <            painCave.isFatal = 1;
1081 <            simError();
1082 <          }    
1083 <        }          
1019 >      for (atom = mol->beginAtom(ai); atom != NULL;
1020 >           atom = mol->nextAtom(ai)) {
1021 >        atom->setSnapshotManager(sman_);
1022 >      }        
1023 >      for (rb = mol->beginRigidBody(rbIter); rb != NULL;
1024 >           rb = mol->nextRigidBody(rbIter)) {
1025 >        rb->setSnapshotManager(sman_);
1026        }
1027 <    }          
1028 <    notifyFortranCutoffPolicy(&cp);
1029 <
1088 <    // Check the Skin Thickness for neighborlists
1089 <    RealType skin;
1090 <    if (simParams_->haveSkinThickness()) {
1091 <      skin = simParams_->getSkinThickness();
1092 <      notifyFortranSkinThickness(&skin);
1093 <    }            
1094 <        
1095 <    // Check if the cutoff was set explicitly:
1096 <    if (simParams_->haveCutoffRadius()) {
1097 <      rcut_ = simParams_->getCutoffRadius();
1098 <      if (simParams_->haveSwitchingRadius()) {
1099 <        rsw_  = simParams_->getSwitchingRadius();
1100 <      } else {
1101 <        if (fInfo_.SIM_uses_Charges |
1102 <            fInfo_.SIM_uses_Dipoles |
1103 <            fInfo_.SIM_uses_RF) {
1104 <          
1105 <          rsw_ = 0.85 * rcut_;
1106 <          sprintf(painCave.errMsg,
1107 <                  "SimCreator Warning: No value was set for the switchingRadius.\n"
1108 <                  "\tOOPSE will use a default value of 85 percent of the cutoffRadius.\n"
1109 <                  "\tswitchingRadius = %f. for this simulation\n", rsw_);
1110 <        painCave.isFatal = 0;
1111 <        simError();
1112 <        } else {
1113 <          rsw_ = rcut_;
1114 <          sprintf(painCave.errMsg,
1115 <                  "SimCreator Warning: No value was set for the switchingRadius.\n"
1116 <                  "\tOOPSE will use the same value as the cutoffRadius.\n"
1117 <                  "\tswitchingRadius = %f. for this simulation\n", rsw_);
1118 <          painCave.isFatal = 0;
1119 <          simError();
1120 <        }
1027 >      for (cg = mol->beginCutoffGroup(cgIter); cg != NULL;
1028 >           cg = mol->nextCutoffGroup(cgIter)) {
1029 >        cg->setSnapshotManager(sman_);
1030        }
1031 <
1032 <      if (simParams_->haveElectrostaticSummationMethod()) {
1033 <        std::string myMethod = simParams_->getElectrostaticSummationMethod();
1125 <        toUpper(myMethod);
1126 <        
1127 <        if (myMethod == "SHIFTED_POTENTIAL") {
1128 <          ljsp_ = 1;
1129 <        } else if (myMethod == "SHIFTED_FORCE") {
1130 <          ljsf_ = 1;
1131 <        }
1031 >      for (bond = mol->beginBond(bondIter); bond != NULL;
1032 >           bond = mol->nextBond(bondIter)) {
1033 >        bond->setSnapshotManager(sman_);
1034        }
1035 <
1036 <      notifyFortranCutoffs(&rcut_, &rsw_, &ljsp_, &ljsf_);
1037 <      
1136 <    } else {
1137 <      
1138 <      // For electrostatic atoms, we'll assume a large safe value:
1139 <      if (fInfo_.SIM_uses_Charges | fInfo_.SIM_uses_Dipoles | fInfo_.SIM_uses_RF) {
1140 <        sprintf(painCave.errMsg,
1141 <                "SimCreator Warning: No value was set for the cutoffRadius.\n"
1142 <                "\tOOPSE will use a default value of 15.0 angstroms"
1143 <                "\tfor the cutoffRadius.\n");
1144 <        painCave.isFatal = 0;
1145 <        simError();
1146 <        rcut_ = 15.0;
1147 <      
1148 <        if (simParams_->haveElectrostaticSummationMethod()) {
1149 <          std::string myMethod = simParams_->getElectrostaticSummationMethod();
1150 <          toUpper(myMethod);
1151 <      
1152 <      // For the time being, we're tethering the LJ shifted behavior to the
1153 <      // electrostaticSummationMethod keyword options
1154 <          if (myMethod == "SHIFTED_POTENTIAL") {
1155 <            ljsp_ = 1;
1156 <          } else if (myMethod == "SHIFTED_FORCE") {
1157 <            ljsf_ = 1;
1158 <          }
1159 <          if (myMethod == "SHIFTED_POTENTIAL" || myMethod == "SHIFTED_FORCE") {
1160 <            if (simParams_->haveSwitchingRadius()){
1161 <              sprintf(painCave.errMsg,
1162 <                      "SimInfo Warning: A value was set for the switchingRadius\n"
1163 <                      "\teven though the electrostaticSummationMethod was\n"
1164 <                      "\tset to %s\n", myMethod.c_str());
1165 <              painCave.isFatal = 1;
1166 <              simError();            
1167 <            }
1168 <          }
1169 <        }
1170 <      
1171 <        if (simParams_->haveSwitchingRadius()){
1172 <          rsw_ = simParams_->getSwitchingRadius();
1173 <        } else {        
1174 <          sprintf(painCave.errMsg,
1175 <                  "SimCreator Warning: No value was set for switchingRadius.\n"
1176 <                  "\tOOPSE will use a default value of\n"
1177 <                  "\t0.85 * cutoffRadius for the switchingRadius\n");
1178 <          painCave.isFatal = 0;
1179 <          simError();
1180 <          rsw_ = 0.85 * rcut_;
1181 <        }
1182 <
1183 <        notifyFortranCutoffs(&rcut_, &rsw_, &ljsp_, &ljsf_);
1184 <
1185 <      } else {
1186 <        // We didn't set rcut explicitly, and we don't have electrostatic atoms, so
1187 <        // We'll punt and let fortran figure out the cutoffs later.
1188 <        
1189 <        notifyFortranYouAreOnYourOwn();
1190 <
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    }
1049  
1195  void SimInfo::setupElectrostaticSummationMethod( int isError ) {    
1196    
1197    int errorOut;
1198    int esm =  NONE;
1199    int sm = UNDAMPED;
1200    RealType alphaVal;
1201    RealType dielectric;
1202    
1203    errorOut = isError;
1050  
1051 <    if (simParams_->haveElectrostaticSummationMethod()) {
1206 <      std::string myMethod = simParams_->getElectrostaticSummationMethod();
1207 <      toUpper(myMethod);
1208 <      if (myMethod == "NONE") {
1209 <        esm = NONE;
1210 <      } else {
1211 <        if (myMethod == "SWITCHING_FUNCTION") {
1212 <          esm = SWITCHING_FUNCTION;
1213 <        } else {
1214 <          if (myMethod == "SHIFTED_POTENTIAL") {
1215 <            esm = SHIFTED_POTENTIAL;
1216 <          } else {
1217 <            if (myMethod == "SHIFTED_FORCE") {            
1218 <              esm = SHIFTED_FORCE;
1219 <            } else {
1220 <              if (myMethod == "REACTION_FIELD") {
1221 <                esm = REACTION_FIELD;
1222 <                dielectric = simParams_->getDielectric();
1223 <                if (!simParams_->haveDielectric()) {
1224 <                  // throw warning
1225 <                  sprintf( painCave.errMsg,
1226 <                           "SimInfo warning: dielectric was not specified in the input file\n\tfor the reaction field correction method.\n"
1227 <                           "\tA default value of %f will be used for the dielectric.\n", dielectric);
1228 <                  painCave.isFatal = 0;
1229 <                  simError();
1230 <                }
1231 <              } else {
1232 <                // throw error        
1233 <                sprintf( painCave.errMsg,
1234 <                         "SimInfo error: Unknown electrostaticSummationMethod.\n"
1235 <                         "\t(Input file specified %s .)\n"
1236 <                         "\telectrostaticSummationMethod must be one of: \"none\",\n"
1237 <                         "\t\"shifted_potential\", \"shifted_force\", or \n"
1238 <                         "\t\"reaction_field\".\n", myMethod.c_str() );
1239 <                painCave.isFatal = 1;
1240 <                simError();
1241 <              }    
1242 <            }          
1243 <          }
1244 <        }
1245 <      }
1246 <    }
1247 <    
1248 <    if (simParams_->haveElectrostaticScreeningMethod()) {
1249 <      std::string myScreen = simParams_->getElectrostaticScreeningMethod();
1250 <      toUpper(myScreen);
1251 <      if (myScreen == "UNDAMPED") {
1252 <        sm = UNDAMPED;
1253 <      } else {
1254 <        if (myScreen == "DAMPED") {
1255 <          sm = DAMPED;
1256 <          if (!simParams_->haveDampingAlpha()) {
1257 <            // first set a cutoff dependent alpha value
1258 <            // we assume alpha depends linearly with rcut from 0 to 20.5 ang
1259 <            alphaVal = 0.5125 - rcut_* 0.025;
1260 <            // for values rcut > 20.5, alpha is zero
1261 <            if (alphaVal < 0) alphaVal = 0;
1051 >  ostream& operator <<(ostream& o, SimInfo& info) {
1052  
1263            // throw warning
1264            sprintf( painCave.errMsg,
1265                     "SimInfo warning: dampingAlpha was not specified in the input file.\n"
1266                     "\tA default value of %f (1/ang) will be used for the cutoff of\n\t%f (ang).\n", alphaVal, rcut_);
1267            painCave.isFatal = 0;
1268            simError();
1269          } else {
1270            alphaVal = simParams_->getDampingAlpha();
1271          }
1272          
1273        } else {
1274          // throw error        
1275          sprintf( painCave.errMsg,
1276                   "SimInfo error: Unknown electrostaticScreeningMethod.\n"
1277                   "\t(Input file specified %s .)\n"
1278                   "\telectrostaticScreeningMethod must be one of: \"undamped\"\n"
1279                   "or \"damped\".\n", myScreen.c_str() );
1280          painCave.isFatal = 1;
1281          simError();
1282        }
1283      }
1284    }
1285    
1286    // let's pass some summation method variables to fortran
1287    setElectrostaticSummationMethod( &esm );
1288    setFortranElectrostaticMethod( &esm );
1289    setScreeningMethod( &sm );
1290    setDampingAlpha( &alphaVal );
1291    setReactionFieldDielectric( &dielectric );
1292    initFortranFF( &errorOut );
1293  }
1294
1295  void SimInfo::setupSwitchingFunction() {    
1296    int ft = CUBIC;
1297
1298    if (simParams_->haveSwitchingFunctionType()) {
1299      std::string funcType = simParams_->getSwitchingFunctionType();
1300      toUpper(funcType);
1301      if (funcType == "CUBIC") {
1302        ft = CUBIC;
1303      } else {
1304        if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
1305          ft = FIFTH_ORDER_POLY;
1306        } else {
1307          // throw error        
1308          sprintf( painCave.errMsg,
1309                   "SimInfo error: Unknown switchingFunctionType. (Input file specified %s .)\n\tswitchingFunctionType must be one of: \"cubic\" or \"fifth_order_polynomial\".", funcType.c_str() );
1310          painCave.isFatal = 1;
1311          simError();
1312        }          
1313      }
1314    }
1315
1316    // send switching function notification to switcheroo
1317    setFunctionType(&ft);
1318
1319  }
1320
1321  void SimInfo::setupAccumulateBoxDipole() {    
1322
1323    // we only call setAccumulateBoxDipole if the accumulateBoxDipole parameter is true
1324    if ( simParams_->haveAccumulateBoxDipole() )
1325      if ( simParams_->getAccumulateBoxDipole() ) {
1326        setAccumulateBoxDipole();
1327        calcBoxDipole_ = true;
1328      }
1329
1330  }
1331
1332  void SimInfo::addProperty(GenericData* genData) {
1333    properties_.addProperty(genData);  
1334  }
1335
1336  void SimInfo::removeProperty(const std::string& propName) {
1337    properties_.removeProperty(propName);  
1338  }
1339
1340  void SimInfo::clearProperties() {
1341    properties_.clearProperties();
1342  }
1343
1344  std::vector<std::string> SimInfo::getPropertyNames() {
1345    return properties_.getPropertyNames();  
1346  }
1347      
1348  std::vector<GenericData*> SimInfo::getProperties() {
1349    return properties_.getProperties();
1350  }
1351
1352  GenericData* SimInfo::getPropertyByName(const std::string& propName) {
1353    return properties_.getPropertyByName(propName);
1354  }
1355
1356  void SimInfo::setSnapshotManager(SnapshotManager* sman) {
1357    if (sman_ == sman) {
1358      return;
1359    }    
1360    delete sman_;
1361    sman_ = sman;
1362
1363    Molecule* mol;
1364    RigidBody* rb;
1365    Atom* atom;
1366    SimInfo::MoleculeIterator mi;
1367    Molecule::RigidBodyIterator rbIter;
1368    Molecule::AtomIterator atomIter;;
1369
1370    for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
1371        
1372      for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) {
1373        atom->setSnapshotManager(sman_);
1374      }
1375        
1376      for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
1377        rb->setSnapshotManager(sman_);
1378      }
1379    }    
1380    
1381  }
1382
1383  Vector3d SimInfo::getComVel(){
1384    SimInfo::MoleculeIterator i;
1385    Molecule* mol;
1386
1387    Vector3d comVel(0.0);
1388    RealType totalMass = 0.0;
1389    
1390
1391    for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1392      RealType mass = mol->getMass();
1393      totalMass += mass;
1394      comVel += mass * mol->getComVel();
1395    }  
1396
1397 #ifdef IS_MPI
1398    RealType tmpMass = totalMass;
1399    Vector3d tmpComVel(comVel);    
1400    MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1401    MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1402 #endif
1403
1404    comVel /= totalMass;
1405
1406    return comVel;
1407  }
1408
1409  Vector3d SimInfo::getCom(){
1410    SimInfo::MoleculeIterator i;
1411    Molecule* mol;
1412
1413    Vector3d com(0.0);
1414    RealType totalMass = 0.0;
1415    
1416    for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1417      RealType mass = mol->getMass();
1418      totalMass += mass;
1419      com += mass * mol->getCom();
1420    }  
1421
1422 #ifdef IS_MPI
1423    RealType tmpMass = totalMass;
1424    Vector3d tmpCom(com);    
1425    MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1426    MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1427 #endif
1428
1429    com /= totalMass;
1430
1431    return com;
1432
1433  }        
1434
1435  std::ostream& operator <<(std::ostream& o, SimInfo& info) {
1436
1053      return o;
1054    }
1055    
1056 <  
1441 <   /*
1442 <   Returns center of mass and center of mass velocity in one function call.
1443 <   */
1444 <  
1445 <   void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){
1446 <      SimInfo::MoleculeIterator i;
1447 <      Molecule* mol;
1448 <      
1449 <    
1450 <      RealType totalMass = 0.0;
1451 <    
1452 <
1453 <      for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1454 <         RealType mass = mol->getMass();
1455 <         totalMass += mass;
1456 <         com += mass * mol->getCom();
1457 <         comVel += mass * mol->getComVel();          
1458 <      }  
1459 <      
1460 < #ifdef IS_MPI
1461 <      RealType tmpMass = totalMass;
1462 <      Vector3d tmpCom(com);  
1463 <      Vector3d tmpComVel(comVel);
1464 <      MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1465 <      MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1466 <      MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1467 < #endif
1468 <      
1469 <      com /= totalMass;
1470 <      comVel /= totalMass;
1471 <   }        
1472 <  
1473 <   /*
1474 <   Return intertia tensor for entire system and angular momentum Vector.
1475 <
1476 <
1477 <       [  Ixx -Ixy  -Ixz ]
1478 <  J =| -Iyx  Iyy  -Iyz |
1479 <       [ -Izx -Iyz   Izz ]
1480 <    */
1481 <
1482 <   void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){
1483 <      
1484 <
1485 <      RealType xx = 0.0;
1486 <      RealType yy = 0.0;
1487 <      RealType zz = 0.0;
1488 <      RealType xy = 0.0;
1489 <      RealType xz = 0.0;
1490 <      RealType yz = 0.0;
1491 <      Vector3d com(0.0);
1492 <      Vector3d comVel(0.0);
1493 <      
1494 <      getComAll(com, comVel);
1495 <      
1496 <      SimInfo::MoleculeIterator i;
1497 <      Molecule* mol;
1498 <      
1499 <      Vector3d thisq(0.0);
1500 <      Vector3d thisv(0.0);
1501 <
1502 <      RealType thisMass = 0.0;
1503 <    
1504 <      
1505 <      
1506 <  
1507 <      for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1508 <        
1509 <         thisq = mol->getCom()-com;
1510 <         thisv = mol->getComVel()-comVel;
1511 <         thisMass = mol->getMass();
1512 <         // Compute moment of intertia coefficients.
1513 <         xx += thisq[0]*thisq[0]*thisMass;
1514 <         yy += thisq[1]*thisq[1]*thisMass;
1515 <         zz += thisq[2]*thisq[2]*thisMass;
1516 <        
1517 <         // compute products of intertia
1518 <         xy += thisq[0]*thisq[1]*thisMass;
1519 <         xz += thisq[0]*thisq[2]*thisMass;
1520 <         yz += thisq[1]*thisq[2]*thisMass;
1521 <            
1522 <         angularMomentum += cross( thisq, thisv ) * thisMass;
1523 <            
1524 <      }  
1525 <      
1526 <      
1527 <      inertiaTensor(0,0) = yy + zz;
1528 <      inertiaTensor(0,1) = -xy;
1529 <      inertiaTensor(0,2) = -xz;
1530 <      inertiaTensor(1,0) = -xy;
1531 <      inertiaTensor(1,1) = xx + zz;
1532 <      inertiaTensor(1,2) = -yz;
1533 <      inertiaTensor(2,0) = -xz;
1534 <      inertiaTensor(2,1) = -yz;
1535 <      inertiaTensor(2,2) = xx + yy;
1536 <      
1537 < #ifdef IS_MPI
1538 <      Mat3x3d tmpI(inertiaTensor);
1539 <      Vector3d tmpAngMom;
1540 <      MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1541 <      MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1542 < #endif
1543 <              
1544 <      return;
1545 <   }
1546 <
1547 <   //Returns the angular momentum of the system
1548 <   Vector3d SimInfo::getAngularMomentum(){
1549 <      
1550 <      Vector3d com(0.0);
1551 <      Vector3d comVel(0.0);
1552 <      Vector3d angularMomentum(0.0);
1553 <      
1554 <      getComAll(com,comVel);
1555 <      
1556 <      SimInfo::MoleculeIterator i;
1557 <      Molecule* mol;
1558 <      
1559 <      Vector3d thisr(0.0);
1560 <      Vector3d thisp(0.0);
1561 <      
1562 <      RealType thisMass;
1563 <      
1564 <      for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {        
1565 <        thisMass = mol->getMass();
1566 <        thisr = mol->getCom()-com;
1567 <        thisp = (mol->getComVel()-comVel)*thisMass;
1568 <        
1569 <        angularMomentum += cross( thisr, thisp );
1570 <        
1571 <      }  
1572 <      
1573 < #ifdef IS_MPI
1574 <      Vector3d tmpAngMom;
1575 <      MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1576 < #endif
1577 <      
1578 <      return angularMomentum;
1579 <   }
1580 <  
1056 >  
1057    StuntDouble* SimInfo::getIOIndexToIntegrableObject(int index) {
1058 <    return IOIndexToIntegrableObject.at(index);
1058 >    if (index >= int(IOIndexToIntegrableObject.size())) {
1059 >      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    }
1068    
1069 <  void SimInfo::setIOIndexToIntegrableObject(const std::vector<StuntDouble*>& v) {
1069 >  void SimInfo::setIOIndexToIntegrableObject(const vector<StuntDouble*>& v) {
1070      IOIndexToIntegrableObject= v;
1071    }
1072  
1073 <  /* Returns the Volume of the simulation based on a ellipsoid with semi-axes
1074 <     based on the radius of gyration V=4/3*Pi*R_1*R_2*R_3
1075 <     where R_i are related to the principle inertia moments R_i = sqrt(C*I_i/N), this reduces to
1076 <     V = 4/3*Pi*(C/N)^3/2*sqrt(det(I)). See S.E. Baltazar et. al. Comp. Mat. Sci. 37 (2006) 526-536.
1077 <  */
1078 <  void SimInfo::getGyrationalVolume(RealType &volume){
1079 <    Mat3x3d intTensor;
1596 <    RealType det;
1597 <    Vector3d dummyAngMom;
1598 <    RealType sysconstants;
1599 <    RealType geomCnst;
1600 <
1601 <    geomCnst = 3.0/2.0;
1602 <    /* Get the inertial tensor and angular momentum for free*/
1603 <    getInertiaTensor(intTensor,dummyAngMom);
1604 <    
1605 <    det = intTensor.determinant();
1606 <    sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1607 <    volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(det);
1608 <    return;
1073 >  void SimInfo::calcNConstraints() {
1074 > #ifdef IS_MPI
1075 >    MPI_Allreduce(&nConstraints_, &nGlobalConstraints_, 1,  
1076 >                  MPI_INT, MPI_SUM, MPI_COMM_WORLD);
1077 > #else
1078 >    nGlobalConstraints_ =  nConstraints_;
1079 > #endif
1080    }
1081  
1082 <  void SimInfo::getGyrationalVolume(RealType &volume, RealType &detI){
1612 <    Mat3x3d intTensor;
1613 <    Vector3d dummyAngMom;
1614 <    RealType sysconstants;
1615 <    RealType geomCnst;
1082 > }//end namespace OpenMD
1083  
1617    geomCnst = 3.0/2.0;
1618    /* Get the inertial tensor and angular momentum for free*/
1619    getInertiaTensor(intTensor,dummyAngMom);
1620    
1621    detI = intTensor.determinant();
1622    sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1623    volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(detI);
1624    return;
1625  }
1626 /*
1627   void SimInfo::setStuntDoubleFromGlobalIndex(std::vector<StuntDouble*> v) {
1628      assert( v.size() == nAtoms_ + nRigidBodies_);
1629      sdByGlobalIndex_ = v;
1630    }
1631
1632    StuntDouble* SimInfo::getStuntDoubleFromGlobalIndex(int index) {
1633      //assert(index < nAtoms_ + nRigidBodies_);
1634      return sdByGlobalIndex_.at(index);
1635    }  
1636 */  
1637 }//end namespace oopse
1638

Comparing trunk/src/brains/SimInfo.cpp (property svn:keywords):
Revision 1386 by gezelter, Fri Oct 23 18:41:09 2009 UTC vs.
Revision 1987 by gezelter, Thu Apr 17 19:07:31 2014 UTC

# Line 0 | Line 1
1 + Author Id Revision Date

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines