ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 1174
Committed: Wed May 12 20:54:10 2004 UTC (20 years, 4 months ago) by gezelter
File size: 53648 byte(s)
Log Message:
Fixes for compilation under Mac OS X with IBM's xl compilers

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