ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 1066
Committed: Tue Feb 24 16:36:33 2004 UTC (20 years, 4 months ago) by tim
File size: 48189 byte(s)
Log Message:
*** empty log message ***

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