ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 1204
Committed: Thu May 27 19:26:42 2004 UTC (20 years, 3 months ago) by gezelter
File size: 56727 byte(s)
Log Message:
bugfix in simsetup?

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 chrisfen 1187 if (globals->getUseThermInt()) {
954     if (globals->haveThermIntLambda() && globals->haveThermIntK()) {
955     info[i].useThermInt = globals->getUseThermInt();
956     info[i].thermIntLambda = globals->getThermIntLambda();
957     info[i].thermIntK = globals->getThermIntK();
958    
959     Restraints *myRestraint = new Restraints(tot_nmol, info[i].thermIntLambda, info[i].thermIntK);
960     info[i].restraint = myRestraint;
961     }
962     else {
963     sprintf(painCave.errMsg,
964     "SimSetup Error:\n"
965     "\tKeyword useThermInt was set to 'true' but\n"
966     "\tthermodynamicIntegrationLambda (and/or\n"
967     "\tthermodynamicIntegrationK) was not specified.\n"
968     "\tPlease provide a lambda value and k value in your .bass file.\n");
969     painCave.isFatal = 1;
970     simError();
971     }
972 chrisfen 1180 }
973 chrisfen 1187 else if(globals->haveThermIntLambda() || globals->haveThermIntK()){
974     sprintf(painCave.errMsg,
975     "SimSetup Warning: If you want to use Thermodynamic\n"
976     "\tIntegration, set useThermInt to 'true' in your .bass file.\n"
977     "\tThe useThermInt keyword is 'false' by default, so your\n"
978     "\tlambda and/or k values are being ignored.\n");
979     painCave.isFatal = 0;
980     simError();
981     }
982 tim 708 }
983 mmeineke 841
984 tim 722 //setup seed for random number generator
985 tim 708 int seedValue;
986    
987 tim 722 if (globals->haveSeed()){
988 tim 708 seedValue = globals->getSeed();
989 tim 722
990     if(seedValue / 1E9 == 0){
991     sprintf(painCave.errMsg,
992     "Seed for sprng library should contain at least 9 digits\n"
993     "OOPSE will generate a seed for user\n");
994     painCave.isFatal = 0;
995     simError();
996    
997     //using seed generated by system instead of invalid seed set by user
998     #ifndef IS_MPI
999     seedValue = make_sprng_seed();
1000     #else
1001     if (worldRank == 0){
1002     seedValue = make_sprng_seed();
1003     }
1004     MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD);
1005     #endif
1006     }
1007     }//end of if branch of globals->haveSeed()
1008 tim 708 else{
1009 tim 722
1010 tim 708 #ifndef IS_MPI
1011 tim 722 seedValue = make_sprng_seed();
1012 tim 708 #else
1013 tim 722 if (worldRank == 0){
1014     seedValue = make_sprng_seed();
1015 tim 708 }
1016 tim 722 MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD);
1017 tim 708 #endif
1018 tim 722 }//end of globals->haveSeed()
1019 tim 708
1020 tim 722 for (int i = 0; i < nInfo; i++){
1021 tim 708 info[i].setSeed(seedValue);
1022     }
1023 tim 1031
1024 mmeineke 614 #ifdef IS_MPI
1025 chrisfen 999 strcpy(checkPointMsg, "Successfully gathered all information from Bass\n");
1026 mmeineke 614 MPIcheckPoint();
1027     #endif // is_mpi
1028     }
1029    
1030    
1031 tim 722 void SimSetup::finalInfoCheck(void){
1032 mmeineke 614 int index;
1033     int usesDipoles;
1034 tim 1113 int usesCharges;
1035 mmeineke 670 int i;
1036 mmeineke 614
1037 tim 722 for (i = 0; i < nInfo; i++){
1038 mmeineke 670 // check electrostatic parameters
1039 tim 722
1040 mmeineke 670 index = 0;
1041     usesDipoles = 0;
1042 tim 722 while ((index < info[i].n_atoms) && !usesDipoles){
1043 mmeineke 670 usesDipoles = (info[i].atoms[index])->hasDipole();
1044     index++;
1045     }
1046 tim 1113 index = 0;
1047     usesCharges = 0;
1048     while ((index < info[i].n_atoms) && !usesCharges){
1049     usesCharges= (info[i].atoms[index])->hasCharge();
1050     index++;
1051     }
1052 mmeineke 614 #ifdef IS_MPI
1053 mmeineke 670 int myUse = usesDipoles;
1054 tim 722 MPI_Allreduce(&myUse, &usesDipoles, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
1055 mmeineke 614 #endif //is_mpi
1056 tim 722
1057 gezelter 1154 double theRcut, theRsw;
1058 tim 722
1059 gezelter 1163 if (globals->haveRcut()) {
1060     theRcut = globals->getRcut();
1061    
1062     if (globals->haveRsw())
1063     theRsw = globals->getRsw();
1064     else
1065     theRsw = theRcut;
1066    
1067     info[i].setDefaultRcut(theRcut, theRsw);
1068    
1069     } else {
1070    
1071     the_ff->calcRcut();
1072     theRcut = info[i].getRcut();
1073    
1074     if (globals->haveRsw())
1075     theRsw = globals->getRsw();
1076     else
1077     theRsw = theRcut;
1078    
1079     info[i].setDefaultRcut(theRcut, theRsw);
1080     }
1081    
1082 tim 722 if (globals->getUseRF()){
1083 mmeineke 670 info[i].useReactionField = 1;
1084 gezelter 1163
1085 gezelter 1154 if (!globals->haveRcut()){
1086 tim 722 sprintf(painCave.errMsg,
1087 gezelter 1154 "SimSetup Warning: No value was set for the cutoffRadius.\n"
1088 gezelter 965 "\tOOPSE will use a default value of 15.0 angstroms"
1089 gezelter 1154 "\tfor the cutoffRadius.\n");
1090 tim 722 painCave.isFatal = 0;
1091     simError();
1092 gezelter 1154 theRcut = 15.0;
1093 mmeineke 614 }
1094 tim 722 else{
1095 gezelter 1154 theRcut = globals->getRcut();
1096 mmeineke 614 }
1097 tim 722
1098 gezelter 1154 if (!globals->haveRsw()){
1099 tim 722 sprintf(painCave.errMsg,
1100 gezelter 1154 "SimSetup Warning: No value was set for switchingRadius.\n"
1101 gezelter 965 "\tOOPSE will use a default value of\n"
1102 gezelter 1154 "\t0.95 * cutoffRadius for the switchingRadius\n");
1103 tim 722 painCave.isFatal = 0;
1104     simError();
1105 gezelter 1154 theRsw = 0.95 * theRcut;
1106 mmeineke 670 }
1107 tim 722 else{
1108 gezelter 1154 theRsw = globals->getRsw();
1109 mmeineke 670 }
1110 tim 722
1111 gezelter 1154 info[i].setDefaultRcut(theRcut, theRsw);
1112 tim 722
1113     if (!globals->haveDielectric()){
1114     sprintf(painCave.errMsg,
1115 gezelter 965 "SimSetup Error: No Dielectric constant was set.\n"
1116     "\tYou are trying to use Reaction Field without"
1117     "\tsetting a dielectric constant!\n");
1118 tim 722 painCave.isFatal = 1;
1119     simError();
1120     }
1121     info[i].dielectric = globals->getDielectric();
1122     }
1123     else{
1124 tim 1113 if (usesDipoles || usesCharges){
1125 gezelter 1154
1126     if (!globals->haveRcut()){
1127 tim 722 sprintf(painCave.errMsg,
1128 gezelter 1154 "SimSetup Warning: No value was set for the cutoffRadius.\n"
1129 gezelter 965 "\tOOPSE will use a default value of 15.0 angstroms"
1130 gezelter 1154 "\tfor the cutoffRadius.\n");
1131     painCave.isFatal = 0;
1132     simError();
1133     theRcut = 15.0;
1134     }
1135 tim 722 else{
1136 gezelter 1154 theRcut = globals->getRcut();
1137 tim 722 }
1138 gezelter 1154
1139     if (!globals->haveRsw()){
1140 tim 722 sprintf(painCave.errMsg,
1141 gezelter 1154 "SimSetup Warning: No value was set for switchingRadius.\n"
1142 gezelter 965 "\tOOPSE will use a default value of\n"
1143 gezelter 1154 "\t0.95 * cutoffRadius for the switchingRadius\n");
1144 tim 722 painCave.isFatal = 0;
1145     simError();
1146 gezelter 1154 theRsw = 0.95 * theRcut;
1147 tim 722 }
1148     else{
1149 gezelter 1154 theRsw = globals->getRsw();
1150 tim 722 }
1151 gezelter 1154
1152     info[i].setDefaultRcut(theRcut, theRsw);
1153 mmeineke 859
1154 tim 722 }
1155     }
1156 mmeineke 670 }
1157 mmeineke 614 #ifdef IS_MPI
1158 tim 722 strcpy(checkPointMsg, "post processing checks out");
1159 mmeineke 614 MPIcheckPoint();
1160     #endif // is_mpi
1161 gezelter 1163
1162     // clean up the forcefield
1163     the_ff->cleanMe();
1164 mmeineke 614 }
1165 mmeineke 859
1166 tim 722 void SimSetup::initSystemCoords(void){
1167 mmeineke 670 int i;
1168 tim 722
1169 tim 689 char* inName;
1170    
1171 tim 722 (info[0].getConfiguration())->createArrays(info[0].n_atoms);
1172 mmeineke 614
1173 tim 722 for (i = 0; i < info[0].n_atoms; i++)
1174     info[0].atoms[i]->setCoords();
1175    
1176     if (globals->haveInitialConfig()){
1177 mmeineke 670 InitializeFromFile* fileInit;
1178 mmeineke 614 #ifdef IS_MPI // is_mpi
1179 tim 722 if (worldRank == 0){
1180 mmeineke 614 #endif //is_mpi
1181 tim 689 inName = globals->getInitialConfig();
1182 tim 722 fileInit = new InitializeFromFile(inName);
1183 mmeineke 614 #ifdef IS_MPI
1184 tim 722 }
1185     else
1186     fileInit = new InitializeFromFile(NULL);
1187 mmeineke 614 #endif
1188 tim 722 fileInit->readInit(info); // default velocities on
1189    
1190 mmeineke 670 delete fileInit;
1191     }
1192     else{
1193 mmeineke 859
1194 mmeineke 670 // no init from bass
1195 mmeineke 859
1196 tim 722 sprintf(painCave.errMsg,
1197 mmeineke 859 "Cannot intialize a simulation without an initial configuration file.\n");
1198 mmeineke 787 painCave.isFatal = 1;;
1199 mmeineke 670 simError();
1200 mmeineke 859
1201 mmeineke 670 }
1202 tim 722
1203 mmeineke 614 #ifdef IS_MPI
1204 tim 722 strcpy(checkPointMsg, "Successfully read in the initial configuration");
1205 mmeineke 614 MPIcheckPoint();
1206     #endif // is_mpi
1207     }
1208    
1209    
1210 tim 722 void SimSetup::makeOutNames(void){
1211 mmeineke 670 int k;
1212 mmeineke 614
1213 mmeineke 670
1214 tim 722 for (k = 0; k < nInfo; k++){
1215 mmeineke 614 #ifdef IS_MPI
1216 tim 722 if (worldRank == 0){
1217 mmeineke 614 #endif // is_mpi
1218 tim 722
1219     if (globals->haveFinalConfig()){
1220     strcpy(info[k].finalName, globals->getFinalConfig());
1221 mmeineke 614 }
1222     else{
1223 tim 722 strcpy(info[k].finalName, inFileName);
1224     char* endTest;
1225     int nameLength = strlen(info[k].finalName);
1226     endTest = &(info[k].finalName[nameLength - 5]);
1227     if (!strcmp(endTest, ".bass")){
1228     strcpy(endTest, ".eor");
1229     }
1230     else if (!strcmp(endTest, ".BASS")){
1231     strcpy(endTest, ".eor");
1232     }
1233     else{
1234     endTest = &(info[k].finalName[nameLength - 4]);
1235     if (!strcmp(endTest, ".bss")){
1236     strcpy(endTest, ".eor");
1237     }
1238     else if (!strcmp(endTest, ".mdl")){
1239     strcpy(endTest, ".eor");
1240     }
1241     else{
1242     strcat(info[k].finalName, ".eor");
1243     }
1244     }
1245 mmeineke 614 }
1246 tim 722
1247 mmeineke 670 // make the sample and status out names
1248 tim 722
1249     strcpy(info[k].sampleName, inFileName);
1250 mmeineke 670 char* endTest;
1251 tim 722 int nameLength = strlen(info[k].sampleName);
1252 mmeineke 670 endTest = &(info[k].sampleName[nameLength - 5]);
1253 tim 722 if (!strcmp(endTest, ".bass")){
1254     strcpy(endTest, ".dump");
1255 mmeineke 614 }
1256 tim 722 else if (!strcmp(endTest, ".BASS")){
1257     strcpy(endTest, ".dump");
1258 mmeineke 614 }
1259     else{
1260 tim 722 endTest = &(info[k].sampleName[nameLength - 4]);
1261     if (!strcmp(endTest, ".bss")){
1262     strcpy(endTest, ".dump");
1263     }
1264     else if (!strcmp(endTest, ".mdl")){
1265     strcpy(endTest, ".dump");
1266     }
1267     else{
1268     strcat(info[k].sampleName, ".dump");
1269     }
1270 mmeineke 614 }
1271 tim 722
1272     strcpy(info[k].statusName, inFileName);
1273     nameLength = strlen(info[k].statusName);
1274 mmeineke 670 endTest = &(info[k].statusName[nameLength - 5]);
1275 tim 722 if (!strcmp(endTest, ".bass")){
1276     strcpy(endTest, ".stat");
1277 mmeineke 614 }
1278 tim 722 else if (!strcmp(endTest, ".BASS")){
1279     strcpy(endTest, ".stat");
1280 mmeineke 614 }
1281     else{
1282 tim 722 endTest = &(info[k].statusName[nameLength - 4]);
1283     if (!strcmp(endTest, ".bss")){
1284     strcpy(endTest, ".stat");
1285     }
1286     else if (!strcmp(endTest, ".mdl")){
1287     strcpy(endTest, ".stat");
1288     }
1289     else{
1290     strcat(info[k].statusName, ".stat");
1291     }
1292 mmeineke 614 }
1293 tim 722
1294 chrisfen 1180 strcpy(info[k].rawPotName, inFileName);
1295     nameLength = strlen(info[k].rawPotName);
1296     endTest = &(info[k].rawPotName[nameLength - 5]);
1297     if (!strcmp(endTest, ".bass")){
1298     strcpy(endTest, ".raw");
1299     }
1300     else if (!strcmp(endTest, ".BASS")){
1301     strcpy(endTest, ".raw");
1302     }
1303     else{
1304     endTest = &(info[k].rawPotName[nameLength - 4]);
1305     if (!strcmp(endTest, ".bss")){
1306     strcpy(endTest, ".raw");
1307     }
1308     else if (!strcmp(endTest, ".mdl")){
1309     strcpy(endTest, ".raw");
1310     }
1311     else{
1312     strcat(info[k].rawPotName, ".raw");
1313     }
1314     }
1315    
1316 mmeineke 670 #ifdef IS_MPI
1317 tim 722
1318 mmeineke 614 }
1319 mmeineke 670 #endif // is_mpi
1320 mmeineke 614 }
1321     }
1322    
1323    
1324 tim 722 void SimSetup::sysObjectsCreation(void){
1325     int i, k;
1326    
1327 mmeineke 614 // create the forceField
1328 tim 689
1329 mmeineke 614 createFF();
1330    
1331     // extract componentList
1332    
1333     compList();
1334    
1335     // calc the number of atoms, bond, bends, and torsions
1336    
1337     calcSysValues();
1338    
1339     #ifdef IS_MPI
1340     // divide the molecules among the processors
1341 tim 722
1342 mmeineke 614 mpiMolDivide();
1343     #endif //is_mpi
1344 tim 722
1345 mmeineke 614 // create the atom and SRI arrays. Also initialize Molecule Stamp ID's
1346 tim 722
1347 mmeineke 614 makeSysArrays();
1348    
1349 mmeineke 616 // make and initialize the molecules (all but atomic coordinates)
1350 tim 722
1351 mmeineke 616 makeMolecules();
1352 tim 722
1353     for (k = 0; k < nInfo; k++){
1354 mmeineke 670 info[k].identArray = new int[info[k].n_atoms];
1355 tim 722 for (i = 0; i < info[k].n_atoms; i++){
1356 mmeineke 670 info[k].identArray[i] = info[k].atoms[i]->getIdent();
1357     }
1358 mmeineke 616 }
1359 mmeineke 614 }
1360    
1361    
1362 tim 722 void SimSetup::createFF(void){
1363     switch (ffCase){
1364     case FF_DUFF:
1365     the_ff = new DUFF();
1366     break;
1367 mmeineke 614
1368 tim 722 case FF_LJ:
1369     the_ff = new LJFF();
1370     break;
1371 mmeineke 614
1372 tim 722 case FF_EAM:
1373     the_ff = new EAM_FF();
1374     break;
1375 mmeineke 614
1376 chrisfen 999 case FF_H2O:
1377     the_ff = new WATER();
1378     break;
1379    
1380 tim 722 default:
1381     sprintf(painCave.errMsg,
1382     "SimSetup Error. Unrecognized force field in case statement.\n");
1383     painCave.isFatal = 1;
1384     simError();
1385 mmeineke 614 }
1386    
1387     #ifdef IS_MPI
1388 tim 722 strcpy(checkPointMsg, "ForceField creation successful");
1389 mmeineke 614 MPIcheckPoint();
1390     #endif // is_mpi
1391     }
1392    
1393    
1394 tim 722 void SimSetup::compList(void){
1395 mmeineke 616 int i;
1396 mmeineke 670 char* id;
1397     LinkedMolStamp* headStamp = new LinkedMolStamp();
1398     LinkedMolStamp* currentStamp = NULL;
1399 tim 722 comp_stamps = new MoleculeStamp * [n_components];
1400 tim 1157 bool haveCutoffGroups;
1401 tim 722
1402 tim 1157 haveCutoffGroups = false;
1403    
1404 mmeineke 614 // make an array of molecule stamps that match the components used.
1405     // also extract the used stamps out into a separate linked list
1406 tim 722
1407     for (i = 0; i < nInfo; i++){
1408 mmeineke 670 info[i].nComponents = n_components;
1409     info[i].componentsNmol = components_nmol;
1410     info[i].compStamps = comp_stamps;
1411     info[i].headStamp = headStamp;
1412     }
1413 mmeineke 614
1414    
1415 tim 722 for (i = 0; i < n_components; i++){
1416 mmeineke 614 id = the_components[i]->getType();
1417     comp_stamps[i] = NULL;
1418 tim 722
1419 mmeineke 614 // check to make sure the component isn't already in the list
1420    
1421 tim 722 comp_stamps[i] = headStamp->match(id);
1422     if (comp_stamps[i] == NULL){
1423 mmeineke 614 // extract the component from the list;
1424 tim 722
1425     currentStamp = stamps->extractMolStamp(id);
1426     if (currentStamp == NULL){
1427     sprintf(painCave.errMsg,
1428     "SimSetup error: Component \"%s\" was not found in the "
1429     "list of declared molecules\n",
1430     id);
1431     painCave.isFatal = 1;
1432     simError();
1433 mmeineke 614 }
1434 tim 722
1435     headStamp->add(currentStamp);
1436     comp_stamps[i] = headStamp->match(id);
1437 mmeineke 614 }
1438 tim 1157
1439     if(comp_stamps[i]->getNCutoffGroups() > 0)
1440     haveCutoffGroups = true;
1441 mmeineke 614 }
1442 tim 1157
1443     for (i = 0; i < nInfo; i++)
1444     info[i].haveCutoffGroups = haveCutoffGroups;
1445 mmeineke 614
1446     #ifdef IS_MPI
1447 tim 722 strcpy(checkPointMsg, "Component stamps successfully extracted\n");
1448 mmeineke 614 MPIcheckPoint();
1449     #endif // is_mpi
1450 tim 722 }
1451 mmeineke 614
1452 tim 722 void SimSetup::calcSysValues(void){
1453 gezelter 1203 int i, j;
1454     int ncutgroups, atomsingroups, ngroupsinstamp;
1455 mmeineke 614
1456 tim 722 int* molMembershipArray;
1457 gezelter 1203 CutoffGroupStamp* cg;
1458 mmeineke 614
1459     tot_atoms = 0;
1460     tot_bonds = 0;
1461     tot_bends = 0;
1462     tot_torsions = 0;
1463 gezelter 1097 tot_rigid = 0;
1464 gezelter 1203 tot_groups = 0;
1465 tim 722 for (i = 0; i < n_components; i++){
1466     tot_atoms += components_nmol[i] * comp_stamps[i]->getNAtoms();
1467     tot_bonds += components_nmol[i] * comp_stamps[i]->getNBonds();
1468     tot_bends += components_nmol[i] * comp_stamps[i]->getNBends();
1469 mmeineke 614 tot_torsions += components_nmol[i] * comp_stamps[i]->getNTorsions();
1470 gezelter 1097 tot_rigid += components_nmol[i] * comp_stamps[i]->getNRigidBodies();
1471 gezelter 1203
1472     ncutgroups = comp_stamps[i]->getNCutoffGroups();
1473     atomsingroups = 0;
1474     for (j=0; j < ncutgroups; j++) {
1475     cg = comp_stamps[i]->getCutoffGroup(j);
1476     atomsingroups += cg->getNMembers();
1477     }
1478     ngroupsinstamp = comp_stamps[i]->getNAtoms() - atomsingroups + ncutgroups;
1479     tot_groups += components_nmol[i] * ngroupsinstamp;
1480 mmeineke 614 }
1481 gezelter 1097
1482 mmeineke 614 tot_SRI = tot_bonds + tot_bends + tot_torsions;
1483 mmeineke 670 molMembershipArray = new int[tot_atoms];
1484 tim 722
1485     for (i = 0; i < nInfo; i++){
1486 mmeineke 670 info[i].n_atoms = tot_atoms;
1487     info[i].n_bonds = tot_bonds;
1488     info[i].n_bends = tot_bends;
1489     info[i].n_torsions = tot_torsions;
1490     info[i].n_SRI = tot_SRI;
1491     info[i].n_mol = tot_nmol;
1492 gezelter 1203 info[i].ngroup = tot_groups;
1493 mmeineke 670 info[i].molMembershipArray = molMembershipArray;
1494 tim 722 }
1495 mmeineke 614 }
1496    
1497     #ifdef IS_MPI
1498    
1499 tim 722 void SimSetup::mpiMolDivide(void){
1500 mmeineke 616 int i, j, k;
1501 mmeineke 614 int localMol, allMol;
1502     int local_atoms, local_bonds, local_bends, local_torsions, local_SRI;
1503 gezelter 1204 int local_rigid, local_groups;
1504 tim 1108 vector<int> globalMolIndex;
1505 gezelter 1204 int ncutgroups, atomsingroups, ngroupsinstamp;
1506     CutoffGroupStamp* cg;
1507 mmeineke 614
1508 tim 722 mpiSim = new mpiSimulation(info);
1509    
1510 tim 1108 mpiSim->divideLabor();
1511     globalAtomIndex = mpiSim->getGlobalAtomIndex();
1512 tim 1129 //globalMolIndex = mpiSim->getGlobalMolIndex();
1513 mmeineke 614
1514     // set up the local variables
1515 tim 722
1516 mmeineke 614 mol2proc = mpiSim->getMolToProcMap();
1517     molCompType = mpiSim->getMolComponentType();
1518 tim 722
1519 mmeineke 614 allMol = 0;
1520     localMol = 0;
1521     local_atoms = 0;
1522     local_bonds = 0;
1523     local_bends = 0;
1524     local_torsions = 0;
1525 gezelter 1097 local_rigid = 0;
1526 gezelter 1204 local_groups = 0;
1527 tim 1108 globalAtomCounter = 0;
1528 mmeineke 614
1529 tim 722 for (i = 0; i < n_components; i++){
1530     for (j = 0; j < components_nmol[i]; j++){
1531     if (mol2proc[allMol] == worldRank){
1532     local_atoms += comp_stamps[i]->getNAtoms();
1533     local_bonds += comp_stamps[i]->getNBonds();
1534     local_bends += comp_stamps[i]->getNBends();
1535     local_torsions += comp_stamps[i]->getNTorsions();
1536 gezelter 1097 local_rigid += comp_stamps[i]->getNRigidBodies();
1537 gezelter 1204
1538     ncutgroups = comp_stamps[i]->getNCutoffGroups();
1539     atomsingroups = 0;
1540     for (k=0; k < ncutgroups; k++) {
1541     cg = comp_stamps[i]->getCutoffGroup(k);
1542     atomsingroups += cg->getNMembers();
1543     }
1544     ngroupsinstamp = comp_stamps[i]->getNAtoms() - atomsingroups +
1545     ncutgroups;
1546     local_groups += ngroupsinstamp;
1547    
1548 tim 722 localMol++;
1549 mmeineke 614 }
1550 tim 722 for (k = 0; k < comp_stamps[i]->getNAtoms(); k++){
1551 tim 1108 info[0].molMembershipArray[globalAtomCounter] = allMol;
1552     globalAtomCounter++;
1553 mmeineke 614 }
1554    
1555 tim 722 allMol++;
1556 mmeineke 614 }
1557     }
1558     local_SRI = local_bonds + local_bends + local_torsions;
1559 tim 722
1560 gezelter 1203 info[0].n_atoms = mpiSim->getNAtomsLocal();
1561 gezelter 1097
1562 tim 722 if (local_atoms != info[0].n_atoms){
1563     sprintf(painCave.errMsg,
1564 gezelter 965 "SimSetup error: mpiSim's localAtom (%d) and SimSetup's\n"
1565     "\tlocalAtom (%d) are not equal.\n",
1566 tim 722 info[0].n_atoms, local_atoms);
1567 mmeineke 614 painCave.isFatal = 1;
1568     simError();
1569     }
1570    
1571 gezelter 1204 info[0].ngroup = mpiSim->getNGroupsLocal();
1572     if (local_groups != info[0].ngroup){
1573     sprintf(painCave.errMsg,
1574     "SimSetup error: mpiSim's localGroups (%d) and SimSetup's\n"
1575     "\tlocalGroups (%d) are not equal.\n",
1576     info[0].ngroup, local_groups);
1577     painCave.isFatal = 1;
1578     simError();
1579     }
1580    
1581 mmeineke 670 info[0].n_bonds = local_bonds;
1582     info[0].n_bends = local_bends;
1583     info[0].n_torsions = local_torsions;
1584     info[0].n_SRI = local_SRI;
1585     info[0].n_mol = localMol;
1586 mmeineke 614
1587 tim 722 strcpy(checkPointMsg, "Passed nlocal consistency check.");
1588 mmeineke 614 MPIcheckPoint();
1589     }
1590 tim 722
1591 mmeineke 614 #endif // is_mpi
1592    
1593    
1594 tim 722 void SimSetup::makeSysArrays(void){
1595 mmeineke 787
1596     #ifndef IS_MPI
1597     int k, j;
1598     #endif // is_mpi
1599     int i, l;
1600 mmeineke 614
1601 mmeineke 670 Atom** the_atoms;
1602     Molecule* the_molecules;
1603 mmeineke 616
1604 tim 722 for (l = 0; l < nInfo; l++){
1605 mmeineke 670 // create the atom and short range interaction arrays
1606 tim 722
1607     the_atoms = new Atom * [info[l].n_atoms];
1608 mmeineke 670 the_molecules = new Molecule[info[l].n_mol];
1609     int molIndex;
1610 mmeineke 614
1611 mmeineke 670 // initialize the molecule's stampID's
1612 tim 722
1613 mmeineke 670 #ifdef IS_MPI
1614 tim 722
1615    
1616 mmeineke 670 molIndex = 0;
1617 gezelter 1203 for (i = 0; i < mpiSim->getNMolGlobal(); i++){
1618 tim 722 if (mol2proc[i] == worldRank){
1619     the_molecules[molIndex].setStampID(molCompType[i]);
1620     the_molecules[molIndex].setMyIndex(molIndex);
1621     the_molecules[molIndex].setGlobalIndex(i);
1622     molIndex++;
1623 mmeineke 670 }
1624 mmeineke 614 }
1625 tim 722
1626 mmeineke 614 #else // is_mpi
1627 tim 722
1628 mmeineke 670 molIndex = 0;
1629 tim 1108 globalAtomCounter = 0;
1630 tim 722 for (i = 0; i < n_components; i++){
1631     for (j = 0; j < components_nmol[i]; j++){
1632     the_molecules[molIndex].setStampID(i);
1633     the_molecules[molIndex].setMyIndex(molIndex);
1634     the_molecules[molIndex].setGlobalIndex(molIndex);
1635     for (k = 0; k < comp_stamps[i]->getNAtoms(); k++){
1636 tim 1108 info[l].molMembershipArray[globalAtomCounter] = molIndex;
1637     globalAtomCounter++;
1638 tim 722 }
1639     molIndex++;
1640 mmeineke 614 }
1641     }
1642 tim 722
1643    
1644 mmeineke 614 #endif // is_mpi
1645    
1646 gezelter 1097 info[l].globalExcludes = new int;
1647     info[l].globalExcludes[0] = 0;
1648    
1649 mmeineke 670 // set the arrays into the SimInfo object
1650 mmeineke 614
1651 mmeineke 670 info[l].atoms = the_atoms;
1652     info[l].molecules = the_molecules;
1653     info[l].nGlobalExcludes = 0;
1654 tim 1157
1655 tim 722 the_ff->setSimInfo(info);
1656 mmeineke 670 }
1657 mmeineke 614 }
1658 mmeineke 616
1659 tim 722 void SimSetup::makeIntegrator(void){
1660 mmeineke 670 int k;
1661    
1662 mmeineke 782 NVE<RealIntegrator>* myNVE = NULL;
1663 tim 722 NVT<RealIntegrator>* myNVT = NULL;
1664 mmeineke 778 NPTi<NPT<RealIntegrator> >* myNPTi = NULL;
1665 mmeineke 780 NPTf<NPT<RealIntegrator> >* myNPTf = NULL;
1666 mmeineke 812 NPTxyz<NPT<RealIntegrator> >* myNPTxyz = NULL;
1667 tim 733
1668 tim 722 for (k = 0; k < nInfo; k++){
1669     switch (ensembleCase){
1670     case NVE_ENS:
1671     if (globals->haveZconstraints()){
1672     setupZConstraint(info[k]);
1673 mmeineke 782 myNVE = new ZConstraint<NVE<RealIntegrator> >(&(info[k]), the_ff);
1674 tim 722 }
1675 mmeineke 782 else{
1676     myNVE = new NVE<RealIntegrator>(&(info[k]), the_ff);
1677     }
1678    
1679     info->the_integrator = myNVE;
1680 tim 722 break;
1681 tim 676
1682 tim 722 case NVT_ENS:
1683     if (globals->haveZconstraints()){
1684     setupZConstraint(info[k]);
1685     myNVT = new ZConstraint<NVT<RealIntegrator> >(&(info[k]), the_ff);
1686     }
1687     else
1688     myNVT = new NVT<RealIntegrator>(&(info[k]), the_ff);
1689    
1690 tim 701 myNVT->setTargetTemp(globals->getTargetTemp());
1691 tim 722
1692     if (globals->haveTauThermostat())
1693 tim 701 myNVT->setTauThermostat(globals->getTauThermostat());
1694 tim 722 else{
1695     sprintf(painCave.errMsg,
1696     "SimSetup error: If you use the NVT\n"
1697 gezelter 965 "\tensemble, you must set tauThermostat.\n");
1698 tim 701 painCave.isFatal = 1;
1699     simError();
1700     }
1701 mmeineke 782
1702     info->the_integrator = myNVT;
1703 tim 701 break;
1704 tim 676
1705 tim 722 case NPTi_ENS:
1706     if (globals->haveZconstraints()){
1707     setupZConstraint(info[k]);
1708 mmeineke 778 myNPTi = new ZConstraint<NPTi<NPT <RealIntegrator> > >(&(info[k]), the_ff);
1709 tim 722 }
1710     else
1711 mmeineke 778 myNPTi = new NPTi<NPT<RealIntegrator> >(&(info[k]), the_ff);
1712 tim 722
1713     myNPTi->setTargetTemp(globals->getTargetTemp());
1714    
1715     if (globals->haveTargetPressure())
1716     myNPTi->setTargetPressure(globals->getTargetPressure());
1717     else{
1718     sprintf(painCave.errMsg,
1719     "SimSetup error: If you use a constant pressure\n"
1720 gezelter 965 "\tensemble, you must set targetPressure in the BASS file.\n");
1721 tim 722 painCave.isFatal = 1;
1722     simError();
1723     }
1724    
1725     if (globals->haveTauThermostat())
1726     myNPTi->setTauThermostat(globals->getTauThermostat());
1727     else{
1728     sprintf(painCave.errMsg,
1729     "SimSetup error: If you use an NPT\n"
1730 gezelter 965 "\tensemble, you must set tauThermostat.\n");
1731 tim 722 painCave.isFatal = 1;
1732     simError();
1733     }
1734    
1735     if (globals->haveTauBarostat())
1736     myNPTi->setTauBarostat(globals->getTauBarostat());
1737     else{
1738     sprintf(painCave.errMsg,
1739 tim 701 "SimSetup error: If you use an NPT\n"
1740 gezelter 965 "\tensemble, you must set tauBarostat.\n");
1741 tim 722 painCave.isFatal = 1;
1742     simError();
1743     }
1744 mmeineke 782
1745     info->the_integrator = myNPTi;
1746 tim 722 break;
1747 tim 676
1748 tim 722 case NPTf_ENS:
1749     if (globals->haveZconstraints()){
1750     setupZConstraint(info[k]);
1751 mmeineke 780 myNPTf = new ZConstraint<NPTf<NPT <RealIntegrator> > >(&(info[k]), the_ff);
1752 tim 722 }
1753     else
1754 mmeineke 780 myNPTf = new NPTf<NPT <RealIntegrator> >(&(info[k]), the_ff);
1755 tim 722
1756     myNPTf->setTargetTemp(globals->getTargetTemp());
1757    
1758     if (globals->haveTargetPressure())
1759     myNPTf->setTargetPressure(globals->getTargetPressure());
1760     else{
1761     sprintf(painCave.errMsg,
1762 tim 701 "SimSetup error: If you use a constant pressure\n"
1763 gezelter 965 "\tensemble, you must set targetPressure in the BASS file.\n");
1764 tim 722 painCave.isFatal = 1;
1765     simError();
1766     }
1767    
1768     if (globals->haveTauThermostat())
1769     myNPTf->setTauThermostat(globals->getTauThermostat());
1770 mmeineke 843
1771 tim 722 else{
1772     sprintf(painCave.errMsg,
1773 tim 701 "SimSetup error: If you use an NPT\n"
1774 gezelter 965 "\tensemble, you must set tauThermostat.\n");
1775 tim 722 painCave.isFatal = 1;
1776     simError();
1777     }
1778    
1779     if (globals->haveTauBarostat())
1780     myNPTf->setTauBarostat(globals->getTauBarostat());
1781 mmeineke 843
1782 tim 722 else{
1783     sprintf(painCave.errMsg,
1784     "SimSetup error: If you use an NPT\n"
1785 gezelter 965 "\tensemble, you must set tauBarostat.\n");
1786 tim 722 painCave.isFatal = 1;
1787     simError();
1788     }
1789 mmeineke 782
1790     info->the_integrator = myNPTf;
1791 tim 722 break;
1792 tim 676
1793 mmeineke 812 case NPTxyz_ENS:
1794     if (globals->haveZconstraints()){
1795     setupZConstraint(info[k]);
1796     myNPTxyz = new ZConstraint<NPTxyz<NPT <RealIntegrator> > >(&(info[k]), the_ff);
1797     }
1798     else
1799     myNPTxyz = new NPTxyz<NPT <RealIntegrator> >(&(info[k]), the_ff);
1800    
1801     myNPTxyz->setTargetTemp(globals->getTargetTemp());
1802    
1803     if (globals->haveTargetPressure())
1804     myNPTxyz->setTargetPressure(globals->getTargetPressure());
1805     else{
1806     sprintf(painCave.errMsg,
1807     "SimSetup error: If you use a constant pressure\n"
1808 gezelter 965 "\tensemble, you must set targetPressure in the BASS file.\n");
1809 mmeineke 812 painCave.isFatal = 1;
1810     simError();
1811     }
1812    
1813     if (globals->haveTauThermostat())
1814     myNPTxyz->setTauThermostat(globals->getTauThermostat());
1815     else{
1816     sprintf(painCave.errMsg,
1817     "SimSetup error: If you use an NPT\n"
1818 gezelter 965 "\tensemble, you must set tauThermostat.\n");
1819 mmeineke 812 painCave.isFatal = 1;
1820     simError();
1821     }
1822    
1823     if (globals->haveTauBarostat())
1824     myNPTxyz->setTauBarostat(globals->getTauBarostat());
1825     else{
1826     sprintf(painCave.errMsg,
1827     "SimSetup error: If you use an NPT\n"
1828 gezelter 965 "\tensemble, you must set tauBarostat.\n");
1829 mmeineke 812 painCave.isFatal = 1;
1830     simError();
1831     }
1832    
1833     info->the_integrator = myNPTxyz;
1834     break;
1835    
1836 tim 722 default:
1837     sprintf(painCave.errMsg,
1838     "SimSetup Error. Unrecognized ensemble in case statement.\n");
1839 tim 701 painCave.isFatal = 1;
1840     simError();
1841 tim 660 }
1842 mmeineke 616 }
1843     }
1844    
1845 tim 722 void SimSetup::initFortran(void){
1846     info[0].refreshSim();
1847 mmeineke 616
1848 tim 722 if (!strcmp(info[0].mixingRule, "standard")){
1849     the_ff->initForceField(LB_MIXING_RULE);
1850 mmeineke 616 }
1851 tim 722 else if (!strcmp(info[0].mixingRule, "explicit")){
1852     the_ff->initForceField(EXPLICIT_MIXING_RULE);
1853 mmeineke 616 }
1854     else{
1855 tim 722 sprintf(painCave.errMsg, "SimSetup Error: unknown mixing rule -> \"%s\"\n",
1856     info[0].mixingRule);
1857 mmeineke 616 painCave.isFatal = 1;
1858     simError();
1859     }
1860    
1861    
1862     #ifdef IS_MPI
1863 tim 722 strcpy(checkPointMsg, "Successfully intialized the mixingRule for Fortran.");
1864 mmeineke 616 MPIcheckPoint();
1865     #endif // is_mpi
1866     }
1867 tim 660
1868 tim 722 void SimSetup::setupZConstraint(SimInfo& theInfo){
1869 tim 701 int nZConstraints;
1870     ZconStamp** zconStamp;
1871 tim 682
1872 tim 722 if (globals->haveZconstraintTime()){
1873 tim 701 //add sample time of z-constraint into SimInfo's property list
1874     DoubleData* zconsTimeProp = new DoubleData();
1875     zconsTimeProp->setID(ZCONSTIME_ID);
1876     zconsTimeProp->setData(globals->getZconsTime());
1877     theInfo.addProperty(zconsTimeProp);
1878     }
1879     else{
1880 tim 722 sprintf(painCave.errMsg,
1881 gezelter 965 "ZConstraint error: If you use a ZConstraint,\n"
1882     "\tyou must set zconsTime.\n");
1883 tim 701 painCave.isFatal = 1;
1884 tim 722 simError();
1885 tim 701 }
1886 tim 682
1887 tim 701 //push zconsTol into siminfo, if user does not specify
1888     //value for zconsTol, a default value will be used
1889     DoubleData* zconsTol = new DoubleData();
1890     zconsTol->setID(ZCONSTOL_ID);
1891 tim 722 if (globals->haveZconsTol()){
1892 tim 701 zconsTol->setData(globals->getZconsTol());
1893     }
1894     else{
1895 tim 722 double defaultZConsTol = 0.01;
1896     sprintf(painCave.errMsg,
1897 gezelter 965 "ZConstraint Warning: Tolerance for z-constraint method is not specified.\n"
1898     "\tOOPSE will use a default value of %f.\n"
1899     "\tTo set the tolerance, use the zconsTol variable.\n",
1900 tim 722 defaultZConsTol);
1901 tim 701 painCave.isFatal = 0;
1902     simError();
1903 tim 699
1904 tim 701 zconsTol->setData(defaultZConsTol);
1905     }
1906     theInfo.addProperty(zconsTol);
1907 tim 699
1908 tim 738 //set Force Subtraction Policy
1909 tim 722 StringData* zconsForcePolicy = new StringData();
1910 tim 701 zconsForcePolicy->setID(ZCONSFORCEPOLICY_ID);
1911 tim 722
1912     if (globals->haveZconsForcePolicy()){
1913 tim 701 zconsForcePolicy->setData(globals->getZconsForcePolicy());
1914 tim 722 }
1915 tim 701 else{
1916 tim 722 sprintf(painCave.errMsg,
1917 gezelter 965 "ZConstraint Warning: No force subtraction policy was set.\n"
1918     "\tOOPSE will use PolicyByMass.\n"
1919     "\tTo set the policy, use the zconsForcePolicy variable.\n");
1920 tim 722 painCave.isFatal = 0;
1921     simError();
1922 tim 736 zconsForcePolicy->setData("BYMASS");
1923 tim 701 }
1924 tim 722
1925     theInfo.addProperty(zconsForcePolicy);
1926    
1927 tim 1091 //set zcons gap
1928     DoubleData* zconsGap = new DoubleData();
1929     zconsGap->setID(ZCONSGAP_ID);
1930    
1931     if (globals->haveZConsGap()){
1932     zconsGap->setData(globals->getZconsGap());
1933     theInfo.addProperty(zconsGap);
1934     }
1935    
1936     //set zcons fixtime
1937     DoubleData* zconsFixtime = new DoubleData();
1938     zconsFixtime->setID(ZCONSFIXTIME_ID);
1939    
1940     if (globals->haveZConsFixTime()){
1941     zconsFixtime->setData(globals->getZconsFixtime());
1942     theInfo.addProperty(zconsFixtime);
1943     }
1944    
1945 tim 1093 //set zconsUsingSMD
1946     IntData* zconsUsingSMD = new IntData();
1947     zconsUsingSMD->setID(ZCONSUSINGSMD_ID);
1948 tim 1091
1949 tim 1093 if (globals->haveZConsUsingSMD()){
1950     zconsUsingSMD->setData(globals->getZconsUsingSMD());
1951     theInfo.addProperty(zconsUsingSMD);
1952     }
1953    
1954 tim 701 //Determine the name of ouput file and add it into SimInfo's property list
1955     //Be careful, do not use inFileName, since it is a pointer which
1956     //point to a string at master node, and slave nodes do not contain that string
1957 tim 722
1958 tim 701 string zconsOutput(theInfo.finalName);
1959 tim 722
1960 tim 701 zconsOutput = zconsOutput.substr(0, zconsOutput.rfind(".")) + ".fz";
1961 tim 722
1962 tim 701 StringData* zconsFilename = new StringData();
1963     zconsFilename->setID(ZCONSFILENAME_ID);
1964     zconsFilename->setData(zconsOutput);
1965 tim 722
1966 tim 701 theInfo.addProperty(zconsFilename);
1967 tim 722
1968 tim 701 //setup index, pos and other parameters of z-constraint molecules
1969     nZConstraints = globals->getNzConstraints();
1970     theInfo.nZconstraints = nZConstraints;
1971    
1972     zconStamp = globals->getZconStamp();
1973     ZConsParaItem tempParaItem;
1974    
1975     ZConsParaData* zconsParaData = new ZConsParaData();
1976     zconsParaData->setID(ZCONSPARADATA_ID);
1977 tim 722
1978     for (int i = 0; i < nZConstraints; i++){
1979 tim 699 tempParaItem.havingZPos = zconStamp[i]->haveZpos();
1980     tempParaItem.zPos = zconStamp[i]->getZpos();
1981     tempParaItem.zconsIndex = zconStamp[i]->getMolIndex();
1982     tempParaItem.kRatio = zconStamp[i]->getKratio();
1983 tim 1093 tempParaItem.havingCantVel = zconStamp[i]->haveCantVel();
1984     tempParaItem.cantVel = zconStamp[i]->getCantVel();
1985 tim 699 zconsParaData->addItem(tempParaItem);
1986 tim 701 }
1987 tim 699
1988 tim 736 //check the uniqueness of index
1989     if(!zconsParaData->isIndexUnique()){
1990     sprintf(painCave.errMsg,
1991 gezelter 965 "ZConstraint Error: molIndex is not unique!\n");
1992 tim 736 painCave.isFatal = 1;
1993     simError();
1994     }
1995    
1996 tim 701 //sort the parameters by index of molecules
1997     zconsParaData->sortByIndex();
1998 tim 736
1999 tim 701 //push data into siminfo, therefore, we can retrieve later
2000     theInfo.addProperty(zconsParaData);
2001 tim 660 }
2002 tim 1031
2003     void SimSetup::makeMinimizer(){
2004 tim 1032
2005 tim 1064 OOPSEMinimizer* myOOPSEMinimizer;
2006 tim 1031 MinimizerParameterSet* param;
2007 tim 1064 char minimizerName[100];
2008 tim 1031
2009     for (int i = 0; i < nInfo; i++){
2010 tim 1064
2011 tim 1031 //prepare parameter set for minimizer
2012     param = new MinimizerParameterSet();
2013     param->setDefaultParameter();
2014    
2015     if (globals->haveMinimizer()){
2016     param->setFTol(globals->getMinFTol());
2017     }
2018    
2019     if (globals->haveMinGTol()){
2020     param->setGTol(globals->getMinGTol());
2021     }
2022    
2023     if (globals->haveMinMaxIter()){
2024     param->setMaxIteration(globals->getMinMaxIter());
2025     }
2026    
2027     if (globals->haveMinWriteFrq()){
2028     param->setMaxIteration(globals->getMinMaxIter());
2029     }
2030    
2031     if (globals->haveMinWriteFrq()){
2032     param->setWriteFrq(globals->getMinWriteFrq());
2033     }
2034    
2035 tim 1064 if (globals->haveMinStepSize()){
2036     param->setStepSize(globals->getMinStepSize());
2037 tim 1031 }
2038    
2039     if (globals->haveMinLSMaxIter()){
2040     param->setLineSearchMaxIteration(globals->getMinLSMaxIter());
2041     }
2042    
2043     if (globals->haveMinLSTol()){
2044     param->setLineSearchTol(globals->getMinLSTol());
2045     }
2046    
2047 tim 1066 strcpy(minimizerName, globals->getMinimizer());
2048 tim 1064
2049     if (!strcasecmp(minimizerName, "CG")){
2050     myOOPSEMinimizer = new PRCGMinimizer(&(info[i]), the_ff, param);
2051     }
2052     else if (!strcasecmp(minimizerName, "SD")){
2053     //myOOPSEMinimizer = MinimizerFactory.creatMinimizer("", &(info[i]), the_ff, param);
2054     myOOPSEMinimizer = new SDMinimizer(&(info[i]), the_ff, param);
2055     }
2056     else{
2057 tim 1066 sprintf(painCave.errMsg,
2058     "SimSetup error: Unrecognized Minimizer, use Conjugate Gradient \n");
2059     painCave.isFatal = 0;
2060     simError();
2061    
2062     myOOPSEMinimizer = new PRCGMinimizer(&(info[i]), the_ff, param);
2063     }
2064 tim 1064 info[i].the_integrator = myOOPSEMinimizer;
2065    
2066 tim 1031 //store the minimizer into simInfo
2067 tim 1064 info[i].the_minimizer = myOOPSEMinimizer;
2068 tim 1031 info[i].has_minimizer = true;
2069     }
2070 tim 1032
2071 tim 1031 }