ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/brains/SimSetup.cpp
Revision: 1650
Committed: Tue Oct 26 22:24:52 2004 UTC (19 years, 8 months ago) by gezelter
File size: 56378 byte(s)
Log Message:
forcefield refactoring for shapes

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