ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 1163
Committed: Wed May 12 14:30:12 2004 UTC (20 years, 4 months ago) by gezelter
File size: 52981 byte(s)
Log Message:
bug fixes for cutoffGroups

File Contents

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