ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 1211
Committed: Tue Jun 1 15:57:30 2004 UTC (20 years, 3 months ago) by tim
File size: 59726 byte(s)
Log Message:
cutoff group in progress

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