ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 812
Committed: Wed Oct 22 21:17:32 2003 UTC (20 years, 8 months ago) by mmeineke
File size: 43816 byte(s)
Log Message:
added a new NPT integrator, NPTxyz. It scales the x, y, and z direction sepeartely. no box skew allowed.

File Contents

# User Rev Content
1 tim 658 #include <algorithm>
2 mmeineke 377 #include <cstdlib>
3     #include <iostream>
4     #include <cmath>
5 tim 658 #include <string>
6 tim 722 #include <sprng.h>
7 mmeineke 377
8     #include "SimSetup.hpp"
9 tim 689 #include "ReadWrite.hpp"
10 mmeineke 377 #include "parse_me.h"
11     #include "Integrator.hpp"
12     #include "simError.h"
13    
14     #ifdef IS_MPI
15     #include "mpiBASS.h"
16     #include "mpiSimulation.hpp"
17     #endif
18    
19 mmeineke 557 // some defines for ensemble and Forcefield cases
20    
21 tim 660 #define NVE_ENS 0
22     #define NVT_ENS 1
23     #define NPTi_ENS 2
24     #define NPTf_ENS 3
25 mmeineke 812 #define NPTxyz_ENS 4
26 mmeineke 557
27 mmeineke 812
28 mmeineke 557 #define FF_DUFF 0
29     #define FF_LJ 1
30 chuckv 653 #define FF_EAM 2
31 mmeineke 557
32 tim 658 using namespace std;
33    
34 mmeineke 377 SimSetup::SimSetup(){
35 mmeineke 656 isInfoArray = 0;
36     nInfo = 1;
37 tim 722
38 mmeineke 377 stamps = new MakeStamps();
39     globals = new Globals();
40 tim 722
41    
42 mmeineke 377 #ifdef IS_MPI
43 tim 722 strcpy(checkPointMsg, "SimSetup creation successful");
44 mmeineke 377 MPIcheckPoint();
45     #endif // IS_MPI
46     }
47    
48     SimSetup::~SimSetup(){
49     delete stamps;
50     delete globals;
51     }
52    
53 tim 722 void SimSetup::setSimInfo(SimInfo* the_info, int theNinfo){
54     info = the_info;
55     nInfo = theNinfo;
56     isInfoArray = 1;
57 mmeineke 670 }
58 mmeineke 656
59    
60 tim 722 void SimSetup::parseFile(char* fileName){
61 mmeineke 377 #ifdef IS_MPI
62 tim 722 if (worldRank == 0){
63 mmeineke 377 #endif // is_mpi
64 tim 722
65 mmeineke 377 inFileName = fileName;
66 tim 722 set_interface_stamps(stamps, globals);
67    
68 mmeineke 377 #ifdef IS_MPI
69     mpiEventInit();
70     #endif
71    
72 tim 722 yacc_BASS(fileName);
73 mmeineke 377
74     #ifdef IS_MPI
75     throwMPIEvent(NULL);
76     }
77 tim 722 else{
78     receiveParse();
79     }
80 mmeineke 377 #endif
81    
82     }
83    
84     #ifdef IS_MPI
85     void SimSetup::receiveParse(void){
86 tim 722 set_interface_stamps(stamps, globals);
87     mpiEventInit();
88     MPIcheckPoint();
89     mpiEventLoop();
90 mmeineke 377 }
91    
92     #endif // is_mpi
93    
94 mmeineke 670 void SimSetup::createSim(void){
95 mmeineke 377
96 mmeineke 614 // gather all of the information from the Bass file
97 tim 689
98 mmeineke 614 gatherInfo();
99 mmeineke 377
100 mmeineke 614 // creation of complex system objects
101 mmeineke 377
102 mmeineke 614 sysObjectsCreation();
103 mmeineke 377
104 mmeineke 614 // check on the post processing info
105 tim 689
106 mmeineke 614 finalInfoCheck();
107 mmeineke 377
108 mmeineke 614 // initialize the system coordinates
109 mmeineke 377
110 tim 722 if (!isInfoArray){
111     initSystemCoords();
112 mmeineke 811
113     if( !(globals->getUseInitTime()) )
114     info[0].currentTime = 0.0;
115 tim 722 }
116 mmeineke 377
117 mmeineke 614 // make the output filenames
118 mmeineke 377
119 mmeineke 614 makeOutNames();
120 tim 722
121 mmeineke 555 // make the integrator
122 tim 722
123 mmeineke 616 makeIntegrator();
124 tim 722
125 chuckv 432 #ifdef IS_MPI
126     mpiSim->mpiRefresh();
127     #endif
128 mmeineke 377
129 chuckv 432 // initialize the Fortran
130 mmeineke 377
131 mmeineke 616 initFortran();
132 mmeineke 377 }
133    
134 mmeineke 407
135 tim 722 void SimSetup::makeMolecules(void){
136 mmeineke 787 int k;
137 mmeineke 412 int i, j, exI, exJ, tempEx, stampID, atomOffset, excludeOffset;
138 mmeineke 616 molInit molInfo;
139 mmeineke 407 DirectionalAtom* dAtom;
140 mmeineke 412 LinkedAssign* extras;
141     LinkedAssign* current_extra;
142 mmeineke 407 AtomStamp* currentAtom;
143     BondStamp* currentBond;
144     BendStamp* currentBend;
145     TorsionStamp* currentTorsion;
146 mmeineke 427
147     bond_pair* theBonds;
148     bend_set* theBends;
149     torsion_set* theTorsions;
150    
151 tim 722
152 mmeineke 407 //init the forceField paramters
153    
154     the_ff->readParams();
155    
156 tim 722
157 mmeineke 427 // init the atoms
158 mmeineke 407
159 mmeineke 427 double ux, uy, uz, u, uSqr;
160 mmeineke 407
161 tim 722 for (k = 0; k < nInfo; k++){
162     the_ff->setSimInfo(&(info[k]));
163    
164 mmeineke 670 atomOffset = 0;
165     excludeOffset = 0;
166 tim 722 for (i = 0; i < info[k].n_mol; i++){
167 mmeineke 670 stampID = info[k].molecules[i].getStampID();
168 mmeineke 412
169 tim 722 molInfo.nAtoms = comp_stamps[stampID]->getNAtoms();
170     molInfo.nBonds = comp_stamps[stampID]->getNBonds();
171     molInfo.nBends = comp_stamps[stampID]->getNBends();
172 mmeineke 670 molInfo.nTorsions = comp_stamps[stampID]->getNTorsions();
173     molInfo.nExcludes = molInfo.nBonds + molInfo.nBends + molInfo.nTorsions;
174 tim 722
175 mmeineke 670 molInfo.myAtoms = &(info[k].atoms[atomOffset]);
176     molInfo.myExcludes = &(info[k].excludes[excludeOffset]);
177 tim 722 molInfo.myBonds = new Bond * [molInfo.nBonds];
178     molInfo.myBends = new Bend * [molInfo.nBends];
179     molInfo.myTorsions = new Torsion * [molInfo.nTorsions];
180 mmeineke 407
181 mmeineke 670 theBonds = new bond_pair[molInfo.nBonds];
182     theBends = new bend_set[molInfo.nBends];
183     theTorsions = new torsion_set[molInfo.nTorsions];
184 tim 722
185 mmeineke 670 // make the Atoms
186 tim 722
187     for (j = 0; j < molInfo.nAtoms; j++){
188     currentAtom = comp_stamps[stampID]->getAtom(j);
189     if (currentAtom->haveOrientation()){
190     dAtom = new DirectionalAtom((j + atomOffset),
191     info[k].getConfiguration());
192     info[k].n_oriented++;
193     molInfo.myAtoms[j] = dAtom;
194    
195     ux = currentAtom->getOrntX();
196     uy = currentAtom->getOrntY();
197     uz = currentAtom->getOrntZ();
198    
199     uSqr = (ux * ux) + (uy * uy) + (uz * uz);
200    
201     u = sqrt(uSqr);
202     ux = ux / u;
203     uy = uy / u;
204     uz = uz / u;
205    
206     dAtom->setSUx(ux);
207     dAtom->setSUy(uy);
208     dAtom->setSUz(uz);
209     }
210     else{
211     molInfo.myAtoms[j] = new GeneralAtom((j + atomOffset),
212     info[k].getConfiguration());
213     }
214     molInfo.myAtoms[j]->setType(currentAtom->getType());
215    
216 mmeineke 407 #ifdef IS_MPI
217 tim 722
218     molInfo.myAtoms[j]->setGlobalIndex(globalIndex[j + atomOffset]);
219    
220 mmeineke 407 #endif // is_mpi
221 mmeineke 670 }
222 tim 722
223     // make the bonds
224     for (j = 0; j < molInfo.nBonds; j++){
225     currentBond = comp_stamps[stampID]->getBond(j);
226     theBonds[j].a = currentBond->getA() + atomOffset;
227     theBonds[j].b = currentBond->getB() + atomOffset;
228    
229     exI = theBonds[j].a;
230     exJ = theBonds[j].b;
231    
232     // exclude_I must always be the smaller of the pair
233     if (exI > exJ){
234     tempEx = exI;
235     exI = exJ;
236     exJ = tempEx;
237     }
238 mmeineke 670 #ifdef IS_MPI
239 tim 722 tempEx = exI;
240     exI = info[k].atoms[tempEx]->getGlobalIndex() + 1;
241     tempEx = exJ;
242     exJ = info[k].atoms[tempEx]->getGlobalIndex() + 1;
243    
244     info[k].excludes[j + excludeOffset]->setPair(exI, exJ);
245 mmeineke 412 #else // isn't MPI
246 tim 722
247     info[k].excludes[j + excludeOffset]->setPair((exI + 1), (exJ + 1));
248 mmeineke 412 #endif //is_mpi
249 mmeineke 670 }
250     excludeOffset += molInfo.nBonds;
251 tim 722
252 mmeineke 670 //make the bends
253 tim 722 for (j = 0; j < molInfo.nBends; j++){
254     currentBend = comp_stamps[stampID]->getBend(j);
255     theBends[j].a = currentBend->getA() + atomOffset;
256     theBends[j].b = currentBend->getB() + atomOffset;
257     theBends[j].c = currentBend->getC() + atomOffset;
258    
259     if (currentBend->haveExtras()){
260     extras = currentBend->getExtras();
261     current_extra = extras;
262    
263     while (current_extra != NULL){
264     if (!strcmp(current_extra->getlhs(), "ghostVectorSource")){
265     switch (current_extra->getType()){
266     case 0:
267     theBends[j].ghost = current_extra->getInt() + atomOffset;
268     theBends[j].isGhost = 1;
269     break;
270    
271     case 1:
272     theBends[j].ghost = (int) current_extra->getDouble() +
273     atomOffset;
274     theBends[j].isGhost = 1;
275     break;
276    
277     default:
278     sprintf(painCave.errMsg,
279     "SimSetup Error: ghostVectorSource was neither a "
280     "double nor an int.\n"
281     "-->Bend[%d] in %s\n",
282     j, comp_stamps[stampID]->getID());
283     painCave.isFatal = 1;
284     simError();
285     }
286     }
287     else{
288     sprintf(painCave.errMsg,
289     "SimSetup Error: unhandled bend assignment:\n"
290     " -->%s in Bend[%d] in %s\n",
291     current_extra->getlhs(), j, comp_stamps[stampID]->getID());
292     painCave.isFatal = 1;
293     simError();
294     }
295    
296     current_extra = current_extra->getNext();
297     }
298 tim 701 }
299 tim 722
300     if (!theBends[j].isGhost){
301     exI = theBends[j].a;
302     exJ = theBends[j].c;
303     }
304     else{
305     exI = theBends[j].a;
306     exJ = theBends[j].b;
307     }
308    
309     // exclude_I must always be the smaller of the pair
310     if (exI > exJ){
311     tempEx = exI;
312     exI = exJ;
313     exJ = tempEx;
314     }
315 mmeineke 670 #ifdef IS_MPI
316 tim 722 tempEx = exI;
317     exI = info[k].atoms[tempEx]->getGlobalIndex() + 1;
318     tempEx = exJ;
319     exJ = info[k].atoms[tempEx]->getGlobalIndex() + 1;
320    
321     info[k].excludes[j + excludeOffset]->setPair(exI, exJ);
322 mmeineke 670 #else // isn't MPI
323 tim 722 info[k].excludes[j + excludeOffset]->setPair((exI + 1), (exJ + 1));
324 mmeineke 670 #endif //is_mpi
325 mmeineke 412 }
326 mmeineke 670 excludeOffset += molInfo.nBends;
327 tim 722
328     for (j = 0; j < molInfo.nTorsions; j++){
329     currentTorsion = comp_stamps[stampID]->getTorsion(j);
330     theTorsions[j].a = currentTorsion->getA() + atomOffset;
331     theTorsions[j].b = currentTorsion->getB() + atomOffset;
332     theTorsions[j].c = currentTorsion->getC() + atomOffset;
333     theTorsions[j].d = currentTorsion->getD() + atomOffset;
334    
335     exI = theTorsions[j].a;
336     exJ = theTorsions[j].d;
337    
338     // exclude_I must always be the smaller of the pair
339     if (exI > exJ){
340     tempEx = exI;
341     exI = exJ;
342     exJ = tempEx;
343     }
344 mmeineke 670 #ifdef IS_MPI
345 tim 722 tempEx = exI;
346     exI = info[k].atoms[tempEx]->getGlobalIndex() + 1;
347     tempEx = exJ;
348     exJ = info[k].atoms[tempEx]->getGlobalIndex() + 1;
349    
350     info[k].excludes[j + excludeOffset]->setPair(exI, exJ);
351 mmeineke 407 #else // isn't MPI
352 tim 722 info[k].excludes[j + excludeOffset]->setPair((exI + 1), (exJ + 1));
353 mmeineke 412 #endif //is_mpi
354 mmeineke 670 }
355     excludeOffset += molInfo.nTorsions;
356 tim 722
357    
358 mmeineke 670 // send the arrays off to the forceField for init.
359 mmeineke 407
360 tim 722 the_ff->initializeAtoms(molInfo.nAtoms, molInfo.myAtoms);
361     the_ff->initializeBonds(molInfo.nBonds, molInfo.myBonds, theBonds);
362     the_ff->initializeBends(molInfo.nBends, molInfo.myBends, theBends);
363     the_ff->initializeTorsions(molInfo.nTorsions, molInfo.myTorsions,
364     theTorsions);
365    
366    
367     info[k].molecules[i].initialize(molInfo);
368    
369    
370 mmeineke 670 atomOffset += molInfo.nAtoms;
371     delete[] theBonds;
372     delete[] theBends;
373     delete[] theTorsions;
374 mmeineke 412 }
375 mmeineke 414 }
376 tim 722
377 chuckv 434 #ifdef IS_MPI
378 tim 722 sprintf(checkPointMsg, "all molecules initialized succesfully");
379 chuckv 434 MPIcheckPoint();
380     #endif // is_mpi
381 tim 722
382 mmeineke 670 // clean up the forcefield
383 chuckv 434
384 mmeineke 420 the_ff->calcRcut();
385 mmeineke 414 the_ff->cleanMe();
386     }
387 mmeineke 407
388 tim 722 void SimSetup::initFromBass(void){
389 mmeineke 377 int i, j, k;
390     int n_cells;
391     double cellx, celly, cellz;
392     double temp1, temp2, temp3;
393     int n_per_extra;
394     int n_extra;
395     int have_extra, done;
396    
397 mmeineke 670 double vel[3];
398     vel[0] = 0.0;
399     vel[1] = 0.0;
400     vel[2] = 0.0;
401    
402 tim 722 temp1 = (double) tot_nmol / 4.0;
403     temp2 = pow(temp1, (1.0 / 3.0));
404     temp3 = ceil(temp2);
405 mmeineke 377
406 tim 722 have_extra = 0;
407     if (temp2 < temp3){
408     // we have a non-complete lattice
409     have_extra = 1;
410 mmeineke 377
411 tim 722 n_cells = (int) temp3 - 1;
412 mmeineke 670 cellx = info[0].boxL[0] / temp3;
413     celly = info[0].boxL[1] / temp3;
414     cellz = info[0].boxL[2] / temp3;
415 tim 722 n_extra = tot_nmol - (4 * n_cells * n_cells * n_cells);
416     temp1 = ((double) n_extra) / (pow(temp3, 3.0) - pow(n_cells, 3.0));
417     n_per_extra = (int) ceil(temp1);
418 mmeineke 377
419 tim 722 if (n_per_extra > 4){
420     sprintf(painCave.errMsg,
421     "SimSetup error. There has been an error in constructing"
422     " the non-complete lattice.\n");
423 mmeineke 377 painCave.isFatal = 1;
424     simError();
425     }
426     }
427     else{
428 tim 722 n_cells = (int) temp3;
429 mmeineke 670 cellx = info[0].boxL[0] / temp3;
430     celly = info[0].boxL[1] / temp3;
431     cellz = info[0].boxL[2] / temp3;
432 mmeineke 377 }
433    
434     current_mol = 0;
435     current_comp_mol = 0;
436     current_comp = 0;
437     current_atom_ndx = 0;
438    
439 tim 722 for (i = 0; i < n_cells ; i++){
440     for (j = 0; j < n_cells; j++){
441     for (k = 0; k < n_cells; k++){
442     makeElement(i * cellx, j * celly, k * cellz);
443 mmeineke 377
444 tim 722 makeElement(i * cellx + 0.5 * cellx, j * celly + 0.5 * celly, k * cellz);
445 mmeineke 377
446 tim 722 makeElement(i * cellx, j * celly + 0.5 * celly, k * cellz + 0.5 * cellz);
447 mmeineke 377
448 tim 722 makeElement(i * cellx + 0.5 * cellx, j * celly, k * cellz + 0.5 * cellz);
449 mmeineke 377 }
450     }
451     }
452    
453 tim 722 if (have_extra){
454 mmeineke 377 done = 0;
455    
456     int start_ndx;
457 tim 722 for (i = 0; i < (n_cells + 1) && !done; i++){
458     for (j = 0; j < (n_cells + 1) && !done; j++){
459     if (i < n_cells){
460     if (j < n_cells){
461     start_ndx = n_cells;
462     }
463     else
464     start_ndx = 0;
465     }
466     else
467     start_ndx = 0;
468 mmeineke 377
469 tim 722 for (k = start_ndx; k < (n_cells + 1) && !done; k++){
470     makeElement(i * cellx, j * celly, k * cellz);
471     done = (current_mol >= tot_nmol);
472 mmeineke 377
473 tim 722 if (!done && n_per_extra > 1){
474     makeElement(i * cellx + 0.5 * cellx, j * celly + 0.5 * celly,
475     k * cellz);
476     done = (current_mol >= tot_nmol);
477     }
478 mmeineke 377
479 tim 722 if (!done && n_per_extra > 2){
480     makeElement(i * cellx, j * celly + 0.5 * celly,
481     k * cellz + 0.5 * cellz);
482     done = (current_mol >= tot_nmol);
483     }
484 mmeineke 377
485 tim 722 if (!done && n_per_extra > 3){
486     makeElement(i * cellx + 0.5 * cellx, j * celly,
487     k * cellz + 0.5 * cellz);
488     done = (current_mol >= tot_nmol);
489     }
490     }
491 mmeineke 377 }
492     }
493     }
494    
495 tim 722 for (i = 0; i < info[0].n_atoms; i++){
496     info[0].atoms[i]->setVel(vel);
497 mmeineke 377 }
498     }
499    
500 tim 722 void SimSetup::makeElement(double x, double y, double z){
501 mmeineke 377 int k;
502     AtomStamp* current_atom;
503     DirectionalAtom* dAtom;
504     double rotMat[3][3];
505 mmeineke 670 double pos[3];
506 mmeineke 377
507 tim 722 for (k = 0; k < comp_stamps[current_comp]->getNAtoms(); k++){
508     current_atom = comp_stamps[current_comp]->getAtom(k);
509     if (!current_atom->havePosition()){
510     sprintf(painCave.errMsg,
511     "SimSetup:initFromBass error.\n"
512     "\tComponent %s, atom %s does not have a position specified.\n"
513     "\tThe initialization routine is unable to give a start"
514     " position.\n",
515     comp_stamps[current_comp]->getID(), current_atom->getType());
516 mmeineke 377 painCave.isFatal = 1;
517     simError();
518     }
519 tim 722
520 mmeineke 670 pos[0] = x + current_atom->getPosX();
521     pos[1] = y + current_atom->getPosY();
522     pos[2] = z + current_atom->getPosZ();
523 mmeineke 377
524 tim 722 info[0].atoms[current_atom_ndx]->setPos(pos);
525 mmeineke 377
526 tim 722 if (info[0].atoms[current_atom_ndx]->isDirectional()){
527     dAtom = (DirectionalAtom *) info[0].atoms[current_atom_ndx];
528 mmeineke 377
529     rotMat[0][0] = 1.0;
530     rotMat[0][1] = 0.0;
531     rotMat[0][2] = 0.0;
532    
533     rotMat[1][0] = 0.0;
534     rotMat[1][1] = 1.0;
535     rotMat[1][2] = 0.0;
536    
537     rotMat[2][0] = 0.0;
538     rotMat[2][1] = 0.0;
539     rotMat[2][2] = 1.0;
540    
541 tim 722 dAtom->setA(rotMat);
542 mmeineke 377 }
543    
544     current_atom_ndx++;
545     }
546    
547     current_mol++;
548     current_comp_mol++;
549    
550 tim 722 if (current_comp_mol >= components_nmol[current_comp]){
551 mmeineke 377 current_comp_mol = 0;
552     current_comp++;
553     }
554     }
555 mmeineke 614
556    
557 tim 722 void SimSetup::gatherInfo(void){
558 mmeineke 787 int i;
559 mmeineke 614
560     ensembleCase = -1;
561     ffCase = -1;
562    
563 mmeineke 670 // set the easy ones first
564 mmeineke 614
565 tim 722 for (i = 0; i < nInfo; i++){
566 mmeineke 670 info[i].target_temp = globals->getTargetTemp();
567     info[i].dt = globals->getDt();
568     info[i].run_time = globals->getRunTime();
569     }
570 mmeineke 616 n_components = globals->getNComponents();
571 mmeineke 614
572    
573     // get the forceField
574    
575 tim 722 strcpy(force_field, globals->getForceField());
576 mmeineke 614
577 tim 722 if (!strcasecmp(force_field, "DUFF")){
578     ffCase = FF_DUFF;
579     }
580     else if (!strcasecmp(force_field, "LJ")){
581     ffCase = FF_LJ;
582     }
583     else if (!strcasecmp(force_field, "EAM")){
584     ffCase = FF_EAM;
585     }
586 mmeineke 614 else{
587 tim 722 sprintf(painCave.errMsg, "SimSetup Error. Unrecognized force field -> %s\n",
588     force_field);
589     painCave.isFatal = 1;
590     simError();
591 mmeineke 614 }
592    
593 tim 722 // get the ensemble
594 mmeineke 614
595 tim 722 strcpy(ensemble, globals->getEnsemble());
596 mmeineke 614
597 tim 722 if (!strcasecmp(ensemble, "NVE")){
598     ensembleCase = NVE_ENS;
599     }
600     else if (!strcasecmp(ensemble, "NVT")){
601     ensembleCase = NVT_ENS;
602     }
603     else if (!strcasecmp(ensemble, "NPTi") || !strcasecmp(ensemble, "NPT")){
604 mmeineke 614 ensembleCase = NPTi_ENS;
605 tim 722 }
606     else if (!strcasecmp(ensemble, "NPTf")){
607     ensembleCase = NPTf_ENS;
608     }
609 mmeineke 812 else if (!strcasecmp(ensemble, "NPTxyz")){
610     ensembleCase = NPTxyz_ENS;
611     }
612 mmeineke 614 else{
613 tim 722 sprintf(painCave.errMsg,
614     "SimSetup Warning. Unrecognized Ensemble -> %s, "
615     "reverting to NVE for this simulation.\n",
616     ensemble);
617     painCave.isFatal = 0;
618     simError();
619     strcpy(ensemble, "NVE");
620     ensembleCase = NVE_ENS;
621 mmeineke 614 }
622    
623 tim 722 for (i = 0; i < nInfo; i++){
624     strcpy(info[i].ensemble, ensemble);
625    
626 mmeineke 670 // get the mixing rule
627 mmeineke 614
628 tim 722 strcpy(info[i].mixingRule, globals->getMixingRule());
629 mmeineke 670 info[i].usePBC = globals->getPBC();
630     }
631 tim 722
632 mmeineke 614 // get the components and calculate the tot_nMol and indvidual n_mol
633 tim 722
634 mmeineke 616 the_components = globals->getComponents();
635 mmeineke 614 components_nmol = new int[n_components];
636    
637    
638 tim 722 if (!globals->haveNMol()){
639 mmeineke 614 // we don't have the total number of molecules, so we assume it is
640     // given in each component
641    
642     tot_nmol = 0;
643 tim 722 for (i = 0; i < n_components; i++){
644     if (!the_components[i]->haveNMol()){
645     // we have a problem
646     sprintf(painCave.errMsg,
647     "SimSetup Error. No global NMol or component NMol"
648     " given. Cannot calculate the number of atoms.\n");
649     painCave.isFatal = 1;
650     simError();
651 mmeineke 614 }
652    
653     tot_nmol += the_components[i]->getNMol();
654     components_nmol[i] = the_components[i]->getNMol();
655     }
656     }
657     else{
658 tim 722 sprintf(painCave.errMsg,
659     "SimSetup error.\n"
660     "\tSorry, the ability to specify total"
661     " nMols and then give molfractions in the components\n"
662     "\tis not currently supported."
663     " Please give nMol in the components.\n");
664 mmeineke 614 painCave.isFatal = 1;
665     simError();
666     }
667    
668     // set the status, sample, and thermal kick times
669    
670 tim 722 for (i = 0; i < nInfo; i++){
671     if (globals->haveSampleTime()){
672 mmeineke 670 info[i].sampleTime = globals->getSampleTime();
673     info[i].statusTime = info[i].sampleTime;
674     info[i].thermalTime = info[i].sampleTime;
675     }
676     else{
677     info[i].sampleTime = globals->getRunTime();
678     info[i].statusTime = info[i].sampleTime;
679     info[i].thermalTime = info[i].sampleTime;
680     }
681 tim 722
682     if (globals->haveStatusTime()){
683 mmeineke 670 info[i].statusTime = globals->getStatusTime();
684     }
685 tim 722
686     if (globals->haveThermalTime()){
687 mmeineke 670 info[i].thermalTime = globals->getThermalTime();
688     }
689 mmeineke 614
690 mmeineke 746 info[i].resetIntegrator = 0;
691     if( globals->haveResetTime() ){
692     info[i].resetTime = globals->getResetTime();
693     info[i].resetIntegrator = 1;
694     }
695    
696 mmeineke 670 // check for the temperature set flag
697 mmeineke 614
698 tim 722 if (globals->haveTempSet())
699     info[i].setTemp = globals->getTempSet();
700    
701 mmeineke 670 // get some of the tricky things that may still be in the globals
702 tim 722
703 mmeineke 670 double boxVector[3];
704 tim 722 if (globals->haveBox()){
705 mmeineke 670 boxVector[0] = globals->getBox();
706     boxVector[1] = globals->getBox();
707     boxVector[2] = globals->getBox();
708 tim 722
709     info[i].setBox(boxVector);
710 mmeineke 670 }
711 tim 722 else if (globals->haveDensity()){
712 mmeineke 670 double vol;
713 tim 722 vol = (double) tot_nmol / globals->getDensity();
714     boxVector[0] = pow(vol, (1.0 / 3.0));
715 mmeineke 670 boxVector[1] = boxVector[0];
716     boxVector[2] = boxVector[0];
717 tim 722
718     info[i].setBox(boxVector);
719     }
720 mmeineke 670 else{
721 tim 722 if (!globals->haveBoxX()){
722     sprintf(painCave.errMsg,
723     "SimSetup error, no periodic BoxX size given.\n");
724     painCave.isFatal = 1;
725     simError();
726 mmeineke 670 }
727     boxVector[0] = globals->getBoxX();
728 tim 722
729     if (!globals->haveBoxY()){
730     sprintf(painCave.errMsg,
731     "SimSetup error, no periodic BoxY size given.\n");
732     painCave.isFatal = 1;
733     simError();
734 mmeineke 670 }
735     boxVector[1] = globals->getBoxY();
736 tim 722
737     if (!globals->haveBoxZ()){
738     sprintf(painCave.errMsg,
739     "SimSetup error, no periodic BoxZ size given.\n");
740     painCave.isFatal = 1;
741     simError();
742 mmeineke 670 }
743     boxVector[2] = globals->getBoxZ();
744 tim 722
745     info[i].setBox(boxVector);
746 mmeineke 614 }
747 tim 708 }
748 mmeineke 614
749 tim 722 //setup seed for random number generator
750 tim 708 int seedValue;
751    
752 tim 722 if (globals->haveSeed()){
753 tim 708 seedValue = globals->getSeed();
754 tim 722
755     if(seedValue / 1E9 == 0){
756     sprintf(painCave.errMsg,
757     "Seed for sprng library should contain at least 9 digits\n"
758     "OOPSE will generate a seed for user\n");
759     painCave.isFatal = 0;
760     simError();
761    
762     //using seed generated by system instead of invalid seed set by user
763     #ifndef IS_MPI
764     seedValue = make_sprng_seed();
765     #else
766     if (worldRank == 0){
767     seedValue = make_sprng_seed();
768     }
769     MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD);
770     #endif
771     }
772     }//end of if branch of globals->haveSeed()
773 tim 708 else{
774 tim 722
775 tim 708 #ifndef IS_MPI
776 tim 722 seedValue = make_sprng_seed();
777 tim 708 #else
778 tim 722 if (worldRank == 0){
779     seedValue = make_sprng_seed();
780 tim 708 }
781 tim 722 MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD);
782 tim 708 #endif
783 tim 722 }//end of globals->haveSeed()
784 tim 708
785 tim 722 for (int i = 0; i < nInfo; i++){
786 tim 708 info[i].setSeed(seedValue);
787     }
788    
789 mmeineke 614 #ifdef IS_MPI
790 tim 722 strcpy(checkPointMsg, "Succesfully gathered all information from Bass\n");
791 mmeineke 614 MPIcheckPoint();
792     #endif // is_mpi
793     }
794    
795    
796 tim 722 void SimSetup::finalInfoCheck(void){
797 mmeineke 614 int index;
798     int usesDipoles;
799 mmeineke 670 int i;
800 mmeineke 614
801 tim 722 for (i = 0; i < nInfo; i++){
802 mmeineke 670 // check electrostatic parameters
803 tim 722
804 mmeineke 670 index = 0;
805     usesDipoles = 0;
806 tim 722 while ((index < info[i].n_atoms) && !usesDipoles){
807 mmeineke 670 usesDipoles = (info[i].atoms[index])->hasDipole();
808     index++;
809     }
810 tim 722
811 mmeineke 614 #ifdef IS_MPI
812 mmeineke 670 int myUse = usesDipoles;
813 tim 722 MPI_Allreduce(&myUse, &usesDipoles, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
814 mmeineke 614 #endif //is_mpi
815 tim 722
816 mmeineke 670 double theEcr, theEst;
817 tim 722
818     if (globals->getUseRF()){
819 mmeineke 670 info[i].useReactionField = 1;
820 tim 722
821     if (!globals->haveECR()){
822     sprintf(painCave.errMsg,
823     "SimSetup Warning: using default value of 1/2 the smallest "
824     "box length for the electrostaticCutoffRadius.\n"
825     "I hope you have a very fast processor!\n");
826     painCave.isFatal = 0;
827     simError();
828     double smallest;
829     smallest = info[i].boxL[0];
830     if (info[i].boxL[1] <= smallest)
831     smallest = info[i].boxL[1];
832     if (info[i].boxL[2] <= smallest)
833     smallest = info[i].boxL[2];
834     theEcr = 0.5 * smallest;
835 mmeineke 614 }
836 tim 722 else{
837     theEcr = globals->getECR();
838 mmeineke 614 }
839 tim 722
840     if (!globals->haveEST()){
841     sprintf(painCave.errMsg,
842     "SimSetup Warning: using default value of 0.05 * the "
843     "electrostaticCutoffRadius for the electrostaticSkinThickness\n");
844     painCave.isFatal = 0;
845     simError();
846     theEst = 0.05 * theEcr;
847 mmeineke 670 }
848 tim 722 else{
849     theEst = globals->getEST();
850 mmeineke 670 }
851 tim 722
852     info[i].setEcr(theEcr, theEst);
853    
854     if (!globals->haveDielectric()){
855     sprintf(painCave.errMsg,
856     "SimSetup Error: You are trying to use Reaction Field without"
857     "setting a dielectric constant!\n");
858     painCave.isFatal = 1;
859     simError();
860     }
861     info[i].dielectric = globals->getDielectric();
862     }
863     else{
864     if (usesDipoles){
865     if (!globals->haveECR()){
866     sprintf(painCave.errMsg,
867     "SimSetup Warning: using default value of 1/2 the smallest "
868     "box length for the electrostaticCutoffRadius.\n"
869     "I hope you have a very fast processor!\n");
870     painCave.isFatal = 0;
871     simError();
872     double smallest;
873     smallest = info[i].boxL[0];
874     if (info[i].boxL[1] <= smallest)
875     smallest = info[i].boxL[1];
876     if (info[i].boxL[2] <= smallest)
877     smallest = info[i].boxL[2];
878     theEcr = 0.5 * smallest;
879     }
880     else{
881     theEcr = globals->getECR();
882     }
883    
884     if (!globals->haveEST()){
885     sprintf(painCave.errMsg,
886     "SimSetup Warning: using default value of 0.05 * the "
887     "electrostaticCutoffRadius for the "
888     "electrostaticSkinThickness\n");
889     painCave.isFatal = 0;
890     simError();
891     theEst = 0.05 * theEcr;
892     }
893     else{
894     theEst = globals->getEST();
895     }
896    
897     info[i].setEcr(theEcr, theEst);
898     }
899     }
900 mmeineke 670 }
901 mmeineke 626
902 mmeineke 614 #ifdef IS_MPI
903 tim 722 strcpy(checkPointMsg, "post processing checks out");
904 mmeineke 614 MPIcheckPoint();
905     #endif // is_mpi
906     }
907    
908 tim 722 void SimSetup::initSystemCoords(void){
909 mmeineke 670 int i;
910 tim 722
911 tim 689 char* inName;
912    
913 tim 722 (info[0].getConfiguration())->createArrays(info[0].n_atoms);
914 mmeineke 614
915 tim 722 for (i = 0; i < info[0].n_atoms; i++)
916     info[0].atoms[i]->setCoords();
917    
918     if (globals->haveInitialConfig()){
919 mmeineke 670 InitializeFromFile* fileInit;
920 mmeineke 614 #ifdef IS_MPI // is_mpi
921 tim 722 if (worldRank == 0){
922 mmeineke 614 #endif //is_mpi
923 tim 689 inName = globals->getInitialConfig();
924 tim 722 fileInit = new InitializeFromFile(inName);
925 mmeineke 614 #ifdef IS_MPI
926 tim 722 }
927     else
928     fileInit = new InitializeFromFile(NULL);
929 mmeineke 614 #endif
930 tim 722 fileInit->readInit(info); // default velocities on
931    
932 mmeineke 670 delete fileInit;
933     }
934     else{
935 mmeineke 614 #ifdef IS_MPI
936 tim 722
937 mmeineke 670 // no init from bass
938 tim 722
939     sprintf(painCave.errMsg,
940     "Cannot intialize a parallel simulation without an initial configuration file.\n");
941 mmeineke 787 painCave.isFatal = 1;;
942 mmeineke 670 simError();
943 tim 722
944 mmeineke 614 #else
945 tim 722
946 mmeineke 670 initFromBass();
947 tim 722
948    
949 mmeineke 614 #endif
950 mmeineke 670 }
951 tim 722
952 mmeineke 614 #ifdef IS_MPI
953 tim 722 strcpy(checkPointMsg, "Successfully read in the initial configuration");
954 mmeineke 614 MPIcheckPoint();
955     #endif // is_mpi
956     }
957    
958    
959 tim 722 void SimSetup::makeOutNames(void){
960 mmeineke 670 int k;
961 mmeineke 614
962 mmeineke 670
963 tim 722 for (k = 0; k < nInfo; k++){
964 mmeineke 614 #ifdef IS_MPI
965 tim 722 if (worldRank == 0){
966 mmeineke 614 #endif // is_mpi
967 tim 722
968     if (globals->haveFinalConfig()){
969     strcpy(info[k].finalName, globals->getFinalConfig());
970 mmeineke 614 }
971     else{
972 tim 722 strcpy(info[k].finalName, inFileName);
973     char* endTest;
974     int nameLength = strlen(info[k].finalName);
975     endTest = &(info[k].finalName[nameLength - 5]);
976     if (!strcmp(endTest, ".bass")){
977     strcpy(endTest, ".eor");
978     }
979     else if (!strcmp(endTest, ".BASS")){
980     strcpy(endTest, ".eor");
981     }
982     else{
983     endTest = &(info[k].finalName[nameLength - 4]);
984     if (!strcmp(endTest, ".bss")){
985     strcpy(endTest, ".eor");
986     }
987     else if (!strcmp(endTest, ".mdl")){
988     strcpy(endTest, ".eor");
989     }
990     else{
991     strcat(info[k].finalName, ".eor");
992     }
993     }
994 mmeineke 614 }
995 tim 722
996 mmeineke 670 // make the sample and status out names
997 tim 722
998     strcpy(info[k].sampleName, inFileName);
999 mmeineke 670 char* endTest;
1000 tim 722 int nameLength = strlen(info[k].sampleName);
1001 mmeineke 670 endTest = &(info[k].sampleName[nameLength - 5]);
1002 tim 722 if (!strcmp(endTest, ".bass")){
1003     strcpy(endTest, ".dump");
1004 mmeineke 614 }
1005 tim 722 else if (!strcmp(endTest, ".BASS")){
1006     strcpy(endTest, ".dump");
1007 mmeineke 614 }
1008     else{
1009 tim 722 endTest = &(info[k].sampleName[nameLength - 4]);
1010     if (!strcmp(endTest, ".bss")){
1011     strcpy(endTest, ".dump");
1012     }
1013     else if (!strcmp(endTest, ".mdl")){
1014     strcpy(endTest, ".dump");
1015     }
1016     else{
1017     strcat(info[k].sampleName, ".dump");
1018     }
1019 mmeineke 614 }
1020 tim 722
1021     strcpy(info[k].statusName, inFileName);
1022     nameLength = strlen(info[k].statusName);
1023 mmeineke 670 endTest = &(info[k].statusName[nameLength - 5]);
1024 tim 722 if (!strcmp(endTest, ".bass")){
1025     strcpy(endTest, ".stat");
1026 mmeineke 614 }
1027 tim 722 else if (!strcmp(endTest, ".BASS")){
1028     strcpy(endTest, ".stat");
1029 mmeineke 614 }
1030     else{
1031 tim 722 endTest = &(info[k].statusName[nameLength - 4]);
1032     if (!strcmp(endTest, ".bss")){
1033     strcpy(endTest, ".stat");
1034     }
1035     else if (!strcmp(endTest, ".mdl")){
1036     strcpy(endTest, ".stat");
1037     }
1038     else{
1039     strcat(info[k].statusName, ".stat");
1040     }
1041 mmeineke 614 }
1042 tim 722
1043 mmeineke 670 #ifdef IS_MPI
1044 tim 722
1045 mmeineke 614 }
1046 mmeineke 670 #endif // is_mpi
1047 mmeineke 614 }
1048     }
1049    
1050    
1051 tim 722 void SimSetup::sysObjectsCreation(void){
1052     int i, k;
1053    
1054 mmeineke 614 // create the forceField
1055 tim 689
1056 mmeineke 614 createFF();
1057    
1058     // extract componentList
1059    
1060     compList();
1061    
1062     // calc the number of atoms, bond, bends, and torsions
1063    
1064     calcSysValues();
1065    
1066     #ifdef IS_MPI
1067     // divide the molecules among the processors
1068 tim 722
1069 mmeineke 614 mpiMolDivide();
1070     #endif //is_mpi
1071 tim 722
1072 mmeineke 614 // create the atom and SRI arrays. Also initialize Molecule Stamp ID's
1073 tim 722
1074 mmeineke 614 makeSysArrays();
1075    
1076 mmeineke 616 // make and initialize the molecules (all but atomic coordinates)
1077 tim 722
1078 mmeineke 616 makeMolecules();
1079 tim 722
1080     for (k = 0; k < nInfo; k++){
1081 mmeineke 670 info[k].identArray = new int[info[k].n_atoms];
1082 tim 722 for (i = 0; i < info[k].n_atoms; i++){
1083 mmeineke 670 info[k].identArray[i] = info[k].atoms[i]->getIdent();
1084     }
1085 mmeineke 616 }
1086 mmeineke 614 }
1087    
1088    
1089 tim 722 void SimSetup::createFF(void){
1090     switch (ffCase){
1091     case FF_DUFF:
1092     the_ff = new DUFF();
1093     break;
1094 mmeineke 614
1095 tim 722 case FF_LJ:
1096     the_ff = new LJFF();
1097     break;
1098 mmeineke 614
1099 tim 722 case FF_EAM:
1100     the_ff = new EAM_FF();
1101     break;
1102 mmeineke 614
1103 tim 722 default:
1104     sprintf(painCave.errMsg,
1105     "SimSetup Error. Unrecognized force field in case statement.\n");
1106     painCave.isFatal = 1;
1107     simError();
1108 mmeineke 614 }
1109    
1110     #ifdef IS_MPI
1111 tim 722 strcpy(checkPointMsg, "ForceField creation successful");
1112 mmeineke 614 MPIcheckPoint();
1113     #endif // is_mpi
1114     }
1115    
1116    
1117 tim 722 void SimSetup::compList(void){
1118 mmeineke 616 int i;
1119 mmeineke 670 char* id;
1120     LinkedMolStamp* headStamp = new LinkedMolStamp();
1121     LinkedMolStamp* currentStamp = NULL;
1122 tim 722 comp_stamps = new MoleculeStamp * [n_components];
1123    
1124 mmeineke 614 // make an array of molecule stamps that match the components used.
1125     // also extract the used stamps out into a separate linked list
1126 tim 722
1127     for (i = 0; i < nInfo; i++){
1128 mmeineke 670 info[i].nComponents = n_components;
1129     info[i].componentsNmol = components_nmol;
1130     info[i].compStamps = comp_stamps;
1131     info[i].headStamp = headStamp;
1132     }
1133 mmeineke 614
1134    
1135 tim 722 for (i = 0; i < n_components; i++){
1136 mmeineke 614 id = the_components[i]->getType();
1137     comp_stamps[i] = NULL;
1138 tim 722
1139 mmeineke 614 // check to make sure the component isn't already in the list
1140    
1141 tim 722 comp_stamps[i] = headStamp->match(id);
1142     if (comp_stamps[i] == NULL){
1143 mmeineke 614 // extract the component from the list;
1144 tim 722
1145     currentStamp = stamps->extractMolStamp(id);
1146     if (currentStamp == NULL){
1147     sprintf(painCave.errMsg,
1148     "SimSetup error: Component \"%s\" was not found in the "
1149     "list of declared molecules\n",
1150     id);
1151     painCave.isFatal = 1;
1152     simError();
1153 mmeineke 614 }
1154 tim 722
1155     headStamp->add(currentStamp);
1156     comp_stamps[i] = headStamp->match(id);
1157 mmeineke 614 }
1158     }
1159    
1160     #ifdef IS_MPI
1161 tim 722 strcpy(checkPointMsg, "Component stamps successfully extracted\n");
1162 mmeineke 614 MPIcheckPoint();
1163     #endif // is_mpi
1164 tim 722 }
1165 mmeineke 614
1166 tim 722 void SimSetup::calcSysValues(void){
1167 mmeineke 787 int i;
1168 mmeineke 614
1169 tim 722 int* molMembershipArray;
1170 mmeineke 614
1171     tot_atoms = 0;
1172     tot_bonds = 0;
1173     tot_bends = 0;
1174     tot_torsions = 0;
1175 tim 722 for (i = 0; i < n_components; i++){
1176     tot_atoms += components_nmol[i] * comp_stamps[i]->getNAtoms();
1177     tot_bonds += components_nmol[i] * comp_stamps[i]->getNBonds();
1178     tot_bends += components_nmol[i] * comp_stamps[i]->getNBends();
1179 mmeineke 614 tot_torsions += components_nmol[i] * comp_stamps[i]->getNTorsions();
1180     }
1181 tim 722
1182 mmeineke 614 tot_SRI = tot_bonds + tot_bends + tot_torsions;
1183 mmeineke 670 molMembershipArray = new int[tot_atoms];
1184 tim 722
1185     for (i = 0; i < nInfo; i++){
1186 mmeineke 670 info[i].n_atoms = tot_atoms;
1187     info[i].n_bonds = tot_bonds;
1188     info[i].n_bends = tot_bends;
1189     info[i].n_torsions = tot_torsions;
1190     info[i].n_SRI = tot_SRI;
1191     info[i].n_mol = tot_nmol;
1192 tim 722
1193 mmeineke 670 info[i].molMembershipArray = molMembershipArray;
1194 tim 722 }
1195 mmeineke 614 }
1196    
1197     #ifdef IS_MPI
1198    
1199 tim 722 void SimSetup::mpiMolDivide(void){
1200 mmeineke 616 int i, j, k;
1201 mmeineke 614 int localMol, allMol;
1202     int local_atoms, local_bonds, local_bends, local_torsions, local_SRI;
1203    
1204 tim 722 mpiSim = new mpiSimulation(info);
1205    
1206 mmeineke 614 globalIndex = mpiSim->divideLabor();
1207    
1208     // set up the local variables
1209 tim 722
1210 mmeineke 614 mol2proc = mpiSim->getMolToProcMap();
1211     molCompType = mpiSim->getMolComponentType();
1212 tim 722
1213 mmeineke 614 allMol = 0;
1214     localMol = 0;
1215     local_atoms = 0;
1216     local_bonds = 0;
1217     local_bends = 0;
1218     local_torsions = 0;
1219     globalAtomIndex = 0;
1220    
1221    
1222 tim 722 for (i = 0; i < n_components; i++){
1223     for (j = 0; j < components_nmol[i]; j++){
1224     if (mol2proc[allMol] == worldRank){
1225     local_atoms += comp_stamps[i]->getNAtoms();
1226     local_bonds += comp_stamps[i]->getNBonds();
1227     local_bends += comp_stamps[i]->getNBends();
1228     local_torsions += comp_stamps[i]->getNTorsions();
1229     localMol++;
1230 mmeineke 614 }
1231 tim 722 for (k = 0; k < comp_stamps[i]->getNAtoms(); k++){
1232 mmeineke 670 info[0].molMembershipArray[globalAtomIndex] = allMol;
1233 mmeineke 614 globalAtomIndex++;
1234     }
1235    
1236 tim 722 allMol++;
1237 mmeineke 614 }
1238     }
1239     local_SRI = local_bonds + local_bends + local_torsions;
1240 tim 722
1241 mmeineke 670 info[0].n_atoms = mpiSim->getMyNlocal();
1242 tim 722
1243     if (local_atoms != info[0].n_atoms){
1244     sprintf(painCave.errMsg,
1245     "SimSetup error: mpiSim's localAtom (%d) and SimSetup's"
1246     " localAtom (%d) are not equal.\n",
1247     info[0].n_atoms, local_atoms);
1248 mmeineke 614 painCave.isFatal = 1;
1249     simError();
1250     }
1251    
1252 mmeineke 670 info[0].n_bonds = local_bonds;
1253     info[0].n_bends = local_bends;
1254     info[0].n_torsions = local_torsions;
1255     info[0].n_SRI = local_SRI;
1256     info[0].n_mol = localMol;
1257 mmeineke 614
1258 tim 722 strcpy(checkPointMsg, "Passed nlocal consistency check.");
1259 mmeineke 614 MPIcheckPoint();
1260     }
1261 tim 722
1262 mmeineke 614 #endif // is_mpi
1263    
1264    
1265 tim 722 void SimSetup::makeSysArrays(void){
1266 mmeineke 787
1267     #ifndef IS_MPI
1268     int k, j;
1269     #endif // is_mpi
1270     int i, l;
1271 mmeineke 614
1272 mmeineke 670 Atom** the_atoms;
1273     Molecule* the_molecules;
1274     Exclude** the_excludes;
1275 mmeineke 616
1276 tim 722
1277     for (l = 0; l < nInfo; l++){
1278 mmeineke 670 // create the atom and short range interaction arrays
1279 tim 722
1280     the_atoms = new Atom * [info[l].n_atoms];
1281 mmeineke 670 the_molecules = new Molecule[info[l].n_mol];
1282     int molIndex;
1283 mmeineke 614
1284 mmeineke 670 // initialize the molecule's stampID's
1285 tim 722
1286 mmeineke 670 #ifdef IS_MPI
1287 tim 722
1288    
1289 mmeineke 670 molIndex = 0;
1290 tim 722 for (i = 0; i < mpiSim->getTotNmol(); i++){
1291     if (mol2proc[i] == worldRank){
1292     the_molecules[molIndex].setStampID(molCompType[i]);
1293     the_molecules[molIndex].setMyIndex(molIndex);
1294     the_molecules[molIndex].setGlobalIndex(i);
1295     molIndex++;
1296 mmeineke 670 }
1297 mmeineke 614 }
1298 tim 722
1299 mmeineke 614 #else // is_mpi
1300 tim 722
1301 mmeineke 670 molIndex = 0;
1302     globalAtomIndex = 0;
1303 tim 722 for (i = 0; i < n_components; i++){
1304     for (j = 0; j < components_nmol[i]; j++){
1305     the_molecules[molIndex].setStampID(i);
1306     the_molecules[molIndex].setMyIndex(molIndex);
1307     the_molecules[molIndex].setGlobalIndex(molIndex);
1308     for (k = 0; k < comp_stamps[i]->getNAtoms(); k++){
1309     info[l].molMembershipArray[globalAtomIndex] = molIndex;
1310     globalAtomIndex++;
1311     }
1312     molIndex++;
1313 mmeineke 614 }
1314     }
1315 tim 722
1316    
1317 mmeineke 614 #endif // is_mpi
1318    
1319    
1320 tim 722 if (info[l].n_SRI){
1321 mmeineke 670 Exclude::createArray(info[l].n_SRI);
1322 tim 722 the_excludes = new Exclude * [info[l].n_SRI];
1323     for (int ex = 0; ex < info[l].n_SRI; ex++){
1324     the_excludes[ex] = new Exclude(ex);
1325 mmeineke 670 }
1326     info[l].globalExcludes = new int;
1327     info[l].n_exclude = info[l].n_SRI;
1328     }
1329     else{
1330 tim 722 Exclude::createArray(1);
1331     the_excludes = new Exclude * ;
1332 mmeineke 670 the_excludes[0] = new Exclude(0);
1333 tim 722 the_excludes[0]->setPair(0, 0);
1334 mmeineke 670 info[l].globalExcludes = new int;
1335     info[l].globalExcludes[0] = 0;
1336     info[l].n_exclude = 0;
1337     }
1338 mmeineke 614
1339 mmeineke 670 // set the arrays into the SimInfo object
1340 mmeineke 614
1341 mmeineke 670 info[l].atoms = the_atoms;
1342     info[l].molecules = the_molecules;
1343     info[l].nGlobalExcludes = 0;
1344     info[l].excludes = the_excludes;
1345 mmeineke 614
1346 tim 722 the_ff->setSimInfo(info);
1347 mmeineke 670 }
1348 mmeineke 614 }
1349 mmeineke 616
1350 tim 722 void SimSetup::makeIntegrator(void){
1351 mmeineke 670 int k;
1352    
1353 mmeineke 782 NVE<RealIntegrator>* myNVE = NULL;
1354 tim 722 NVT<RealIntegrator>* myNVT = NULL;
1355 mmeineke 778 NPTi<NPT<RealIntegrator> >* myNPTi = NULL;
1356 mmeineke 780 NPTf<NPT<RealIntegrator> >* myNPTf = NULL;
1357 mmeineke 812 NPTxyz<NPT<RealIntegrator> >* myNPTxyz = NULL;
1358 tim 733
1359 tim 722 for (k = 0; k < nInfo; k++){
1360     switch (ensembleCase){
1361     case NVE_ENS:
1362     if (globals->haveZconstraints()){
1363     setupZConstraint(info[k]);
1364 mmeineke 782 myNVE = new ZConstraint<NVE<RealIntegrator> >(&(info[k]), the_ff);
1365 tim 722 }
1366 mmeineke 782 else{
1367     myNVE = new NVE<RealIntegrator>(&(info[k]), the_ff);
1368     }
1369    
1370     info->the_integrator = myNVE;
1371 tim 722 break;
1372 tim 676
1373 tim 722 case NVT_ENS:
1374     if (globals->haveZconstraints()){
1375     setupZConstraint(info[k]);
1376     myNVT = new ZConstraint<NVT<RealIntegrator> >(&(info[k]), the_ff);
1377     }
1378     else
1379     myNVT = new NVT<RealIntegrator>(&(info[k]), the_ff);
1380    
1381 tim 701 myNVT->setTargetTemp(globals->getTargetTemp());
1382 tim 722
1383     if (globals->haveTauThermostat())
1384 tim 701 myNVT->setTauThermostat(globals->getTauThermostat());
1385 tim 722 else{
1386     sprintf(painCave.errMsg,
1387     "SimSetup error: If you use the NVT\n"
1388     " ensemble, you must set tauThermostat.\n");
1389 tim 701 painCave.isFatal = 1;
1390     simError();
1391     }
1392 mmeineke 782
1393     info->the_integrator = myNVT;
1394 tim 701 break;
1395 tim 676
1396 tim 722 case NPTi_ENS:
1397     if (globals->haveZconstraints()){
1398     setupZConstraint(info[k]);
1399 mmeineke 778 myNPTi = new ZConstraint<NPTi<NPT <RealIntegrator> > >(&(info[k]), the_ff);
1400 tim 722 }
1401     else
1402 mmeineke 778 myNPTi = new NPTi<NPT<RealIntegrator> >(&(info[k]), the_ff);
1403 tim 722
1404     myNPTi->setTargetTemp(globals->getTargetTemp());
1405    
1406     if (globals->haveTargetPressure())
1407     myNPTi->setTargetPressure(globals->getTargetPressure());
1408     else{
1409     sprintf(painCave.errMsg,
1410     "SimSetup error: If you use a constant pressure\n"
1411     " ensemble, you must set targetPressure in the BASS file.\n");
1412     painCave.isFatal = 1;
1413     simError();
1414     }
1415    
1416     if (globals->haveTauThermostat())
1417     myNPTi->setTauThermostat(globals->getTauThermostat());
1418     else{
1419     sprintf(painCave.errMsg,
1420     "SimSetup error: If you use an NPT\n"
1421 tim 701 " ensemble, you must set tauThermostat.\n");
1422 tim 722 painCave.isFatal = 1;
1423     simError();
1424     }
1425    
1426     if (globals->haveTauBarostat())
1427     myNPTi->setTauBarostat(globals->getTauBarostat());
1428     else{
1429     sprintf(painCave.errMsg,
1430 tim 701 "SimSetup error: If you use an NPT\n"
1431     " ensemble, you must set tauBarostat.\n");
1432 tim 722 painCave.isFatal = 1;
1433     simError();
1434     }
1435 mmeineke 782
1436     info->the_integrator = myNPTi;
1437 tim 722 break;
1438 tim 676
1439 tim 722 case NPTf_ENS:
1440     if (globals->haveZconstraints()){
1441     setupZConstraint(info[k]);
1442 mmeineke 780 myNPTf = new ZConstraint<NPTf<NPT <RealIntegrator> > >(&(info[k]), the_ff);
1443 tim 722 }
1444     else
1445 mmeineke 780 myNPTf = new NPTf<NPT <RealIntegrator> >(&(info[k]), the_ff);
1446 tim 722
1447     myNPTf->setTargetTemp(globals->getTargetTemp());
1448    
1449     if (globals->haveTargetPressure())
1450     myNPTf->setTargetPressure(globals->getTargetPressure());
1451     else{
1452     sprintf(painCave.errMsg,
1453 tim 701 "SimSetup error: If you use a constant pressure\n"
1454     " ensemble, you must set targetPressure in the BASS file.\n");
1455 tim 722 painCave.isFatal = 1;
1456     simError();
1457     }
1458    
1459     if (globals->haveTauThermostat())
1460     myNPTf->setTauThermostat(globals->getTauThermostat());
1461     else{
1462     sprintf(painCave.errMsg,
1463 tim 701 "SimSetup error: If you use an NPT\n"
1464 tim 722 " ensemble, you must set tauThermostat.\n");
1465     painCave.isFatal = 1;
1466     simError();
1467     }
1468    
1469     if (globals->haveTauBarostat())
1470     myNPTf->setTauBarostat(globals->getTauBarostat());
1471     else{
1472     sprintf(painCave.errMsg,
1473     "SimSetup error: If you use an NPT\n"
1474 tim 701 " ensemble, you must set tauBarostat.\n");
1475 tim 722 painCave.isFatal = 1;
1476     simError();
1477     }
1478 mmeineke 782
1479     info->the_integrator = myNPTf;
1480 tim 722 break;
1481 tim 676
1482 mmeineke 812 case NPTxyz_ENS:
1483     if (globals->haveZconstraints()){
1484     setupZConstraint(info[k]);
1485     myNPTxyz = new ZConstraint<NPTxyz<NPT <RealIntegrator> > >(&(info[k]), the_ff);
1486     }
1487     else
1488     myNPTxyz = new NPTxyz<NPT <RealIntegrator> >(&(info[k]), the_ff);
1489    
1490     myNPTxyz->setTargetTemp(globals->getTargetTemp());
1491    
1492     if (globals->haveTargetPressure())
1493     myNPTxyz->setTargetPressure(globals->getTargetPressure());
1494     else{
1495     sprintf(painCave.errMsg,
1496     "SimSetup error: If you use a constant pressure\n"
1497     " ensemble, you must set targetPressure in the BASS file.\n");
1498     painCave.isFatal = 1;
1499     simError();
1500     }
1501    
1502     if (globals->haveTauThermostat())
1503     myNPTxyz->setTauThermostat(globals->getTauThermostat());
1504     else{
1505     sprintf(painCave.errMsg,
1506     "SimSetup error: If you use an NPT\n"
1507     " ensemble, you must set tauThermostat.\n");
1508     painCave.isFatal = 1;
1509     simError();
1510     }
1511    
1512     if (globals->haveTauBarostat())
1513     myNPTxyz->setTauBarostat(globals->getTauBarostat());
1514     else{
1515     sprintf(painCave.errMsg,
1516     "SimSetup error: If you use an NPT\n"
1517     " ensemble, you must set tauBarostat.\n");
1518     painCave.isFatal = 1;
1519     simError();
1520     }
1521    
1522     info->the_integrator = myNPTxyz;
1523     break;
1524    
1525 tim 722 default:
1526     sprintf(painCave.errMsg,
1527     "SimSetup Error. Unrecognized ensemble in case statement.\n");
1528 tim 701 painCave.isFatal = 1;
1529     simError();
1530 tim 660 }
1531 mmeineke 616 }
1532     }
1533    
1534 tim 722 void SimSetup::initFortran(void){
1535     info[0].refreshSim();
1536 mmeineke 616
1537 tim 722 if (!strcmp(info[0].mixingRule, "standard")){
1538     the_ff->initForceField(LB_MIXING_RULE);
1539 mmeineke 616 }
1540 tim 722 else if (!strcmp(info[0].mixingRule, "explicit")){
1541     the_ff->initForceField(EXPLICIT_MIXING_RULE);
1542 mmeineke 616 }
1543     else{
1544 tim 722 sprintf(painCave.errMsg, "SimSetup Error: unknown mixing rule -> \"%s\"\n",
1545     info[0].mixingRule);
1546 mmeineke 616 painCave.isFatal = 1;
1547     simError();
1548     }
1549    
1550    
1551     #ifdef IS_MPI
1552 tim 722 strcpy(checkPointMsg, "Successfully intialized the mixingRule for Fortran.");
1553 mmeineke 616 MPIcheckPoint();
1554     #endif // is_mpi
1555     }
1556 tim 660
1557 tim 722 void SimSetup::setupZConstraint(SimInfo& theInfo){
1558 tim 701 int nZConstraints;
1559     ZconStamp** zconStamp;
1560 tim 682
1561 tim 722 if (globals->haveZconstraintTime()){
1562 tim 701 //add sample time of z-constraint into SimInfo's property list
1563     DoubleData* zconsTimeProp = new DoubleData();
1564     zconsTimeProp->setID(ZCONSTIME_ID);
1565     zconsTimeProp->setData(globals->getZconsTime());
1566     theInfo.addProperty(zconsTimeProp);
1567     }
1568     else{
1569 tim 722 sprintf(painCave.errMsg,
1570     "ZConstraint error: If you use an ZConstraint\n"
1571     " , you must set sample time.\n");
1572 tim 701 painCave.isFatal = 1;
1573 tim 722 simError();
1574 tim 701 }
1575 tim 682
1576 tim 701 //push zconsTol into siminfo, if user does not specify
1577     //value for zconsTol, a default value will be used
1578     DoubleData* zconsTol = new DoubleData();
1579     zconsTol->setID(ZCONSTOL_ID);
1580 tim 722 if (globals->haveZconsTol()){
1581 tim 701 zconsTol->setData(globals->getZconsTol());
1582     }
1583     else{
1584 tim 722 double defaultZConsTol = 0.01;
1585     sprintf(painCave.errMsg,
1586     "ZConstraint Waring: Tolerance for z-constraint methodl is not specified\n"
1587     " , default value %f is used.\n",
1588     defaultZConsTol);
1589 tim 701 painCave.isFatal = 0;
1590     simError();
1591 tim 699
1592 tim 701 zconsTol->setData(defaultZConsTol);
1593     }
1594     theInfo.addProperty(zconsTol);
1595 tim 699
1596 tim 738 //set Force Subtraction Policy
1597 tim 722 StringData* zconsForcePolicy = new StringData();
1598 tim 701 zconsForcePolicy->setID(ZCONSFORCEPOLICY_ID);
1599 tim 722
1600     if (globals->haveZconsForcePolicy()){
1601 tim 701 zconsForcePolicy->setData(globals->getZconsForcePolicy());
1602 tim 722 }
1603 tim 701 else{
1604 tim 722 sprintf(painCave.errMsg,
1605 tim 738 "ZConstraint Warning: User does not set force Subtraction policy, "
1606 tim 736 "PolicyByMass is used\n");
1607 tim 722 painCave.isFatal = 0;
1608     simError();
1609 tim 736 zconsForcePolicy->setData("BYMASS");
1610 tim 701 }
1611 tim 722
1612     theInfo.addProperty(zconsForcePolicy);
1613    
1614 tim 701 //Determine the name of ouput file and add it into SimInfo's property list
1615     //Be careful, do not use inFileName, since it is a pointer which
1616     //point to a string at master node, and slave nodes do not contain that string
1617 tim 722
1618 tim 701 string zconsOutput(theInfo.finalName);
1619 tim 722
1620 tim 701 zconsOutput = zconsOutput.substr(0, zconsOutput.rfind(".")) + ".fz";
1621 tim 722
1622 tim 701 StringData* zconsFilename = new StringData();
1623     zconsFilename->setID(ZCONSFILENAME_ID);
1624     zconsFilename->setData(zconsOutput);
1625 tim 722
1626 tim 701 theInfo.addProperty(zconsFilename);
1627 tim 722
1628 tim 701 //setup index, pos and other parameters of z-constraint molecules
1629     nZConstraints = globals->getNzConstraints();
1630     theInfo.nZconstraints = nZConstraints;
1631    
1632     zconStamp = globals->getZconStamp();
1633     ZConsParaItem tempParaItem;
1634    
1635     ZConsParaData* zconsParaData = new ZConsParaData();
1636     zconsParaData->setID(ZCONSPARADATA_ID);
1637 tim 722
1638     for (int i = 0; i < nZConstraints; i++){
1639 tim 699 tempParaItem.havingZPos = zconStamp[i]->haveZpos();
1640     tempParaItem.zPos = zconStamp[i]->getZpos();
1641     tempParaItem.zconsIndex = zconStamp[i]->getMolIndex();
1642     tempParaItem.kRatio = zconStamp[i]->getKratio();
1643    
1644     zconsParaData->addItem(tempParaItem);
1645 tim 701 }
1646 tim 699
1647 tim 736 //check the uniqueness of index
1648     if(!zconsParaData->isIndexUnique()){
1649     sprintf(painCave.errMsg,
1650     "ZConstraint Error: molIndex is not unique\n");
1651     painCave.isFatal = 1;
1652     simError();
1653     }
1654    
1655 tim 701 //sort the parameters by index of molecules
1656     zconsParaData->sortByIndex();
1657 tim 736
1658 tim 701 //push data into siminfo, therefore, we can retrieve later
1659     theInfo.addProperty(zconsParaData);
1660 tim 660 }