ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 778
Committed: Fri Sep 19 20:00:27 2003 UTC (20 years, 9 months ago) by mmeineke
File size: 45061 byte(s)
Log Message:
added NPT base class. NPTi is up to date. NPTf is not.

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