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 3 by tim, Fri Sep 24 16:27:58 2004 UTC vs.
branches/development/src/brains/SimInfo.cpp (file contents), Revision 1769 by gezelter, Mon Jul 9 14:15:52 2012 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]  Kuang & Gezelter,  J. Chem. Phys. 133, 164101 (2010).
40 > * [5]  Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 > */
42 >
43 > /**
44 > * @file SimInfo.cpp
45 > * @author    tlin
46 > * @date  11/02/2004
47 > * @version 1.0
48 > */
49  
50 < #include <iostream>
51 < using namespace std;
50 > #include <algorithm>
51 > #include <set>
52 > #include <map>
53  
54   #include "brains/SimInfo.hpp"
55 < #define __C
56 < #include "brains/fSimulation.h"
55 > #include "math/Vector3.hpp"
56 > #include "primitives/Molecule.hpp"
57 > #include "primitives/StuntDouble.hpp"
58 > #include "utils/MemoryUtils.hpp"
59   #include "utils/simError.h"
60 <
61 < #include "UseTheForce/fortranWrappers.hpp"
62 <
63 < #include "math/MatVec3.h"
16 <
60 > #include "selection/SelectionManager.hpp"
61 > #include "io/ForceFieldOptions.hpp"
62 > #include "brains/ForceField.hpp"
63 > #include "nonbonded/SwitchingFunction.hpp"
64   #ifdef IS_MPI
65 < #include "brains/mpiSimulation.hpp"
65 > #include <mpi.h>
66   #endif
67  
68 < inline double roundMe( double x ){
69 <  return ( x >= 0 ) ? floor( x + 0.5 ) : ceil( x - 0.5 );
70 < }
71 <          
72 < inline double min( double a, double b ){
73 <  return (a < b ) ? a : b;
74 < }
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), nGlobalFluctuatingCharges_(0),
76 >    nAtoms_(0), nBonds_(0),  nBends_(0), nTorsions_(0), nInversions_(0),
77 >    nRigidBodies_(0), nIntegrableObjects_(0), nCutoffGroups_(0),
78 >    nConstraints_(0), nFluctuatingCharges_(0), sman_(NULL), topologyDone_(false),
79 >    calcBoxDipole_(false), useAtomicVirial_(true) {    
80 >    
81 >    MoleculeStamp* molStamp;
82 >    int nMolWithSameStamp;
83 >    int nCutoffAtoms = 0; // number of atoms belong to cutoff groups
84 >    int nGroups = 0;       //total cutoff groups defined in meta-data file
85 >    CutoffGroupStamp* cgStamp;    
86 >    RigidBodyStamp* rbStamp;
87 >    int nRigidAtoms = 0;
88 >    
89 >    vector<Component*> components = simParams->getComponents();
90 >    
91 >    for (vector<Component*>::iterator i = components.begin(); i !=components.end(); ++i) {
92 >      molStamp = (*i)->getMoleculeStamp();
93 >      nMolWithSameStamp = (*i)->getNMol();
94 >      
95 >      addMoleculeStamp(molStamp, nMolWithSameStamp);
96 >      
97 >      //calculate atoms in molecules
98 >      nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp;  
99 >      
100 >      //calculate atoms in cutoff groups
101 >      int nAtomsInGroups = 0;
102 >      int nCutoffGroupsInStamp = molStamp->getNCutoffGroups();
103 >      
104 >      for (int j=0; j < nCutoffGroupsInStamp; j++) {
105 >        cgStamp = molStamp->getCutoffGroupStamp(j);
106 >        nAtomsInGroups += cgStamp->getNMembers();
107 >      }
108 >      
109 >      nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
110 >      
111 >      nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;            
112 >      
113 >      //calculate atoms in rigid bodies
114 >      int nAtomsInRigidBodies = 0;
115 >      int nRigidBodiesInStamp = molStamp->getNRigidBodies();
116 >      
117 >      for (int j=0; j < nRigidBodiesInStamp; j++) {
118 >        rbStamp = molStamp->getRigidBodyStamp(j);
119 >        nAtomsInRigidBodies += rbStamp->getNMembers();
120 >      }
121 >      
122 >      nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
123 >      nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;            
124 >      
125 >    }
126 >    
127 >    //every free atom (atom does not belong to cutoff groups) is a cutoff
128 >    //group therefore the total number of cutoff groups in the system is
129 >    //equal to the total number of atoms minus number of atoms belong to
130 >    //cutoff group defined in meta-data file plus the number of cutoff
131 >    //groups defined in meta-data file
132  
133 < SimInfo* currentInfo;
133 >    nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
134 >    
135 >    //every free atom (atom does not belong to rigid bodies) is an
136 >    //integrable object therefore the total number of integrable objects
137 >    //in the system is equal to the total number of atoms minus number of
138 >    //atoms belong to rigid body defined in meta-data file plus the number
139 >    //of rigid bodies defined in meta-data file
140 >    nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms
141 >      + nGlobalRigidBodies_;
142 >    
143 >    nGlobalMols_ = molStampIds_.size();
144 >    molToProcMap_.resize(nGlobalMols_);
145 >  }
146 >  
147 >  SimInfo::~SimInfo() {
148 >    map<int, Molecule*>::iterator i;
149 >    for (i = molecules_.begin(); i != molecules_.end(); ++i) {
150 >      delete i->second;
151 >    }
152 >    molecules_.clear();
153 >      
154 >    delete sman_;
155 >    delete simParams_;
156 >    delete forceField_;
157 >  }
158  
31 SimInfo::SimInfo(){
159  
160 <  n_constraints = 0;
161 <  nZconstraints = 0;
162 <  n_oriented = 0;
163 <  n_dipoles = 0;
164 <  ndf = 0;
165 <  ndfRaw = 0;
166 <  nZconstraints = 0;
167 <  the_integrator = NULL;
168 <  setTemp = 0;
169 <  thermalTime = 0.0;
170 <  currentTime = 0.0;
171 <  rCut = 0.0;
172 <  rSw = 0.0;
173 <
174 <  haveRcut = 0;
175 <  haveRsw = 0;
176 <  boxIsInit = 0;
160 >  bool SimInfo::addMolecule(Molecule* mol) {
161 >    MoleculeIterator i;
162 >    
163 >    i = molecules_.find(mol->getGlobalIndex());
164 >    if (i == molecules_.end() ) {
165 >      
166 >      molecules_.insert(make_pair(mol->getGlobalIndex(), mol));
167 >      
168 >      nAtoms_ += mol->getNAtoms();
169 >      nBonds_ += mol->getNBonds();
170 >      nBends_ += mol->getNBends();
171 >      nTorsions_ += mol->getNTorsions();
172 >      nInversions_ += mol->getNInversions();
173 >      nRigidBodies_ += mol->getNRigidBodies();
174 >      nIntegrableObjects_ += mol->getNIntegrableObjects();
175 >      nCutoffGroups_ += mol->getNCutoffGroups();
176 >      nConstraints_ += mol->getNConstraintPairs();
177 >      
178 >      addInteractionPairs(mol);
179 >      
180 >      return true;
181 >    } else {
182 >      return false;
183 >    }
184 >  }
185    
186 <  resetTime = 1e99;
186 >  bool SimInfo::removeMolecule(Molecule* mol) {
187 >    MoleculeIterator i;
188 >    i = molecules_.find(mol->getGlobalIndex());
189  
190 <  orthoRhombic = 0;
54 <  orthoTolerance = 1E-6;
55 <  useInitXSstate = true;
190 >    if (i != molecules_.end() ) {
191  
192 <  usePBC = 0;
193 <  useLJ = 0;
194 <  useSticky = 0;
195 <  useCharges = 0;
196 <  useDipoles = 0;
197 <  useReactionField = 0;
198 <  useGB = 0;
199 <  useEAM = 0;
200 <  useSolidThermInt = 0;
201 <  useLiquidThermInt = 0;
192 >      assert(mol == i->second);
193 >        
194 >      nAtoms_ -= mol->getNAtoms();
195 >      nBonds_ -= mol->getNBonds();
196 >      nBends_ -= mol->getNBends();
197 >      nTorsions_ -= mol->getNTorsions();
198 >      nInversions_ -= mol->getNInversions();
199 >      nRigidBodies_ -= mol->getNRigidBodies();
200 >      nIntegrableObjects_ -= mol->getNIntegrableObjects();
201 >      nCutoffGroups_ -= mol->getNCutoffGroups();
202 >      nConstraints_ -= mol->getNConstraintPairs();
203  
204 <  haveCutoffGroups = false;
204 >      removeInteractionPairs(mol);
205 >      molecules_.erase(mol->getGlobalIndex());
206  
207 <  excludes = Exclude::Instance();
207 >      delete mol;
208 >        
209 >      return true;
210 >    } else {
211 >      return false;
212 >    }
213 >  }    
214  
215 <  myConfiguration = new SimState();
215 >        
216 >  Molecule* SimInfo::beginMolecule(MoleculeIterator& i) {
217 >    i = molecules_.begin();
218 >    return i == molecules_.end() ? NULL : i->second;
219 >  }    
220  
221 <  has_minimizer = false;
222 <  the_minimizer =NULL;
221 >  Molecule* SimInfo::nextMolecule(MoleculeIterator& i) {
222 >    ++i;
223 >    return i == molecules_.end() ? NULL : i->second;    
224 >  }
225  
77  ngroup = 0;
226  
227 <  wrapMeSimInfo( this );
228 < }
227 >  void SimInfo::calcNdf() {
228 >    int ndf_local, nfq_local;
229 >    MoleculeIterator i;
230 >    vector<StuntDouble*>::iterator j;
231 >    vector<Atom*>::iterator k;
232  
233 +    Molecule* mol;
234 +    StuntDouble* sd;
235 +    Atom* atom;
236  
237 < SimInfo::~SimInfo(){
237 >    ndf_local = 0;
238 >    nfq_local = 0;
239 >    
240 >    for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
241  
242 <  delete myConfiguration;
242 >      for (sd = mol->beginIntegrableObject(j); sd != NULL;
243 >           sd = mol->nextIntegrableObject(j)) {
244  
245 <  map<string, GenericData*>::iterator i;
88 <  
89 <  for(i = properties.begin(); i != properties.end(); i++)
90 <    delete (*i).second;
245 >        ndf_local += 3;
246  
247 < }
247 >        if (sd->isDirectional()) {
248 >          if (sd->isLinear()) {
249 >            ndf_local += 2;
250 >          } else {
251 >            ndf_local += 3;
252 >          }
253 >        }
254 >      }
255  
256 < void SimInfo::setBox(double newBox[3]) {
257 <  
258 <  int i, j;
259 <  double tempMat[3][3];
256 >      for (atom = mol->beginFluctuatingCharge(k); atom != NULL;
257 >           atom = mol->nextFluctuatingCharge(k)) {
258 >        if (atom->isFluctuatingCharge()) {
259 >          nfq_local++;
260 >        }
261 >      }
262 >    }
263 >    
264 >    ndfLocal_ = ndf_local;
265  
266 <  for(i=0; i<3; i++)
267 <    for (j=0; j<3; j++) tempMat[i][j] = 0.0;;
266 >    // n_constraints is local, so subtract them on each processor
267 >    ndf_local -= nConstraints_;
268  
269 <  tempMat[0][0] = newBox[0];
270 <  tempMat[1][1] = newBox[1];
271 <  tempMat[2][2] = newBox[2];
269 > #ifdef IS_MPI
270 >    MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
271 >    MPI_Allreduce(&nfq_local,&nGlobalFluctuatingCharges_,1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
272 > #else
273 >    ndf_ = ndf_local;
274 >    nGlobalFluctuatingCharges_ = nfq_local;
275 > #endif
276  
277 <  setBoxM( tempMat );
277 >    // nZconstraints_ is global, as are the 3 COM translations for the
278 >    // entire system:
279 >    ndf_ = ndf_ - 3 - nZconstraint_;
280  
281 < }
281 >  }
282  
283 < void SimInfo::setBoxM( double theBox[3][3] ){
283 >  int SimInfo::getFdf() {
284 > #ifdef IS_MPI
285 >    MPI_Allreduce(&fdf_local,&fdf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
286 > #else
287 >    fdf_ = fdf_local;
288 > #endif
289 >    return fdf_;
290 >  }
291    
292 <  int i, j;
293 <  double FortranHmat[9]; // to preserve compatibility with Fortran the
294 <                         // ordering in the array is as follows:
295 <                         // [ 0 3 6 ]
296 <                         // [ 1 4 7 ]
297 <                         // [ 2 5 8 ]
298 <  double FortranHmatInv[9]; // the inverted Hmat (for Fortran);
299 <
300 <  if( !boxIsInit ) boxIsInit = 1;
301 <
302 <  for(i=0; i < 3; i++)
303 <    for (j=0; j < 3; j++) Hmat[i][j] = theBox[i][j];
304 <  
305 <  calcBoxL();
126 <  calcHmatInv();
127 <
128 <  for(i=0; i < 3; i++) {
129 <    for (j=0; j < 3; j++) {
130 <      FortranHmat[3*j + i] = Hmat[i][j];
131 <      FortranHmatInv[3*j + i] = HmatInv[i][j];
292 >  unsigned int SimInfo::getNLocalCutoffGroups(){
293 >    int nLocalCutoffAtoms = 0;
294 >    Molecule* mol;
295 >    MoleculeIterator mi;
296 >    CutoffGroup* cg;
297 >    Molecule::CutoffGroupIterator ci;
298 >    
299 >    for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) {
300 >      
301 >      for (cg = mol->beginCutoffGroup(ci); cg != NULL;
302 >           cg = mol->nextCutoffGroup(ci)) {
303 >        nLocalCutoffAtoms += cg->getNumAtom();
304 >        
305 >      }        
306      }
307 +    
308 +    return nAtoms_ - nLocalCutoffAtoms + nCutoffGroups_;
309    }
310 +    
311 +  void SimInfo::calcNdfRaw() {
312 +    int ndfRaw_local;
313  
314 <  setFortranBoxSize(FortranHmat, FortranHmatInv, &orthoRhombic);
315 <
316 < }
317 <
314 >    MoleculeIterator i;
315 >    vector<StuntDouble*>::iterator j;
316 >    Molecule* mol;
317 >    StuntDouble* sd;
318  
319 < void SimInfo::getBoxM (double theBox[3][3]) {
319 >    // Raw degrees of freedom that we have to set
320 >    ndfRaw_local = 0;
321 >    
322 >    for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
323  
324 <  int i, j;
325 <  for(i=0; i<3; i++)
144 <    for (j=0; j<3; j++) theBox[i][j] = Hmat[i][j];
145 < }
324 >      for (sd = mol->beginIntegrableObject(j); sd != NULL;
325 >           sd = mol->nextIntegrableObject(j)) {
326  
327 +        ndfRaw_local += 3;
328  
329 < void SimInfo::scaleBox(double scale) {
330 <  double theBox[3][3];
331 <  int i, j;
332 <
333 <  // cerr << "Scaling box by " << scale << "\n";
334 <
335 <  for(i=0; i<3; i++)
336 <    for (j=0; j<3; j++) theBox[i][j] = Hmat[i][j]*scale;
156 <
157 <  setBoxM(theBox);
158 <
159 < }
160 <
161 < void SimInfo::calcHmatInv( void ) {
162 <  
163 <  int oldOrtho;
164 <  int i,j;
165 <  double smallDiag;
166 <  double tol;
167 <  double sanity[3][3];
168 <
169 <  invertMat3( Hmat, HmatInv );
170 <
171 <  // check to see if Hmat is orthorhombic
172 <  
173 <  oldOrtho = orthoRhombic;
174 <
175 <  smallDiag = fabs(Hmat[0][0]);
176 <  if(smallDiag > fabs(Hmat[1][1])) smallDiag = fabs(Hmat[1][1]);
177 <  if(smallDiag > fabs(Hmat[2][2])) smallDiag = fabs(Hmat[2][2]);
178 <  tol = smallDiag * orthoTolerance;
179 <
180 <  orthoRhombic = 1;
181 <  
182 <  for (i = 0; i < 3; i++ ) {
183 <    for (j = 0 ; j < 3; j++) {
184 <      if (i != j) {
185 <        if (orthoRhombic) {
186 <          if ( fabs(Hmat[i][j]) >= tol) orthoRhombic = 0;
187 <        }        
329 >        if (sd->isDirectional()) {
330 >          if (sd->isLinear()) {
331 >            ndfRaw_local += 2;
332 >          } else {
333 >            ndfRaw_local += 3;
334 >          }
335 >        }
336 >            
337        }
338      }
190  }
191
192  if( oldOrtho != orthoRhombic ){
339      
340 <    if( orthoRhombic ) {
341 <      sprintf( painCave.errMsg,
342 <               "OOPSE is switching from the default Non-Orthorhombic\n"
343 <               "\tto the faster Orthorhombic periodic boundary computations.\n"
344 <               "\tThis is usually a good thing, but if you wan't the\n"
199 <               "\tNon-Orthorhombic computations, make the orthoBoxTolerance\n"
200 <               "\tvariable ( currently set to %G ) smaller.\n",
201 <               orthoTolerance);
202 <      painCave.severity = OOPSE_INFO;
203 <      simError();
204 <    }
205 <    else {
206 <      sprintf( painCave.errMsg,
207 <               "OOPSE is switching from the faster Orthorhombic to the more\n"
208 <               "\tflexible Non-Orthorhombic periodic boundary computations.\n"
209 <               "\tThis is usually because the box has deformed under\n"
210 <               "\tNPTf integration. If you wan't to live on the edge with\n"
211 <               "\tthe Orthorhombic computations, make the orthoBoxTolerance\n"
212 <               "\tvariable ( currently set to %G ) larger.\n",
213 <               orthoTolerance);
214 <      painCave.severity = OOPSE_WARNING;
215 <      simError();
216 <    }
340 > #ifdef IS_MPI
341 >    MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
342 > #else
343 >    ndfRaw_ = ndfRaw_local;
344 > #endif
345    }
218 }
346  
347 < void SimInfo::calcBoxL( void ){
347 >  void SimInfo::calcNdfTrans() {
348 >    int ndfTrans_local;
349  
350 <  double dx, dy, dz, dsq;
350 >    ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
351  
224  // boxVol = Determinant of Hmat
352  
353 <  boxVol = matDet3( Hmat );
353 > #ifdef IS_MPI
354 >    MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
355 > #else
356 >    ndfTrans_ = ndfTrans_local;
357 > #endif
358  
359 <  // boxLx
360 <  
361 <  dx = Hmat[0][0]; dy = Hmat[1][0]; dz = Hmat[2][0];
231 <  dsq = dx*dx + dy*dy + dz*dz;
232 <  boxL[0] = sqrt( dsq );
233 <  //maxCutoff = 0.5 * boxL[0];
359 >    ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
360 >
361 >  }
362  
363 <  // boxLy
364 <  
365 <  dx = Hmat[0][1]; dy = Hmat[1][1]; dz = Hmat[2][1];
366 <  dsq = dx*dx + dy*dy + dz*dz;
367 <  boxL[1] = sqrt( dsq );
368 <  //if( (0.5 * boxL[1]) < maxCutoff ) maxCutoff = 0.5 * boxL[1];
363 >  void SimInfo::addInteractionPairs(Molecule* mol) {
364 >    ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
365 >    vector<Bond*>::iterator bondIter;
366 >    vector<Bend*>::iterator bendIter;
367 >    vector<Torsion*>::iterator torsionIter;
368 >    vector<Inversion*>::iterator inversionIter;
369 >    Bond* bond;
370 >    Bend* bend;
371 >    Torsion* torsion;
372 >    Inversion* inversion;
373 >    int a;
374 >    int b;
375 >    int c;
376 >    int d;
377  
378 +    // atomGroups can be used to add special interaction maps between
379 +    // groups of atoms that are in two separate rigid bodies.
380 +    // However, most site-site interactions between two rigid bodies
381 +    // are probably not special, just the ones between the physically
382 +    // bonded atoms.  Interactions *within* a single rigid body should
383 +    // always be excluded.  These are done at the bottom of this
384 +    // function.
385  
386 <  // boxLz
387 <  
388 <  dx = Hmat[0][2]; dy = Hmat[1][2]; dz = Hmat[2][2];
389 <  dsq = dx*dx + dy*dy + dz*dz;
390 <  boxL[2] = sqrt( dsq );
391 <  //if( (0.5 * boxL[2]) < maxCutoff ) maxCutoff = 0.5 * boxL[2];
386 >    map<int, set<int> > atomGroups;
387 >    Molecule::RigidBodyIterator rbIter;
388 >    RigidBody* rb;
389 >    Molecule::IntegrableObjectIterator ii;
390 >    StuntDouble* sd;
391 >    
392 >    for (sd = mol->beginIntegrableObject(ii); sd != NULL;
393 >         sd = mol->nextIntegrableObject(ii)) {
394 >      
395 >      if (sd->isRigidBody()) {
396 >        rb = static_cast<RigidBody*>(sd);
397 >        vector<Atom*> atoms = rb->getAtoms();
398 >        set<int> rigidAtoms;
399 >        for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
400 >          rigidAtoms.insert(atoms[i]->getGlobalIndex());
401 >        }
402 >        for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
403 >          atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
404 >        }      
405 >      } else {
406 >        set<int> oneAtomSet;
407 >        oneAtomSet.insert(sd->getGlobalIndex());
408 >        atomGroups.insert(map<int, set<int> >::value_type(sd->getGlobalIndex(), oneAtomSet));        
409 >      }
410 >    }  
411 >          
412 >    for (bond= mol->beginBond(bondIter); bond != NULL;
413 >         bond = mol->nextBond(bondIter)) {
414  
415 <  //calculate the max cutoff
416 <  maxCutoff =  calcMaxCutOff();
417 <  
418 <  checkCutOffs();
415 >      a = bond->getAtomA()->getGlobalIndex();
416 >      b = bond->getAtomB()->getGlobalIndex();  
417 >    
418 >      if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
419 >        oneTwoInteractions_.addPair(a, b);
420 >      } else {
421 >        excludedInteractions_.addPair(a, b);
422 >      }
423 >    }
424  
425 < }
425 >    for (bend= mol->beginBend(bendIter); bend != NULL;
426 >         bend = mol->nextBend(bendIter)) {
427  
428 +      a = bend->getAtomA()->getGlobalIndex();
429 +      b = bend->getAtomB()->getGlobalIndex();        
430 +      c = bend->getAtomC()->getGlobalIndex();
431 +      
432 +      if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
433 +        oneTwoInteractions_.addPair(a, b);      
434 +        oneTwoInteractions_.addPair(b, c);
435 +      } else {
436 +        excludedInteractions_.addPair(a, b);
437 +        excludedInteractions_.addPair(b, c);
438 +      }
439  
440 < double SimInfo::calcMaxCutOff(){
440 >      if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
441 >        oneThreeInteractions_.addPair(a, c);      
442 >      } else {
443 >        excludedInteractions_.addPair(a, c);
444 >      }
445 >    }
446  
447 <  double ri[3], rj[3], rk[3];
448 <  double rij[3], rjk[3], rki[3];
262 <  double minDist;
447 >    for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
448 >         torsion = mol->nextTorsion(torsionIter)) {
449  
450 <  ri[0] = Hmat[0][0];
451 <  ri[1] = Hmat[1][0];
452 <  ri[2] = Hmat[2][0];
450 >      a = torsion->getAtomA()->getGlobalIndex();
451 >      b = torsion->getAtomB()->getGlobalIndex();        
452 >      c = torsion->getAtomC()->getGlobalIndex();        
453 >      d = torsion->getAtomD()->getGlobalIndex();      
454  
455 <  rj[0] = Hmat[0][1];
456 <  rj[1] = Hmat[1][1];
457 <  rj[2] = Hmat[2][1];
455 >      if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
456 >        oneTwoInteractions_.addPair(a, b);      
457 >        oneTwoInteractions_.addPair(b, c);
458 >        oneTwoInteractions_.addPair(c, d);
459 >      } else {
460 >        excludedInteractions_.addPair(a, b);
461 >        excludedInteractions_.addPair(b, c);
462 >        excludedInteractions_.addPair(c, d);
463 >      }
464  
465 <  rk[0] = Hmat[0][2];
466 <  rk[1] = Hmat[1][2];
467 <  rk[2] = Hmat[2][2];
468 <    
469 <  crossProduct3(ri, rj, rij);
470 <  distXY = dotProduct3(rk,rij) / norm3(rij);
465 >      if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
466 >        oneThreeInteractions_.addPair(a, c);      
467 >        oneThreeInteractions_.addPair(b, d);      
468 >      } else {
469 >        excludedInteractions_.addPair(a, c);
470 >        excludedInteractions_.addPair(b, d);
471 >      }
472  
473 <  crossProduct3(rj,rk, rjk);
474 <  distYZ = dotProduct3(ri,rjk) / norm3(rjk);
473 >      if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
474 >        oneFourInteractions_.addPair(a, d);      
475 >      } else {
476 >        excludedInteractions_.addPair(a, d);
477 >      }
478 >    }
479  
480 <  crossProduct3(rk,ri, rki);
481 <  distZX = dotProduct3(rj,rki) / norm3(rki);
480 >    for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
481 >         inversion = mol->nextInversion(inversionIter)) {
482  
483 <  minDist = min(min(distXY, distYZ), distZX);
484 <  return minDist/2;
485 <  
486 < }
483 >      a = inversion->getAtomA()->getGlobalIndex();
484 >      b = inversion->getAtomB()->getGlobalIndex();        
485 >      c = inversion->getAtomC()->getGlobalIndex();        
486 >      d = inversion->getAtomD()->getGlobalIndex();        
487  
488 < void SimInfo::wrapVector( double thePos[3] ){
488 >      if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
489 >        oneTwoInteractions_.addPair(a, b);      
490 >        oneTwoInteractions_.addPair(a, c);
491 >        oneTwoInteractions_.addPair(a, d);
492 >      } else {
493 >        excludedInteractions_.addPair(a, b);
494 >        excludedInteractions_.addPair(a, c);
495 >        excludedInteractions_.addPair(a, d);
496 >      }
497  
498 <  int i;
499 <  double scaled[3];
498 >      if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
499 >        oneThreeInteractions_.addPair(b, c);    
500 >        oneThreeInteractions_.addPair(b, d);    
501 >        oneThreeInteractions_.addPair(c, d);      
502 >      } else {
503 >        excludedInteractions_.addPair(b, c);
504 >        excludedInteractions_.addPair(b, d);
505 >        excludedInteractions_.addPair(c, d);
506 >      }
507 >    }
508  
509 <  if( !orthoRhombic ){
510 <    // calc the scaled coordinates.
511 <  
509 >    for (rb = mol->beginRigidBody(rbIter); rb != NULL;
510 >         rb = mol->nextRigidBody(rbIter)) {
511 >      vector<Atom*> atoms = rb->getAtoms();
512 >      for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
513 >        for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
514 >          a = atoms[i]->getGlobalIndex();
515 >          b = atoms[j]->getGlobalIndex();
516 >          excludedInteractions_.addPair(a, b);
517 >        }
518 >      }
519 >    }        
520  
299    matVecMul3(HmatInv, thePos, scaled);
300    
301    for(i=0; i<3; i++)
302      scaled[i] -= roundMe(scaled[i]);
303    
304    // calc the wrapped real coordinates from the wrapped scaled coordinates
305    
306    matVecMul3(Hmat, scaled, thePos);
307
521    }
309  else{
310    // calc the scaled coordinates.
311    
312    for(i=0; i<3; i++)
313      scaled[i] = thePos[i]*HmatInv[i][i];
314    
315    // wrap the scaled coordinates
316    
317    for(i=0; i<3; i++)
318      scaled[i] -= roundMe(scaled[i]);
319    
320    // calc the wrapped real coordinates from the wrapped scaled coordinates
321    
322    for(i=0; i<3; i++)
323      thePos[i] = scaled[i]*Hmat[i][i];
324  }
325    
326 }
522  
523 +  void SimInfo::removeInteractionPairs(Molecule* mol) {
524 +    ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
525 +    vector<Bond*>::iterator bondIter;
526 +    vector<Bend*>::iterator bendIter;
527 +    vector<Torsion*>::iterator torsionIter;
528 +    vector<Inversion*>::iterator inversionIter;
529 +    Bond* bond;
530 +    Bend* bend;
531 +    Torsion* torsion;
532 +    Inversion* inversion;
533 +    int a;
534 +    int b;
535 +    int c;
536 +    int d;
537  
538 < int SimInfo::getNDF(){
539 <  int ndf_local;
540 <
541 <  ndf_local = 0;
542 <  
543 <  for(int i = 0; i < integrableObjects.size(); i++){
544 <    ndf_local += 3;
545 <    if (integrableObjects[i]->isDirectional()) {
546 <      if (integrableObjects[i]->isLinear())
547 <        ndf_local += 2;
548 <      else
549 <        ndf_local += 3;
538 >    map<int, set<int> > atomGroups;
539 >    Molecule::RigidBodyIterator rbIter;
540 >    RigidBody* rb;
541 >    Molecule::IntegrableObjectIterator ii;
542 >    StuntDouble* sd;
543 >    
544 >    for (sd = mol->beginIntegrableObject(ii); sd != NULL;
545 >         sd = mol->nextIntegrableObject(ii)) {
546 >      
547 >      if (sd->isRigidBody()) {
548 >        rb = static_cast<RigidBody*>(sd);
549 >        vector<Atom*> atoms = rb->getAtoms();
550 >        set<int> rigidAtoms;
551 >        for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
552 >          rigidAtoms.insert(atoms[i]->getGlobalIndex());
553 >        }
554 >        for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
555 >          atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
556 >        }      
557 >      } else {
558 >        set<int> oneAtomSet;
559 >        oneAtomSet.insert(sd->getGlobalIndex());
560 >        atomGroups.insert(map<int, set<int> >::value_type(sd->getGlobalIndex(), oneAtomSet));        
561 >      }
562 >    }  
563 >
564 >    for (bond= mol->beginBond(bondIter); bond != NULL;
565 >         bond = mol->nextBond(bondIter)) {
566 >      
567 >      a = bond->getAtomA()->getGlobalIndex();
568 >      b = bond->getAtomB()->getGlobalIndex();  
569 >    
570 >      if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
571 >        oneTwoInteractions_.removePair(a, b);
572 >      } else {
573 >        excludedInteractions_.removePair(a, b);
574 >      }
575      }
342  }
576  
577 <  // n_constraints is local, so subtract them on each processor:
577 >    for (bend= mol->beginBend(bendIter); bend != NULL;
578 >         bend = mol->nextBend(bendIter)) {
579  
580 <  ndf_local -= n_constraints;
580 >      a = bend->getAtomA()->getGlobalIndex();
581 >      b = bend->getAtomB()->getGlobalIndex();        
582 >      c = bend->getAtomC()->getGlobalIndex();
583 >      
584 >      if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
585 >        oneTwoInteractions_.removePair(a, b);      
586 >        oneTwoInteractions_.removePair(b, c);
587 >      } else {
588 >        excludedInteractions_.removePair(a, b);
589 >        excludedInteractions_.removePair(b, c);
590 >      }
591  
592 < #ifdef IS_MPI
593 <  MPI_Allreduce(&ndf_local,&ndf,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
594 < #else
595 <  ndf = ndf_local;
596 < #endif
592 >      if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
593 >        oneThreeInteractions_.removePair(a, c);      
594 >      } else {
595 >        excludedInteractions_.removePair(a, c);
596 >      }
597 >    }
598  
599 <  // nZconstraints is global, as are the 3 COM translations for the
600 <  // entire system:
599 >    for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
600 >         torsion = mol->nextTorsion(torsionIter)) {
601  
602 <  ndf = ndf - 3 - nZconstraints;
602 >      a = torsion->getAtomA()->getGlobalIndex();
603 >      b = torsion->getAtomB()->getGlobalIndex();        
604 >      c = torsion->getAtomC()->getGlobalIndex();        
605 >      d = torsion->getAtomD()->getGlobalIndex();      
606 >  
607 >      if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
608 >        oneTwoInteractions_.removePair(a, b);      
609 >        oneTwoInteractions_.removePair(b, c);
610 >        oneTwoInteractions_.removePair(c, d);
611 >      } else {
612 >        excludedInteractions_.removePair(a, b);
613 >        excludedInteractions_.removePair(b, c);
614 >        excludedInteractions_.removePair(c, d);
615 >      }
616  
617 <  return ndf;
618 < }
617 >      if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
618 >        oneThreeInteractions_.removePair(a, c);      
619 >        oneThreeInteractions_.removePair(b, d);      
620 >      } else {
621 >        excludedInteractions_.removePair(a, c);
622 >        excludedInteractions_.removePair(b, d);
623 >      }
624  
625 < int SimInfo::getNDFraw() {
626 <  int ndfRaw_local;
625 >      if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
626 >        oneFourInteractions_.removePair(a, d);      
627 >      } else {
628 >        excludedInteractions_.removePair(a, d);
629 >      }
630 >    }
631  
632 <  // Raw degrees of freedom that we have to set
633 <  ndfRaw_local = 0;
632 >    for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
633 >         inversion = mol->nextInversion(inversionIter)) {
634  
635 <  for(int i = 0; i < integrableObjects.size(); i++){
636 <    ndfRaw_local += 3;
637 <    if (integrableObjects[i]->isDirectional()) {
638 <       if (integrableObjects[i]->isLinear())
639 <        ndfRaw_local += 2;
640 <      else
641 <        ndfRaw_local += 3;
635 >      a = inversion->getAtomA()->getGlobalIndex();
636 >      b = inversion->getAtomB()->getGlobalIndex();        
637 >      c = inversion->getAtomC()->getGlobalIndex();        
638 >      d = inversion->getAtomD()->getGlobalIndex();        
639 >
640 >      if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
641 >        oneTwoInteractions_.removePair(a, b);      
642 >        oneTwoInteractions_.removePair(a, c);
643 >        oneTwoInteractions_.removePair(a, d);
644 >      } else {
645 >        excludedInteractions_.removePair(a, b);
646 >        excludedInteractions_.removePair(a, c);
647 >        excludedInteractions_.removePair(a, d);
648 >      }
649 >
650 >      if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
651 >        oneThreeInteractions_.removePair(b, c);    
652 >        oneThreeInteractions_.removePair(b, d);    
653 >        oneThreeInteractions_.removePair(c, d);      
654 >      } else {
655 >        excludedInteractions_.removePair(b, c);
656 >        excludedInteractions_.removePair(b, d);
657 >        excludedInteractions_.removePair(c, d);
658 >      }
659      }
660 +
661 +    for (rb = mol->beginRigidBody(rbIter); rb != NULL;
662 +         rb = mol->nextRigidBody(rbIter)) {
663 +      vector<Atom*> atoms = rb->getAtoms();
664 +      for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
665 +        for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
666 +          a = atoms[i]->getGlobalIndex();
667 +          b = atoms[j]->getGlobalIndex();
668 +          excludedInteractions_.removePair(a, b);
669 +        }
670 +      }
671 +    }        
672 +    
673    }
674 +  
675 +  
676 +  void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) {
677 +    int curStampId;
678      
679 +    //index from 0
680 +    curStampId = moleculeStamps_.size();
681 +
682 +    moleculeStamps_.push_back(molStamp);
683 +    molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
684 +  }
685 +
686 +
687 +  /**
688 +   * update
689 +   *
690 +   *  Performs the global checks and variable settings after the
691 +   *  objects have been created.
692 +   *
693 +   */
694 +  void SimInfo::update() {  
695 +    setupSimVariables();
696 +    calcNdf();
697 +    calcNdfRaw();
698 +    calcNdfTrans();
699 +  }
700 +  
701 +  /**
702 +   * getSimulatedAtomTypes
703 +   *
704 +   * Returns an STL set of AtomType* that are actually present in this
705 +   * simulation.  Must query all processors to assemble this information.
706 +   *
707 +   */
708 +  set<AtomType*> SimInfo::getSimulatedAtomTypes() {
709 +    SimInfo::MoleculeIterator mi;
710 +    Molecule* mol;
711 +    Molecule::AtomIterator ai;
712 +    Atom* atom;
713 +    set<AtomType*> atomTypes;
714 +    
715 +    for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
716 +      for(atom = mol->beginAtom(ai); atom != NULL;
717 +          atom = mol->nextAtom(ai)) {
718 +        atomTypes.insert(atom->getAtomType());
719 +      }      
720 +    }    
721 +    
722   #ifdef IS_MPI
379  MPI_Allreduce(&ndfRaw_local,&ndfRaw,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
380 #else
381  ndfRaw = ndfRaw_local;
382 #endif
723  
724 <  return ndfRaw;
725 < }
724 >    // loop over the found atom types on this processor, and add their
725 >    // numerical idents to a vector:
726 >    
727 >    vector<int> foundTypes;
728 >    set<AtomType*>::iterator i;
729 >    for (i = atomTypes.begin(); i != atomTypes.end(); ++i)
730 >      foundTypes.push_back( (*i)->getIdent() );
731  
732 < int SimInfo::getNDFtranslational() {
733 <  int ndfTrans_local;
732 >    // count_local holds the number of found types on this processor
733 >    int count_local = foundTypes.size();
734  
735 <  ndfTrans_local = 3 * integrableObjects.size() - n_constraints;
735 >    int nproc = MPI::COMM_WORLD.Get_size();
736  
737 +    // we need arrays to hold the counts and displacement vectors for
738 +    // all processors
739 +    vector<int> counts(nproc, 0);
740 +    vector<int> disps(nproc, 0);
741  
742 < #ifdef IS_MPI
743 <  MPI_Allreduce(&ndfTrans_local,&ndfTrans,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
744 < #else
745 <  ndfTrans = ndfTrans_local;
746 < #endif
742 >    // fill the counts array
743 >    MPI::COMM_WORLD.Allgather(&count_local, 1, MPI::INT, &counts[0],
744 >                              1, MPI::INT);
745 >  
746 >    // use the processor counts to compute the displacement array
747 >    disps[0] = 0;    
748 >    int totalCount = counts[0];
749 >    for (int iproc = 1; iproc < nproc; iproc++) {
750 >      disps[iproc] = disps[iproc-1] + counts[iproc-1];
751 >      totalCount += counts[iproc];
752 >    }
753  
754 <  ndfTrans = ndfTrans - 3 - nZconstraints;
754 >    // we need a (possibly redundant) set of all found types:
755 >    vector<int> ftGlobal(totalCount);
756 >    
757 >    // now spray out the foundTypes to all the other processors:    
758 >    MPI::COMM_WORLD.Allgatherv(&foundTypes[0], count_local, MPI::INT,
759 >                               &ftGlobal[0], &counts[0], &disps[0],
760 >                               MPI::INT);
761  
762 <  return ndfTrans;
402 < }
762 >    vector<int>::iterator j;
763  
764 < int SimInfo::getTotIntegrableObjects() {
765 <  int nObjs_local;
766 <  int nObjs;
764 >    // foundIdents is a stl set, so inserting an already found ident
765 >    // will have no effect.
766 >    set<int> foundIdents;
767  
768 <  nObjs_local =  integrableObjects.size();
768 >    for (j = ftGlobal.begin(); j != ftGlobal.end(); ++j)
769 >      foundIdents.insert((*j));
770 >    
771 >    // now iterate over the foundIdents and get the actual atom types
772 >    // that correspond to these:
773 >    set<int>::iterator it;
774 >    for (it = foundIdents.begin(); it != foundIdents.end(); ++it)
775 >      atomTypes.insert( forceField_->getAtomType((*it)) );
776 >
777 > #endif
778  
779 +    return atomTypes;        
780 +  }
781  
782 +  void SimInfo::setupSimVariables() {
783 +    useAtomicVirial_ = simParams_->getUseAtomicVirial();
784 +    // we only call setAccumulateBoxDipole if the accumulateBoxDipole parameter is true
785 +    calcBoxDipole_ = false;
786 +    if ( simParams_->haveAccumulateBoxDipole() )
787 +      if ( simParams_->getAccumulateBoxDipole() ) {
788 +        calcBoxDipole_ = true;      
789 +      }
790 +    
791 +    set<AtomType*>::iterator i;
792 +    set<AtomType*> atomTypes;
793 +    atomTypes = getSimulatedAtomTypes();    
794 +    bool usesElectrostatic = false;
795 +    bool usesMetallic = false;
796 +    bool usesDirectional = false;
797 +    bool usesFluctuatingCharges =  false;
798 +    //loop over all of the atom types
799 +    for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
800 +      usesElectrostatic |= (*i)->isElectrostatic();
801 +      usesMetallic |= (*i)->isMetal();
802 +      usesDirectional |= (*i)->isDirectional();
803 +      usesFluctuatingCharges |= (*i)->isFluctuatingCharge();
804 +    }
805 +
806   #ifdef IS_MPI
807 <  MPI_Allreduce(&nObjs_local,&nObjs,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
807 >    bool temp;
808 >    temp = usesDirectional;
809 >    MPI::COMM_WORLD.Allreduce(&temp, &usesDirectionalAtoms_, 1, MPI::BOOL,
810 >                              MPI::LOR);
811 >        
812 >    temp = usesMetallic;
813 >    MPI::COMM_WORLD.Allreduce(&temp, &usesMetallicAtoms_, 1, MPI::BOOL,
814 >                              MPI::LOR);
815 >    
816 >    temp = usesElectrostatic;
817 >    MPI::COMM_WORLD.Allreduce(&temp, &usesElectrostaticAtoms_, 1, MPI::BOOL,
818 >                              MPI::LOR);
819 >
820 >    temp = usesFluctuatingCharges;
821 >    MPI::COMM_WORLD.Allreduce(&temp, &usesFluctuatingCharges_, 1, MPI::BOOL,
822 >                              MPI::LOR);
823   #else
824 <  nObjs = nObjs_local;
824 >
825 >    usesDirectionalAtoms_ = usesDirectional;
826 >    usesMetallicAtoms_ = usesMetallic;
827 >    usesElectrostaticAtoms_ = usesElectrostatic;
828 >    usesFluctuatingCharges_ = usesFluctuatingCharges;
829 >
830   #endif
831 +    
832 +    requiresPrepair_ = usesMetallicAtoms_ ? true : false;
833 +    requiresSkipCorrection_ = usesElectrostaticAtoms_ ? true : false;
834 +    requiresSelfCorrection_ = usesElectrostaticAtoms_ ? true : false;    
835 +  }
836  
837  
838 <  return nObjs;
839 < }
838 >  vector<int> SimInfo::getGlobalAtomIndices() {
839 >    SimInfo::MoleculeIterator mi;
840 >    Molecule* mol;
841 >    Molecule::AtomIterator ai;
842 >    Atom* atom;
843  
844 < void SimInfo::refreshSim(){
844 >    vector<int> GlobalAtomIndices(getNAtoms(), 0);
845 >    
846 >    for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) {
847 >      
848 >      for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
849 >        GlobalAtomIndices[atom->getLocalIndex()] = atom->getGlobalIndex();
850 >      }
851 >    }
852 >    return GlobalAtomIndices;
853 >  }
854  
423  simtype fInfo;
424  int isError;
425  int n_global;
426  int* excl;
855  
856 <  fInfo.dielect = 0.0;
856 >  vector<int> SimInfo::getGlobalGroupIndices() {
857 >    SimInfo::MoleculeIterator mi;
858 >    Molecule* mol;
859 >    Molecule::CutoffGroupIterator ci;
860 >    CutoffGroup* cg;
861  
862 <  if( useDipoles ){
863 <    if( useReactionField )fInfo.dielect = dielectric;
862 >    vector<int> GlobalGroupIndices;
863 >    
864 >    for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) {
865 >      
866 >      //local index of cutoff group is trivial, it only depends on the
867 >      //order of travesing
868 >      for (cg = mol->beginCutoffGroup(ci); cg != NULL;
869 >           cg = mol->nextCutoffGroup(ci)) {
870 >        GlobalGroupIndices.push_back(cg->getGlobalIndex());
871 >      }        
872 >    }
873 >    return GlobalGroupIndices;
874    }
875  
434  fInfo.SIM_uses_PBC = usePBC;
435  //fInfo.SIM_uses_LJ = 0;
436  fInfo.SIM_uses_LJ = useLJ;
437  fInfo.SIM_uses_sticky = useSticky;
438  //fInfo.SIM_uses_sticky = 0;
439  fInfo.SIM_uses_charges = useCharges;
440  fInfo.SIM_uses_dipoles = useDipoles;
441  //fInfo.SIM_uses_dipoles = 0;
442  fInfo.SIM_uses_RF = useReactionField;
443  //fInfo.SIM_uses_RF = 0;
444  fInfo.SIM_uses_GB = useGB;
445  fInfo.SIM_uses_EAM = useEAM;
876  
877 <  n_exclude = excludes->getSize();
878 <  excl = excludes->getFortranArray();
449 <  
450 < #ifdef IS_MPI
451 <  n_global = mpiSim->getNAtomsGlobal();
452 < #else
453 <  n_global = n_atoms;
454 < #endif
455 <  
456 <  isError = 0;
457 <  
458 <  getFortranGroupArrays(this, FglobalGroupMembership, mfact);
459 <  //it may not be a good idea to pass the address of first element in vector
460 <  //since c++ standard does not require vector to be stored continuously in meomory
461 <  //Most of the compilers will organize the memory of vector continuously
462 <  setFsimulation( &fInfo, &n_global, &n_atoms, identArray, &n_exclude, excl,
463 <                  &nGlobalExcludes, globalExcludes, molMembershipArray,
464 <                  &mfact[0], &ngroup, &FglobalGroupMembership[0], &isError);
877 >  void SimInfo::prepareTopology() {
878 >    int nExclude, nOneTwo, nOneThree, nOneFour;
879  
880 <  if( isError ){
880 >    //calculate mass ratio of cutoff group
881 >    SimInfo::MoleculeIterator mi;
882 >    Molecule* mol;
883 >    Molecule::CutoffGroupIterator ci;
884 >    CutoffGroup* cg;
885 >    Molecule::AtomIterator ai;
886 >    Atom* atom;
887 >    RealType totalMass;
888 >
889 >    /**
890 >     * The mass factor is the relative mass of an atom to the total
891 >     * mass of the cutoff group it belongs to.  By default, all atoms
892 >     * are their own cutoff groups, and therefore have mass factors of
893 >     * 1.  We need some special handling for massless atoms, which
894 >     * will be treated as carrying the entire mass of the cutoff
895 >     * group.
896 >     */
897 >    massFactors_.clear();
898 >    massFactors_.resize(getNAtoms(), 1.0);
899      
900 <    sprintf( painCave.errMsg,
901 <             "There was an error setting the simulation information in fortran.\n" );
902 <    painCave.isFatal = 1;
471 <    painCave.severity = OOPSE_ERROR;
472 <    simError();
473 <  }
474 <  
475 < #ifdef IS_MPI
476 <  sprintf( checkPointMsg,
477 <           "succesfully sent the simulation information to fortran.\n");
478 <  MPIcheckPoint();
479 < #endif // is_mpi
480 <  
481 <  this->ndf = this->getNDF();
482 <  this->ndfRaw = this->getNDFraw();
483 <  this->ndfTrans = this->getNDFtranslational();
484 < }
900 >    for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {        
901 >      for (cg = mol->beginCutoffGroup(ci); cg != NULL;
902 >           cg = mol->nextCutoffGroup(ci)) {
903  
904 < void SimInfo::setDefaultRcut( double theRcut ){
905 <  
906 <  haveRcut = 1;
907 <  rCut = theRcut;
908 <  rList = rCut + 1.0;
909 <  
910 <  notifyFortranCutOffs( &rCut, &rSw, &rList );
911 < }
904 >        totalMass = cg->getMass();
905 >        for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
906 >          // Check for massless groups - set mfact to 1 if true
907 >          if (totalMass != 0)
908 >            massFactors_[atom->getLocalIndex()] = atom->getMass()/totalMass;
909 >          else
910 >            massFactors_[atom->getLocalIndex()] = 1.0;
911 >        }
912 >      }      
913 >    }
914  
915 < void SimInfo::setDefaultRcut( double theRcut, double theRsw ){
915 >    // Build the identArray_
916  
917 <  rSw = theRsw;
918 <  setDefaultRcut( theRcut );
919 < }
917 >    identArray_.clear();
918 >    identArray_.reserve(getNAtoms());    
919 >    for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {        
920 >      for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
921 >        identArray_.push_back(atom->getIdent());
922 >      }
923 >    }    
924 >    
925 >    //scan topology
926  
927 +    nExclude = excludedInteractions_.getSize();
928 +    nOneTwo = oneTwoInteractions_.getSize();
929 +    nOneThree = oneThreeInteractions_.getSize();
930 +    nOneFour = oneFourInteractions_.getSize();
931  
932 < void SimInfo::checkCutOffs( void ){
933 <  
934 <  if( boxIsInit ){
935 <    
936 <    //we need to check cutOffs against the box
937 <    
508 <    if( rCut > maxCutoff ){
509 <      sprintf( painCave.errMsg,
510 <               "cutoffRadius is too large for the current periodic box.\n"
511 <               "\tCurrent Value of cutoffRadius = %G at time %G\n "
512 <               "\tThis is larger than half of at least one of the\n"
513 <               "\tperiodic box vectors.  Right now, the Box matrix is:\n"
514 <               "\n"
515 <               "\t[ %G %G %G ]\n"
516 <               "\t[ %G %G %G ]\n"
517 <               "\t[ %G %G %G ]\n",
518 <               rCut, currentTime,
519 <               Hmat[0][0], Hmat[0][1], Hmat[0][2],
520 <               Hmat[1][0], Hmat[1][1], Hmat[1][2],
521 <               Hmat[2][0], Hmat[2][1], Hmat[2][2]);
522 <      painCave.severity = OOPSE_ERROR;
523 <      painCave.isFatal = 1;
524 <      simError();
525 <    }    
526 <  } else {
527 <    // initialize this stuff before using it, OK?
528 <    sprintf( painCave.errMsg,
529 <             "Trying to check cutoffs without a box.\n"
530 <             "\tOOPSE should have better programmers than that.\n" );
531 <    painCave.severity = OOPSE_ERROR;
532 <    painCave.isFatal = 1;
533 <    simError();      
932 >    int* excludeList = excludedInteractions_.getPairList();
933 >    int* oneTwoList = oneTwoInteractions_.getPairList();
934 >    int* oneThreeList = oneThreeInteractions_.getPairList();
935 >    int* oneFourList = oneFourInteractions_.getPairList();
936 >
937 >    topologyDone_ = true;
938    }
535  
536 }
939  
940 < void SimInfo::addProperty(GenericData* prop){
940 >  void SimInfo::addProperty(GenericData* genData) {
941 >    properties_.addProperty(genData);  
942 >  }
943  
944 <  map<string, GenericData*>::iterator result;
945 <  result = properties.find(prop->getID());
542 <  
543 <  //we can't simply use  properties[prop->getID()] = prop,
544 <  //it will cause memory leak if we already contain a propery which has the same name of prop
545 <  
546 <  if(result != properties.end()){
547 <    
548 <    delete (*result).second;
549 <    (*result).second = prop;
550 <      
944 >  void SimInfo::removeProperty(const string& propName) {
945 >    properties_.removeProperty(propName);  
946    }
552  else{
947  
948 <    properties[prop->getID()] = prop;
948 >  void SimInfo::clearProperties() {
949 >    properties_.clearProperties();
950 >  }
951  
952 +  vector<string> SimInfo::getPropertyNames() {
953 +    return properties_.getPropertyNames();  
954    }
955 <    
956 < }
955 >      
956 >  vector<GenericData*> SimInfo::getProperties() {
957 >    return properties_.getProperties();
958 >  }
959  
960 < GenericData* SimInfo::getProperty(const string& propName){
960 >  GenericData* SimInfo::getPropertyByName(const string& propName) {
961 >    return properties_.getPropertyByName(propName);
962 >  }
963 >
964 >  void SimInfo::setSnapshotManager(SnapshotManager* sman) {
965 >    if (sman_ == sman) {
966 >      return;
967 >    }    
968 >    delete sman_;
969 >    sman_ = sman;
970 >
971 >    Molecule* mol;
972 >    RigidBody* rb;
973 >    Atom* atom;
974 >    CutoffGroup* cg;
975 >    SimInfo::MoleculeIterator mi;
976 >    Molecule::RigidBodyIterator rbIter;
977 >    Molecule::AtomIterator atomIter;
978 >    Molecule::CutoffGroupIterator cgIter;
979  
980 <  map<string, GenericData*>::iterator result;
981 <  
982 <  //string lowerCaseName = ();
983 <  
984 <  result = properties.find(propName);
985 <  
986 <  if(result != properties.end())
987 <    return (*result).second;  
988 <  else  
571 <    return NULL;  
572 < }
980 >    for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
981 >        
982 >      for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) {
983 >        atom->setSnapshotManager(sman_);
984 >      }
985 >        
986 >      for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
987 >        rb->setSnapshotManager(sman_);
988 >      }
989  
990 +      for (cg = mol->beginCutoffGroup(cgIter); cg != NULL; cg = mol->nextCutoffGroup(cgIter)) {
991 +        cg->setSnapshotManager(sman_);
992 +      }
993 +    }    
994 +    
995 +  }
996  
997 < void SimInfo::getFortranGroupArrays(SimInfo* info,
998 <                                    vector<int>& FglobalGroupMembership,
999 <                                    vector<double>& mfact){
997 >
998 >  ostream& operator <<(ostream& o, SimInfo& info) {
999 >
1000 >    return o;
1001 >  }
1002 >  
1003    
1004 <  Molecule* myMols;
1005 <  Atom** myAtoms;
1006 <  int numAtom;
582 <  double mtot;
583 <  int numMol;
584 <  int numCutoffGroups;
585 <  CutoffGroup* myCutoffGroup;
586 <  vector<CutoffGroup*>::iterator iterCutoff;
587 <  Atom* cutoffAtom;
588 <  vector<Atom*>::iterator iterAtom;
589 <  int atomIndex;
590 <  double totalMass;
1004 >  StuntDouble* SimInfo::getIOIndexToIntegrableObject(int index) {
1005 >    return IOIndexToIntegrableObject.at(index);
1006 >  }
1007    
1008 <  mfact.clear();
1009 <  FglobalGroupMembership.clear();
1010 <  
1008 >  void SimInfo::setIOIndexToIntegrableObject(const vector<StuntDouble*>& v) {
1009 >    IOIndexToIntegrableObject= v;
1010 >  }
1011 > /*
1012 >   void SimInfo::setStuntDoubleFromGlobalIndex(vector<StuntDouble*> v) {
1013 >      assert( v.size() == nAtoms_ + nRigidBodies_);
1014 >      sdByGlobalIndex_ = v;
1015 >    }
1016  
1017 <  // Fix the silly fortran indexing problem
1017 >    StuntDouble* SimInfo::getStuntDoubleFromGlobalIndex(int index) {
1018 >      //assert(index < nAtoms_ + nRigidBodies_);
1019 >      return sdByGlobalIndex_.at(index);
1020 >    }  
1021 > */  
1022 >  int SimInfo::getNGlobalConstraints() {
1023 >    int nGlobalConstraints;
1024   #ifdef IS_MPI
1025 <  numAtom = mpiSim->getNAtomsGlobal();
1025 >    MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM,
1026 >                  MPI_COMM_WORLD);    
1027   #else
1028 <  numAtom = n_atoms;
1028 >    nGlobalConstraints =  nConstraints_;
1029   #endif
1030 <  for (int i = 0; i < numAtom; i++)
603 <    FglobalGroupMembership.push_back(globalGroupMembership[i] + 1);
604 <  
605 <
606 <  myMols = info->molecules;
607 <  numMol = info->n_mol;
608 <  for(int i  = 0; i < numMol; i++){
609 <    numCutoffGroups = myMols[i].getNCutoffGroups();
610 <    for(myCutoffGroup =myMols[i].beginCutoffGroup(iterCutoff);
611 <        myCutoffGroup != NULL;
612 <        myCutoffGroup =myMols[i].nextCutoffGroup(iterCutoff)){
613 <
614 <      totalMass = myCutoffGroup->getMass();
615 <      
616 <      for(cutoffAtom = myCutoffGroup->beginAtom(iterAtom);
617 <          cutoffAtom != NULL;
618 <          cutoffAtom = myCutoffGroup->nextAtom(iterAtom)){
619 <        mfact.push_back(cutoffAtom->getMass()/totalMass);
620 <      }  
621 <    }
1030 >    return nGlobalConstraints;
1031    }
1032  
1033 < }
1033 > }//end namespace OpenMD
1034 >

Comparing:
trunk/src/brains/SimInfo.cpp (property svn:keywords), Revision 3 by tim, Fri Sep 24 16:27:58 2004 UTC vs.
branches/development/src/brains/SimInfo.cpp (property svn:keywords), Revision 1769 by gezelter, Mon Jul 9 14:15:52 2012 UTC

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

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines