ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 1104
Committed: Tue Apr 13 16:26:03 2004 UTC (20 years, 5 months ago) by gezelter
File size: 50934 byte(s)
Log Message:
Now molecules can keep track of their own IntegrableObjects (and
RigidBodies).  Also a bug-fix so that SimInfo can keep track of
RigidBodies (which was done incorrectly before).

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