ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 1214
Committed: Tue Jun 1 18:42:58 2004 UTC (20 years, 3 months ago) by gezelter
File size: 62504 byte(s)
Log Message:
Cutoff Groups for MPI

File Contents

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