ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 999
Committed: Fri Jan 30 15:01:09 2004 UTC (20 years, 5 months ago) by chrisfen
File size: 46117 byte(s)
Log Message:
Substantial changes. OOPSE now has a working WATER.cpp forcefield and parser.
This involved changes to WATER.cpp and ForceFields amoung other files. One important
note: a hardwiring of LJ_rcut was made in calc_LJ_FF.F90. This will be removed on
the next commit...

File Contents

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