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

Comparing:
trunk/src/brains/SimInfo.cpp (file contents), Revision 143 by chrisfen, Fri Oct 22 22:54:01 2004 UTC vs.
branches/development/src/brains/SimInfo.cpp (file contents), Revision 1547 by gezelter, Mon Apr 11 18:44:16 2011 UTC

# Line 1 | Line 1
1 < #include <stdlib.h>
2 < #include <string.h>
3 < #include <math.h>
1 > /*
2 > * 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 > * 1. Redistributions of source code must retain the above copyright
10 > *    notice, this list of conditions and the following disclaimer.
11 > *
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.
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 > *
32 > * SUPPORT OPEN SCIENCE!  If you use OpenMD or its source code in your
33 > * research, please cite the appropriate papers when you publish your
34 > * work.  Good starting points are:
35 > *                                                                      
36 > * [1]  Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).            
37 > * [2]  Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).          
38 > * [3]  Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).          
39 > * [4]  Vardeman & Gezelter, in progress (2009).                        
40 > */
41 >
42 > /**
43 > * @file SimInfo.cpp
44 > * @author    tlin
45 > * @date  11/02/2004
46 > * @version 1.0
47 > */
48  
49 < #include <iostream>
50 < using namespace std;
49 > #include <algorithm>
50 > #include <set>
51 > #include <map>
52  
53   #include "brains/SimInfo.hpp"
54 < #define __C
55 < #include "brains/fSimulation.h"
54 > #include "math/Vector3.hpp"
55 > #include "primitives/Molecule.hpp"
56 > #include "primitives/StuntDouble.hpp"
57 > #include "UseTheForce/DarkSide/neighborLists_interface.h"
58 > #include "UseTheForce/doForces_interface.h"
59 > #include "utils/MemoryUtils.hpp"
60   #include "utils/simError.h"
61 < #include "UseTheForce/DarkSide/simulation_interface.h"
62 < #include "UseTheForce/notifyCutoffs_interface.h"
61 > #include "selection/SelectionManager.hpp"
62 > #include "io/ForceFieldOptions.hpp"
63 > #include "UseTheForce/ForceField.hpp"
64 > #include "nonbonded/SwitchingFunction.hpp"
65  
15 //#include "UseTheForce/fortranWrappers.hpp"
16
17 #include "math/MatVec3.h"
18
66   #ifdef IS_MPI
67 < #include "brains/mpiSimulation.hpp"
68 < #endif
67 > #include "UseTheForce/mpiComponentPlan.h"
68 > #include "UseTheForce/DarkSide/simParallel_interface.h"
69 > #endif
70  
71 < inline double roundMe( double x ){
72 <  return ( x >= 0 ) ? floor( x + 0.5 ) : ceil( x - 0.5 );
25 < }
26 <          
27 < inline double min( double a, double b ){
28 <  return (a < b ) ? a : b;
29 < }
30 <
31 < SimInfo* currentInfo;
32 <
33 < SimInfo::SimInfo(){
34 <
35 <  n_constraints = 0;
36 <  nZconstraints = 0;
37 <  n_oriented = 0;
38 <  n_dipoles = 0;
39 <  ndf = 0;
40 <  ndfRaw = 0;
41 <  nZconstraints = 0;
42 <  the_integrator = NULL;
43 <  setTemp = 0;
44 <  thermalTime = 0.0;
45 <  currentTime = 0.0;
46 <  rCut = 0.0;
47 <  rSw = 0.0;
48 <
49 <  haveRcut = 0;
50 <  haveRsw = 0;
51 <  boxIsInit = 0;
71 > using namespace std;
72 > namespace OpenMD {
73    
74 <  resetTime = 1e99;
74 >  SimInfo::SimInfo(ForceField* ff, Globals* simParams) :
75 >    forceField_(ff), simParams_(simParams),
76 >    ndf_(0), fdf_local(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0),
77 >    nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0),
78 >    nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
79 >    nAtoms_(0), nBonds_(0),  nBends_(0), nTorsions_(0), nInversions_(0),
80 >    nRigidBodies_(0), nIntegrableObjects_(0), nCutoffGroups_(0),
81 >    nConstraints_(0), sman_(NULL), fortranInitialized_(false),
82 >    calcBoxDipole_(false), useAtomicVirial_(true) {    
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(); i !=components.end(); ++i) {
95 >      molStamp = (*i)->getMoleculeStamp();
96 >      nMolWithSameStamp = (*i)->getNMol();
97 >      
98 >      addMoleculeStamp(molStamp, nMolWithSameStamp);
99 >      
100 >      //calculate atoms in molecules
101 >      nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp;  
102 >      
103 >      //calculate atoms in cutoff groups
104 >      int nAtomsInGroups = 0;
105 >      int nCutoffGroupsInStamp = molStamp->getNCutoffGroups();
106 >      
107 >      for (int j=0; j < nCutoffGroupsInStamp; j++) {
108 >        cgStamp = molStamp->getCutoffGroupStamp(j);
109 >        nAtomsInGroups += cgStamp->getNMembers();
110 >      }
111 >      
112 >      nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
113 >      
114 >      nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;            
115 >      
116 >      //calculate atoms in rigid bodies
117 >      int nAtomsInRigidBodies = 0;
118 >      int nRigidBodiesInStamp = molStamp->getNRigidBodies();
119 >      
120 >      for (int j=0; j < nRigidBodiesInStamp; j++) {
121 >        rbStamp = molStamp->getRigidBodyStamp(j);
122 >        nAtomsInRigidBodies += rbStamp->getNMembers();
123 >      }
124 >      
125 >      nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
126 >      nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;            
127 >      
128 >    }
129 >    
130 >    //every free atom (atom does not belong to cutoff groups) is a cutoff
131 >    //group therefore the total number of cutoff groups in the system is
132 >    //equal to the total number of atoms minus number of atoms belong to
133 >    //cutoff group defined in meta-data file plus the number of cutoff
134 >    //groups defined in meta-data file
135 >    std::cerr << "nGA = " << nGlobalAtoms_ << "\n";
136 >    std::cerr << "nCA = " << nCutoffAtoms << "\n";
137 >    std::cerr << "nG = " << nGroups << "\n";
138  
139 <  orthoRhombic = 0;
56 <  orthoTolerance = 1E-6;
57 <  useInitXSstate = true;
139 >    nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
140  
141 <  usePBC = 0;
142 <  useDirectionalAtoms = 0;
143 <  useLennardJones = 0;
144 <  useElectrostatics = 0;
145 <  useCharges = 0;
146 <  useDipoles = 0;
147 <  useSticky = 0;
148 <  useGayBerne = 0;
149 <  useEAM = 0;
150 <  useShapes = 0;
151 <  useFLARB = 0;
152 <
153 <  useSolidThermInt = 0;
154 <  useLiquidThermInt = 0;
155 <
156 <  haveCutoffGroups = false;
157 <
158 <  excludes = Exclude::Instance();
159 <
160 <  myConfiguration = new SimState();
161 <
162 <  has_minimizer = false;
163 <  the_minimizer =NULL;
164 <
165 <  ngroup = 0;
84 <
85 < }
86 <
87 <
88 < SimInfo::~SimInfo(){
89 <
90 <  delete myConfiguration;
91 <
92 <  map<string, GenericData*>::iterator i;
93 <  
94 <  for(i = properties.begin(); i != properties.end(); i++)
95 <    delete (*i).second;
96 <
97 < }
98 <
99 < void SimInfo::setBox(double newBox[3]) {
100 <  
101 <  int i, j;
102 <  double tempMat[3][3];
103 <
104 <  for(i=0; i<3; i++)
105 <    for (j=0; j<3; j++) tempMat[i][j] = 0.0;;
106 <
107 <  tempMat[0][0] = newBox[0];
108 <  tempMat[1][1] = newBox[1];
109 <  tempMat[2][2] = newBox[2];
110 <
111 <  setBoxM( tempMat );
112 <
113 < }
114 <
115 < void SimInfo::setBoxM( double theBox[3][3] ){
116 <  
117 <  int i, j;
118 <  double FortranHmat[9]; // to preserve compatibility with Fortran the
119 <                         // ordering in the array is as follows:
120 <                         // [ 0 3 6 ]
121 <                         // [ 1 4 7 ]
122 <                         // [ 2 5 8 ]
123 <  double FortranHmatInv[9]; // the inverted Hmat (for Fortran);
124 <
125 <  if( !boxIsInit ) boxIsInit = 1;
126 <
127 <  for(i=0; i < 3; i++)
128 <    for (j=0; j < 3; j++) Hmat[i][j] = theBox[i][j];
129 <  
130 <  calcBoxL();
131 <  calcHmatInv();
141 >    std::cerr << "nGCG = " << nGlobalCutoffGroups_ << "\n";
142 >    
143 >    //every free atom (atom does not belong to rigid bodies) is an
144 >    //integrable object therefore the total number of integrable objects
145 >    //in the system is equal to the total number of atoms minus number of
146 >    //atoms belong to rigid body defined in meta-data file plus the number
147 >    //of rigid bodies defined in meta-data file
148 >    nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms
149 >      + nGlobalRigidBodies_;
150 >    
151 >    nGlobalMols_ = molStampIds_.size();
152 >    molToProcMap_.resize(nGlobalMols_);
153 >  }
154 >  
155 >  SimInfo::~SimInfo() {
156 >    map<int, Molecule*>::iterator i;
157 >    for (i = molecules_.begin(); i != molecules_.end(); ++i) {
158 >      delete i->second;
159 >    }
160 >    molecules_.clear();
161 >      
162 >    delete sman_;
163 >    delete simParams_;
164 >    delete forceField_;
165 >  }
166  
167 <  for(i=0; i < 3; i++) {
168 <    for (j=0; j < 3; j++) {
169 <      FortranHmat[3*j + i] = Hmat[i][j];
170 <      FortranHmatInv[3*j + i] = HmatInv[i][j];
167 >
168 >  bool SimInfo::addMolecule(Molecule* mol) {
169 >    MoleculeIterator i;
170 >    
171 >    i = molecules_.find(mol->getGlobalIndex());
172 >    if (i == molecules_.end() ) {
173 >      
174 >      molecules_.insert(make_pair(mol->getGlobalIndex(), mol));
175 >      
176 >      nAtoms_ += mol->getNAtoms();
177 >      nBonds_ += mol->getNBonds();
178 >      nBends_ += mol->getNBends();
179 >      nTorsions_ += mol->getNTorsions();
180 >      nInversions_ += mol->getNInversions();
181 >      nRigidBodies_ += mol->getNRigidBodies();
182 >      nIntegrableObjects_ += mol->getNIntegrableObjects();
183 >      nCutoffGroups_ += mol->getNCutoffGroups();
184 >      nConstraints_ += mol->getNConstraintPairs();
185 >      
186 >      addInteractionPairs(mol);
187 >      
188 >      return true;
189 >    } else {
190 >      return false;
191      }
192    }
193 +  
194 +  bool SimInfo::removeMolecule(Molecule* mol) {
195 +    MoleculeIterator i;
196 +    i = molecules_.find(mol->getGlobalIndex());
197  
198 <  setFortranBox(FortranHmat, FortranHmatInv, &orthoRhombic);
141 <
142 < }
143 <
198 >    if (i != molecules_.end() ) {
199  
200 < void SimInfo::getBoxM (double theBox[3][3]) {
200 >      assert(mol == i->second);
201 >        
202 >      nAtoms_ -= mol->getNAtoms();
203 >      nBonds_ -= mol->getNBonds();
204 >      nBends_ -= mol->getNBends();
205 >      nTorsions_ -= mol->getNTorsions();
206 >      nInversions_ -= mol->getNInversions();
207 >      nRigidBodies_ -= mol->getNRigidBodies();
208 >      nIntegrableObjects_ -= mol->getNIntegrableObjects();
209 >      nCutoffGroups_ -= mol->getNCutoffGroups();
210 >      nConstraints_ -= mol->getNConstraintPairs();
211  
212 <  int i, j;
213 <  for(i=0; i<3; i++)
149 <    for (j=0; j<3; j++) theBox[i][j] = Hmat[i][j];
150 < }
212 >      removeInteractionPairs(mol);
213 >      molecules_.erase(mol->getGlobalIndex());
214  
215 +      delete mol;
216 +        
217 +      return true;
218 +    } else {
219 +      return false;
220 +    }
221 +  }    
222  
223 < void SimInfo::scaleBox(double scale) {
224 <  double theBox[3][3];
225 <  int i, j;
223 >        
224 >  Molecule* SimInfo::beginMolecule(MoleculeIterator& i) {
225 >    i = molecules_.begin();
226 >    return i == molecules_.end() ? NULL : i->second;
227 >  }    
228  
229 <  // cerr << "Scaling box by " << scale << "\n";
229 >  Molecule* SimInfo::nextMolecule(MoleculeIterator& i) {
230 >    ++i;
231 >    return i == molecules_.end() ? NULL : i->second;    
232 >  }
233  
159  for(i=0; i<3; i++)
160    for (j=0; j<3; j++) theBox[i][j] = Hmat[i][j]*scale;
234  
235 <  setBoxM(theBox);
235 >  void SimInfo::calcNdf() {
236 >    int ndf_local;
237 >    MoleculeIterator i;
238 >    vector<StuntDouble*>::iterator j;
239 >    Molecule* mol;
240 >    StuntDouble* integrableObject;
241  
242 < }
242 >    ndf_local = 0;
243 >    
244 >    for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
245 >      for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
246 >           integrableObject = mol->nextIntegrableObject(j)) {
247  
248 < void SimInfo::calcHmatInv( void ) {
167 <  
168 <  int oldOrtho;
169 <  int i,j;
170 <  double smallDiag;
171 <  double tol;
172 <  double sanity[3][3];
248 >        ndf_local += 3;
249  
250 <  invertMat3( Hmat, HmatInv );
250 >        if (integrableObject->isDirectional()) {
251 >          if (integrableObject->isLinear()) {
252 >            ndf_local += 2;
253 >          } else {
254 >            ndf_local += 3;
255 >          }
256 >        }
257 >            
258 >      }
259 >    }
260 >    
261 >    // n_constraints is local, so subtract them on each processor
262 >    ndf_local -= nConstraints_;
263  
264 <  // check to see if Hmat is orthorhombic
265 <  
266 <  oldOrtho = orthoRhombic;
267 <
268 <  smallDiag = fabs(Hmat[0][0]);
269 <  if(smallDiag > fabs(Hmat[1][1])) smallDiag = fabs(Hmat[1][1]);
270 <  if(smallDiag > fabs(Hmat[2][2])) smallDiag = fabs(Hmat[2][2]);
271 <  tol = smallDiag * orthoTolerance;
272 <
273 <  orthoRhombic = 1;
274 <  
275 <  for (i = 0; i < 3; i++ ) {
276 <    for (j = 0 ; j < 3; j++) {
277 <      if (i != j) {
278 <        if (orthoRhombic) {
279 <          if ( fabs(Hmat[i][j]) >= tol) orthoRhombic = 0;
280 <        }        
264 > #ifdef IS_MPI
265 >    MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
266 > #else
267 >    ndf_ = ndf_local;
268 > #endif
269 >
270 >    // nZconstraints_ is global, as are the 3 COM translations for the
271 >    // entire system:
272 >    ndf_ = ndf_ - 3 - nZconstraint_;
273 >
274 >  }
275 >
276 >  int SimInfo::getFdf() {
277 > #ifdef IS_MPI
278 >    MPI_Allreduce(&fdf_local,&fdf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
279 > #else
280 >    fdf_ = fdf_local;
281 > #endif
282 >    return fdf_;
283 >  }
284 >    
285 >  void SimInfo::calcNdfRaw() {
286 >    int ndfRaw_local;
287 >
288 >    MoleculeIterator i;
289 >    vector<StuntDouble*>::iterator j;
290 >    Molecule* mol;
291 >    StuntDouble* integrableObject;
292 >
293 >    // Raw degrees of freedom that we have to set
294 >    ndfRaw_local = 0;
295 >    
296 >    for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
297 >      for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
298 >           integrableObject = mol->nextIntegrableObject(j)) {
299 >
300 >        ndfRaw_local += 3;
301 >
302 >        if (integrableObject->isDirectional()) {
303 >          if (integrableObject->isLinear()) {
304 >            ndfRaw_local += 2;
305 >          } else {
306 >            ndfRaw_local += 3;
307 >          }
308 >        }
309 >            
310 >      }
311 >    }
312 >    
313 > #ifdef IS_MPI
314 >    MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
315 > #else
316 >    ndfRaw_ = ndfRaw_local;
317 > #endif
318 >  }
319 >
320 >  void SimInfo::calcNdfTrans() {
321 >    int ndfTrans_local;
322 >
323 >    ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
324 >
325 >
326 > #ifdef IS_MPI
327 >    MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
328 > #else
329 >    ndfTrans_ = ndfTrans_local;
330 > #endif
331 >
332 >    ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
333 >
334 >  }
335 >
336 >  void SimInfo::addInteractionPairs(Molecule* mol) {
337 >    ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
338 >    vector<Bond*>::iterator bondIter;
339 >    vector<Bend*>::iterator bendIter;
340 >    vector<Torsion*>::iterator torsionIter;
341 >    vector<Inversion*>::iterator inversionIter;
342 >    Bond* bond;
343 >    Bend* bend;
344 >    Torsion* torsion;
345 >    Inversion* inversion;
346 >    int a;
347 >    int b;
348 >    int c;
349 >    int d;
350 >
351 >    // atomGroups can be used to add special interaction maps between
352 >    // groups of atoms that are in two separate rigid bodies.
353 >    // However, most site-site interactions between two rigid bodies
354 >    // are probably not special, just the ones between the physically
355 >    // bonded atoms.  Interactions *within* a single rigid body should
356 >    // always be excluded.  These are done at the bottom of this
357 >    // function.
358 >
359 >    map<int, set<int> > atomGroups;
360 >    Molecule::RigidBodyIterator rbIter;
361 >    RigidBody* rb;
362 >    Molecule::IntegrableObjectIterator ii;
363 >    StuntDouble* integrableObject;
364 >    
365 >    for (integrableObject = mol->beginIntegrableObject(ii);
366 >         integrableObject != NULL;
367 >         integrableObject = mol->nextIntegrableObject(ii)) {
368 >      
369 >      if (integrableObject->isRigidBody()) {
370 >        rb = static_cast<RigidBody*>(integrableObject);
371 >        vector<Atom*> atoms = rb->getAtoms();
372 >        set<int> rigidAtoms;
373 >        for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
374 >          rigidAtoms.insert(atoms[i]->getGlobalIndex());
375 >        }
376 >        for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
377 >          atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
378 >        }      
379 >      } else {
380 >        set<int> oneAtomSet;
381 >        oneAtomSet.insert(integrableObject->getGlobalIndex());
382 >        atomGroups.insert(map<int, set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));        
383        }
384 <    }
385 <  }
384 >    }  
385 >          
386 >    for (bond= mol->beginBond(bondIter); bond != NULL;
387 >         bond = mol->nextBond(bondIter)) {
388  
389 <  if( oldOrtho != orthoRhombic ){
389 >      a = bond->getAtomA()->getGlobalIndex();
390 >      b = bond->getAtomB()->getGlobalIndex();  
391      
392 <    if( orthoRhombic ) {
393 <      sprintf( painCave.errMsg,
394 <               "OOPSE is switching from the default Non-Orthorhombic\n"
395 <               "\tto the faster Orthorhombic periodic boundary computations.\n"
396 <               "\tThis is usually a good thing, but if you wan't the\n"
204 <               "\tNon-Orthorhombic computations, make the orthoBoxTolerance\n"
205 <               "\tvariable ( currently set to %G ) smaller.\n",
206 <               orthoTolerance);
207 <      painCave.severity = OOPSE_INFO;
208 <      simError();
392 >      if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
393 >        oneTwoInteractions_.addPair(a, b);
394 >      } else {
395 >        excludedInteractions_.addPair(a, b);
396 >      }
397      }
210    else {
211      sprintf( painCave.errMsg,
212               "OOPSE is switching from the faster Orthorhombic to the more\n"
213               "\tflexible Non-Orthorhombic periodic boundary computations.\n"
214               "\tThis is usually because the box has deformed under\n"
215               "\tNPTf integration. If you wan't to live on the edge with\n"
216               "\tthe Orthorhombic computations, make the orthoBoxTolerance\n"
217               "\tvariable ( currently set to %G ) larger.\n",
218               orthoTolerance);
219      painCave.severity = OOPSE_WARNING;
220      simError();
221    }
222  }
223 }
398  
399 < void SimInfo::calcBoxL( void ){
399 >    for (bend= mol->beginBend(bendIter); bend != NULL;
400 >         bend = mol->nextBend(bendIter)) {
401  
402 <  double dx, dy, dz, dsq;
402 >      a = bend->getAtomA()->getGlobalIndex();
403 >      b = bend->getAtomB()->getGlobalIndex();        
404 >      c = bend->getAtomC()->getGlobalIndex();
405 >      
406 >      if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
407 >        oneTwoInteractions_.addPair(a, b);      
408 >        oneTwoInteractions_.addPair(b, c);
409 >      } else {
410 >        excludedInteractions_.addPair(a, b);
411 >        excludedInteractions_.addPair(b, c);
412 >      }
413  
414 <  // boxVol = Determinant of Hmat
414 >      if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
415 >        oneThreeInteractions_.addPair(a, c);      
416 >      } else {
417 >        excludedInteractions_.addPair(a, c);
418 >      }
419 >    }
420  
421 <  boxVol = matDet3( Hmat );
421 >    for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
422 >         torsion = mol->nextTorsion(torsionIter)) {
423  
424 <  // boxLx
425 <  
426 <  dx = Hmat[0][0]; dy = Hmat[1][0]; dz = Hmat[2][0];
427 <  dsq = dx*dx + dy*dy + dz*dz;
237 <  boxL[0] = sqrt( dsq );
238 <  //maxCutoff = 0.5 * boxL[0];
424 >      a = torsion->getAtomA()->getGlobalIndex();
425 >      b = torsion->getAtomB()->getGlobalIndex();        
426 >      c = torsion->getAtomC()->getGlobalIndex();        
427 >      d = torsion->getAtomD()->getGlobalIndex();      
428  
429 <  // boxLy
430 <  
431 <  dx = Hmat[0][1]; dy = Hmat[1][1]; dz = Hmat[2][1];
432 <  dsq = dx*dx + dy*dy + dz*dz;
433 <  boxL[1] = sqrt( dsq );
434 <  //if( (0.5 * boxL[1]) < maxCutoff ) maxCutoff = 0.5 * boxL[1];
429 >      if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
430 >        oneTwoInteractions_.addPair(a, b);      
431 >        oneTwoInteractions_.addPair(b, c);
432 >        oneTwoInteractions_.addPair(c, d);
433 >      } else {
434 >        excludedInteractions_.addPair(a, b);
435 >        excludedInteractions_.addPair(b, c);
436 >        excludedInteractions_.addPair(c, d);
437 >      }
438  
439 +      if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
440 +        oneThreeInteractions_.addPair(a, c);      
441 +        oneThreeInteractions_.addPair(b, d);      
442 +      } else {
443 +        excludedInteractions_.addPair(a, c);
444 +        excludedInteractions_.addPair(b, d);
445 +      }
446  
447 <  // boxLz
448 <  
449 <  dx = Hmat[0][2]; dy = Hmat[1][2]; dz = Hmat[2][2];
450 <  dsq = dx*dx + dy*dy + dz*dz;
451 <  boxL[2] = sqrt( dsq );
452 <  //if( (0.5 * boxL[2]) < maxCutoff ) maxCutoff = 0.5 * boxL[2];
447 >      if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
448 >        oneFourInteractions_.addPair(a, d);      
449 >      } else {
450 >        excludedInteractions_.addPair(a, d);
451 >      }
452 >    }
453  
454 <  //calculate the max cutoff
455 <  maxCutoff =  calcMaxCutOff();
257 <  
258 <  checkCutOffs();
454 >    for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
455 >         inversion = mol->nextInversion(inversionIter)) {
456  
457 < }
457 >      a = inversion->getAtomA()->getGlobalIndex();
458 >      b = inversion->getAtomB()->getGlobalIndex();        
459 >      c = inversion->getAtomC()->getGlobalIndex();        
460 >      d = inversion->getAtomD()->getGlobalIndex();        
461  
462 +      if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
463 +        oneTwoInteractions_.addPair(a, b);      
464 +        oneTwoInteractions_.addPair(a, c);
465 +        oneTwoInteractions_.addPair(a, d);
466 +      } else {
467 +        excludedInteractions_.addPair(a, b);
468 +        excludedInteractions_.addPair(a, c);
469 +        excludedInteractions_.addPair(a, d);
470 +      }
471  
472 < double SimInfo::calcMaxCutOff(){
472 >      if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
473 >        oneThreeInteractions_.addPair(b, c);    
474 >        oneThreeInteractions_.addPair(b, d);    
475 >        oneThreeInteractions_.addPair(c, d);      
476 >      } else {
477 >        excludedInteractions_.addPair(b, c);
478 >        excludedInteractions_.addPair(b, d);
479 >        excludedInteractions_.addPair(c, d);
480 >      }
481 >    }
482  
483 <  double ri[3], rj[3], rk[3];
484 <  double rij[3], rjk[3], rki[3];
485 <  double minDist;
483 >    for (rb = mol->beginRigidBody(rbIter); rb != NULL;
484 >         rb = mol->nextRigidBody(rbIter)) {
485 >      vector<Atom*> atoms = rb->getAtoms();
486 >      for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
487 >        for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
488 >          a = atoms[i]->getGlobalIndex();
489 >          b = atoms[j]->getGlobalIndex();
490 >          excludedInteractions_.addPair(a, b);
491 >        }
492 >      }
493 >    }        
494  
495 <  ri[0] = Hmat[0][0];
270 <  ri[1] = Hmat[1][0];
271 <  ri[2] = Hmat[2][0];
495 >  }
496  
497 <  rj[0] = Hmat[0][1];
498 <  rj[1] = Hmat[1][1];
499 <  rj[2] = Hmat[2][1];
497 >  void SimInfo::removeInteractionPairs(Molecule* mol) {
498 >    ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
499 >    vector<Bond*>::iterator bondIter;
500 >    vector<Bend*>::iterator bendIter;
501 >    vector<Torsion*>::iterator torsionIter;
502 >    vector<Inversion*>::iterator inversionIter;
503 >    Bond* bond;
504 >    Bend* bend;
505 >    Torsion* torsion;
506 >    Inversion* inversion;
507 >    int a;
508 >    int b;
509 >    int c;
510 >    int d;
511  
512 <  rk[0] = Hmat[0][2];
513 <  rk[1] = Hmat[1][2];
514 <  rk[2] = Hmat[2][2];
512 >    map<int, set<int> > atomGroups;
513 >    Molecule::RigidBodyIterator rbIter;
514 >    RigidBody* rb;
515 >    Molecule::IntegrableObjectIterator ii;
516 >    StuntDouble* integrableObject;
517      
518 <  crossProduct3(ri, rj, rij);
519 <  distXY = dotProduct3(rk,rij) / norm3(rij);
518 >    for (integrableObject = mol->beginIntegrableObject(ii);
519 >         integrableObject != NULL;
520 >         integrableObject = mol->nextIntegrableObject(ii)) {
521 >      
522 >      if (integrableObject->isRigidBody()) {
523 >        rb = static_cast<RigidBody*>(integrableObject);
524 >        vector<Atom*> atoms = rb->getAtoms();
525 >        set<int> rigidAtoms;
526 >        for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
527 >          rigidAtoms.insert(atoms[i]->getGlobalIndex());
528 >        }
529 >        for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
530 >          atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
531 >        }      
532 >      } else {
533 >        set<int> oneAtomSet;
534 >        oneAtomSet.insert(integrableObject->getGlobalIndex());
535 >        atomGroups.insert(map<int, set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));        
536 >      }
537 >    }  
538  
539 <  crossProduct3(rj,rk, rjk);
540 <  distYZ = dotProduct3(ri,rjk) / norm3(rjk);
539 >    for (bond= mol->beginBond(bondIter); bond != NULL;
540 >         bond = mol->nextBond(bondIter)) {
541 >      
542 >      a = bond->getAtomA()->getGlobalIndex();
543 >      b = bond->getAtomB()->getGlobalIndex();  
544 >    
545 >      if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
546 >        oneTwoInteractions_.removePair(a, b);
547 >      } else {
548 >        excludedInteractions_.removePair(a, b);
549 >      }
550 >    }
551  
552 <  crossProduct3(rk,ri, rki);
553 <  distZX = dotProduct3(rj,rki) / norm3(rki);
552 >    for (bend= mol->beginBend(bendIter); bend != NULL;
553 >         bend = mol->nextBend(bendIter)) {
554  
555 <  minDist = min(min(distXY, distYZ), distZX);
556 <  return minDist/2;
557 <  
558 < }
555 >      a = bend->getAtomA()->getGlobalIndex();
556 >      b = bend->getAtomB()->getGlobalIndex();        
557 >      c = bend->getAtomC()->getGlobalIndex();
558 >      
559 >      if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
560 >        oneTwoInteractions_.removePair(a, b);      
561 >        oneTwoInteractions_.removePair(b, c);
562 >      } else {
563 >        excludedInteractions_.removePair(a, b);
564 >        excludedInteractions_.removePair(b, c);
565 >      }
566  
567 < void SimInfo::wrapVector( double thePos[3] ){
567 >      if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
568 >        oneThreeInteractions_.removePair(a, c);      
569 >      } else {
570 >        excludedInteractions_.removePair(a, c);
571 >      }
572 >    }
573  
574 <  int i;
575 <  double scaled[3];
574 >    for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
575 >         torsion = mol->nextTorsion(torsionIter)) {
576  
577 <  if( !orthoRhombic ){
578 <    // calc the scaled coordinates.
577 >      a = torsion->getAtomA()->getGlobalIndex();
578 >      b = torsion->getAtomB()->getGlobalIndex();        
579 >      c = torsion->getAtomC()->getGlobalIndex();        
580 >      d = torsion->getAtomD()->getGlobalIndex();      
581    
582 +      if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
583 +        oneTwoInteractions_.removePair(a, b);      
584 +        oneTwoInteractions_.removePair(b, c);
585 +        oneTwoInteractions_.removePair(c, d);
586 +      } else {
587 +        excludedInteractions_.removePair(a, b);
588 +        excludedInteractions_.removePair(b, c);
589 +        excludedInteractions_.removePair(c, d);
590 +      }
591  
592 <    matVecMul3(HmatInv, thePos, scaled);
593 <    
594 <    for(i=0; i<3; i++)
595 <      scaled[i] -= roundMe(scaled[i]);
596 <    
597 <    // calc the wrapped real coordinates from the wrapped scaled coordinates
598 <    
311 <    matVecMul3(Hmat, scaled, thePos);
592 >      if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
593 >        oneThreeInteractions_.removePair(a, c);      
594 >        oneThreeInteractions_.removePair(b, d);      
595 >      } else {
596 >        excludedInteractions_.removePair(a, c);
597 >        excludedInteractions_.removePair(b, d);
598 >      }
599  
600 <  }
601 <  else{
602 <    // calc the scaled coordinates.
603 <    
604 <    for(i=0; i<3; i++)
605 <      scaled[i] = thePos[i]*HmatInv[i][i];
319 <    
320 <    // wrap the scaled coordinates
321 <    
322 <    for(i=0; i<3; i++)
323 <      scaled[i] -= roundMe(scaled[i]);
324 <    
325 <    // calc the wrapped real coordinates from the wrapped scaled coordinates
326 <    
327 <    for(i=0; i<3; i++)
328 <      thePos[i] = scaled[i]*Hmat[i][i];
329 <  }
330 <    
331 < }
600 >      if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
601 >        oneFourInteractions_.removePair(a, d);      
602 >      } else {
603 >        excludedInteractions_.removePair(a, d);
604 >      }
605 >    }
606  
607 +    for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
608 +         inversion = mol->nextInversion(inversionIter)) {
609  
610 < int SimInfo::getNDF(){
611 <  int ndf_local;
610 >      a = inversion->getAtomA()->getGlobalIndex();
611 >      b = inversion->getAtomB()->getGlobalIndex();        
612 >      c = inversion->getAtomC()->getGlobalIndex();        
613 >      d = inversion->getAtomD()->getGlobalIndex();        
614  
615 <  ndf_local = 0;
616 <  
617 <  for(int i = 0; i < integrableObjects.size(); i++){
618 <    ndf_local += 3;
619 <    if (integrableObjects[i]->isDirectional()) {
620 <      if (integrableObjects[i]->isLinear())
621 <        ndf_local += 2;
622 <      else
623 <        ndf_local += 3;
615 >      if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
616 >        oneTwoInteractions_.removePair(a, b);      
617 >        oneTwoInteractions_.removePair(a, c);
618 >        oneTwoInteractions_.removePair(a, d);
619 >      } else {
620 >        excludedInteractions_.removePair(a, b);
621 >        excludedInteractions_.removePair(a, c);
622 >        excludedInteractions_.removePair(a, d);
623 >      }
624 >
625 >      if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
626 >        oneThreeInteractions_.removePair(b, c);    
627 >        oneThreeInteractions_.removePair(b, d);    
628 >        oneThreeInteractions_.removePair(c, d);      
629 >      } else {
630 >        excludedInteractions_.removePair(b, c);
631 >        excludedInteractions_.removePair(b, d);
632 >        excludedInteractions_.removePair(c, d);
633 >      }
634      }
635 +
636 +    for (rb = mol->beginRigidBody(rbIter); rb != NULL;
637 +         rb = mol->nextRigidBody(rbIter)) {
638 +      vector<Atom*> atoms = rb->getAtoms();
639 +      for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
640 +        for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
641 +          a = atoms[i]->getGlobalIndex();
642 +          b = atoms[j]->getGlobalIndex();
643 +          excludedInteractions_.removePair(a, b);
644 +        }
645 +      }
646 +    }        
647 +    
648    }
649 +  
650 +  
651 +  void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) {
652 +    int curStampId;
653 +    
654 +    //index from 0
655 +    curStampId = moleculeStamps_.size();
656  
657 <  // n_constraints is local, so subtract them on each processor:
657 >    moleculeStamps_.push_back(molStamp);
658 >    molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
659 >  }
660  
661 <  ndf_local -= n_constraints;
661 >
662 >  /**
663 >   * update
664 >   *
665 >   *  Performs the global checks and variable settings after the
666 >   *  objects have been created.
667 >   *
668 >   */
669 >  void SimInfo::update() {  
670 >    setupSimVariables();
671 >    calcNdf();
672 >    calcNdfRaw();
673 >    calcNdfTrans();
674 >  }
675 >  
676 >  /**
677 >   * getSimulatedAtomTypes
678 >   *
679 >   * Returns an STL set of AtomType* that are actually present in this
680 >   * simulation.  Must query all processors to assemble this information.
681 >   *
682 >   */
683 >  set<AtomType*> SimInfo::getSimulatedAtomTypes() {
684 >    SimInfo::MoleculeIterator mi;
685 >    Molecule* mol;
686 >    Molecule::AtomIterator ai;
687 >    Atom* atom;
688 >    set<AtomType*> atomTypes;
689 >    
690 >    for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {      
691 >      for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
692 >        atomTypes.insert(atom->getAtomType());
693 >      }      
694 >    }    
695  
696   #ifdef IS_MPI
354  MPI_Allreduce(&ndf_local,&ndf,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
355 #else
356  ndf = ndf_local;
357 #endif
697  
698 <  // nZconstraints is global, as are the 3 COM translations for the
699 <  // entire system:
698 >    // loop over the found atom types on this processor, and add their
699 >    // numerical idents to a vector:
700  
701 <  ndf = ndf - 3 - nZconstraints;
701 >    vector<int> foundTypes;
702 >    set<AtomType*>::iterator i;
703 >    for (i = atomTypes.begin(); i != atomTypes.end(); ++i)
704 >      foundTypes.push_back( (*i)->getIdent() );
705  
706 <  return ndf;
707 < }
706 >    // count_local holds the number of found types on this processor
707 >    int count_local = foundTypes.size();
708  
709 < int SimInfo::getNDFraw() {
710 <  int ndfRaw_local;
709 >    // count holds the total number of found types on all processors
710 >    // (some will be redundant with the ones found locally):
711 >    int count;
712 >    MPI::COMM_WORLD.Allreduce(&count_local, &count, 1, MPI::INT, MPI::SUM);
713  
714 <  // Raw degrees of freedom that we have to set
715 <  ndfRaw_local = 0;
714 >    // create a vector to hold the globally found types, and resize it:
715 >    vector<int> ftGlobal;
716 >    ftGlobal.resize(count);
717 >    vector<int> counts;
718  
719 <  for(int i = 0; i < integrableObjects.size(); i++){
720 <    ndfRaw_local += 3;
721 <    if (integrableObjects[i]->isDirectional()) {
722 <       if (integrableObjects[i]->isLinear())
723 <        ndfRaw_local += 2;
724 <      else
379 <        ndfRaw_local += 3;
380 <    }
381 <  }
719 >    int nproc = MPI::COMM_WORLD.Get_size();
720 >    counts.resize(nproc);
721 >    vector<int> disps;
722 >    disps.resize(nproc);
723 >
724 >    // now spray out the foundTypes to all the other processors:
725      
726 < #ifdef IS_MPI
727 <  MPI_Allreduce(&ndfRaw_local,&ndfRaw,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
728 < #else
729 <  ndfRaw = ndfRaw_local;
726 >    MPI::COMM_WORLD.Allgatherv(&foundTypes[0], count_local, MPI::INT,
727 >                               &ftGlobal[0], &counts[0], &disps[0], MPI::INT);
728 >
729 >    // foundIdents is a stl set, so inserting an already found ident
730 >    // will have no effect.
731 >    set<int> foundIdents;
732 >    vector<int>::iterator j;
733 >    for (j = ftGlobal.begin(); j != ftGlobal.end(); ++j)
734 >      foundIdents.insert((*j));
735 >    
736 >    // now iterate over the foundIdents and get the actual atom types
737 >    // that correspond to these:
738 >    set<int>::iterator it;
739 >    for (it = foundIdents.begin(); it != foundIdents.end(); ++it)
740 >      atomTypes.insert( forceField_->getAtomType((*it)) );
741 >
742   #endif
743 +    
744 +    return atomTypes;        
745 +  }
746  
747 <  return ndfRaw;
748 < }
747 >  void SimInfo::setupSimVariables() {
748 >    useAtomicVirial_ = simParams_->getUseAtomicVirial();
749 >    // we only call setAccumulateBoxDipole if the accumulateBoxDipole parameter is true
750 >    calcBoxDipole_ = false;
751 >    if ( simParams_->haveAccumulateBoxDipole() )
752 >      if ( simParams_->getAccumulateBoxDipole() ) {
753 >        calcBoxDipole_ = true;      
754 >      }
755  
756 < int SimInfo::getNDFtranslational() {
757 <  int ndfTrans_local;
756 >    set<AtomType*>::iterator i;
757 >    set<AtomType*> atomTypes;
758 >    atomTypes = getSimulatedAtomTypes();    
759 >    int usesElectrostatic = 0;
760 >    int usesMetallic = 0;
761 >    int usesDirectional = 0;
762 >    //loop over all of the atom types
763 >    for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
764 >      usesElectrostatic |= (*i)->isElectrostatic();
765 >      usesMetallic |= (*i)->isMetal();
766 >      usesDirectional |= (*i)->isDirectional();
767 >    }
768  
769 <  ndfTrans_local = 3 * integrableObjects.size() - n_constraints;
769 > #ifdef IS_MPI    
770 >    int temp;
771 >    temp = usesDirectional;
772 >    MPI_Allreduce(&temp, &usesDirectionalAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
773  
774 +    temp = usesMetallic;
775 +    MPI_Allreduce(&temp, &usesMetallicAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
776  
777 < #ifdef IS_MPI
778 <  MPI_Allreduce(&ndfTrans_local,&ndfTrans,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
400 < #else
401 <  ndfTrans = ndfTrans_local;
777 >    temp = usesElectrostatic;
778 >    MPI_Allreduce(&temp, &usesElectrostaticAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
779   #endif
780 +    fInfo_.SIM_uses_PBC = usesPeriodicBoundaries_;    
781 +    fInfo_.SIM_uses_DirectionalAtoms = usesDirectionalAtoms_;
782 +    fInfo_.SIM_uses_MetallicAtoms = usesMetallicAtoms_;
783 +    fInfo_.SIM_requires_SkipCorrection = usesElectrostaticAtoms_;
784 +    fInfo_.SIM_requires_SelfCorrection = usesElectrostaticAtoms_;
785 +    fInfo_.SIM_uses_AtomicVirial = usesAtomicVirial_;
786 +  }
787  
788 <  ndfTrans = ndfTrans - 3 - nZconstraints;
788 >  void SimInfo::setupFortran() {
789 >    int isError;
790 >    int nExclude, nOneTwo, nOneThree, nOneFour;
791 >    vector<int> fortranGlobalGroupMembership;
792 >    
793 >    isError = 0;
794  
795 <  return ndfTrans;
796 < }
795 >    //globalGroupMembership_ is filled by SimCreator    
796 >    for (int i = 0; i < nGlobalAtoms_; i++) {
797 >      fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
798 >    }
799  
800 < int SimInfo::getTotIntegrableObjects() {
801 <  int nObjs_local;
802 <  int nObjs;
800 >    //calculate mass ratio of cutoff group
801 >    vector<RealType> mfact;
802 >    SimInfo::MoleculeIterator mi;
803 >    Molecule* mol;
804 >    Molecule::CutoffGroupIterator ci;
805 >    CutoffGroup* cg;
806 >    Molecule::AtomIterator ai;
807 >    Atom* atom;
808 >    RealType totalMass;
809  
810 <  nObjs_local =  integrableObjects.size();
810 >    //to avoid memory reallocation, reserve enough space for mfact
811 >    mfact.reserve(getNCutoffGroups());
812 >    
813 >    for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {        
814 >      for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
815  
816 +        totalMass = cg->getMass();
817 +        for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
818 +          // Check for massless groups - set mfact to 1 if true
819 +          if (totalMass != 0)
820 +            mfact.push_back(atom->getMass()/totalMass);
821 +          else
822 +            mfact.push_back( 1.0 );
823 +        }
824 +      }      
825 +    }
826  
827 < #ifdef IS_MPI
417 <  MPI_Allreduce(&nObjs_local,&nObjs,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
418 < #else
419 <  nObjs = nObjs_local;
420 < #endif
827 >    // Build the identArray_
828  
829 +    identArray_.clear();
830 +    identArray_.reserve(getNAtoms());    
831 +    for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {        
832 +      for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
833 +        identArray_.push_back(atom->getIdent());
834 +      }
835 +    }    
836  
837 <  return nObjs;
838 < }
837 >    //fill molMembershipArray
838 >    //molMembershipArray is filled by SimCreator    
839 >    vector<int> molMembershipArray(nGlobalAtoms_);
840 >    for (int i = 0; i < nGlobalAtoms_; i++) {
841 >      molMembershipArray[i] = globalMolMembership_[i] + 1;
842 >    }
843 >    
844 >    //setup fortran simulation
845  
846 < void SimInfo::refreshSim(){
846 >    nExclude = excludedInteractions_.getSize();
847 >    nOneTwo = oneTwoInteractions_.getSize();
848 >    nOneThree = oneThreeInteractions_.getSize();
849 >    nOneFour = oneFourInteractions_.getSize();
850  
851 <  simtype fInfo;
852 <  int isError;
853 <  int n_global;
854 <  int* excl;
851 >    int* excludeList = excludedInteractions_.getPairList();
852 >    int* oneTwoList = oneTwoInteractions_.getPairList();
853 >    int* oneThreeList = oneThreeInteractions_.getPairList();
854 >    int* oneFourList = oneFourInteractions_.getPairList();
855  
856 <  fInfo.dielect = 0.0;
856 >    setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray_[0],
857 >                   &nExclude, excludeList,
858 >                   &nOneTwo, oneTwoList,
859 >                   &nOneThree, oneThreeList,
860 >                   &nOneFour, oneFourList,
861 >                   &molMembershipArray[0], &mfact[0], &nCutoffGroups_,
862 >                   &fortranGlobalGroupMembership[0], &isError);
863 >    
864 >    if( isError ){
865 >      
866 >      sprintf( painCave.errMsg,
867 >               "There was an error setting the simulation information in fortran.\n" );
868 >      painCave.isFatal = 1;
869 >      painCave.severity = OPENMD_ERROR;
870 >      simError();
871 >    }
872 >    
873 >    
874 >    sprintf( checkPointMsg,
875 >             "succesfully sent the simulation information to fortran.\n");
876 >    
877 >    errorCheckPoint();
878 >    
879 >    // Setup number of neighbors in neighbor list if present
880 >    if (simParams_->haveNeighborListNeighbors()) {
881 >      int nlistNeighbors = simParams_->getNeighborListNeighbors();
882 >      setNeighbors(&nlistNeighbors);
883 >    }
884 >  
885 > #ifdef IS_MPI    
886 >    //SimInfo is responsible for creating localToGlobalAtomIndex and
887 >    //localToGlobalGroupIndex
888 >    vector<int> localToGlobalAtomIndex(getNAtoms(), 0);
889 >    vector<int> localToGlobalCutoffGroupIndex;
890 >    mpiSimData parallelData;
891  
892 <  if( useDipoles ){
436 <    if( useReactionField )fInfo.dielect = dielectric;
437 <  }
892 >    for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) {
893  
894 <  fInfo.SIM_uses_PBC = usePBC;
894 >      //local index(index in DataStorge) of atom is important
895 >      for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
896 >        localToGlobalAtomIndex[atom->getLocalIndex()] = atom->getGlobalIndex() + 1;
897 >      }
898  
899 <  if (useSticky || useDipoles || useGayBerne || useShapes) {
900 <    useDirectionalAtoms = 1;
901 <    fInfo.SIM_uses_DirectionalAtoms = useDirectionalAtoms;
902 <  }
899 >      //local index of cutoff group is trivial, it only depends on the order of travesing
900 >      for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
901 >        localToGlobalCutoffGroupIndex.push_back(cg->getGlobalIndex() + 1);
902 >      }        
903 >        
904 >    }
905  
906 <  fInfo.SIM_uses_LennardJones = useLennardJones;
906 >    //fill up mpiSimData struct
907 >    parallelData.nMolGlobal = getNGlobalMolecules();
908 >    parallelData.nMolLocal = getNMolecules();
909 >    parallelData.nAtomsGlobal = getNGlobalAtoms();
910 >    parallelData.nAtomsLocal = getNAtoms();
911 >    parallelData.nGroupsGlobal = getNGlobalCutoffGroups();
912 >    parallelData.nGroupsLocal = getNCutoffGroups();
913 >    parallelData.myNode = worldRank;
914 >    MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors));
915  
916 <  if (useCharges || useDipoles) {
917 <    useElectrostatics = 1;
918 <    fInfo.SIM_uses_Electrostatics = useElectrostatics;
919 <  }
916 >    //pass mpiSimData struct and index arrays to fortran
917 >    setFsimParallel(&parallelData, &(parallelData.nAtomsLocal),
918 >                    &localToGlobalAtomIndex[0],  &(parallelData.nGroupsLocal),
919 >                    &localToGlobalCutoffGroupIndex[0], &isError);
920  
921 <  fInfo.SIM_uses_Charges = useCharges;
922 <  fInfo.SIM_uses_Dipoles = useDipoles;
923 <  fInfo.SIM_uses_Sticky = useSticky;
924 <  fInfo.SIM_uses_GayBerne = useGayBerne;
925 <  fInfo.SIM_uses_EAM = useEAM;
926 <  fInfo.SIM_uses_Shapes = useShapes;
459 <  fInfo.SIM_uses_FLARB = useFLARB;
460 <  fInfo.SIM_uses_RF = useReactionField;
921 >    if (isError) {
922 >      sprintf(painCave.errMsg,
923 >              "mpiRefresh errror: fortran didn't like something we gave it.\n");
924 >      painCave.isFatal = 1;
925 >      simError();
926 >    }
927  
928 <  n_exclude = excludes->getSize();
929 <  excl = excludes->getFortranArray();
464 <  
465 < #ifdef IS_MPI
466 <  n_global = mpiSim->getNAtomsGlobal();
467 < #else
468 <  n_global = n_atoms;
928 >    sprintf(checkPointMsg, " mpiRefresh successful.\n");
929 >    errorCheckPoint();
930   #endif
470  
471  isError = 0;
472  
473  getFortranGroupArrays(this, FglobalGroupMembership, mfact);
474  //it may not be a good idea to pass the address of first element in vector
475  //since c++ standard does not require vector to be stored continuously in meomory
476  //Most of the compilers will organize the memory of vector continuously
477  setFortranSim( &fInfo, &n_global, &n_atoms, identArray, &n_exclude, excl,
478                  &nGlobalExcludes, globalExcludes, molMembershipArray,
479                  &mfact[0], &ngroup, &FglobalGroupMembership[0], &isError);
931  
932 <  if( isError ){
933 <    
934 <    sprintf( painCave.errMsg,
935 <             "There was an error setting the simulation information in fortran.\n" );
936 <    painCave.isFatal = 1;
937 <    painCave.severity = OOPSE_ERROR;
938 <    simError();
932 >    initFortranFF(&isError);
933 >    if (isError) {
934 >      sprintf(painCave.errMsg,
935 >              "initFortranFF errror: fortran didn't like something we gave it.\n");
936 >      painCave.isFatal = 1;
937 >      simError();
938 >    }
939 >    fortranInitialized_ = true;
940    }
489  
490 #ifdef IS_MPI
491  sprintf( checkPointMsg,
492           "succesfully sent the simulation information to fortran.\n");
493  MPIcheckPoint();
494 #endif // is_mpi
495  
496  this->ndf = this->getNDF();
497  this->ndfRaw = this->getNDFraw();
498  this->ndfTrans = this->getNDFtranslational();
499 }
941  
942 < void SimInfo::setDefaultRcut( double theRcut ){
943 <  
944 <  haveRcut = 1;
504 <  rCut = theRcut;
505 <  rList = rCut + 1.0;
506 <  
507 <  notifyFortranCutoffs( &rCut, &rSw, &rList );
508 < }
942 >  void SimInfo::addProperty(GenericData* genData) {
943 >    properties_.addProperty(genData);  
944 >  }
945  
946 < void SimInfo::setDefaultRcut( double theRcut, double theRsw ){
946 >  void SimInfo::removeProperty(const string& propName) {
947 >    properties_.removeProperty(propName);  
948 >  }
949  
950 <  rSw = theRsw;
951 <  setDefaultRcut( theRcut );
952 < }
950 >  void SimInfo::clearProperties() {
951 >    properties_.clearProperties();
952 >  }
953  
954 +  vector<string> SimInfo::getPropertyNames() {
955 +    return properties_.getPropertyNames();  
956 +  }
957 +      
958 +  vector<GenericData*> SimInfo::getProperties() {
959 +    return properties_.getProperties();
960 +  }
961  
962 < void SimInfo::checkCutOffs( void ){
963 <  
964 <  if( boxIsInit ){
965 <    
966 <    //we need to check cutOffs against the box
967 <    
968 <    if( rCut > maxCutoff ){
969 <      sprintf( painCave.errMsg,
970 <               "cutoffRadius is too large for the current periodic box.\n"
971 <               "\tCurrent Value of cutoffRadius = %G at time %G\n "
972 <               "\tThis is larger than half of at least one of the\n"
973 <               "\tperiodic box vectors.  Right now, the Box matrix is:\n"
974 <               "\n"
975 <               "\t[ %G %G %G ]\n"
976 <               "\t[ %G %G %G ]\n"
977 <               "\t[ %G %G %G ]\n",
978 <               rCut, currentTime,
979 <               Hmat[0][0], Hmat[0][1], Hmat[0][2],
980 <               Hmat[1][0], Hmat[1][1], Hmat[1][2],
981 <               Hmat[2][0], Hmat[2][1], Hmat[2][2]);
982 <      painCave.severity = OOPSE_ERROR;
983 <      painCave.isFatal = 1;
984 <      simError();
962 >  GenericData* SimInfo::getPropertyByName(const string& propName) {
963 >    return properties_.getPropertyByName(propName);
964 >  }
965 >
966 >  void SimInfo::setSnapshotManager(SnapshotManager* sman) {
967 >    if (sman_ == sman) {
968 >      return;
969 >    }    
970 >    delete sman_;
971 >    sman_ = sman;
972 >
973 >    Molecule* mol;
974 >    RigidBody* rb;
975 >    Atom* atom;
976 >    CutoffGroup* cg;
977 >    SimInfo::MoleculeIterator mi;
978 >    Molecule::RigidBodyIterator rbIter;
979 >    Molecule::AtomIterator atomIter;
980 >    Molecule::CutoffGroupIterator cgIter;
981 >
982 >    for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
983 >        
984 >      for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) {
985 >        atom->setSnapshotManager(sman_);
986 >      }
987 >        
988 >      for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
989 >        rb->setSnapshotManager(sman_);
990 >      }
991 >
992 >      for (cg = mol->beginCutoffGroup(cgIter); cg != NULL; cg = mol->nextCutoffGroup(cgIter)) {
993 >        cg->setSnapshotManager(sman_);
994 >      }
995      }    
996 <  } else {
542 <    // initialize this stuff before using it, OK?
543 <    sprintf( painCave.errMsg,
544 <             "Trying to check cutoffs without a box.\n"
545 <             "\tOOPSE should have better programmers than that.\n" );
546 <    painCave.severity = OOPSE_ERROR;
547 <    painCave.isFatal = 1;
548 <    simError();      
996 >    
997    }
550  
551 }
998  
999 < void SimInfo::addProperty(GenericData* prop){
999 >  Vector3d SimInfo::getComVel(){
1000 >    SimInfo::MoleculeIterator i;
1001 >    Molecule* mol;
1002  
1003 <  map<string, GenericData*>::iterator result;
1004 <  result = properties.find(prop->getID());
557 <  
558 <  //we can't simply use  properties[prop->getID()] = prop,
559 <  //it will cause memory leak if we already contain a propery which has the same name of prop
560 <  
561 <  if(result != properties.end()){
1003 >    Vector3d comVel(0.0);
1004 >    RealType totalMass = 0.0;
1005      
1006 <    delete (*result).second;
1007 <    (*result).second = prop;
1008 <      
1006 >
1007 >    for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1008 >      RealType mass = mol->getMass();
1009 >      totalMass += mass;
1010 >      comVel += mass * mol->getComVel();
1011 >    }  
1012 >
1013 > #ifdef IS_MPI
1014 >    RealType tmpMass = totalMass;
1015 >    Vector3d tmpComVel(comVel);    
1016 >    MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1017 >    MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1018 > #endif
1019 >
1020 >    comVel /= totalMass;
1021 >
1022 >    return comVel;
1023    }
567  else{
1024  
1025 <    properties[prop->getID()] = prop;
1025 >  Vector3d SimInfo::getCom(){
1026 >    SimInfo::MoleculeIterator i;
1027 >    Molecule* mol;
1028  
1029 +    Vector3d com(0.0);
1030 +    RealType totalMass = 0.0;
1031 +    
1032 +    for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1033 +      RealType mass = mol->getMass();
1034 +      totalMass += mass;
1035 +      com += mass * mol->getCom();
1036 +    }  
1037 +
1038 + #ifdef IS_MPI
1039 +    RealType tmpMass = totalMass;
1040 +    Vector3d tmpCom(com);    
1041 +    MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1042 +    MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1043 + #endif
1044 +
1045 +    com /= totalMass;
1046 +
1047 +    return com;
1048 +
1049 +  }        
1050 +
1051 +  ostream& operator <<(ostream& o, SimInfo& info) {
1052 +
1053 +    return o;
1054    }
1055 +  
1056 +  
1057 +   /*
1058 +   Returns center of mass and center of mass velocity in one function call.
1059 +   */
1060 +  
1061 +   void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){
1062 +      SimInfo::MoleculeIterator i;
1063 +      Molecule* mol;
1064 +      
1065      
1066 < }
1066 >      RealType totalMass = 0.0;
1067 >    
1068  
1069 < GenericData* SimInfo::getProperty(const string& propName){
1070 <
1071 <  map<string, GenericData*>::iterator result;
1072 <  
1073 <  //string lowerCaseName = ();
1074 <  
1075 <  result = properties.find(propName);
1076 <  
1077 <  if(result != properties.end())
1078 <    return (*result).second;  
1079 <  else  
1080 <    return NULL;  
1081 < }
1069 >      for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1070 >         RealType mass = mol->getMass();
1071 >         totalMass += mass;
1072 >         com += mass * mol->getCom();
1073 >         comVel += mass * mol->getComVel();          
1074 >      }  
1075 >      
1076 > #ifdef IS_MPI
1077 >      RealType tmpMass = totalMass;
1078 >      Vector3d tmpCom(com);  
1079 >      Vector3d tmpComVel(comVel);
1080 >      MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1081 >      MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1082 >      MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1083 > #endif
1084 >      
1085 >      com /= totalMass;
1086 >      comVel /= totalMass;
1087 >   }        
1088 >  
1089 >   /*
1090 >   Return intertia tensor for entire system and angular momentum Vector.
1091  
1092  
1093 < void SimInfo::getFortranGroupArrays(SimInfo* info,
1094 <                                    vector<int>& FglobalGroupMembership,
1095 <                                    vector<double>& mfact){
1096 <  
594 <  Molecule* myMols;
595 <  Atom** myAtoms;
596 <  int numAtom;
597 <  double mtot;
598 <  int numMol;
599 <  int numCutoffGroups;
600 <  CutoffGroup* myCutoffGroup;
601 <  vector<CutoffGroup*>::iterator iterCutoff;
602 <  Atom* cutoffAtom;
603 <  vector<Atom*>::iterator iterAtom;
604 <  int atomIndex;
605 <  double totalMass;
606 <  
607 <  mfact.clear();
608 <  FglobalGroupMembership.clear();
609 <  
1093 >       [  Ixx -Ixy  -Ixz ]
1094 >    J =| -Iyx  Iyy  -Iyz |
1095 >       [ -Izx -Iyz   Izz ]
1096 >    */
1097  
1098 <  // Fix the silly fortran indexing problem
1098 >   void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){
1099 >      
1100 >
1101 >      RealType xx = 0.0;
1102 >      RealType yy = 0.0;
1103 >      RealType zz = 0.0;
1104 >      RealType xy = 0.0;
1105 >      RealType xz = 0.0;
1106 >      RealType yz = 0.0;
1107 >      Vector3d com(0.0);
1108 >      Vector3d comVel(0.0);
1109 >      
1110 >      getComAll(com, comVel);
1111 >      
1112 >      SimInfo::MoleculeIterator i;
1113 >      Molecule* mol;
1114 >      
1115 >      Vector3d thisq(0.0);
1116 >      Vector3d thisv(0.0);
1117 >
1118 >      RealType thisMass = 0.0;
1119 >    
1120 >      
1121 >      
1122 >  
1123 >      for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1124 >        
1125 >         thisq = mol->getCom()-com;
1126 >         thisv = mol->getComVel()-comVel;
1127 >         thisMass = mol->getMass();
1128 >         // Compute moment of intertia coefficients.
1129 >         xx += thisq[0]*thisq[0]*thisMass;
1130 >         yy += thisq[1]*thisq[1]*thisMass;
1131 >         zz += thisq[2]*thisq[2]*thisMass;
1132 >        
1133 >         // compute products of intertia
1134 >         xy += thisq[0]*thisq[1]*thisMass;
1135 >         xz += thisq[0]*thisq[2]*thisMass;
1136 >         yz += thisq[1]*thisq[2]*thisMass;
1137 >            
1138 >         angularMomentum += cross( thisq, thisv ) * thisMass;
1139 >            
1140 >      }  
1141 >      
1142 >      
1143 >      inertiaTensor(0,0) = yy + zz;
1144 >      inertiaTensor(0,1) = -xy;
1145 >      inertiaTensor(0,2) = -xz;
1146 >      inertiaTensor(1,0) = -xy;
1147 >      inertiaTensor(1,1) = xx + zz;
1148 >      inertiaTensor(1,2) = -yz;
1149 >      inertiaTensor(2,0) = -xz;
1150 >      inertiaTensor(2,1) = -yz;
1151 >      inertiaTensor(2,2) = xx + yy;
1152 >      
1153   #ifdef IS_MPI
1154 <  numAtom = mpiSim->getNAtomsGlobal();
1155 < #else
1156 <  numAtom = n_atoms;
1154 >      Mat3x3d tmpI(inertiaTensor);
1155 >      Vector3d tmpAngMom;
1156 >      MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1157 >      MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1158   #endif
1159 <  for (int i = 0; i < numAtom; i++)
1160 <    FglobalGroupMembership.push_back(globalGroupMembership[i] + 1);
1161 <  
1159 >              
1160 >      return;
1161 >   }
1162  
1163 <  myMols = info->molecules;
1164 <  numMol = info->n_mol;
623 <  for(int i  = 0; i < numMol; i++){
624 <    numCutoffGroups = myMols[i].getNCutoffGroups();
625 <    for(myCutoffGroup =myMols[i].beginCutoffGroup(iterCutoff);
626 <        myCutoffGroup != NULL;
627 <        myCutoffGroup =myMols[i].nextCutoffGroup(iterCutoff)){
628 <
629 <      totalMass = myCutoffGroup->getMass();
1163 >   //Returns the angular momentum of the system
1164 >   Vector3d SimInfo::getAngularMomentum(){
1165        
1166 <      for(cutoffAtom = myCutoffGroup->beginAtom(iterAtom);
1167 <          cutoffAtom != NULL;
1168 <          cutoffAtom = myCutoffGroup->nextAtom(iterAtom)){
1169 <        mfact.push_back(cutoffAtom->getMass()/totalMass);
1166 >      Vector3d com(0.0);
1167 >      Vector3d comVel(0.0);
1168 >      Vector3d angularMomentum(0.0);
1169 >      
1170 >      getComAll(com,comVel);
1171 >      
1172 >      SimInfo::MoleculeIterator i;
1173 >      Molecule* mol;
1174 >      
1175 >      Vector3d thisr(0.0);
1176 >      Vector3d thisp(0.0);
1177 >      
1178 >      RealType thisMass;
1179 >      
1180 >      for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {        
1181 >        thisMass = mol->getMass();
1182 >        thisr = mol->getCom()-com;
1183 >        thisp = (mol->getComVel()-comVel)*thisMass;
1184 >        
1185 >        angularMomentum += cross( thisr, thisp );
1186 >        
1187        }  
1188 +      
1189 + #ifdef IS_MPI
1190 +      Vector3d tmpAngMom;
1191 +      MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1192 + #endif
1193 +      
1194 +      return angularMomentum;
1195 +   }
1196 +  
1197 +  StuntDouble* SimInfo::getIOIndexToIntegrableObject(int index) {
1198 +    return IOIndexToIntegrableObject.at(index);
1199 +  }
1200 +  
1201 +  void SimInfo::setIOIndexToIntegrableObject(const vector<StuntDouble*>& v) {
1202 +    IOIndexToIntegrableObject= v;
1203 +  }
1204 +
1205 +  /* Returns the Volume of the simulation based on a ellipsoid with semi-axes
1206 +     based on the radius of gyration V=4/3*Pi*R_1*R_2*R_3
1207 +     where R_i are related to the principle inertia moments R_i = sqrt(C*I_i/N), this reduces to
1208 +     V = 4/3*Pi*(C/N)^3/2*sqrt(det(I)). See S.E. Baltazar et. al. Comp. Mat. Sci. 37 (2006) 526-536.
1209 +  */
1210 +  void SimInfo::getGyrationalVolume(RealType &volume){
1211 +    Mat3x3d intTensor;
1212 +    RealType det;
1213 +    Vector3d dummyAngMom;
1214 +    RealType sysconstants;
1215 +    RealType geomCnst;
1216 +
1217 +    geomCnst = 3.0/2.0;
1218 +    /* Get the inertial tensor and angular momentum for free*/
1219 +    getInertiaTensor(intTensor,dummyAngMom);
1220 +    
1221 +    det = intTensor.determinant();
1222 +    sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1223 +    volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(det);
1224 +    return;
1225 +  }
1226 +
1227 +  void SimInfo::getGyrationalVolume(RealType &volume, RealType &detI){
1228 +    Mat3x3d intTensor;
1229 +    Vector3d dummyAngMom;
1230 +    RealType sysconstants;
1231 +    RealType geomCnst;
1232 +
1233 +    geomCnst = 3.0/2.0;
1234 +    /* Get the inertial tensor and angular momentum for free*/
1235 +    getInertiaTensor(intTensor,dummyAngMom);
1236 +    
1237 +    detI = intTensor.determinant();
1238 +    sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1239 +    volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(detI);
1240 +    return;
1241 +  }
1242 + /*
1243 +   void SimInfo::setStuntDoubleFromGlobalIndex(vector<StuntDouble*> v) {
1244 +      assert( v.size() == nAtoms_ + nRigidBodies_);
1245 +      sdByGlobalIndex_ = v;
1246      }
1247 +
1248 +    StuntDouble* SimInfo::getStuntDoubleFromGlobalIndex(int index) {
1249 +      //assert(index < nAtoms_ + nRigidBodies_);
1250 +      return sdByGlobalIndex_.at(index);
1251 +    }  
1252 + */  
1253 +  int SimInfo::getNGlobalConstraints() {
1254 +    int nGlobalConstraints;
1255 + #ifdef IS_MPI
1256 +    MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM,
1257 +                  MPI_COMM_WORLD);    
1258 + #else
1259 +    nGlobalConstraints =  nConstraints_;
1260 + #endif
1261 +    return nGlobalConstraints;
1262    }
1263  
1264 < }
1264 > }//end namespace OpenMD
1265 >

Comparing:
trunk/src/brains/SimInfo.cpp (property svn:keywords), Revision 143 by chrisfen, Fri Oct 22 22:54:01 2004 UTC vs.
branches/development/src/brains/SimInfo.cpp (property svn:keywords), Revision 1547 by gezelter, Mon Apr 11 18:44:16 2011 UTC

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

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines