ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
(Generate patch)

Comparing trunk/OOPSE/libmdtools/SimSetup.cpp (file contents):
Revision 572 by mmeineke, Wed Jul 2 21:26:55 2003 UTC vs.
Revision 1174 by gezelter, Wed May 12 20:54:10 2004 UTC

# Line 1 | Line 1
1 < #include <cstdlib>
1 > #include <algorithm>
2 > #include <stdlib.h>
3   #include <iostream>
4 < #include <cmath>
5 <
4 > #include <math.h>
5 > #include <string>
6 > #include <sprng.h>
7   #include "SimSetup.hpp"
8 + #include "ReadWrite.hpp"
9   #include "parse_me.h"
10   #include "Integrator.hpp"
11   #include "simError.h"
12 + #include "RigidBody.hpp"
13 + //#include "ConjugateMinimizer.hpp"
14 + #include "OOPSEMinimizer.hpp"
15  
16   #ifdef IS_MPI
17   #include "mpiBASS.h"
# Line 14 | Line 20
20  
21   // some defines for ensemble and Forcefield  cases
22  
23 < #define NVE_ENS 0
24 < #define NVT_ENS 1
25 < #define NPT_ENS 2
23 > #define NVE_ENS        0
24 > #define NVT_ENS        1
25 > #define NPTi_ENS       2
26 > #define NPTf_ENS       3
27 > #define NPTxyz_ENS     4
28  
21 #define FF_DUFF 0
22 #define FF_LJ   1
29  
30 + #define FF_DUFF  0
31 + #define FF_LJ    1
32 + #define FF_EAM   2
33 + #define FF_H2O   3
34  
35 + using namespace std;
36 +
37 + /**
38 + * Check whether dividend is divisble by divisor or not
39 + */
40 + bool isDivisible(double dividend, double divisor){
41 +  double tolerance = 0.000001;
42 +  double quotient;
43 +  double diff;
44 +  int intQuotient;
45 +  
46 +  quotient = dividend / divisor;
47 +
48 +  if (quotient < 0)
49 +    quotient = -quotient;
50 +
51 +  intQuotient = int (quotient + tolerance);
52 +
53 +  diff = fabs(fabs(dividend) - intQuotient  * fabs(divisor));
54 +
55 +  if (diff <= tolerance)
56 +    return true;
57 +  else
58 +    return false;  
59 + }
60 +
61   SimSetup::SimSetup(){
62 +  
63 +  initSuspend = false;
64 +  isInfoArray = 0;
65 +  nInfo = 1;
66 +
67    stamps = new MakeStamps();
68    globals = new Globals();
69 <  
69 >
70 >
71   #ifdef IS_MPI
72 <  strcpy( checkPointMsg, "SimSetup creation successful" );
72 >  strcpy(checkPointMsg, "SimSetup creation successful");
73    MPIcheckPoint();
74   #endif // IS_MPI
75   }
# Line 37 | Line 79 | SimSetup::~SimSetup(){
79    delete globals;
80   }
81  
82 < void SimSetup::parseFile( char* fileName ){
82 > void SimSetup::setSimInfo(SimInfo* the_info, int theNinfo){
83 >  info = the_info;
84 >  nInfo = theNinfo;
85 >  isInfoArray = 1;
86 >  initSuspend = true;
87 > }
88  
89 +
90 + void SimSetup::parseFile(char* fileName){
91   #ifdef IS_MPI
92 <  if( worldRank == 0 ){
92 >  if (worldRank == 0){
93   #endif // is_mpi
94 <    
94 >
95      inFileName = fileName;
96 <    set_interface_stamps( stamps, globals );
97 <    
96 >    set_interface_stamps(stamps, globals);
97 >
98   #ifdef IS_MPI
99      mpiEventInit();
100   #endif
101  
102 <    yacc_BASS( fileName );
102 >    yacc_BASS(fileName);
103  
104   #ifdef IS_MPI
105      throwMPIEvent(NULL);
106    }
107 <  else receiveParse();
107 >  else{
108 >    receiveParse();
109 >  }
110   #endif
111  
112   }
113  
114   #ifdef IS_MPI
115   void SimSetup::receiveParse(void){
116 <
117 <    set_interface_stamps( stamps, globals );
118 <    mpiEventInit();
119 <    MPIcheckPoint();
69 <    mpiEventLoop();
70 <
116 >  set_interface_stamps(stamps, globals);
117 >  mpiEventInit();
118 >  MPIcheckPoint();
119 >  mpiEventLoop();
120   }
121  
122   #endif // is_mpi
123  
124 < void SimSetup::createSim( void ){
124 > void SimSetup::createSim(void){
125  
126 <  MakeStamps *the_stamps;
78 <  Globals* the_globals;
79 <  int i, j, k, globalAtomIndex;
80 <  
81 <  int ensembleCase;
82 <  int ffCase;
83 <  
84 <  ensembleCase = -1;
85 <  ffCase = -1;
126 >  // gather all of the information from the Bass file
127  
128 <  // get the stamps and globals;
88 <  the_stamps = stamps;
89 <  the_globals = globals;
128 >  gatherInfo();
129  
130 <  // set the easy ones first
92 <  simnfo->target_temp = the_globals->getTargetTemp();
93 <  simnfo->dt = the_globals->getDt();
94 <  simnfo->run_time = the_globals->getRunTime();
130 >  // creation of complex system objects
131  
132 <  // get the ones we know are there, yet still may need some work.
97 <  n_components = the_globals->getNComponents();
98 <  strcpy( force_field, the_globals->getForceField() );
132 >  sysObjectsCreation();
133  
134 <  if( !strcasecmp( force_field, "DUFF" )) ffCase = FF_DUFF;
101 <  else if( !strcasecmp( force_field, "LJ" )) ffCase = FF_LJ;
102 <  else{
103 <    sprintf( painCave.errMsg,
104 <             "SimSetup Error. Unrecognized force field -> %s\n",
105 <             force_field );
106 <    painCave.isFatal = 1;
107 <    simError();
108 <  }
134 >  // check on the post processing info
135  
136 <  // get the ensemble:
111 <  strcpy( ensemble, the_globals->getEnsemble() );
136 >  finalInfoCheck();
137  
138 <  if( !strcasecmp( ensemble, "NVE" ))      ensembleCase = NVE_ENS;
139 <  else if( !strcasecmp( ensemble, "NVT" )) ensembleCase = NVT_ENS;
140 <  else if( !strcasecmp( ensemble, "NPT" )) ensembleCase = NPT_ENS;
141 <  else{
142 <    sprintf( painCave.errMsg,
143 <             "SimSetup Warning. Unrecognized Ensemble -> %s, "
144 <             "reverting to NVE for this simulation.\n",
120 <             ensemble );
121 <    painCave.isFatal = 0;
122 <    simError();
123 <    strcpy( ensemble, "NVE" );
124 <    ensembleCase = NVE_ENS;
138 >  // initialize the system coordinates
139 >
140 >  if ( !initSuspend ){
141 >    initSystemCoords();
142 >
143 >    if( !(globals->getUseInitTime()) )
144 >      info[0].currentTime = 0.0;
145    }  
126  strcpy( simnfo->ensemble, ensemble );
146  
147 +  // make the output filenames
148  
149 < //   if( !strcasecmp( ensemble, "NPT" ) ) {
150 < //     the_extendedsystem = new ExtendedSystem( simnfo );
151 < //     the_extendedsystem->setTargetTemp(the_globals->getTargetTemp());
152 < //     if (the_globals->haveTargetPressure())
153 < //       the_extendedsystem->setTargetPressure(the_globals->getTargetPressure());
134 < //     else {
135 < //       sprintf( painCave.errMsg,
136 < //                "SimSetup error: If you use the constant pressure\n"
137 < //                "    ensemble, you must set targetPressure.\n"
138 < //                "    This was found in the BASS file.\n");
139 < //       painCave.isFatal = 1;
140 < //       simError();
141 < //     }
149 >  makeOutNames();
150 >  
151 > #ifdef IS_MPI
152 >  mpiSim->mpiRefresh();
153 > #endif
154  
155 < //     if (the_globals->haveTauThermostat())
144 < //       the_extendedsystem->setTauThermostat(the_globals->getTauThermostat());
145 < //     else if (the_globals->haveQmass())
146 < //       the_extendedsystem->setQmass(the_globals->getQmass());
147 < //     else {
148 < //       sprintf( painCave.errMsg,
149 < //                "SimSetup error: If you use one of the constant temperature\n"
150 < //                "    ensembles, you must set either tauThermostat or qMass.\n"
151 < //                "    Neither of these was found in the BASS file.\n");
152 < //       painCave.isFatal = 1;
153 < //       simError();
154 < //     }
155 >  // initialize the Fortran
156  
157 < //     if (the_globals->haveTauBarostat())
157 < //       the_extendedsystem->setTauBarostat(the_globals->getTauBarostat());
158 < //     else {
159 < //       sprintf( painCave.errMsg,
160 < //                "SimSetup error: If you use the constant pressure\n"
161 < //                "    ensemble, you must set tauBarostat.\n"
162 < //                "    This was found in the BASS file.\n");
163 < //       painCave.isFatal = 1;
164 < //       simError();
165 < //     }
157 >  initFortran();
158  
159 < //   } else if ( !strcasecmp( ensemble, "NVT") ) {
160 < //     the_extendedsystem = new ExtendedSystem( simnfo );
161 < //     the_extendedsystem->setTargetTemp(the_globals->getTargetTemp());
159 >  if (globals->haveMinimizer())
160 >    // make minimizer
161 >    makeMinimizer();
162 >  else
163 >    // make the integrator
164 >    makeIntegrator();
165  
166 < //     if (the_globals->haveTauThermostat())
172 < //       the_extendedsystem->setTauThermostat(the_globals->getTauThermostat());
173 < //     else if (the_globals->haveQmass())
174 < //       the_extendedsystem->setQmass(the_globals->getQmass());
175 < //     else {
176 < //       sprintf( painCave.errMsg,
177 < //                "SimSetup error: If you use one of the constant temperature\n"
178 < //                "    ensembles, you must set either tauThermostat or qMass.\n"
179 < //                "    Neither of these was found in the BASS file.\n");
180 < //       painCave.isFatal = 1;
181 < //       simError();
182 < //     }
166 > }
167  
184  strcpy( simnfo->mixingRule, the_globals->getMixingRule() );
185  simnfo->usePBC = the_globals->getPBC();
186          
187  int usesDipoles = 0;
188  switch( ffCase ){
168  
169 <  case FF_DUFF:
170 <    the_ff = new DUFF();
171 <    usesDipoles = 1;
172 <    break;
169 > void SimSetup::makeMolecules(void){
170 >  int i, j, k;
171 >  int exI, exJ, exK, exL, slI, slJ;
172 >  int tempI, tempJ, tempK, tempL;
173 >  int molI;
174 >  int stampID, atomOffset, rbOffset;
175 >  molInit molInfo;
176 >  DirectionalAtom* dAtom;
177 >  RigidBody* myRB;
178 >  StuntDouble* mySD;
179 >  LinkedAssign* extras;
180 >  LinkedAssign* current_extra;
181 >  AtomStamp* currentAtom;
182 >  BondStamp* currentBond;
183 >  BendStamp* currentBend;
184 >  TorsionStamp* currentTorsion;
185 >  RigidBodyStamp* currentRigidBody;
186 >  CutoffGroupStamp* currentCutoffGroup;
187 >  CutoffGroup* myCutoffGroup;
188 >  int nCutoffGroups;// number of cutoff group of a molecule defined in mdl file
189 >  set<int> cutoffAtomSet; //atoms belong to  cutoffgroup defined at mdl file
190  
191 <  case FF_LJ:
192 <    the_ff = new LJFF();
193 <    break;
191 >  bond_pair* theBonds;
192 >  bend_set* theBends;
193 >  torsion_set* theTorsions;
194  
195 <  default:
200 <    sprintf( painCave.errMsg,
201 <             "SimSetup Error. Unrecognized force field in case statement.\n");
202 <    painCave.isFatal = 1;
203 <    simError();
204 <  }
195 >  set<int> skipList;
196  
197 < #ifdef IS_MPI
198 <  strcpy( checkPointMsg, "ForceField creation successful" );
199 <  MPIcheckPoint();
209 < #endif // is_mpi
197 >  double phi, theta, psi;
198 >  char* molName;
199 >  char rbName[100];
200  
201 <  // get the components and calculate the tot_nMol and indvidual n_mol
212 <  the_components = the_globals->getComponents();
213 <  components_nmol = new int[n_components];
214 <  comp_stamps = new MoleculeStamp*[n_components];
201 >  //init the forceField paramters
202  
203 <  if( !the_globals->haveNMol() ){
217 <    // we don't have the total number of molecules, so we assume it is
218 <    // given in each component
203 >  the_ff->readParams();
204  
205 <    tot_nmol = 0;
221 <    for( i=0; i<n_components; i++ ){
205 >  // init the atoms
206  
207 <      if( !the_components[i]->haveNMol() ){
208 <        // we have a problem
209 <        sprintf( painCave.errMsg,
210 <                 "SimSetup Error. No global NMol or component NMol"
211 <                 " given. Cannot calculate the number of atoms.\n" );
212 <        painCave.isFatal = 1;
213 <        simError();
214 <      }
207 >  int nMembers, nNew, rb1, rb2;
208 >
209 >  for (k = 0; k < nInfo; k++){
210 >    the_ff->setSimInfo(&(info[k]));
211 >
212 >    atomOffset = 0;
213 >
214 >    for (i = 0; i < info[k].n_mol; i++){
215 >      stampID = info[k].molecules[i].getStampID();
216 >      molName = comp_stamps[stampID]->getID();
217  
218 <      tot_nmol += the_components[i]->getNMol();
219 <      components_nmol[i] = the_components[i]->getNMol();
220 <    }
221 <  }
222 <  else{
237 <    sprintf( painCave.errMsg,
238 <             "SimSetup error.\n"
239 <             "\tSorry, the ability to specify total"
240 <             " nMols and then give molfractions in the components\n"
241 <             "\tis not currently supported."
242 <             " Please give nMol in the components.\n" );
243 <    painCave.isFatal = 1;
244 <    simError();
245 <    
246 <    
247 <    //     tot_nmol = the_globals->getNMol();
248 <    
249 <    //   //we have the total number of molecules, now we check for molfractions
250 <    //     for( i=0; i<n_components; i++ ){
251 <    
252 <    //       if( !the_components[i]->haveMolFraction() ){
253 <    
254 <    //  if( !the_components[i]->haveNMol() ){
255 <    //    //we have a problem
256 <    //    std::cerr << "SimSetup error. Neither molFraction nor "
257 <    //              << " nMol was given in component
258 <    
259 <  }
218 >      molInfo.nAtoms = comp_stamps[stampID]->getNAtoms();
219 >      molInfo.nBonds = comp_stamps[stampID]->getNBonds();
220 >      molInfo.nBends = comp_stamps[stampID]->getNBends();
221 >      molInfo.nTorsions = comp_stamps[stampID]->getNTorsions();
222 >      molInfo.nRigidBodies = comp_stamps[stampID]->getNRigidBodies();
223  
224 < #ifdef IS_MPI
225 <  strcpy( checkPointMsg, "Have the number of components" );
226 <  MPIcheckPoint();
264 < #endif // is_mpi
224 >      nCutoffGroups = comp_stamps[stampID]->getNCutoffGroups();
225 >      
226 >      molInfo.myAtoms = &(info[k].atoms[atomOffset]);
227  
228 <  // make an array of molecule stamps that match the components used.
229 <  // also extract the used stamps out into a separate linked list
228 >      if (molInfo.nBonds > 0)
229 >        molInfo.myBonds = new Bond*[molInfo.nBonds];
230 >      else
231 >        molInfo.myBonds = NULL;
232  
233 <  simnfo->nComponents = n_components;
234 <  simnfo->componentsNmol = components_nmol;
235 <  simnfo->compStamps = comp_stamps;
236 <  simnfo->headStamp = new LinkedMolStamp();
273 <  
274 <  char* id;
275 <  LinkedMolStamp* headStamp = simnfo->headStamp;
276 <  LinkedMolStamp* currentStamp = NULL;
277 <  for( i=0; i<n_components; i++ ){
233 >      if (molInfo.nBends > 0)
234 >        molInfo.myBends = new Bend*[molInfo.nBends];
235 >      else
236 >        molInfo.myBends = NULL;
237  
238 <    id = the_components[i]->getType();
239 <    comp_stamps[i] = NULL;
240 <    
241 <    // check to make sure the component isn't already in the list
238 >      if (molInfo.nTorsions > 0)
239 >        molInfo.myTorsions = new Torsion *[molInfo.nTorsions];
240 >      else
241 >        molInfo.myTorsions = NULL;
242  
243 <    comp_stamps[i] = headStamp->match( id );
244 <    if( comp_stamps[i] == NULL ){
243 >      theBonds = new bond_pair[molInfo.nBonds];
244 >      theBends = new bend_set[molInfo.nBends];
245 >      theTorsions = new torsion_set[molInfo.nTorsions];
246        
247 <      // extract the component from the list;
288 <      
289 <      currentStamp = the_stamps->extractMolStamp( id );
290 <      if( currentStamp == NULL ){
291 <        sprintf( painCave.errMsg,
292 <                 "SimSetup error: Component \"%s\" was not found in the "
293 <                 "list of declared molecules\n",
294 <                 id );
295 <        painCave.isFatal = 1;
296 <        simError();
297 <      }
298 <      
299 <      headStamp->add( currentStamp );
300 <      comp_stamps[i] = headStamp->match( id );
301 <    }
302 <  }
247 >      // make the Atoms
248  
249 < #ifdef IS_MPI
250 <  strcpy( checkPointMsg, "Component stamps successfully extracted\n" );
306 <  MPIcheckPoint();
307 < #endif // is_mpi
308 <  
249 >      for (j = 0; j < molInfo.nAtoms; j++){
250 >        currentAtom = comp_stamps[stampID]->getAtom(j);
251  
252 +        if (currentAtom->haveOrientation()){
253 +          dAtom = new DirectionalAtom((j + atomOffset),
254 +                                      info[k].getConfiguration());
255 +          info[k].n_oriented++;
256 +          molInfo.myAtoms[j] = dAtom;
257  
258 +          // Directional Atoms have standard unit vectors which are oriented
259 +          // in space using the three Euler angles.  We assume the standard
260 +          // unit vector was originally along the z axis below.
261  
262 <  // caclulate the number of atoms, bonds, bends and torsions
262 >          phi = currentAtom->getEulerPhi() * M_PI / 180.0;
263 >          theta = currentAtom->getEulerTheta() * M_PI / 180.0;
264 >          psi = currentAtom->getEulerPsi()* M_PI / 180.0;
265  
266 <  tot_atoms = 0;
267 <  tot_bonds = 0;
268 <  tot_bends = 0;
269 <  tot_torsions = 0;
318 <  for( i=0; i<n_components; i++ ){
319 <    
320 <    tot_atoms +=    components_nmol[i] * comp_stamps[i]->getNAtoms();
321 <    tot_bonds +=    components_nmol[i] * comp_stamps[i]->getNBonds();
322 <    tot_bends +=    components_nmol[i] * comp_stamps[i]->getNBends();
323 <    tot_torsions += components_nmol[i] * comp_stamps[i]->getNTorsions();
324 <  }
266 >          dAtom->setUnitFrameFromEuler(phi, theta, psi);
267 >            
268 >        }
269 >        else{
270  
271 <  tot_SRI = tot_bonds + tot_bends + tot_torsions;
271 >          molInfo.myAtoms[j] = new Atom((j + atomOffset), info[k].getConfiguration());
272  
273 <  simnfo->n_atoms = tot_atoms;
329 <  simnfo->n_bonds = tot_bonds;
330 <  simnfo->n_bends = tot_bends;
331 <  simnfo->n_torsions = tot_torsions;
332 <  simnfo->n_SRI = tot_SRI;
333 <  simnfo->n_mol = tot_nmol;
334 <  
335 <  simnfo->molMembershipArray = new int[tot_atoms];
273 >        }
274  
275 +        molInfo.myAtoms[j]->setType(currentAtom->getType());
276   #ifdef IS_MPI
277  
278 <  // divide the molecules among processors here.
340 <  
341 <  mpiSim = new mpiSimulation( simnfo );
342 <  
343 <  globalIndex = mpiSim->divideLabor();
278 >        molInfo.myAtoms[j]->setGlobalIndex(globalAtomIndex[j + atomOffset]);
279  
280 <  // set up the local variables
281 <  
347 <  int localMol, allMol;
348 <  int local_atoms, local_bonds, local_bends, local_torsions, local_SRI;
280 > #endif // is_mpi
281 >      }
282  
283 <  int* mol2proc = mpiSim->getMolToProcMap();
284 <  int* molCompType = mpiSim->getMolComponentType();
285 <  
286 <  allMol = 0;
287 <  localMol = 0;
355 <  local_atoms = 0;
356 <  local_bonds = 0;
357 <  local_bends = 0;
358 <  local_torsions = 0;
359 <  globalAtomIndex = 0;
283 >      // make the bonds
284 >      for (j = 0; j < molInfo.nBonds; j++){
285 >        currentBond = comp_stamps[stampID]->getBond(j);
286 >        theBonds[j].a = currentBond->getA() + atomOffset;
287 >        theBonds[j].b = currentBond->getB() + atomOffset;
288  
289 +        tempI = theBonds[j].a;
290 +        tempJ = theBonds[j].b;
291  
292 <  for( i=0; i<n_components; i++ ){
292 > #ifdef IS_MPI
293 >        exI = info[k].atoms[tempI]->getGlobalIndex() + 1;
294 >        exJ = info[k].atoms[tempJ]->getGlobalIndex() + 1;
295 > #else
296 >        exI = tempI + 1;
297 >        exJ = tempJ + 1;
298 > #endif
299  
300 <    for( j=0; j<components_nmol[i]; j++ ){
365 <      
366 <      if( mol2proc[allMol] == worldRank ){
367 <        
368 <        local_atoms +=    comp_stamps[i]->getNAtoms();
369 <        local_bonds +=    comp_stamps[i]->getNBonds();
370 <        local_bends +=    comp_stamps[i]->getNBends();
371 <        local_torsions += comp_stamps[i]->getNTorsions();
372 <        localMol++;
373 <      }      
374 <      for (k = 0; k < comp_stamps[i]->getNAtoms(); k++) {
375 <        simnfo->molMembershipArray[globalAtomIndex] = allMol;
376 <        globalAtomIndex++;
300 >        info[k].excludes->addPair(exI, exJ);
301        }
302  
303 <      allMol++;      
304 <    }
305 <  }
306 <  local_SRI = local_bonds + local_bends + local_torsions;
307 <  
308 <  simnfo->n_atoms = mpiSim->getMyNlocal();  
385 <  
386 <  if( local_atoms != simnfo->n_atoms ){
387 <    sprintf( painCave.errMsg,
388 <             "SimSetup error: mpiSim's localAtom (%d) and SimSetup's"
389 <             " localAtom (%d) are not equal.\n",
390 <             simnfo->n_atoms,
391 <             local_atoms );
392 <    painCave.isFatal = 1;
393 <    simError();
394 <  }
303 >      //make the bends
304 >      for (j = 0; j < molInfo.nBends; j++){
305 >        currentBend = comp_stamps[stampID]->getBend(j);
306 >        theBends[j].a = currentBend->getA() + atomOffset;
307 >        theBends[j].b = currentBend->getB() + atomOffset;
308 >        theBends[j].c = currentBend->getC() + atomOffset;
309  
310 <  simnfo->n_bonds = local_bonds;
311 <  simnfo->n_bends = local_bends;
312 <  simnfo->n_torsions = local_torsions;
399 <  simnfo->n_SRI = local_SRI;
400 <  simnfo->n_mol = localMol;
310 >        if (currentBend->haveExtras()){
311 >          extras = currentBend->getExtras();
312 >          current_extra = extras;
313  
314 <  strcpy( checkPointMsg, "Passed nlocal consistency check." );
315 <  MPIcheckPoint();
316 <  
317 <  
318 < #endif // is_mpi
319 <  
314 >          while (current_extra != NULL){
315 >            if (!strcmp(current_extra->getlhs(), "ghostVectorSource")){
316 >              switch (current_extra->getType()){
317 >                case 0:
318 >                  theBends[j].ghost = current_extra->getInt() + atomOffset;
319 >                  theBends[j].isGhost = 1;
320 >                  break;
321  
322 <  // create the atom and short range interaction arrays
322 >                case 1:
323 >                  theBends[j].ghost = (int) current_extra->getDouble() +
324 >                                      atomOffset;
325 >                  theBends[j].isGhost = 1;
326 >                  break;
327  
328 <  Atom::createArrays(simnfo->n_atoms);
329 <  the_atoms = new Atom*[simnfo->n_atoms];
330 <  the_molecules = new Molecule[simnfo->n_mol];
331 <  int molIndex;
332 <
333 <  // initialize the molecule's stampID's
328 >                default:
329 >                  sprintf(painCave.errMsg,
330 >                          "SimSetup Error: ghostVectorSource was neither a "
331 >                          "double nor an int.\n"
332 >                          "-->Bend[%d] in %s\n",
333 >                          j, comp_stamps[stampID]->getID());
334 >                  painCave.isFatal = 1;
335 >                  simError();
336 >              }
337 >            }
338 >            else{
339 >              sprintf(painCave.errMsg,
340 >                      "SimSetup Error: unhandled bend assignment:\n"
341 >                      "    -->%s in Bend[%d] in %s\n",
342 >                      current_extra->getlhs(), j, comp_stamps[stampID]->getID());
343 >              painCave.isFatal = 1;
344 >              simError();
345 >            }
346  
347 +            current_extra = current_extra->getNext();
348 +          }
349 +        }
350 +
351 +        if (theBends[j].isGhost) {
352 +          
353 +          tempI = theBends[j].a;
354 +          tempJ = theBends[j].b;
355 +          
356   #ifdef IS_MPI
357 <  
357 >          exI = info[k].atoms[tempI]->getGlobalIndex() + 1;
358 >          exJ = info[k].atoms[tempJ]->getGlobalIndex() + 1;
359 > #else
360 >          exI = tempI + 1;
361 >          exJ = tempJ + 1;
362 > #endif          
363 >          info[k].excludes->addPair(exI, exJ);
364  
365 <  molIndex = 0;
422 <  for(i=0; i<mpiSim->getTotNmol(); i++){
423 <    
424 <    if(mol2proc[i] == worldRank ){
425 <      the_molecules[molIndex].setStampID( molCompType[i] );
426 <      the_molecules[molIndex].setMyIndex( molIndex );
427 <      the_molecules[molIndex].setGlobalIndex( i );
428 <      molIndex++;
429 <    }
430 <  }
365 >        } else {
366  
367 < #else // is_mpi
368 <  
369 <  molIndex = 0;
370 <  globalAtomIndex = 0;
371 <  for(i=0; i<n_components; i++){
372 <    for(j=0; j<components_nmol[i]; j++ ){
373 <      the_molecules[molIndex].setStampID( i );
374 <      the_molecules[molIndex].setMyIndex( molIndex );
375 <      the_molecules[molIndex].setGlobalIndex( molIndex );
376 <      for (k = 0; k < comp_stamps[i]->getNAtoms(); k++) {
377 <        simnfo->molMembershipArray[globalAtomIndex] = molIndex;
378 <        globalAtomIndex++;
367 >          tempI = theBends[j].a;
368 >          tempJ = theBends[j].b;
369 >          tempK = theBends[j].c;
370 >          
371 > #ifdef IS_MPI
372 >          exI = info[k].atoms[tempI]->getGlobalIndex() + 1;
373 >          exJ = info[k].atoms[tempJ]->getGlobalIndex() + 1;
374 >          exK = info[k].atoms[tempK]->getGlobalIndex() + 1;
375 > #else
376 >          exI = tempI + 1;
377 >          exJ = tempJ + 1;
378 >          exK = tempK + 1;
379 > #endif
380 >          
381 >          info[k].excludes->addPair(exI, exK);
382 >          info[k].excludes->addPair(exI, exJ);
383 >          info[k].excludes->addPair(exJ, exK);
384 >        }
385        }
445      molIndex++;
446    }
447  }
448    
386  
387 < #endif // is_mpi
387 >      for (j = 0; j < molInfo.nTorsions; j++){
388 >        currentTorsion = comp_stamps[stampID]->getTorsion(j);
389 >        theTorsions[j].a = currentTorsion->getA() + atomOffset;
390 >        theTorsions[j].b = currentTorsion->getB() + atomOffset;
391 >        theTorsions[j].c = currentTorsion->getC() + atomOffset;
392 >        theTorsions[j].d = currentTorsion->getD() + atomOffset;
393  
394 +        tempI = theTorsions[j].a;      
395 +        tempJ = theTorsions[j].b;
396 +        tempK = theTorsions[j].c;
397 +        tempL = theTorsions[j].d;
398  
399 <  if( simnfo->n_SRI ){
400 <    
401 <    Exclude::createArray(simnfo->n_SRI);
402 <    the_excludes = new Exclude*[simnfo->n_SRI];
403 <    for( int ex=0; ex<simnfo->n_SRI; ex++) the_excludes[ex] = new Exclude(ex);
404 <    simnfo->globalExcludes = new int;
405 <    simnfo->n_exclude = simnfo->n_SRI;
406 <  }
407 <  else{
408 <    
409 <    Exclude::createArray( 1 );
464 <    the_excludes = new Exclude*;
465 <    the_excludes[0] = new Exclude(0);
466 <    the_excludes[0]->setPair( 0,0 );
467 <    simnfo->globalExcludes = new int;
468 <    simnfo->globalExcludes[0] = 0;
469 <    simnfo->n_exclude = 0;
470 <  }
399 > #ifdef IS_MPI
400 >        exI = info[k].atoms[tempI]->getGlobalIndex() + 1;
401 >        exJ = info[k].atoms[tempJ]->getGlobalIndex() + 1;
402 >        exK = info[k].atoms[tempK]->getGlobalIndex() + 1;
403 >        exL = info[k].atoms[tempL]->getGlobalIndex() + 1;
404 > #else
405 >        exI = tempI + 1;
406 >        exJ = tempJ + 1;
407 >        exK = tempK + 1;
408 >        exL = tempL + 1;
409 > #endif
410  
411 <  // set the arrays into the SimInfo object
411 >        info[k].excludes->addPair(exI, exJ);
412 >        info[k].excludes->addPair(exI, exK);
413 >        info[k].excludes->addPair(exI, exL);        
414 >        info[k].excludes->addPair(exJ, exK);
415 >        info[k].excludes->addPair(exJ, exL);
416 >        info[k].excludes->addPair(exK, exL);
417 >      }
418  
419 <  simnfo->atoms = the_atoms;
420 <  simnfo->molecules = the_molecules;
421 <  simnfo->nGlobalExcludes = 0;
422 <  simnfo->excludes = the_excludes;
419 >      
420 >      molInfo.myRigidBodies.clear();
421 >      
422 >      for (j = 0; j < molInfo.nRigidBodies; j++){
423  
424 +        currentRigidBody = comp_stamps[stampID]->getRigidBody(j);
425 +        nMembers = currentRigidBody->getNMembers();
426  
427 <  // get some of the tricky things that may still be in the globals
427 >        // Create the Rigid Body:
428  
429 <  double boxVector[3];
483 <  if( the_globals->haveBox() ){
484 <    boxVector[0] = the_globals->getBox();
485 <    boxVector[1] = the_globals->getBox();
486 <    boxVector[2] = the_globals->getBox();
487 <    
488 <    simnfo->setBox( boxVector );
489 <  }
490 <  else if( the_globals->haveDensity() ){
429 >        myRB = new RigidBody();
430  
431 <    double vol;
432 <    vol = (double)tot_nmol / the_globals->getDensity();
433 <     boxVector[0] = pow( vol, ( 1.0 / 3.0 ) );
434 <     boxVector[1] = boxVector[0];
496 <     boxVector[2] = boxVector[0];
431 >        sprintf(rbName,"%s_RB_%d", molName, j);
432 >        myRB->setType(rbName);
433 >        
434 >        for (rb1 = 0; rb1 < nMembers; rb1++) {
435  
436 <    simnfo->setBox( boxVector );
437 <  }
500 <  else{
501 <    if( !the_globals->haveBoxX() ){
502 <      sprintf( painCave.errMsg,
503 <               "SimSetup error, no periodic BoxX size given.\n" );
504 <      painCave.isFatal = 1;
505 <      simError();
506 <    }
507 <    boxVector[0] = the_globals->getBoxX();
436 >          // molI is atom numbering inside this molecule
437 >          molI = currentRigidBody->getMember(rb1);    
438  
439 <    if( !the_globals->haveBoxY() ){
440 <      sprintf( painCave.errMsg,
511 <               "SimSetup error, no periodic BoxY size given.\n" );
512 <      painCave.isFatal = 1;
513 <      simError();
514 <    }
515 <    boxVector[1] = the_globals->getBoxY();
439 >          // tempI is atom numbering on local processor
440 >          tempI = molI + atomOffset;
441  
442 <    if( !the_globals->haveBoxZ() ){
443 <      sprintf( painCave.errMsg,
444 <               "SimSetup error, no periodic BoxZ size given.\n" );
520 <      painCave.isFatal = 1;
521 <      simError();
522 <    }
523 <    boxVector[2] = the_globals->getBoxZ();
442 >          // currentAtom is the AtomStamp (which we need for
443 >          // rigid body reference positions)
444 >          currentAtom = comp_stamps[stampID]->getAtom(molI);
445  
446 <    simnfo->setBox( boxVector );
447 <  }
446 >          // When we add to the rigid body, add the atom itself and
447 >          // the stamp info:
448  
449 +          myRB->addAtom(info[k].atoms[tempI], currentAtom);
450 +          
451 +          // Add this atom to the Skip List for the integrators
452   #ifdef IS_MPI
453 <  strcpy( checkPointMsg, "Box size set up" );
454 <  MPIcheckPoint();
455 < #endif // is_mpi
453 >          slI = info[k].atoms[tempI]->getGlobalIndex();
454 > #else
455 >          slI = tempI;
456 > #endif
457 >          skipList.insert(slI);
458 >          
459 >        }
460 >        
461 >        for(rb1 = 0; rb1 < nMembers - 1; rb1++) {
462 >          for(rb2 = rb1+1; rb2 < nMembers; rb2++) {
463 >            
464 >            tempI = currentRigidBody->getMember(rb1);
465 >            tempJ = currentRigidBody->getMember(rb2);
466 >            
467 >            // Some explanation is required here.
468 >            // Fortran indexing starts at 1, while c indexing starts at 0
469 >            // Also, in parallel computations, the GlobalIndex is
470 >            // used for the exclude list:
471 >            
472 > #ifdef IS_MPI
473 >            exI = molInfo.myAtoms[tempI]->getGlobalIndex() + 1;
474 >            exJ = molInfo.myAtoms[tempJ]->getGlobalIndex() + 1;
475 > #else
476 >            exI = molInfo.myAtoms[tempI]->getIndex() + 1;
477 >            exJ = molInfo.myAtoms[tempJ]->getIndex() + 1;
478 > #endif
479 >            
480 >            info[k].excludes->addPair(exI, exJ);
481 >            
482 >          }
483 >        }
484  
485 +        molInfo.myRigidBodies.push_back(myRB);
486 +        info[k].rigidBodies.push_back(myRB);
487 +      }
488 +      
489  
490 <  // initialize the arrays
490 >      //create cutoff group for molecule
491  
492 <  the_ff->setSimInfo( simnfo );
492 >      cutoffAtomSet.clear();
493 >      molInfo.myCutoffGroups.clear();
494 >      
495 >      for (j = 0; j < nCutoffGroups; j++){
496  
497 <  makeMolecules();
498 <  simnfo->identArray = new int[simnfo->n_atoms];
540 <  for(i=0; i<simnfo->n_atoms; i++){
541 <    simnfo->identArray[i] = the_atoms[i]->getIdent();
542 <  }
543 <  
544 <  if (the_globals->getUseRF() ) {
545 <    simnfo->useReactionField = 1;
546 <  
547 <    if( !the_globals->haveECR() ){
548 <      sprintf( painCave.errMsg,
549 <               "SimSetup Warning: using default value of 1/2 the smallest "
550 <               "box length for the electrostaticCutoffRadius.\n"
551 <               "I hope you have a very fast processor!\n");
552 <      painCave.isFatal = 0;
553 <      simError();
554 <      double smallest;
555 <      smallest = simnfo->boxLx;
556 <      if (simnfo->boxLy <= smallest) smallest = simnfo->boxLy;
557 <      if (simnfo->boxLz <= smallest) smallest = simnfo->boxLz;
558 <      simnfo->ecr = 0.5 * smallest;
559 <    } else {
560 <      simnfo->ecr        = the_globals->getECR();
561 <    }
497 >        currentCutoffGroup = comp_stamps[stampID]->getCutoffGroup(j);
498 >        nMembers = currentCutoffGroup->getNMembers();
499  
500 <    if( !the_globals->haveEST() ){
501 <      sprintf( painCave.errMsg,
502 <               "SimSetup Warning: using default value of 0.05 * the "
503 <               "electrostaticCutoffRadius for the electrostaticSkinThickness\n"
504 <               );
505 <      painCave.isFatal = 0;
506 <      simError();
507 <      simnfo->est = 0.05 * simnfo->ecr;
508 <    } else {
509 <      simnfo->est        = the_globals->getEST();
510 <    }
511 <    
512 <    if(!the_globals->haveDielectric() ){
513 <      sprintf( painCave.errMsg,
514 <               "SimSetup Error: You are trying to use Reaction Field without"
515 <               "setting a dielectric constant!\n"
516 <               );
517 <      painCave.isFatal = 1;
518 <      simError();
519 <    }
520 <    simnfo->dielectric = the_globals->getDielectric();  
521 <  } else {
522 <    if (usesDipoles) {
523 <      
524 <      if( !the_globals->haveECR() ){
525 <        sprintf( painCave.errMsg,
526 <                 "SimSetup Warning: using default value of 1/2 the smallest "
527 <                 "box length for the electrostaticCutoffRadius.\n"
591 <                 "I hope you have a very fast processor!\n");
592 <        painCave.isFatal = 0;
593 <        simError();
594 <        double smallest;
595 <        smallest = simnfo->boxLx;
596 <        if (simnfo->boxLy <= smallest) smallest = simnfo->boxLy;
597 <        if (simnfo->boxLz <= smallest) smallest = simnfo->boxLz;
598 <        simnfo->ecr = 0.5 * smallest;
599 <      } else {
600 <        simnfo->ecr        = the_globals->getECR();
500 >        myCutoffGroup = new CutoffGroup();
501 >        
502 >        for (int cg = 0; cg < nMembers; cg++) {
503 >
504 >          // molI is atom numbering inside this molecule
505 >          molI = currentCutoffGroup->getMember(cg);    
506 >
507 >          // tempI is atom numbering on local processor
508 >          tempI = molI + atomOffset;
509 >          
510 >          myCutoffGroup->addAtom(info[k].atoms[tempI]);          
511 >
512 >          cutoffAtomSet.insert(tempI);
513 >        }
514 >
515 >        molInfo.myCutoffGroups.push_back(myCutoffGroup);
516 >      }//end for (j = 0; j < molInfo.nCutoffGroups; j++)
517 >
518 >      //creat a cutoff group for every atom  in current molecule which does not belong to cutoffgroup defined at mdl file
519 >
520 >      for(j = 0; j < molInfo.nAtoms; j++){
521 >
522 >        if(cutoffAtomSet.find(molInfo.myAtoms[j]->getIndex()) == cutoffAtomSet.end()){
523 >          myCutoffGroup = new CutoffGroup();
524 >          myCutoffGroup->addAtom(molInfo.myAtoms[j]);
525 >          molInfo.myCutoffGroups.push_back(myCutoffGroup);
526 >        }
527 >          
528        }
529 +
530 +              
531 +
532 +
533 +      // After this is all set up, scan through the atoms to
534 +      // see if they can be added to the integrableObjects:
535 +
536 +      molInfo.myIntegrableObjects.clear();
537        
603      if( !the_globals->haveEST() ){
604        sprintf( painCave.errMsg,
605                 "SimSetup Warning: using default value of 5%% of the "
606                 "electrostaticCutoffRadius for the "
607                 "electrostaticSkinThickness\n"
608                 );
609        painCave.isFatal = 0;
610        simError();
611        simnfo->est = 0.05 * simnfo->ecr;
612      } else {
613        simnfo->est        = the_globals->getEST();
614      }
615    }
616  }  
538  
539 < #ifdef IS_MPI
619 <  strcpy( checkPointMsg, "electrostatic parameters check out" );
620 <  MPIcheckPoint();
621 < #endif // is_mpi
539 >      for (j = 0; j < molInfo.nAtoms; j++){
540  
623 if( the_globals->haveInitialConfig() ){
624
625     InitializeFromFile* fileInit;
626 #ifdef IS_MPI // is_mpi
627     if( worldRank == 0 ){
628 #endif //is_mpi
629   fileInit = new InitializeFromFile( the_globals->getInitialConfig() );
541   #ifdef IS_MPI
542 <     }else fileInit = new InitializeFromFile( NULL );
542 >        slJ = molInfo.myAtoms[j]->getGlobalIndex();
543 > #else
544 >        slJ = j+atomOffset;
545   #endif
633   fileInit->read_xyz( simnfo ); // default velocities on
546  
547 <   delete fileInit;
636 < }
637 < else{
547 >        // if they aren't on the skip list, then they can be integrated
548  
549 < #ifdef IS_MPI
549 >        if (skipList.find(slJ) == skipList.end()) {
550 >          mySD = (StuntDouble *) molInfo.myAtoms[j];
551 >          info[k].integrableObjects.push_back(mySD);
552 >          molInfo.myIntegrableObjects.push_back(mySD);
553 >        }
554 >      }
555  
556 <  // no init from bass
642 <  
643 <  sprintf( painCave.errMsg,
644 <           "Cannot intialize a parallel simulation without an initial configuration file.\n" );
645 <  painCave.isFatal;
646 <  simError();
647 <  
648 < #else
556 >      // all rigid bodies are integrated:
557  
558 <  initFromBass();
558 >      for (j = 0; j < molInfo.nRigidBodies; j++) {
559 >        mySD = (StuntDouble *) molInfo.myRigidBodies[j];
560 >        info[k].integrableObjects.push_back(mySD);      
561 >        molInfo.myIntegrableObjects.push_back(mySD);
562 >      }
563 >    
564 >      
565 >      // send the arrays off to the forceField for init.
566 >      
567 >      the_ff->initializeAtoms(molInfo.nAtoms, molInfo.myAtoms);
568 >      the_ff->initializeBonds(molInfo.nBonds, molInfo.myBonds, theBonds);
569 >      the_ff->initializeBends(molInfo.nBends, molInfo.myBends, theBends);
570 >      the_ff->initializeTorsions(molInfo.nTorsions, molInfo.myTorsions,
571 >                                 theTorsions);
572  
573 +      info[k].molecules[i].initialize(molInfo);
574  
653 #endif
654 }
575  
576 +      atomOffset += molInfo.nAtoms;
577 +      delete[] theBonds;
578 +      delete[] theBends;
579 +      delete[] theTorsions;
580 +    }    
581 +  }
582 +
583   #ifdef IS_MPI
584 <  strcpy( checkPointMsg, "Successfully read in the initial configuration" );
584 >  sprintf(checkPointMsg, "all molecules initialized succesfully");
585    MPIcheckPoint();
586   #endif // is_mpi
587  
588 + }
589  
590 <  
591 <
592 <  
590 > void SimSetup::initFromBass(void){
591 >  int i, j, k;
592 >  int n_cells;
593 >  double cellx, celly, cellz;
594 >  double temp1, temp2, temp3;
595 >  int n_per_extra;
596 >  int n_extra;
597 >  int have_extra, done;
598  
599 <  
600 < #ifdef IS_MPI
601 <  if( worldRank == 0 ){
602 < #endif // is_mpi
603 <    
604 <    if( the_globals->haveFinalConfig() ){
605 <      strcpy( simnfo->finalName, the_globals->getFinalConfig() );
599 >  double vel[3];
600 >  vel[0] = 0.0;
601 >  vel[1] = 0.0;
602 >  vel[2] = 0.0;
603 >
604 >  temp1 = (double) tot_nmol / 4.0;
605 >  temp2 = pow(temp1, (1.0 / 3.0));
606 >  temp3 = ceil(temp2);
607 >
608 >  have_extra = 0;
609 >  if (temp2 < temp3){
610 >    // we have a non-complete lattice
611 >    have_extra = 1;
612 >
613 >    n_cells = (int) temp3 - 1;
614 >    cellx = info[0].boxL[0] / temp3;
615 >    celly = info[0].boxL[1] / temp3;
616 >    cellz = info[0].boxL[2] / temp3;
617 >    n_extra = tot_nmol - (4 * n_cells * n_cells * n_cells);
618 >    temp1 = ((double) n_extra) / (pow(temp3, 3.0) - pow(n_cells, 3.0));
619 >    n_per_extra = (int) ceil(temp1);
620 >
621 >    if (n_per_extra > 4){
622 >      sprintf(painCave.errMsg,
623 >              "SimSetup error. There has been an error in constructing"
624 >              " the non-complete lattice.\n");
625 >      painCave.isFatal = 1;
626 >      simError();
627      }
674    else{
675      strcpy( simnfo->finalName, inFileName );
676      char* endTest;
677      int nameLength = strlen( simnfo->finalName );
678      endTest = &(simnfo->finalName[nameLength - 5]);
679      if( !strcmp( endTest, ".bass" ) ){
680        strcpy( endTest, ".eor" );
681      }
682      else if( !strcmp( endTest, ".BASS" ) ){
683        strcpy( endTest, ".eor" );
684      }
685      else{
686        endTest = &(simnfo->finalName[nameLength - 4]);
687        if( !strcmp( endTest, ".bss" ) ){
688          strcpy( endTest, ".eor" );
689        }
690        else if( !strcmp( endTest, ".mdl" ) ){
691          strcpy( endTest, ".eor" );
692        }
693        else{
694          strcat( simnfo->finalName, ".eor" );
695        }
696      }
697    }
698    
699    // make the sample and status out names
700    
701    strcpy( simnfo->sampleName, inFileName );
702    char* endTest;
703    int nameLength = strlen( simnfo->sampleName );
704    endTest = &(simnfo->sampleName[nameLength - 5]);
705    if( !strcmp( endTest, ".bass" ) ){
706      strcpy( endTest, ".dump" );
707    }
708    else if( !strcmp( endTest, ".BASS" ) ){
709      strcpy( endTest, ".dump" );
710    }
711    else{
712      endTest = &(simnfo->sampleName[nameLength - 4]);
713      if( !strcmp( endTest, ".bss" ) ){
714        strcpy( endTest, ".dump" );
715      }
716      else if( !strcmp( endTest, ".mdl" ) ){
717        strcpy( endTest, ".dump" );
718      }
719      else{
720        strcat( simnfo->sampleName, ".dump" );
721      }
722    }
723    
724    strcpy( simnfo->statusName, inFileName );
725    nameLength = strlen( simnfo->statusName );
726    endTest = &(simnfo->statusName[nameLength - 5]);
727    if( !strcmp( endTest, ".bass" ) ){
728      strcpy( endTest, ".stat" );
729    }
730    else if( !strcmp( endTest, ".BASS" ) ){
731      strcpy( endTest, ".stat" );
732    }
733    else{
734      endTest = &(simnfo->statusName[nameLength - 4]);
735      if( !strcmp( endTest, ".bss" ) ){
736        strcpy( endTest, ".stat" );
737      }
738      else if( !strcmp( endTest, ".mdl" ) ){
739        strcpy( endTest, ".stat" );
740      }
741      else{
742        strcat( simnfo->statusName, ".stat" );
743      }
744    }
745    
746 #ifdef IS_MPI
628    }
748 #endif // is_mpi
749  
750  // set the status, sample, and themal kick times
751  
752  if( the_globals->haveSampleTime() ){
753    simnfo->sampleTime = the_globals->getSampleTime();
754    simnfo->statusTime = simnfo->sampleTime;
755    simnfo->thermalTime = simnfo->sampleTime;
756  }
629    else{
630 <    simnfo->sampleTime = the_globals->getRunTime();
631 <    simnfo->statusTime = simnfo->sampleTime;
632 <    simnfo->thermalTime = simnfo->sampleTime;
630 >    n_cells = (int) temp3;
631 >    cellx = info[0].boxL[0] / temp3;
632 >    celly = info[0].boxL[1] / temp3;
633 >    cellz = info[0].boxL[2] / temp3;
634    }
635  
636 <  if( the_globals->haveStatusTime() ){
637 <    simnfo->statusTime = the_globals->getStatusTime();
638 <  }
636 >  current_mol = 0;
637 >  current_comp_mol = 0;
638 >  current_comp = 0;
639 >  current_atom_ndx = 0;
640  
641 <  if( the_globals->haveThermalTime() ){
642 <    simnfo->thermalTime = the_globals->getThermalTime();
641 >  for (i = 0; i < n_cells ; i++){
642 >    for (j = 0; j < n_cells; j++){
643 >      for (k = 0; k < n_cells; k++){
644 >        makeElement(i * cellx, j * celly, k * cellz);
645 >
646 >        makeElement(i * cellx + 0.5 * cellx, j * celly + 0.5 * celly, k * cellz);
647 >
648 >        makeElement(i * cellx, j * celly + 0.5 * celly, k * cellz + 0.5 * cellz);
649 >
650 >        makeElement(i * cellx + 0.5 * cellx, j * celly, k * cellz + 0.5 * cellz);
651 >      }
652 >    }
653    }
654  
655 <  // check for the temperature set flag
655 >  if (have_extra){
656 >    done = 0;
657  
658 <  if( the_globals->haveTempSet() ) simnfo->setTemp = the_globals->getTempSet();
658 >    int start_ndx;
659 >    for (i = 0; i < (n_cells + 1) && !done; i++){
660 >      for (j = 0; j < (n_cells + 1) && !done; j++){
661 >        if (i < n_cells){
662 >          if (j < n_cells){
663 >            start_ndx = n_cells;
664 >          }
665 >          else
666 >            start_ndx = 0;
667 >        }
668 >        else
669 >          start_ndx = 0;
670  
671 +        for (k = start_ndx; k < (n_cells + 1) && !done; k++){
672 +          makeElement(i * cellx, j * celly, k * cellz);
673 +          done = (current_mol >= tot_nmol);
674  
675 <  // make the integrator
676 <  
677 <  
678 <  NVT* myNVT = NULL;
679 <  switch( ensembleCase ){
675 >          if (!done && n_per_extra > 1){
676 >            makeElement(i * cellx + 0.5 * cellx, j * celly + 0.5 * celly,
677 >                        k * cellz);
678 >            done = (current_mol >= tot_nmol);
679 >          }
680  
681 <  case NVE_ENS:
682 <    new NVE( simnfo, the_ff );
683 <    break;
681 >          if (!done && n_per_extra > 2){
682 >            makeElement(i * cellx, j * celly + 0.5 * celly,
683 >                        k * cellz + 0.5 * cellz);
684 >            done = (current_mol >= tot_nmol);
685 >          }
686  
687 <  case NVT_ENS:
688 <    myNVT = new NVT( simnfo, the_ff );
689 <    myNVT->setTargetTemp(the_globals->getTargetTemp());
687 >          if (!done && n_per_extra > 3){
688 >            makeElement(i * cellx + 0.5 * cellx, j * celly,
689 >                        k * cellz + 0.5 * cellz);
690 >            done = (current_mol >= tot_nmol);
691 >          }
692 >        }
693 >      }
694 >    }
695 >  }
696  
697 <    if (the_globals->haveTauThermostat())
698 <      myNVT->setTauThermostat(the_globals->getTauThermostat());
699 < //     else if (the_globals->haveQmass())
700 < //       myNVT->setQmass(the_globals->getQmass());
701 <    else {
702 <      sprintf( painCave.errMsg,
703 <               "SimSetup error: If you use the NVT\n"
704 <               "    ensemble, you must set either tauThermostat or qMass.\n"
705 <               "    Neither of these was found in the BASS file.\n");
697 >  for (i = 0; i < info[0].n_atoms; i++){
698 >    info[0].atoms[i]->setVel(vel);
699 >  }
700 > }
701 >
702 > void SimSetup::makeElement(double x, double y, double z){
703 >  int k;
704 >  AtomStamp* current_atom;
705 >  DirectionalAtom* dAtom;
706 >  double rotMat[3][3];
707 >  double pos[3];
708 >
709 >  for (k = 0; k < comp_stamps[current_comp]->getNAtoms(); k++){
710 >    current_atom = comp_stamps[current_comp]->getAtom(k);
711 >    if (!current_atom->havePosition()){
712 >      sprintf(painCave.errMsg,
713 >              "SimSetup:initFromBass error.\n"
714 >              "\tComponent %s, atom %s does not have a position specified.\n"
715 >              "\tThe initialization routine is unable to give a start"
716 >              " position.\n",
717 >              comp_stamps[current_comp]->getID(), current_atom->getType());
718        painCave.isFatal = 1;
719        simError();
720      }
802    break;
721  
722 <  default:
723 <    sprintf( painCave.errMsg,
724 <             "SimSetup Error. Unrecognized ensemble in case statement.\n");
807 <    painCave.isFatal = 1;
808 <    simError();
809 <  }
722 >    pos[0] = x + current_atom->getPosX();
723 >    pos[1] = y + current_atom->getPosY();
724 >    pos[2] = z + current_atom->getPosZ();
725  
726 +    info[0].atoms[current_atom_ndx]->setPos(pos);
727  
728 < #ifdef IS_MPI
729 <  mpiSim->mpiRefresh();
814 < #endif
728 >    if (info[0].atoms[current_atom_ndx]->isDirectional()){
729 >      dAtom = (DirectionalAtom *) info[0].atoms[current_atom_ndx];
730  
731 <  // initialize the Fortran
731 >      rotMat[0][0] = 1.0;
732 >      rotMat[0][1] = 0.0;
733 >      rotMat[0][2] = 0.0;
734  
735 +      rotMat[1][0] = 0.0;
736 +      rotMat[1][1] = 1.0;
737 +      rotMat[1][2] = 0.0;
738  
739 <  simnfo->refreshSim();
740 <  
741 <  if( !strcmp( simnfo->mixingRule, "standard") ){
742 <    the_ff->initForceField( LB_MIXING_RULE );
739 >      rotMat[2][0] = 0.0;
740 >      rotMat[2][1] = 0.0;
741 >      rotMat[2][2] = 1.0;
742 >
743 >      dAtom->setA(rotMat);
744 >    }
745 >
746 >    current_atom_ndx++;
747    }
824  else if( !strcmp( simnfo->mixingRule, "explicit") ){
825    the_ff->initForceField( EXPLICIT_MIXING_RULE );
826  }
827  else{
828    sprintf( painCave.errMsg,
829             "SimSetup Error: unknown mixing rule -> \"%s\"\n",
830             simnfo->mixingRule );
831    painCave.isFatal = 1;
832    simError();
833  }
748  
749 +  current_mol++;
750 +  current_comp_mol++;
751  
752 < #ifdef IS_MPI
753 <  strcpy( checkPointMsg,
754 <          "Successfully intialized the mixingRule for Fortran." );
755 <  MPIcheckPoint();
840 < #endif // is_mpi
752 >  if (current_comp_mol >= components_nmol[current_comp]){
753 >    current_comp_mol = 0;
754 >    current_comp++;
755 >  }
756   }
757  
758  
759 < void SimSetup::makeMolecules( void ){
759 > void SimSetup::gatherInfo(void){
760 >  int i;
761  
762 <  int i, j, exI, exJ, tempEx, stampID, atomOffset, excludeOffset;
763 <  molInit info;
848 <  DirectionalAtom* dAtom;
849 <  LinkedAssign* extras;
850 <  LinkedAssign* current_extra;
851 <  AtomStamp* currentAtom;
852 <  BondStamp* currentBond;
853 <  BendStamp* currentBend;
854 <  TorsionStamp* currentTorsion;
762 >  ensembleCase = -1;
763 >  ffCase = -1;
764  
765 <  bond_pair* theBonds;
857 <  bend_set* theBends;
858 <  torsion_set* theTorsions;
859 <
860 <  
861 <  //init the forceField paramters
765 >  // set the easy ones first
766  
767 <  the_ff->readParams();
767 >  for (i = 0; i < nInfo; i++){
768 >    info[i].target_temp = globals->getTargetTemp();
769 >    info[i].dt = globals->getDt();
770 >    info[i].run_time = globals->getRunTime();
771 >  }
772 >  n_components = globals->getNComponents();
773  
865  
866  // init the atoms
774  
775 <  double ux, uy, uz, u, uSqr;
869 <  
870 <  atomOffset = 0;
871 <  excludeOffset = 0;
872 <  for(i=0; i<simnfo->n_mol; i++){
873 <    
874 <    stampID = the_molecules[i].getStampID();
775 >  // get the forceField
776  
777 <    info.nAtoms    = comp_stamps[stampID]->getNAtoms();
877 <    info.nBonds    = comp_stamps[stampID]->getNBonds();
878 <    info.nBends    = comp_stamps[stampID]->getNBends();
879 <    info.nTorsions = comp_stamps[stampID]->getNTorsions();
880 <    info.nExcludes = info.nBonds + info.nBends + info.nTorsions;
777 >  strcpy(force_field, globals->getForceField());
778  
779 <    info.myAtoms = &the_atoms[atomOffset];
780 <    info.myExcludes = &the_excludes[excludeOffset];
781 <    info.myBonds = new Bond*[info.nBonds];
782 <    info.myBends = new Bend*[info.nBends];
783 <    info.myTorsions = new Torsion*[info.nTorsions];
779 >  if (!strcasecmp(force_field, "DUFF")){
780 >    ffCase = FF_DUFF;
781 >  }
782 >  else if (!strcasecmp(force_field, "LJ")){
783 >    ffCase = FF_LJ;
784 >  }
785 >  else if (!strcasecmp(force_field, "EAM")){
786 >    ffCase = FF_EAM;
787 >  }
788 >  else if (!strcasecmp(force_field, "WATER")){
789 >    ffCase = FF_H2O;
790 >  }
791 >  else{
792 >    sprintf(painCave.errMsg, "SimSetup Error. Unrecognized force field -> %s\n",
793 >            force_field);
794 >         painCave.isFatal = 1;
795 >         simError();
796 >  }
797  
798 <    theBonds = new bond_pair[info.nBonds];
799 <    theBends = new bend_set[info.nBends];
800 <    theTorsions = new torsion_set[info.nTorsions];
798 >    // get the ensemble
799 >
800 >  strcpy(ensemble, globals->getEnsemble());
801 >
802 >  if (!strcasecmp(ensemble, "NVE")){
803 >    ensembleCase = NVE_ENS;
804 >  }
805 >  else if (!strcasecmp(ensemble, "NVT")){
806 >    ensembleCase = NVT_ENS;
807 >  }
808 >  else if (!strcasecmp(ensemble, "NPTi") || !strcasecmp(ensemble, "NPT")){
809 >    ensembleCase = NPTi_ENS;
810 >  }
811 >  else if (!strcasecmp(ensemble, "NPTf")){
812 >    ensembleCase = NPTf_ENS;
813 >  }
814 >  else if (!strcasecmp(ensemble, "NPTxyz")){
815 >    ensembleCase = NPTxyz_ENS;
816 >  }
817 >  else{
818 >    sprintf(painCave.errMsg,
819 >            "SimSetup Warning. Unrecognized Ensemble -> %s \n"
820 >            "\treverting to NVE for this simulation.\n",
821 >            ensemble);
822 >         painCave.isFatal = 0;
823 >         simError();
824 >         strcpy(ensemble, "NVE");
825 >         ensembleCase = NVE_ENS;
826 >  }  
827 >
828 >  for (i = 0; i < nInfo; i++){
829 >    strcpy(info[i].ensemble, ensemble);
830 >
831 >    // get the mixing rule
832 >
833 >    strcpy(info[i].mixingRule, globals->getMixingRule());
834 >    info[i].usePBC = globals->getPBC();
835 >  }
836 >
837 >  // get the components and calculate the tot_nMol and indvidual n_mol
838 >
839 >  the_components = globals->getComponents();
840 >  components_nmol = new int[n_components];
841 >
842 >
843 >  if (!globals->haveNMol()){
844 >    // we don't have the total number of molecules, so we assume it is
845 >    // given in each component
846 >
847 >    tot_nmol = 0;
848 >    for (i = 0; i < n_components; i++){
849 >      if (!the_components[i]->haveNMol()){
850 >        // we have a problem
851 >        sprintf(painCave.errMsg,
852 >                "SimSetup Error. No global NMol or component NMol given.\n"
853 >                "\tCannot calculate the number of atoms.\n");
854 >        painCave.isFatal = 1;
855 >        simError();
856 >      }
857 >
858 >      tot_nmol += the_components[i]->getNMol();
859 >      components_nmol[i] = the_components[i]->getNMol();
860 >    }
861 >  }
862 >  else{
863 >    sprintf(painCave.errMsg,
864 >            "SimSetup error.\n"
865 >            "\tSorry, the ability to specify total"
866 >            " nMols and then give molfractions in the components\n"
867 >            "\tis not currently supported."
868 >            " Please give nMol in the components.\n");
869 >    painCave.isFatal = 1;
870 >    simError();
871 >  }
872 >
873 >  //check whether sample time, status time, thermal time and reset time are divisble by dt
874 >  if (globals->haveSampleTime() && !isDivisible(globals->getSampleTime(), globals->getDt())){
875 >    sprintf(painCave.errMsg,
876 >            "Sample time is not divisible by dt.\n"
877 >            "\tThis will result in samples that are not uniformly\n"
878 >            "\tdistributed in time.  If this is a problem, change\n"
879 >            "\tyour sampleTime variable.\n");
880 >    painCave.isFatal = 0;
881 >    simError();    
882 >  }
883 >
884 >  if (globals->haveStatusTime() && !isDivisible(globals->getStatusTime(), globals->getDt())){
885 >    sprintf(painCave.errMsg,
886 >            "Status time is not divisible by dt.\n"
887 >            "\tThis will result in status reports that are not uniformly\n"
888 >            "\tdistributed in time.  If this is a problem, change \n"
889 >            "\tyour statusTime variable.\n");
890 >    painCave.isFatal = 0;
891 >    simError();    
892 >  }
893 >
894 >  if (globals->haveThermalTime() && !isDivisible(globals->getThermalTime(), globals->getDt())){
895 >    sprintf(painCave.errMsg,
896 >            "Thermal time is not divisible by dt.\n"
897 >            "\tThis will result in thermalizations that are not uniformly\n"
898 >            "\tdistributed in time.  If this is a problem, change \n"
899 >            "\tyour thermalTime variable.\n");
900 >    painCave.isFatal = 0;
901 >    simError();    
902 >  }  
903 >
904 >  if (globals->haveResetTime() && !isDivisible(globals->getResetTime(), globals->getDt())){
905 >    sprintf(painCave.errMsg,
906 >            "Reset time is not divisible by dt.\n"
907 >            "\tThis will result in integrator resets that are not uniformly\n"
908 >            "\tdistributed in time.  If this is a problem, change\n"
909 >            "\tyour resetTime variable.\n");
910 >    painCave.isFatal = 0;
911 >    simError();    
912 >  }
913 >
914 >  // set the status, sample, and thermal kick times
915 >
916 >  for (i = 0; i < nInfo; i++){
917 >    if (globals->haveSampleTime()){
918 >      info[i].sampleTime = globals->getSampleTime();
919 >      info[i].statusTime = info[i].sampleTime;
920 >    }
921 >    else{
922 >      info[i].sampleTime = globals->getRunTime();
923 >      info[i].statusTime = info[i].sampleTime;
924 >    }
925 >
926 >    if (globals->haveStatusTime()){
927 >      info[i].statusTime = globals->getStatusTime();
928 >    }
929 >
930 >    if (globals->haveThermalTime()){
931 >      info[i].thermalTime = globals->getThermalTime();
932 >    } else {
933 >      info[i].thermalTime = globals->getRunTime();
934 >    }
935 >
936 >    info[i].resetIntegrator = 0;
937 >    if( globals->haveResetTime() ){
938 >      info[i].resetTime = globals->getResetTime();
939 >      info[i].resetIntegrator = 1;
940 >    }
941 >
942 >    // check for the temperature set flag
943      
944 <    // make the Atoms
944 >    if (globals->haveTempSet())
945 >      info[i].setTemp = globals->getTempSet();
946 >
947 >    // check for the extended State init
948 >
949 >    info[i].useInitXSstate = globals->getUseInitXSstate();
950 >    info[i].orthoTolerance = globals->getOrthoBoxTolerance();
951      
952 <    for(j=0; j<info.nAtoms; j++){
953 <      
954 <      currentAtom = comp_stamps[stampID]->getAtom( j );
955 <      if( currentAtom->haveOrientation() ){
956 <        
957 <        dAtom = new DirectionalAtom(j + atomOffset);
958 <        simnfo->n_oriented++;
959 <        info.myAtoms[j] = dAtom;
960 <        
961 <        ux = currentAtom->getOrntX();
962 <        uy = currentAtom->getOrntY();
963 <        uz = currentAtom->getOrntZ();
964 <        
965 <        uSqr = (ux * ux) + (uy * uy) + (uz * uz);
966 <        
967 <        u = sqrt( uSqr );
968 <        ux = ux / u;
969 <        uy = uy / u;
970 <        uz = uz / u;
971 <        
972 <        dAtom->setSUx( ux );
973 <        dAtom->setSUy( uy );
974 <        dAtom->setSUz( uz );
975 <      }
976 <      else{
977 <        info.myAtoms[j] = new GeneralAtom(j + atomOffset);
978 <      }
921 <      info.myAtoms[j]->setType( currentAtom->getType() );
952 >  }
953 >  
954 >  //setup seed for random number generator
955 >  int seedValue;
956 >
957 >  if (globals->haveSeed()){
958 >    seedValue = globals->getSeed();
959 >
960 >    if(seedValue / 1E9 == 0){
961 >      sprintf(painCave.errMsg,
962 >              "Seed for sprng library should contain at least 9 digits\n"
963 >              "OOPSE will generate a seed for user\n");
964 >      painCave.isFatal = 0;
965 >      simError();
966 >
967 >      //using seed generated by system instead of invalid seed set by user
968 > #ifndef IS_MPI
969 >      seedValue = make_sprng_seed();
970 > #else
971 >      if (worldRank == 0){
972 >        seedValue = make_sprng_seed();
973 >      }
974 >      MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD);  
975 > #endif      
976 >    }
977 >  }//end of if branch of globals->haveSeed()
978 >  else{
979      
980 + #ifndef IS_MPI
981 +    seedValue = make_sprng_seed();
982 + #else
983 +    if (worldRank == 0){
984 +      seedValue = make_sprng_seed();
985 +    }
986 +    MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD);  
987 + #endif
988 +  }//end of globals->haveSeed()
989 +
990 +  for (int i = 0; i < nInfo; i++){
991 +    info[i].setSeed(seedValue);
992 +  }
993 +  
994   #ifdef IS_MPI
995 <      
996 <      info.myAtoms[j]->setGlobalIndex( globalIndex[j+atomOffset] );
926 <      
995 >  strcpy(checkPointMsg, "Successfully gathered all information from Bass\n");
996 >  MPIcheckPoint();
997   #endif // is_mpi
998 <    }
929 <    
930 <    // make the bonds
931 <    for(j=0; j<info.nBonds; j++){
932 <      
933 <      currentBond = comp_stamps[stampID]->getBond( j );
934 <      theBonds[j].a = currentBond->getA() + atomOffset;
935 <      theBonds[j].b = currentBond->getB() + atomOffset;
998 > }
999  
937      exI = theBonds[j].a;
938      exJ = theBonds[j].b;
1000  
1001 <      // exclude_I must always be the smaller of the pair
1002 <      if( exI > exJ ){
1003 <        tempEx = exI;
1004 <        exI = exJ;
1005 <        exJ = tempEx;
1006 <      }
1001 > void SimSetup::finalInfoCheck(void){
1002 >  int index;
1003 >  int usesDipoles;
1004 >  int usesCharges;
1005 >  int i;
1006 >
1007 >  for (i = 0; i < nInfo; i++){
1008 >    // check electrostatic parameters
1009 >
1010 >    index = 0;
1011 >    usesDipoles = 0;
1012 >    while ((index < info[i].n_atoms) && !usesDipoles){
1013 >      usesDipoles = (info[i].atoms[index])->hasDipole();
1014 >      index++;
1015 >    }
1016 >    index = 0;
1017 >    usesCharges = 0;
1018 >    while ((index < info[i].n_atoms) && !usesCharges){
1019 >      usesCharges= (info[i].atoms[index])->hasCharge();
1020 >      index++;
1021 >    }
1022   #ifdef IS_MPI
1023 <      tempEx = exI;
1024 <      exI = the_atoms[tempEx]->getGlobalIndex() + 1;
1025 <      tempEx = exJ;
1026 <      exJ = the_atoms[tempEx]->getGlobalIndex() + 1;
1023 >    int myUse = usesDipoles;
1024 >    MPI_Allreduce(&myUse, &usesDipoles, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
1025 > #endif //is_mpi
1026 >
1027 >    double theRcut, theRsw;
1028 >
1029 >    if (globals->haveRcut()) {
1030 >      theRcut = globals->getRcut();
1031 >
1032 >      if (globals->haveRsw())
1033 >        theRsw = globals->getRsw();
1034 >      else
1035 >        theRsw = theRcut;
1036        
1037 <      the_excludes[j+excludeOffset]->setPair( exI, exJ );
953 < #else  // isn't MPI
1037 >      info[i].setDefaultRcut(theRcut, theRsw);
1038  
1039 <      the_excludes[j+excludeOffset]->setPair( (exI+1), (exJ+1) );
1040 < #endif  //is_mpi
1039 >    } else {
1040 >      
1041 >      the_ff->calcRcut();
1042 >      theRcut = info[i].getRcut();
1043 >
1044 >      if (globals->haveRsw())
1045 >        theRsw = globals->getRsw();
1046 >      else
1047 >        theRsw = theRcut;
1048 >      
1049 >      info[i].setDefaultRcut(theRcut, theRsw);
1050      }
958    excludeOffset += info.nBonds;
1051  
1052 <    //make the bends
1053 <    for(j=0; j<info.nBends; j++){
1052 >    if (globals->getUseRF()){
1053 >      info[i].useReactionField = 1;
1054        
1055 <      currentBend = comp_stamps[stampID]->getBend( j );
1056 <      theBends[j].a = currentBend->getA() + atomOffset;
1057 <      theBends[j].b = currentBend->getB() + atomOffset;
1058 <      theBends[j].c = currentBend->getC() + atomOffset;
1059 <          
1060 <      if( currentBend->haveExtras() ){
1061 <            
1062 <        extras = currentBend->getExtras();
971 <        current_extra = extras;
972 <            
973 <        while( current_extra != NULL ){
974 <          if( !strcmp( current_extra->getlhs(), "ghostVectorSource" )){
975 <                
976 <            switch( current_extra->getType() ){
977 <              
978 <            case 0:
979 <              theBends[j].ghost =
980 <                current_extra->getInt() + atomOffset;
981 <              theBends[j].isGhost = 1;
982 <              break;
983 <                  
984 <            case 1:
985 <              theBends[j].ghost =
986 <                (int)current_extra->getDouble() + atomOffset;
987 <              theBends[j].isGhost = 1;
988 <              break;
989 <              
990 <            default:
991 <              sprintf( painCave.errMsg,
992 <                       "SimSetup Error: ghostVectorSource was neither a "
993 <                       "double nor an int.\n"
994 <                       "-->Bend[%d] in %s\n",
995 <                       j, comp_stamps[stampID]->getID() );
996 <              painCave.isFatal = 1;
997 <              simError();
998 <            }
999 <          }
1000 <          
1001 <          else{
1002 <            
1003 <            sprintf( painCave.errMsg,
1004 <                     "SimSetup Error: unhandled bend assignment:\n"
1005 <                     "    -->%s in Bend[%d] in %s\n",
1006 <                     current_extra->getlhs(),
1007 <                     j, comp_stamps[stampID]->getID() );
1008 <            painCave.isFatal = 1;
1009 <            simError();
1010 <          }
1011 <          
1012 <          current_extra = current_extra->getNext();
1013 <        }
1055 >      if (!globals->haveRcut()){
1056 >        sprintf(painCave.errMsg,
1057 >                "SimSetup Warning: No value was set for the cutoffRadius.\n"
1058 >                "\tOOPSE will use a default value of 15.0 angstroms"
1059 >                "\tfor the cutoffRadius.\n");
1060 >        painCave.isFatal = 0;
1061 >        simError();
1062 >        theRcut = 15.0;
1063        }
1064 <          
1065 <      if( !theBends[j].isGhost ){
1017 <            
1018 <        exI = theBends[j].a;
1019 <        exJ = theBends[j].c;
1064 >      else{
1065 >        theRcut = globals->getRcut();
1066        }
1067 +
1068 +      if (!globals->haveRsw()){
1069 +        sprintf(painCave.errMsg,
1070 +                "SimSetup Warning: No value was set for switchingRadius.\n"
1071 +                "\tOOPSE will use a default value of\n"
1072 +                "\t0.95 * cutoffRadius for the switchingRadius\n");
1073 +        painCave.isFatal = 0;
1074 +        simError();
1075 +        theRsw = 0.95 * theRcut;
1076 +      }
1077        else{
1078 <        
1023 <        exI = theBends[j].a;
1024 <        exJ = theBends[j].b;
1078 >        theRsw = globals->getRsw();
1079        }
1080 <      
1081 <      // exclude_I must always be the smaller of the pair
1082 <      if( exI > exJ ){
1083 <        tempEx = exI;
1084 <        exI = exJ;
1085 <        exJ = tempEx;
1080 >
1081 >      info[i].setDefaultRcut(theRcut, theRsw);
1082 >
1083 >      if (!globals->haveDielectric()){
1084 >        sprintf(painCave.errMsg,
1085 >                "SimSetup Error: No Dielectric constant was set.\n"
1086 >                "\tYou are trying to use Reaction Field without"
1087 >                "\tsetting a dielectric constant!\n");
1088 >        painCave.isFatal = 1;
1089 >        simError();
1090        }
1091 < #ifdef IS_MPI
1034 <      tempEx = exI;
1035 <      exI = the_atoms[tempEx]->getGlobalIndex() + 1;
1036 <      tempEx = exJ;
1037 <      exJ = the_atoms[tempEx]->getGlobalIndex() + 1;
1038 <      
1039 <      the_excludes[j+excludeOffset]->setPair( exI, exJ );
1040 < #else  // isn't MPI
1041 <      the_excludes[j+excludeOffset]->setPair( (exI+1), (exJ+1) );
1042 < #endif  //is_mpi
1091 >      info[i].dielectric = globals->getDielectric();
1092      }
1093 <    excludeOffset += info.nBends;
1093 >    else{
1094 >      if (usesDipoles || usesCharges){
1095  
1096 <    for(j=0; j<info.nTorsions; j++){
1097 <      
1098 <      currentTorsion = comp_stamps[stampID]->getTorsion( j );
1099 <      theTorsions[j].a = currentTorsion->getA() + atomOffset;
1100 <      theTorsions[j].b = currentTorsion->getB() + atomOffset;
1101 <      theTorsions[j].c = currentTorsion->getC() + atomOffset;
1102 <      theTorsions[j].d = currentTorsion->getD() + atomOffset;
1103 <      
1054 <      exI = theTorsions[j].a;
1055 <      exJ = theTorsions[j].d;
1056 <
1057 <      // exclude_I must always be the smaller of the pair
1058 <      if( exI > exJ ){
1059 <        tempEx = exI;
1060 <        exI = exJ;
1061 <        exJ = tempEx;
1096 >        if (!globals->haveRcut()){
1097 >          sprintf(painCave.errMsg,
1098 >                  "SimSetup Warning: No value was set for the cutoffRadius.\n"
1099 >                  "\tOOPSE will use a default value of 15.0 angstroms"
1100 >                  "\tfor the cutoffRadius.\n");
1101 >          painCave.isFatal = 0;
1102 >          simError();
1103 >          theRcut = 15.0;
1104        }
1105 < #ifdef IS_MPI
1106 <      tempEx = exI;
1107 <      exI = the_atoms[tempEx]->getGlobalIndex() + 1;
1108 <      tempEx = exJ;
1109 <      exJ = the_atoms[tempEx]->getGlobalIndex() + 1;
1110 <      
1111 <      the_excludes[j+excludeOffset]->setPair( exI, exJ );
1112 < #else  // isn't MPI
1113 <      the_excludes[j+excludeOffset]->setPair( (exI+1), (exJ+1) );
1114 < #endif  //is_mpi
1105 >        else{
1106 >          theRcut = globals->getRcut();
1107 >        }
1108 >        
1109 >        if (!globals->haveRsw()){
1110 >          sprintf(painCave.errMsg,
1111 >                  "SimSetup Warning: No value was set for switchingRadius.\n"
1112 >                  "\tOOPSE will use a default value of\n"
1113 >                  "\t0.95 * cutoffRadius for the switchingRadius\n");
1114 >          painCave.isFatal = 0;
1115 >          simError();
1116 >          theRsw = 0.95 * theRcut;
1117 >        }
1118 >        else{
1119 >          theRsw = globals->getRsw();
1120 >        }
1121 >        
1122 >        info[i].setDefaultRcut(theRcut, theRsw);
1123 >        
1124 >      }
1125      }
1126 <    excludeOffset += info.nTorsions;
1126 >  }
1127 > #ifdef IS_MPI
1128 >  strcpy(checkPointMsg, "post processing checks out");
1129 >  MPIcheckPoint();
1130 > #endif // is_mpi
1131  
1132 <    
1133 <    // send the arrays off to the forceField for init.
1132 >  // clean up the forcefield
1133 >  the_ff->cleanMe();
1134 > }
1135 >  
1136 > void SimSetup::initSystemCoords(void){
1137 >  int i;
1138  
1139 <    the_ff->initializeAtoms( info.nAtoms, info.myAtoms );
1080 <    the_ff->initializeBonds( info.nBonds, info.myBonds, theBonds );
1081 <    the_ff->initializeBends( info.nBends, info.myBends, theBends );
1082 <    the_ff->initializeTorsions( info.nTorsions, info.myTorsions, theTorsions );
1139 >  char* inName;
1140  
1141 +  (info[0].getConfiguration())->createArrays(info[0].n_atoms);
1142  
1143 <    the_molecules[i].initialize( info );
1143 >  for (i = 0; i < info[0].n_atoms; i++)
1144 >    info[0].atoms[i]->setCoords();
1145  
1146 +  if (globals->haveInitialConfig()){
1147 +    InitializeFromFile* fileInit;
1148 + #ifdef IS_MPI // is_mpi
1149 +    if (worldRank == 0){
1150 + #endif //is_mpi
1151 +      inName = globals->getInitialConfig();
1152 +      fileInit = new InitializeFromFile(inName);
1153 + #ifdef IS_MPI
1154 +    }
1155 +    else
1156 +      fileInit = new InitializeFromFile(NULL);
1157 + #endif
1158 +    fileInit->readInit(info); // default velocities on
1159  
1160 <    atomOffset += info.nAtoms;
1089 <    delete[] theBonds;
1090 <    delete[] theBends;
1091 <    delete[] theTorsions;
1160 >    delete fileInit;
1161    }
1162 +  else{
1163 +    
1164 +    // no init from bass
1165 +    
1166 +    sprintf(painCave.errMsg,
1167 +            "Cannot intialize a simulation without an initial configuration file.\n");
1168 +    painCave.isFatal = 1;;
1169 +    simError();
1170 +    
1171 +  }
1172  
1173   #ifdef IS_MPI
1174 <  sprintf( checkPointMsg, "all molecules initialized succesfully" );
1174 >  strcpy(checkPointMsg, "Successfully read in the initial configuration");
1175    MPIcheckPoint();
1176   #endif // is_mpi
1177 + }
1178  
1099  // clean up the forcefield
1100  the_ff->calcRcut();
1101  the_ff->cleanMe();
1179  
1180 < }
1180 > void SimSetup::makeOutNames(void){
1181 >  int k;
1182  
1105 void SimSetup::initFromBass( void ){
1183  
1184 <  int i, j, k;
1185 <  int n_cells;
1186 <  double cellx, celly, cellz;
1187 <  double temp1, temp2, temp3;
1111 <  int n_per_extra;
1112 <  int n_extra;
1113 <  int have_extra, done;
1184 >  for (k = 0; k < nInfo; k++){
1185 > #ifdef IS_MPI
1186 >    if (worldRank == 0){
1187 > #endif // is_mpi
1188  
1189 <  temp1 = (double)tot_nmol / 4.0;
1190 <  temp2 = pow( temp1, ( 1.0 / 3.0 ) );
1191 <  temp3 = ceil( temp2 );
1189 >      if (globals->haveFinalConfig()){
1190 >        strcpy(info[k].finalName, globals->getFinalConfig());
1191 >      }
1192 >      else{
1193 >        strcpy(info[k].finalName, inFileName);
1194 >        char* endTest;
1195 >        int nameLength = strlen(info[k].finalName);
1196 >        endTest = &(info[k].finalName[nameLength - 5]);
1197 >        if (!strcmp(endTest, ".bass")){
1198 >          strcpy(endTest, ".eor");
1199 >        }
1200 >        else if (!strcmp(endTest, ".BASS")){
1201 >          strcpy(endTest, ".eor");
1202 >        }
1203 >        else{
1204 >          endTest = &(info[k].finalName[nameLength - 4]);
1205 >          if (!strcmp(endTest, ".bss")){
1206 >            strcpy(endTest, ".eor");
1207 >          }
1208 >          else if (!strcmp(endTest, ".mdl")){
1209 >            strcpy(endTest, ".eor");
1210 >          }
1211 >          else{
1212 >            strcat(info[k].finalName, ".eor");
1213 >          }
1214 >        }
1215 >      }
1216  
1217 <  have_extra =0;
1120 <  if( temp2 < temp3 ){ // we have a non-complete lattice
1121 <    have_extra =1;
1217 >      // make the sample and status out names
1218  
1219 <    n_cells = (int)temp3 - 1;
1220 <    cellx = simnfo->boxLx / temp3;
1221 <    celly = simnfo->boxLy / temp3;
1222 <    cellz = simnfo->boxLz / temp3;
1223 <    n_extra = tot_nmol - ( 4 * n_cells * n_cells * n_cells );
1224 <    temp1 = ((double)n_extra) / ( pow( temp3, 3.0 ) - pow( n_cells, 3.0 ) );
1225 <    n_per_extra = (int)ceil( temp1 );
1219 >      strcpy(info[k].sampleName, inFileName);
1220 >      char* endTest;
1221 >      int nameLength = strlen(info[k].sampleName);
1222 >      endTest = &(info[k].sampleName[nameLength - 5]);
1223 >      if (!strcmp(endTest, ".bass")){
1224 >        strcpy(endTest, ".dump");
1225 >      }
1226 >      else if (!strcmp(endTest, ".BASS")){
1227 >        strcpy(endTest, ".dump");
1228 >      }
1229 >      else{
1230 >        endTest = &(info[k].sampleName[nameLength - 4]);
1231 >        if (!strcmp(endTest, ".bss")){
1232 >          strcpy(endTest, ".dump");
1233 >        }
1234 >        else if (!strcmp(endTest, ".mdl")){
1235 >          strcpy(endTest, ".dump");
1236 >        }
1237 >        else{
1238 >          strcat(info[k].sampleName, ".dump");
1239 >        }
1240 >      }
1241  
1242 <    if( n_per_extra > 4){
1243 <      sprintf( painCave.errMsg,
1244 <               "SimSetup error. There has been an error in constructing"
1245 <               " the non-complete lattice.\n" );
1246 <      painCave.isFatal = 1;
1247 <      simError();
1242 >      strcpy(info[k].statusName, inFileName);
1243 >      nameLength = strlen(info[k].statusName);
1244 >      endTest = &(info[k].statusName[nameLength - 5]);
1245 >      if (!strcmp(endTest, ".bass")){
1246 >        strcpy(endTest, ".stat");
1247 >      }
1248 >      else if (!strcmp(endTest, ".BASS")){
1249 >        strcpy(endTest, ".stat");
1250 >      }
1251 >      else{
1252 >        endTest = &(info[k].statusName[nameLength - 4]);
1253 >        if (!strcmp(endTest, ".bss")){
1254 >          strcpy(endTest, ".stat");
1255 >        }
1256 >        else if (!strcmp(endTest, ".mdl")){
1257 >          strcpy(endTest, ".stat");
1258 >        }
1259 >        else{
1260 >          strcat(info[k].statusName, ".stat");
1261 >        }
1262 >      }
1263 >
1264 > #ifdef IS_MPI
1265 >
1266      }
1267 + #endif // is_mpi
1268    }
1269 <  else{
1270 <    n_cells = (int)temp3;
1271 <    cellx = simnfo->boxLx / temp3;
1272 <    celly = simnfo->boxLy / temp3;
1273 <    cellz = simnfo->boxLz / temp3;
1269 > }
1270 >
1271 >
1272 > void SimSetup::sysObjectsCreation(void){
1273 >  int i, k;
1274 >
1275 >  // create the forceField
1276 >
1277 >  createFF();
1278 >
1279 >  // extract componentList
1280 >
1281 >  compList();
1282 >
1283 >  // calc the number of atoms, bond, bends, and torsions
1284 >
1285 >  calcSysValues();
1286 >
1287 > #ifdef IS_MPI
1288 >  // divide the molecules among the processors
1289 >
1290 >  mpiMolDivide();
1291 > #endif //is_mpi
1292 >
1293 >  // create the atom and SRI arrays. Also initialize Molecule Stamp ID's
1294 >
1295 >  makeSysArrays();
1296 >
1297 >  // make and initialize the molecules (all but atomic coordinates)
1298 >
1299 >  makeMolecules();
1300 >
1301 >  for (k = 0; k < nInfo; k++){
1302 >    info[k].identArray = new int[info[k].n_atoms];
1303 >    for (i = 0; i < info[k].n_atoms; i++){
1304 >      info[k].identArray[i] = info[k].atoms[i]->getIdent();
1305 >    }
1306    }
1307 + }
1308  
1146  current_mol = 0;
1147  current_comp_mol = 0;
1148  current_comp = 0;
1149  current_atom_ndx = 0;
1309  
1310 <  for( i=0; i < n_cells ; i++ ){
1311 <    for( j=0; j < n_cells; j++ ){
1312 <      for( k=0; k < n_cells; k++ ){
1310 > void SimSetup::createFF(void){
1311 >  switch (ffCase){
1312 >    case FF_DUFF:
1313 >      the_ff = new DUFF();
1314 >      break;
1315  
1316 <        makeElement( i * cellx,
1317 <                     j * celly,
1318 <                     k * cellz );
1316 >    case FF_LJ:
1317 >      the_ff = new LJFF();
1318 >      break;
1319  
1320 <        makeElement( i * cellx + 0.5 * cellx,
1321 <                     j * celly + 0.5 * celly,
1322 <                     k * cellz );
1320 >    case FF_EAM:
1321 >      the_ff = new EAM_FF();
1322 >      break;
1323  
1324 <        makeElement( i * cellx,
1325 <                     j * celly + 0.5 * celly,
1326 <                     k * cellz + 0.5 * cellz );
1324 >    case FF_H2O:
1325 >      the_ff = new WATER();
1326 >      break;
1327  
1328 <        makeElement( i * cellx + 0.5 * cellx,
1329 <                     j * celly,
1330 <                     k * cellz + 0.5 * cellz );
1328 >    default:
1329 >      sprintf(painCave.errMsg,
1330 >              "SimSetup Error. Unrecognized force field in case statement.\n");
1331 >      painCave.isFatal = 1;
1332 >      simError();
1333 >  }
1334 >
1335 > #ifdef IS_MPI
1336 >  strcpy(checkPointMsg, "ForceField creation successful");
1337 >  MPIcheckPoint();
1338 > #endif // is_mpi
1339 > }
1340 >
1341 >
1342 > void SimSetup::compList(void){
1343 >  int i;
1344 >  char* id;
1345 >  LinkedMolStamp* headStamp = new LinkedMolStamp();
1346 >  LinkedMolStamp* currentStamp = NULL;
1347 >  comp_stamps = new MoleculeStamp * [n_components];
1348 >  bool haveCutoffGroups;
1349 >
1350 >  haveCutoffGroups = false;
1351 >  
1352 >  // make an array of molecule stamps that match the components used.
1353 >  // also extract the used stamps out into a separate linked list
1354 >
1355 >  for (i = 0; i < nInfo; i++){
1356 >    info[i].nComponents = n_components;
1357 >    info[i].componentsNmol = components_nmol;
1358 >    info[i].compStamps = comp_stamps;
1359 >    info[i].headStamp = headStamp;
1360 >  }
1361 >
1362 >
1363 >  for (i = 0; i < n_components; i++){
1364 >    id = the_components[i]->getType();
1365 >    comp_stamps[i] = NULL;
1366 >
1367 >    // check to make sure the component isn't already in the list
1368 >
1369 >    comp_stamps[i] = headStamp->match(id);
1370 >    if (comp_stamps[i] == NULL){
1371 >      // extract the component from the list;
1372 >
1373 >      currentStamp = stamps->extractMolStamp(id);
1374 >      if (currentStamp == NULL){
1375 >        sprintf(painCave.errMsg,
1376 >                "SimSetup error: Component \"%s\" was not found in the "
1377 >                "list of declared molecules\n",
1378 >                id);
1379 >        painCave.isFatal = 1;
1380 >        simError();
1381        }
1382 +
1383 +      headStamp->add(currentStamp);
1384 +      comp_stamps[i] = headStamp->match(id);
1385      }
1386 +
1387 +    if(comp_stamps[i]->getNCutoffGroups() > 0)
1388 +      haveCutoffGroups = true;    
1389    }
1390 +    
1391 +  for (i = 0; i < nInfo; i++)
1392 +    info[i].haveCutoffGroups = haveCutoffGroups;
1393  
1394 <  if( have_extra ){
1395 <    done = 0;
1394 > #ifdef IS_MPI
1395 >  strcpy(checkPointMsg, "Component stamps successfully extracted\n");
1396 >  MPIcheckPoint();
1397 > #endif // is_mpi
1398 > }
1399  
1400 <    int start_ndx;
1401 <    for( i=0; i < (n_cells+1) && !done; i++ ){
1179 <      for( j=0; j < (n_cells+1) && !done; j++ ){
1400 > void SimSetup::calcSysValues(void){
1401 >  int i;
1402  
1403 <        if( i < n_cells ){
1403 >  int* molMembershipArray;
1404  
1405 <          if( j < n_cells ){
1406 <            start_ndx = n_cells;
1407 <          }
1408 <          else start_ndx = 0;
1409 <        }
1410 <        else start_ndx = 0;
1405 >  tot_atoms = 0;
1406 >  tot_bonds = 0;
1407 >  tot_bends = 0;
1408 >  tot_torsions = 0;
1409 >  tot_rigid = 0;
1410 >  for (i = 0; i < n_components; i++){
1411 >    tot_atoms += components_nmol[i] * comp_stamps[i]->getNAtoms();
1412 >    tot_bonds += components_nmol[i] * comp_stamps[i]->getNBonds();
1413 >    tot_bends += components_nmol[i] * comp_stamps[i]->getNBends();
1414 >    tot_torsions += components_nmol[i] * comp_stamps[i]->getNTorsions();
1415 >    tot_rigid += components_nmol[i] * comp_stamps[i]->getNRigidBodies();
1416 >  }
1417 >  
1418 >  tot_SRI = tot_bonds + tot_bends + tot_torsions;
1419 >  molMembershipArray = new int[tot_atoms];
1420  
1421 <        for( k=start_ndx; k < (n_cells+1) && !done; k++ ){
1421 >  for (i = 0; i < nInfo; i++){
1422 >    info[i].n_atoms = tot_atoms;
1423 >    info[i].n_bonds = tot_bonds;
1424 >    info[i].n_bends = tot_bends;
1425 >    info[i].n_torsions = tot_torsions;
1426 >    info[i].n_SRI = tot_SRI;
1427 >    info[i].n_mol = tot_nmol;
1428  
1429 <          makeElement( i * cellx,
1430 <                       j * celly,
1431 <                       k * cellz );
1195 <          done = ( current_mol >= tot_nmol );
1429 >    info[i].molMembershipArray = molMembershipArray;
1430 >  }
1431 > }
1432  
1433 <          if( !done && n_per_extra > 1 ){
1198 <            makeElement( i * cellx + 0.5 * cellx,
1199 <                         j * celly + 0.5 * celly,
1200 <                         k * cellz );
1201 <            done = ( current_mol >= tot_nmol );
1202 <          }
1433 > #ifdef IS_MPI
1434  
1435 <          if( !done && n_per_extra > 2){
1436 <            makeElement( i * cellx,
1437 <                         j * celly + 0.5 * celly,
1438 <                         k * cellz + 0.5 * cellz );
1439 <            done = ( current_mol >= tot_nmol );
1440 <          }
1435 > void SimSetup::mpiMolDivide(void){
1436 >  int i, j, k;
1437 >  int localMol, allMol;
1438 >  int local_atoms, local_bonds, local_bends, local_torsions, local_SRI;
1439 >  int local_rigid;
1440 >  vector<int> globalMolIndex;
1441  
1442 <          if( !done && n_per_extra > 3){
1443 <            makeElement( i * cellx + 0.5 * cellx,
1444 <                         j * celly,
1445 <                         k * cellz + 0.5 * cellz );
1446 <            done = ( current_mol >= tot_nmol );
1447 <          }
1448 <        }
1442 >  mpiSim = new mpiSimulation(info);
1443 >
1444 >  mpiSim->divideLabor();
1445 >  globalAtomIndex = mpiSim->getGlobalAtomIndex();
1446 >  //globalMolIndex = mpiSim->getGlobalMolIndex();
1447 >
1448 >  // set up the local variables
1449 >
1450 >  mol2proc = mpiSim->getMolToProcMap();
1451 >  molCompType = mpiSim->getMolComponentType();
1452 >
1453 >  allMol = 0;
1454 >  localMol = 0;
1455 >  local_atoms = 0;
1456 >  local_bonds = 0;
1457 >  local_bends = 0;
1458 >  local_torsions = 0;
1459 >  local_rigid = 0;
1460 >  globalAtomCounter = 0;
1461 >
1462 >  for (i = 0; i < n_components; i++){
1463 >    for (j = 0; j < components_nmol[i]; j++){
1464 >      if (mol2proc[allMol] == worldRank){
1465 >        local_atoms += comp_stamps[i]->getNAtoms();
1466 >        local_bonds += comp_stamps[i]->getNBonds();
1467 >        local_bends += comp_stamps[i]->getNBends();
1468 >        local_torsions += comp_stamps[i]->getNTorsions();
1469 >        local_rigid += comp_stamps[i]->getNRigidBodies();
1470 >        localMol++;
1471 >      }      
1472 >      for (k = 0; k < comp_stamps[i]->getNAtoms(); k++){
1473 >        info[0].molMembershipArray[globalAtomCounter] = allMol;
1474 >        globalAtomCounter++;
1475        }
1476 +
1477 +      allMol++;
1478      }
1479    }
1480 +  local_SRI = local_bonds + local_bends + local_torsions;
1481  
1482 +  info[0].n_atoms = mpiSim->getMyNlocal();  
1483 +  
1484  
1485 <  for( i=0; i<simnfo->n_atoms; i++ ){
1486 <    simnfo->atoms[i]->set_vx( 0.0 );
1487 <    simnfo->atoms[i]->set_vy( 0.0 );
1488 <    simnfo->atoms[i]->set_vz( 0.0 );
1485 >  if (local_atoms != info[0].n_atoms){
1486 >    sprintf(painCave.errMsg,
1487 >            "SimSetup error: mpiSim's localAtom (%d) and SimSetup's\n"
1488 >            "\tlocalAtom (%d) are not equal.\n",
1489 >            info[0].n_atoms, local_atoms);
1490 >    painCave.isFatal = 1;
1491 >    simError();
1492    }
1493 +
1494 +  info[0].n_bonds = local_bonds;
1495 +  info[0].n_bends = local_bends;
1496 +  info[0].n_torsions = local_torsions;
1497 +  info[0].n_SRI = local_SRI;
1498 +  info[0].n_mol = localMol;
1499 +
1500 +  strcpy(checkPointMsg, "Passed nlocal consistency check.");
1501 +  MPIcheckPoint();
1502   }
1503  
1504 < void SimSetup::makeElement( double x, double y, double z ){
1504 > #endif // is_mpi
1505  
1232  int k;
1233  AtomStamp* current_atom;
1234  DirectionalAtom* dAtom;
1235  double rotMat[3][3];
1506  
1507 <  for( k=0; k<comp_stamps[current_comp]->getNAtoms(); k++ ){
1507 > void SimSetup::makeSysArrays(void){
1508 >
1509 > #ifndef IS_MPI
1510 >  int k, j;
1511 > #endif // is_mpi
1512 >  int i, l;
1513  
1514 <    current_atom = comp_stamps[current_comp]->getAtom( k );
1515 <    if( !current_atom->havePosition() ){
1241 <      sprintf( painCave.errMsg,
1242 <               "SimSetup:initFromBass error.\n"
1243 <               "\tComponent %s, atom %s does not have a position specified.\n"
1244 <               "\tThe initialization routine is unable to give a start"
1245 <               " position.\n",
1246 <               comp_stamps[current_comp]->getID(),
1247 <               current_atom->getType() );
1248 <      painCave.isFatal = 1;
1249 <      simError();
1250 <    }
1514 >  Atom** the_atoms;
1515 >  Molecule* the_molecules;
1516  
1517 <    the_atoms[current_atom_ndx]->setX( x + current_atom->getPosX() );
1518 <    the_atoms[current_atom_ndx]->setY( y + current_atom->getPosY() );
1254 <    the_atoms[current_atom_ndx]->setZ( z + current_atom->getPosZ() );
1517 >  for (l = 0; l < nInfo; l++){
1518 >    // create the atom and short range interaction arrays
1519  
1520 <    if( the_atoms[current_atom_ndx]->isDirectional() ){
1520 >    the_atoms = new Atom * [info[l].n_atoms];
1521 >    the_molecules = new Molecule[info[l].n_mol];
1522 >    int molIndex;
1523  
1524 <      dAtom = (DirectionalAtom *)the_atoms[current_atom_ndx];
1524 >    // initialize the molecule's stampID's
1525  
1526 <      rotMat[0][0] = 1.0;
1261 <      rotMat[0][1] = 0.0;
1262 <      rotMat[0][2] = 0.0;
1526 > #ifdef IS_MPI
1527  
1264      rotMat[1][0] = 0.0;
1265      rotMat[1][1] = 1.0;
1266      rotMat[1][2] = 0.0;
1528  
1529 <      rotMat[2][0] = 0.0;
1530 <      rotMat[2][1] = 0.0;
1531 <      rotMat[2][2] = 1.0;
1529 >    molIndex = 0;
1530 >    for (i = 0; i < mpiSim->getTotNmol(); i++){
1531 >      if (mol2proc[i] == worldRank){
1532 >        the_molecules[molIndex].setStampID(molCompType[i]);
1533 >        the_molecules[molIndex].setMyIndex(molIndex);
1534 >        the_molecules[molIndex].setGlobalIndex(i);
1535 >        molIndex++;
1536 >      }
1537 >    }
1538  
1539 <      dAtom->setA( rotMat );
1539 > #else // is_mpi
1540 >
1541 >    molIndex = 0;
1542 >    globalAtomCounter = 0;
1543 >    for (i = 0; i < n_components; i++){
1544 >      for (j = 0; j < components_nmol[i]; j++){
1545 >        the_molecules[molIndex].setStampID(i);
1546 >        the_molecules[molIndex].setMyIndex(molIndex);
1547 >        the_molecules[molIndex].setGlobalIndex(molIndex);
1548 >        for (k = 0; k < comp_stamps[i]->getNAtoms(); k++){
1549 >          info[l].molMembershipArray[globalAtomCounter] = molIndex;
1550 >          globalAtomCounter++;
1551 >        }
1552 >        molIndex++;
1553 >      }
1554      }
1555  
1275    current_atom_ndx++;
1276  }
1556  
1557 <  current_mol++;
1279 <  current_comp_mol++;
1557 > #endif // is_mpi
1558  
1559 <  if( current_comp_mol >= components_nmol[current_comp] ){
1559 >    info[l].globalExcludes = new int;
1560 >    info[l].globalExcludes[0] = 0;
1561 >    
1562 >    // set the arrays into the SimInfo object
1563  
1564 <    current_comp_mol = 0;
1565 <    current_comp++;
1564 >    info[l].atoms = the_atoms;
1565 >    info[l].molecules = the_molecules;
1566 >    info[l].nGlobalExcludes = 0;
1567 >    
1568 >    the_ff->setSimInfo(info);
1569 >  }
1570 > }
1571 >
1572 > void SimSetup::makeIntegrator(void){
1573 >  int k;
1574 >
1575 >  NVE<RealIntegrator>* myNVE = NULL;
1576 >  NVT<RealIntegrator>* myNVT = NULL;
1577 >  NPTi<NPT<RealIntegrator> >* myNPTi = NULL;
1578 >  NPTf<NPT<RealIntegrator> >* myNPTf = NULL;
1579 >  NPTxyz<NPT<RealIntegrator> >* myNPTxyz = NULL;
1580 >  
1581 >  for (k = 0; k < nInfo; k++){
1582 >    switch (ensembleCase){
1583 >      case NVE_ENS:
1584 >        if (globals->haveZconstraints()){
1585 >          setupZConstraint(info[k]);
1586 >          myNVE = new ZConstraint<NVE<RealIntegrator> >(&(info[k]), the_ff);
1587 >        }
1588 >        else{
1589 >          myNVE = new NVE<RealIntegrator>(&(info[k]), the_ff);
1590 >        }
1591 >        
1592 >        info->the_integrator = myNVE;
1593 >        break;
1594 >
1595 >      case NVT_ENS:
1596 >        if (globals->haveZconstraints()){
1597 >          setupZConstraint(info[k]);
1598 >          myNVT = new ZConstraint<NVT<RealIntegrator> >(&(info[k]), the_ff);
1599 >        }
1600 >        else
1601 >          myNVT = new NVT<RealIntegrator>(&(info[k]), the_ff);
1602 >
1603 >        myNVT->setTargetTemp(globals->getTargetTemp());
1604 >
1605 >        if (globals->haveTauThermostat())
1606 >          myNVT->setTauThermostat(globals->getTauThermostat());
1607 >        else{
1608 >          sprintf(painCave.errMsg,
1609 >                  "SimSetup error: If you use the NVT\n"
1610 >                  "\tensemble, you must set tauThermostat.\n");
1611 >          painCave.isFatal = 1;
1612 >          simError();
1613 >        }
1614 >
1615 >        info->the_integrator = myNVT;
1616 >        break;
1617 >
1618 >      case NPTi_ENS:
1619 >        if (globals->haveZconstraints()){
1620 >          setupZConstraint(info[k]);
1621 >          myNPTi = new ZConstraint<NPTi<NPT <RealIntegrator> > >(&(info[k]), the_ff);
1622 >        }
1623 >        else
1624 >          myNPTi = new NPTi<NPT<RealIntegrator> >(&(info[k]), the_ff);
1625 >
1626 >        myNPTi->setTargetTemp(globals->getTargetTemp());
1627 >
1628 >        if (globals->haveTargetPressure())
1629 >          myNPTi->setTargetPressure(globals->getTargetPressure());
1630 >        else{
1631 >          sprintf(painCave.errMsg,
1632 >                  "SimSetup error: If you use a constant pressure\n"
1633 >                  "\tensemble, you must set targetPressure in the BASS file.\n");
1634 >          painCave.isFatal = 1;
1635 >          simError();
1636 >        }
1637 >
1638 >        if (globals->haveTauThermostat())
1639 >          myNPTi->setTauThermostat(globals->getTauThermostat());
1640 >        else{
1641 >          sprintf(painCave.errMsg,
1642 >                  "SimSetup error: If you use an NPT\n"
1643 >                  "\tensemble, you must set tauThermostat.\n");
1644 >          painCave.isFatal = 1;
1645 >          simError();
1646 >        }
1647 >
1648 >        if (globals->haveTauBarostat())
1649 >          myNPTi->setTauBarostat(globals->getTauBarostat());
1650 >        else{
1651 >          sprintf(painCave.errMsg,
1652 >                  "SimSetup error: If you use an NPT\n"
1653 >                  "\tensemble, you must set tauBarostat.\n");
1654 >          painCave.isFatal = 1;
1655 >          simError();
1656 >        }
1657 >
1658 >        info->the_integrator = myNPTi;
1659 >        break;
1660 >
1661 >      case NPTf_ENS:
1662 >        if (globals->haveZconstraints()){
1663 >          setupZConstraint(info[k]);
1664 >          myNPTf = new ZConstraint<NPTf<NPT <RealIntegrator> > >(&(info[k]), the_ff);
1665 >        }
1666 >        else
1667 >          myNPTf = new NPTf<NPT <RealIntegrator> >(&(info[k]), the_ff);
1668 >
1669 >        myNPTf->setTargetTemp(globals->getTargetTemp());
1670 >
1671 >        if (globals->haveTargetPressure())
1672 >          myNPTf->setTargetPressure(globals->getTargetPressure());
1673 >        else{
1674 >          sprintf(painCave.errMsg,
1675 >                  "SimSetup error: If you use a constant pressure\n"
1676 >                  "\tensemble, you must set targetPressure in the BASS file.\n");
1677 >          painCave.isFatal = 1;
1678 >          simError();
1679 >        }    
1680 >
1681 >        if (globals->haveTauThermostat())
1682 >          myNPTf->setTauThermostat(globals->getTauThermostat());
1683 >
1684 >        else{
1685 >          sprintf(painCave.errMsg,
1686 >                  "SimSetup error: If you use an NPT\n"
1687 >                  "\tensemble, you must set tauThermostat.\n");
1688 >          painCave.isFatal = 1;
1689 >          simError();
1690 >        }
1691 >
1692 >        if (globals->haveTauBarostat())
1693 >          myNPTf->setTauBarostat(globals->getTauBarostat());
1694 >
1695 >        else{
1696 >          sprintf(painCave.errMsg,
1697 >                  "SimSetup error: If you use an NPT\n"
1698 >                  "\tensemble, you must set tauBarostat.\n");
1699 >          painCave.isFatal = 1;
1700 >          simError();
1701 >        }
1702 >
1703 >        info->the_integrator = myNPTf;
1704 >        break;
1705 >
1706 >      case NPTxyz_ENS:
1707 >        if (globals->haveZconstraints()){
1708 >          setupZConstraint(info[k]);
1709 >          myNPTxyz = new ZConstraint<NPTxyz<NPT <RealIntegrator> > >(&(info[k]), the_ff);
1710 >        }
1711 >        else
1712 >          myNPTxyz = new NPTxyz<NPT <RealIntegrator> >(&(info[k]), the_ff);
1713 >
1714 >        myNPTxyz->setTargetTemp(globals->getTargetTemp());
1715 >
1716 >        if (globals->haveTargetPressure())
1717 >          myNPTxyz->setTargetPressure(globals->getTargetPressure());
1718 >        else{
1719 >          sprintf(painCave.errMsg,
1720 >                  "SimSetup error: If you use a constant pressure\n"
1721 >                  "\tensemble, you must set targetPressure in the BASS file.\n");
1722 >          painCave.isFatal = 1;
1723 >          simError();
1724 >        }    
1725 >
1726 >        if (globals->haveTauThermostat())
1727 >          myNPTxyz->setTauThermostat(globals->getTauThermostat());
1728 >        else{
1729 >          sprintf(painCave.errMsg,
1730 >                  "SimSetup error: If you use an NPT\n"
1731 >                  "\tensemble, you must set tauThermostat.\n");
1732 >          painCave.isFatal = 1;
1733 >          simError();
1734 >        }
1735 >
1736 >        if (globals->haveTauBarostat())
1737 >          myNPTxyz->setTauBarostat(globals->getTauBarostat());
1738 >        else{
1739 >          sprintf(painCave.errMsg,
1740 >                  "SimSetup error: If you use an NPT\n"
1741 >                  "\tensemble, you must set tauBarostat.\n");
1742 >          painCave.isFatal = 1;
1743 >          simError();
1744 >        }
1745 >
1746 >        info->the_integrator = myNPTxyz;
1747 >        break;
1748 >
1749 >      default:
1750 >        sprintf(painCave.errMsg,
1751 >                "SimSetup Error. Unrecognized ensemble in case statement.\n");
1752 >        painCave.isFatal = 1;
1753 >        simError();
1754 >    }
1755    }
1756   }
1757 +
1758 + void SimSetup::initFortran(void){
1759 +  info[0].refreshSim();
1760 +
1761 +  if (!strcmp(info[0].mixingRule, "standard")){
1762 +    the_ff->initForceField(LB_MIXING_RULE);
1763 +  }
1764 +  else if (!strcmp(info[0].mixingRule, "explicit")){
1765 +    the_ff->initForceField(EXPLICIT_MIXING_RULE);
1766 +  }
1767 +  else{
1768 +    sprintf(painCave.errMsg, "SimSetup Error: unknown mixing rule -> \"%s\"\n",
1769 +            info[0].mixingRule);
1770 +    painCave.isFatal = 1;
1771 +    simError();
1772 +  }
1773 +
1774 +
1775 + #ifdef IS_MPI
1776 +  strcpy(checkPointMsg, "Successfully intialized the mixingRule for Fortran.");
1777 +  MPIcheckPoint();
1778 + #endif // is_mpi
1779 + }
1780 +
1781 + void SimSetup::setupZConstraint(SimInfo& theInfo){
1782 +  int nZConstraints;
1783 +  ZconStamp** zconStamp;
1784 +
1785 +  if (globals->haveZconstraintTime()){
1786 +    //add sample time of z-constraint  into SimInfo's property list                    
1787 +    DoubleData* zconsTimeProp = new DoubleData();
1788 +    zconsTimeProp->setID(ZCONSTIME_ID);
1789 +    zconsTimeProp->setData(globals->getZconsTime());
1790 +    theInfo.addProperty(zconsTimeProp);
1791 +  }
1792 +  else{
1793 +    sprintf(painCave.errMsg,
1794 +            "ZConstraint error: If you use a ZConstraint,\n"
1795 +            "\tyou must set zconsTime.\n");
1796 +    painCave.isFatal = 1;
1797 +    simError();
1798 +  }
1799 +
1800 +  //push zconsTol into siminfo, if user does not specify
1801 +  //value for zconsTol, a default value will be used
1802 +  DoubleData* zconsTol = new DoubleData();
1803 +  zconsTol->setID(ZCONSTOL_ID);
1804 +  if (globals->haveZconsTol()){
1805 +    zconsTol->setData(globals->getZconsTol());
1806 +  }
1807 +  else{
1808 +    double defaultZConsTol = 0.01;
1809 +    sprintf(painCave.errMsg,
1810 +            "ZConstraint Warning: Tolerance for z-constraint method is not specified.\n"
1811 +            "\tOOPSE will use a default value of %f.\n"
1812 +            "\tTo set the tolerance, use the zconsTol variable.\n",
1813 +            defaultZConsTol);
1814 +    painCave.isFatal = 0;
1815 +    simError();      
1816 +
1817 +    zconsTol->setData(defaultZConsTol);
1818 +  }
1819 +  theInfo.addProperty(zconsTol);
1820 +
1821 +  //set Force Subtraction Policy
1822 +  StringData* zconsForcePolicy = new StringData();
1823 +  zconsForcePolicy->setID(ZCONSFORCEPOLICY_ID);
1824 +
1825 +  if (globals->haveZconsForcePolicy()){
1826 +    zconsForcePolicy->setData(globals->getZconsForcePolicy());
1827 +  }
1828 +  else{
1829 +    sprintf(painCave.errMsg,
1830 +            "ZConstraint Warning: No force subtraction policy was set.\n"
1831 +            "\tOOPSE will use PolicyByMass.\n"
1832 +            "\tTo set the policy, use the zconsForcePolicy variable.\n");
1833 +    painCave.isFatal = 0;
1834 +    simError();
1835 +    zconsForcePolicy->setData("BYMASS");
1836 +  }
1837 +
1838 +  theInfo.addProperty(zconsForcePolicy);
1839 +
1840 +  //set zcons gap
1841 +  DoubleData* zconsGap = new DoubleData();
1842 +  zconsGap->setID(ZCONSGAP_ID);
1843 +
1844 +  if (globals->haveZConsGap()){
1845 +    zconsGap->setData(globals->getZconsGap());
1846 +    theInfo.addProperty(zconsGap);  
1847 +  }
1848 +
1849 +  //set zcons fixtime
1850 +  DoubleData* zconsFixtime = new DoubleData();
1851 +  zconsFixtime->setID(ZCONSFIXTIME_ID);
1852 +
1853 +  if (globals->haveZConsFixTime()){
1854 +    zconsFixtime->setData(globals->getZconsFixtime());
1855 +    theInfo.addProperty(zconsFixtime);  
1856 +  }
1857 +
1858 +  //set zconsUsingSMD
1859 +  IntData* zconsUsingSMD = new IntData();
1860 +  zconsUsingSMD->setID(ZCONSUSINGSMD_ID);
1861 +
1862 +  if (globals->haveZConsUsingSMD()){
1863 +    zconsUsingSMD->setData(globals->getZconsUsingSMD());
1864 +    theInfo.addProperty(zconsUsingSMD);  
1865 +  }
1866 +
1867 +  //Determine the name of ouput file and add it into SimInfo's property list
1868 +  //Be careful, do not use inFileName, since it is a pointer which
1869 +  //point to a string at master node, and slave nodes do not contain that string
1870 +
1871 +  string zconsOutput(theInfo.finalName);
1872 +
1873 +  zconsOutput = zconsOutput.substr(0, zconsOutput.rfind(".")) + ".fz";
1874 +
1875 +  StringData* zconsFilename = new StringData();
1876 +  zconsFilename->setID(ZCONSFILENAME_ID);
1877 +  zconsFilename->setData(zconsOutput);
1878 +
1879 +  theInfo.addProperty(zconsFilename);
1880 +
1881 +  //setup index, pos and other parameters of z-constraint molecules
1882 +  nZConstraints = globals->getNzConstraints();
1883 +  theInfo.nZconstraints = nZConstraints;
1884 +
1885 +  zconStamp = globals->getZconStamp();
1886 +  ZConsParaItem tempParaItem;
1887 +
1888 +  ZConsParaData* zconsParaData = new ZConsParaData();
1889 +  zconsParaData->setID(ZCONSPARADATA_ID);
1890 +
1891 +  for (int i = 0; i < nZConstraints; i++){
1892 +    tempParaItem.havingZPos = zconStamp[i]->haveZpos();
1893 +    tempParaItem.zPos = zconStamp[i]->getZpos();
1894 +    tempParaItem.zconsIndex = zconStamp[i]->getMolIndex();
1895 +    tempParaItem.kRatio = zconStamp[i]->getKratio();
1896 +    tempParaItem.havingCantVel = zconStamp[i]->haveCantVel();
1897 +    tempParaItem.cantVel = zconStamp[i]->getCantVel();    
1898 +    zconsParaData->addItem(tempParaItem);
1899 +  }
1900 +
1901 +  //check the uniqueness of index  
1902 +  if(!zconsParaData->isIndexUnique()){
1903 +    sprintf(painCave.errMsg,
1904 +            "ZConstraint Error: molIndex is not unique!\n");
1905 +    painCave.isFatal = 1;
1906 +    simError();
1907 +  }
1908 +
1909 +  //sort the parameters by index of molecules
1910 +  zconsParaData->sortByIndex();
1911 +  
1912 +  //push data into siminfo, therefore, we can retrieve later
1913 +  theInfo.addProperty(zconsParaData);
1914 + }
1915 +
1916 + void SimSetup::makeMinimizer(){
1917 +
1918 +  OOPSEMinimizer* myOOPSEMinimizer;
1919 +  MinimizerParameterSet* param;
1920 +  char minimizerName[100];
1921 +  
1922 +  for (int i = 0; i < nInfo; i++){
1923 +    
1924 +    //prepare parameter set for minimizer
1925 +    param = new MinimizerParameterSet();
1926 +    param->setDefaultParameter();
1927 +
1928 +    if (globals->haveMinimizer()){
1929 +      param->setFTol(globals->getMinFTol());
1930 +    }
1931 +
1932 +    if (globals->haveMinGTol()){
1933 +      param->setGTol(globals->getMinGTol());
1934 +    }
1935 +
1936 +    if (globals->haveMinMaxIter()){
1937 +      param->setMaxIteration(globals->getMinMaxIter());
1938 +    }
1939 +
1940 +    if (globals->haveMinWriteFrq()){
1941 +      param->setMaxIteration(globals->getMinMaxIter());
1942 +    }
1943 +
1944 +    if (globals->haveMinWriteFrq()){
1945 +      param->setWriteFrq(globals->getMinWriteFrq());
1946 +    }
1947 +    
1948 +    if (globals->haveMinStepSize()){
1949 +      param->setStepSize(globals->getMinStepSize());
1950 +    }
1951 +
1952 +    if (globals->haveMinLSMaxIter()){
1953 +      param->setLineSearchMaxIteration(globals->getMinLSMaxIter());
1954 +    }    
1955 +
1956 +    if (globals->haveMinLSTol()){
1957 +      param->setLineSearchTol(globals->getMinLSTol());
1958 +    }    
1959 +
1960 +    strcpy(minimizerName, globals->getMinimizer());
1961 +
1962 +    if (!strcasecmp(minimizerName, "CG")){
1963 +      myOOPSEMinimizer = new PRCGMinimizer(&(info[i]), the_ff, param);
1964 +    }
1965 +    else if (!strcasecmp(minimizerName, "SD")){
1966 +    //myOOPSEMinimizer = MinimizerFactory.creatMinimizer("", &(info[i]), the_ff, param);
1967 +      myOOPSEMinimizer = new SDMinimizer(&(info[i]), the_ff, param);
1968 +    }
1969 +    else{
1970 +          sprintf(painCave.errMsg,
1971 +                  "SimSetup error: Unrecognized Minimizer, use Conjugate Gradient \n");
1972 +          painCave.isFatal = 0;
1973 +          simError();
1974 +
1975 +      myOOPSEMinimizer = new PRCGMinimizer(&(info[i]), the_ff, param);          
1976 +    }
1977 +     info[i].the_integrator = myOOPSEMinimizer;
1978 +
1979 +     //store the minimizer into simInfo
1980 +     info[i].the_minimizer = myOOPSEMinimizer;
1981 +     info[i].has_minimizer = true;
1982 +  }
1983 +
1984 + }

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines