ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 1180
Committed: Thu May 20 20:24:07 2004 UTC (20 years, 4 months ago) by chrisfen
File size: 54684 byte(s)
Log Message:
Several additions... Restraints.cpp and .hpp were included for restraining particles in thermodynamic integration.  By including these, changes were made in Integrator, SimInfo, ForceFeilds, SimSetup, StatWriter, and possibly some other files.  Two bass keywords were also added for performing thermodynamic integration: a lambda value one and a k power one.

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