ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 814
Committed: Thu Oct 23 19:57:25 2003 UTC (20 years, 8 months ago) by mmeineke
File size: 43815 byte(s)
Log Message:
added eam ForceField files to the init

fixed an eam mpi parmeter setup bug

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