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

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

Comparing trunk/src/brains/SimInfo.cpp (property svn:keywords):
Revision 124 by chuckv, Wed Oct 20 20:46:20 2004 UTC vs.
Revision 1953 by gezelter, Thu Dec 5 18:19:26 2013 UTC

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

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines