ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 1129
Committed: Thu Apr 22 03:29:30 2004 UTC (20 years, 4 months ago) by tim
File size: 51593 byte(s)
Log Message:
fixed two bugs in MPI version of InitfromFile and one unmatch MPI command in DumpWriter

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