ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/brains/SimSetup.cpp
Revision: 1628
Committed: Thu Oct 21 20:15:31 2004 UTC (19 years, 8 months ago) by gezelter
File size: 56222 byte(s)
Log Message:
Breaky Breaky.   Fixey Fixey.

File Contents

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