ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-1.0/libmdtools/SimSetup.cpp
Revision: 1451
Committed: Mon Aug 9 14:50:35 2004 UTC (19 years, 11 months ago) by chrisfen
File size: 56459 byte(s)
Log Message:
Restraint loop is now more versatile

File Contents

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