ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-1.0/libmdtools/SimSetup.cpp
Revision: 1419
Committed: Tue Jul 27 18:14:16 2004 UTC (20 years, 9 months ago) by gezelter
File size: 56446 byte(s)
Log Message:
BASS eradication project (part 3)

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 "SimSetup.hpp"
8 #include "ReadWrite.hpp"
9 #include "parse_me.h"
10 #include "Integrator.hpp"
11 #include "simError.h"
12 #include "RigidBody.hpp"
13 #include "OOPSEMinimizer.hpp"
14
15 #ifdef IS_MPI
16 #include "mpiBASS.h"
17 #include "mpiSimulation.hpp"
18 #endif
19
20 // some defines for ensemble and Forcefield cases
21
22 #define NVE_ENS 0
23 #define NVT_ENS 1
24 #define NPTi_ENS 2
25 #define NPTf_ENS 3
26 #define NPTxyz_ENS 4
27
28
29 #define FF_DUFF 0
30 #define FF_LJ 1
31 #define FF_EAM 2
32 #define FF_H2O 3
33
34 using namespace std;
35
36 /**
37 * Check whether dividend is divisble by divisor or not
38 */
39 bool isDivisible(double dividend, double divisor){
40 double tolerance = 0.000001;
41 double quotient;
42 double diff;
43 int intQuotient;
44
45 quotient = dividend / divisor;
46
47 if (quotient < 0)
48 quotient = -quotient;
49
50 intQuotient = int (quotient + tolerance);
51
52 diff = fabs(fabs(dividend) - intQuotient * fabs(divisor));
53
54 if (diff <= tolerance)
55 return true;
56 else
57 return false;
58 }
59
60 string getPrefix(const string& str ){
61 string prefix;
62 string suffix;
63 int pos;
64
65 pos = str.rfind(".");
66
67 if (pos >= 0) {
68 prefix = str.substr(0, pos);
69 suffix = str.substr(pos, str.size());
70
71 // leave .bass there in case we've reverted to old habits
72 if (LowerCase(suffix) == ".md" || LowerCase(suffix) == ".bass")
73 return prefix;
74 else
75 return str;
76
77 } else
78 return str;
79 };
80
81
82 SimSetup::SimSetup(){
83
84 initSuspend = false;
85 isInfoArray = 0;
86 nInfo = 1;
87
88 stamps = new MakeStamps();
89 globals = new Globals();
90
91
92 #ifdef IS_MPI
93 strcpy(checkPointMsg, "SimSetup creation successful");
94 MPIcheckPoint();
95 #endif // IS_MPI
96 }
97
98 SimSetup::~SimSetup(){
99 delete stamps;
100 delete globals;
101 }
102
103 void SimSetup::setSimInfo(SimInfo* the_info, int theNinfo){
104 info = the_info;
105 nInfo = theNinfo;
106 isInfoArray = 1;
107 initSuspend = true;
108 }
109
110
111 void SimSetup::parseFile(char* fileName){
112 #ifdef IS_MPI
113 if (worldRank == 0){
114 #endif // is_mpi
115
116 inFileName = fileName;
117 set_interface_stamps(stamps, globals);
118
119 #ifdef IS_MPI
120 mpiEventInit();
121 #endif
122
123 yacc_BASS(fileName);
124
125 #ifdef IS_MPI
126 throwMPIEvent(NULL);
127 }
128 else{
129 receiveParse();
130 }
131 #endif
132
133 }
134
135 #ifdef IS_MPI
136 void SimSetup::receiveParse(void){
137 set_interface_stamps(stamps, globals);
138 mpiEventInit();
139 MPIcheckPoint();
140 mpiEventLoop();
141 }
142
143 #endif // is_mpi
144
145 void SimSetup::createSim(void){
146
147 // gather all of the information from the meta-data file
148
149 gatherInfo();
150
151 // creation of complex system objects
152
153 sysObjectsCreation();
154
155 // check on the post processing info
156
157 finalInfoCheck();
158
159 // initialize the system coordinates
160
161 if ( !initSuspend ){
162 initSystemCoords();
163
164 if( !(globals->getUseInitTime()) )
165 info[0].currentTime = 0.0;
166 }
167
168 // make the output filenames
169
170 makeOutNames();
171
172 #ifdef IS_MPI
173 mpiSim->mpiRefresh();
174 #endif
175
176 // initialize the Fortran
177
178 initFortran();
179
180 if (globals->haveMinimizer())
181 // make minimizer
182 makeMinimizer();
183 else
184 // make the integrator
185 makeIntegrator();
186
187 }
188
189
190 void SimSetup::makeMolecules(void){
191 int i, j, k;
192 int exI, exJ, exK, exL, slI, slJ;
193 int tempI, tempJ, tempK, tempL;
194 int molI, globalID;
195 int stampID, atomOffset, rbOffset, groupOffset;
196 molInit molInfo;
197 DirectionalAtom* dAtom;
198 RigidBody* myRB;
199 StuntDouble* mySD;
200 LinkedAssign* extras;
201 LinkedAssign* current_extra;
202 AtomStamp* currentAtom;
203 BondStamp* currentBond;
204 BendStamp* currentBend;
205 TorsionStamp* currentTorsion;
206 RigidBodyStamp* currentRigidBody;
207 CutoffGroupStamp* currentCutoffGroup;
208 CutoffGroup* myCutoffGroup;
209 int nCutoffGroups;// number of cutoff group of a molecule defined in mdl file
210 set<int> cutoffAtomSet; //atoms belong to cutoffgroup defined at mdl file
211
212 bond_pair* theBonds;
213 bend_set* theBends;
214 torsion_set* theTorsions;
215
216 set<int> skipList;
217
218 double phi, theta, psi;
219 char* molName;
220 char rbName[100];
221
222 int whichRigidBody;
223 int consAtomIndex; //index of constraint atom in rigid body's atom array
224 double bondLength2;
225 //init the forceField paramters
226
227 the_ff->readParams();
228
229 // init the atoms
230
231 int nMembers, nNew, rb1, rb2;
232
233 for (k = 0; k < nInfo; k++){
234 the_ff->setSimInfo(&(info[k]));
235
236 #ifdef IS_MPI
237 info[k].globalGroupMembership = new int[mpiSim->getNAtomsGlobal()];
238 for (i = 0; i < mpiSim->getNAtomsGlobal(); i++)
239 info[k].globalGroupMembership[i] = 0;
240 #else
241 info[k].globalGroupMembership = new int[info[k].n_atoms];
242 for (i = 0; i < info[k].n_atoms; i++)
243 info[k].globalGroupMembership[i] = 0;
244 #endif
245
246 atomOffset = 0;
247 groupOffset = 0;
248
249 for (i = 0; i < info[k].n_mol; i++){
250 stampID = info[k].molecules[i].getStampID();
251 molName = comp_stamps[stampID]->getID();
252
253 molInfo.nAtoms = comp_stamps[stampID]->getNAtoms();
254 molInfo.nBonds = comp_stamps[stampID]->getNBonds();
255 molInfo.nBends = comp_stamps[stampID]->getNBends();
256 molInfo.nTorsions = comp_stamps[stampID]->getNTorsions();
257 molInfo.nRigidBodies = comp_stamps[stampID]->getNRigidBodies();
258
259 nCutoffGroups = comp_stamps[stampID]->getNCutoffGroups();
260
261 molInfo.myAtoms = &(info[k].atoms[atomOffset]);
262
263 if (molInfo.nBonds > 0)
264 molInfo.myBonds = new Bond*[molInfo.nBonds];
265 else
266 molInfo.myBonds = NULL;
267
268 if (molInfo.nBends > 0)
269 molInfo.myBends = new Bend*[molInfo.nBends];
270 else
271 molInfo.myBends = NULL;
272
273 if (molInfo.nTorsions > 0)
274 molInfo.myTorsions = new Torsion *[molInfo.nTorsions];
275 else
276 molInfo.myTorsions = NULL;
277
278 theBonds = new bond_pair[molInfo.nBonds];
279 theBends = new bend_set[molInfo.nBends];
280 theTorsions = new torsion_set[molInfo.nTorsions];
281
282 // make the Atoms
283
284 for (j = 0; j < molInfo.nAtoms; j++){
285 currentAtom = comp_stamps[stampID]->getAtom(j);
286
287 if (currentAtom->haveOrientation()){
288 dAtom = new DirectionalAtom((j + atomOffset),
289 info[k].getConfiguration());
290 info[k].n_oriented++;
291 molInfo.myAtoms[j] = dAtom;
292
293 // Directional Atoms have standard unit vectors which are oriented
294 // in space using the three Euler angles. We assume the standard
295 // unit vector was originally along the z axis below.
296
297 phi = currentAtom->getEulerPhi() * M_PI / 180.0;
298 theta = currentAtom->getEulerTheta() * M_PI / 180.0;
299 psi = currentAtom->getEulerPsi()* M_PI / 180.0;
300
301 dAtom->setUnitFrameFromEuler(phi, theta, psi);
302
303 }
304 else{
305
306 molInfo.myAtoms[j] = new Atom((j + atomOffset), info[k].getConfiguration());
307
308 }
309
310 molInfo.myAtoms[j]->setType(currentAtom->getType());
311 #ifdef IS_MPI
312 molInfo.myAtoms[j]->setGlobalIndex(globalAtomIndex[j + atomOffset]);
313 #endif // is_mpi
314 }
315
316 // make the bonds
317 for (j = 0; j < molInfo.nBonds; j++){
318 currentBond = comp_stamps[stampID]->getBond(j);
319 theBonds[j].a = currentBond->getA() + atomOffset;
320 theBonds[j].b = currentBond->getB() + atomOffset;
321
322 tempI = theBonds[j].a;
323 tempJ = theBonds[j].b;
324
325 #ifdef IS_MPI
326 exI = info[k].atoms[tempI]->getGlobalIndex() + 1;
327 exJ = info[k].atoms[tempJ]->getGlobalIndex() + 1;
328 #else
329 exI = tempI + 1;
330 exJ = tempJ + 1;
331 #endif
332
333 info[k].excludes->addPair(exI, exJ);
334 }
335
336 //make the bends
337 for (j = 0; j < molInfo.nBends; j++){
338 currentBend = comp_stamps[stampID]->getBend(j);
339 theBends[j].a = currentBend->getA() + atomOffset;
340 theBends[j].b = currentBend->getB() + atomOffset;
341 theBends[j].c = currentBend->getC() + atomOffset;
342
343 if (currentBend->haveExtras()){
344 extras = currentBend->getExtras();
345 current_extra = extras;
346
347 while (current_extra != NULL){
348 if (!strcmp(current_extra->getlhs(), "ghostVectorSource")){
349 switch (current_extra->getType()){
350 case 0:
351 theBends[j].ghost = current_extra->getInt() + atomOffset;
352 theBends[j].isGhost = 1;
353 break;
354
355 case 1:
356 theBends[j].ghost = (int) current_extra->getDouble() +
357 atomOffset;
358 theBends[j].isGhost = 1;
359 break;
360
361 default:
362 sprintf(painCave.errMsg,
363 "SimSetup Error: ghostVectorSource was neither a "
364 "double nor an int.\n"
365 "-->Bend[%d] in %s\n",
366 j, comp_stamps[stampID]->getID());
367 painCave.isFatal = 1;
368 simError();
369 }
370 }
371 else{
372 sprintf(painCave.errMsg,
373 "SimSetup Error: unhandled bend assignment:\n"
374 " -->%s in Bend[%d] in %s\n",
375 current_extra->getlhs(), j, comp_stamps[stampID]->getID());
376 painCave.isFatal = 1;
377 simError();
378 }
379
380 current_extra = current_extra->getNext();
381 }
382 }
383
384 if (theBends[j].isGhost) {
385
386 tempI = theBends[j].a;
387 tempJ = theBends[j].b;
388
389 #ifdef IS_MPI
390 exI = info[k].atoms[tempI]->getGlobalIndex() + 1;
391 exJ = info[k].atoms[tempJ]->getGlobalIndex() + 1;
392 #else
393 exI = tempI + 1;
394 exJ = tempJ + 1;
395 #endif
396 info[k].excludes->addPair(exI, exJ);
397
398 } else {
399
400 tempI = theBends[j].a;
401 tempJ = theBends[j].b;
402 tempK = theBends[j].c;
403
404 #ifdef IS_MPI
405 exI = info[k].atoms[tempI]->getGlobalIndex() + 1;
406 exJ = info[k].atoms[tempJ]->getGlobalIndex() + 1;
407 exK = info[k].atoms[tempK]->getGlobalIndex() + 1;
408 #else
409 exI = tempI + 1;
410 exJ = tempJ + 1;
411 exK = tempK + 1;
412 #endif
413
414 info[k].excludes->addPair(exI, exK);
415 info[k].excludes->addPair(exI, exJ);
416 info[k].excludes->addPair(exJ, exK);
417 }
418 }
419
420 for (j = 0; j < molInfo.nTorsions; j++){
421 currentTorsion = comp_stamps[stampID]->getTorsion(j);
422 theTorsions[j].a = currentTorsion->getA() + atomOffset;
423 theTorsions[j].b = currentTorsion->getB() + atomOffset;
424 theTorsions[j].c = currentTorsion->getC() + atomOffset;
425 theTorsions[j].d = currentTorsion->getD() + atomOffset;
426
427 tempI = theTorsions[j].a;
428 tempJ = theTorsions[j].b;
429 tempK = theTorsions[j].c;
430 tempL = theTorsions[j].d;
431
432 #ifdef IS_MPI
433 exI = info[k].atoms[tempI]->getGlobalIndex() + 1;
434 exJ = info[k].atoms[tempJ]->getGlobalIndex() + 1;
435 exK = info[k].atoms[tempK]->getGlobalIndex() + 1;
436 exL = info[k].atoms[tempL]->getGlobalIndex() + 1;
437 #else
438 exI = tempI + 1;
439 exJ = tempJ + 1;
440 exK = tempK + 1;
441 exL = tempL + 1;
442 #endif
443
444 info[k].excludes->addPair(exI, exJ);
445 info[k].excludes->addPair(exI, exK);
446 info[k].excludes->addPair(exI, exL);
447 info[k].excludes->addPair(exJ, exK);
448 info[k].excludes->addPair(exJ, exL);
449 info[k].excludes->addPair(exK, exL);
450 }
451
452
453 molInfo.myRigidBodies.clear();
454
455 for (j = 0; j < molInfo.nRigidBodies; j++){
456
457 currentRigidBody = comp_stamps[stampID]->getRigidBody(j);
458 nMembers = currentRigidBody->getNMembers();
459
460 // Create the Rigid Body:
461
462 myRB = new RigidBody();
463
464 sprintf(rbName,"%s_RB_%d", molName, j);
465 myRB->setType(rbName);
466
467 for (rb1 = 0; rb1 < nMembers; rb1++) {
468
469 // molI is atom numbering inside this molecule
470 molI = currentRigidBody->getMember(rb1);
471
472 // tempI is atom numbering on local processor
473 tempI = molI + atomOffset;
474
475 // currentAtom is the AtomStamp (which we need for
476 // rigid body reference positions)
477 currentAtom = comp_stamps[stampID]->getAtom(molI);
478
479 // When we add to the rigid body, add the atom itself and
480 // the stamp info:
481
482 myRB->addAtom(info[k].atoms[tempI], currentAtom);
483
484 // Add this atom to the Skip List for the integrators
485 #ifdef IS_MPI
486 slI = info[k].atoms[tempI]->getGlobalIndex();
487 #else
488 slI = tempI;
489 #endif
490 skipList.insert(slI);
491
492 }
493
494 for(rb1 = 0; rb1 < nMembers - 1; rb1++) {
495 for(rb2 = rb1+1; rb2 < nMembers; rb2++) {
496
497 tempI = currentRigidBody->getMember(rb1);
498 tempJ = currentRigidBody->getMember(rb2);
499
500 // Some explanation is required here.
501 // Fortran indexing starts at 1, while c indexing starts at 0
502 // Also, in parallel computations, the GlobalIndex is
503 // used for the exclude list:
504
505 #ifdef IS_MPI
506 exI = molInfo.myAtoms[tempI]->getGlobalIndex() + 1;
507 exJ = molInfo.myAtoms[tempJ]->getGlobalIndex() + 1;
508 #else
509 exI = molInfo.myAtoms[tempI]->getIndex() + 1;
510 exJ = molInfo.myAtoms[tempJ]->getIndex() + 1;
511 #endif
512
513 info[k].excludes->addPair(exI, exJ);
514
515 }
516 }
517
518 molInfo.myRigidBodies.push_back(myRB);
519 info[k].rigidBodies.push_back(myRB);
520 }
521
522
523 //create cutoff group for molecule
524
525 cutoffAtomSet.clear();
526 molInfo.myCutoffGroups.clear();
527
528 for (j = 0; j < nCutoffGroups; j++){
529
530 currentCutoffGroup = comp_stamps[stampID]->getCutoffGroup(j);
531 nMembers = currentCutoffGroup->getNMembers();
532
533 myCutoffGroup = new CutoffGroup();
534
535 #ifdef IS_MPI
536 myCutoffGroup->setGlobalIndex(globalGroupIndex[groupOffset]);
537 #else
538 myCutoffGroup->setGlobalIndex(groupOffset);
539 #endif
540
541 for (int cg = 0; cg < nMembers; cg++) {
542
543 // molI is atom numbering inside this molecule
544 molI = currentCutoffGroup->getMember(cg);
545
546 // tempI is atom numbering on local processor
547 tempI = molI + atomOffset;
548
549 #ifdef IS_MPI
550 globalID = info[k].atoms[tempI]->getGlobalIndex();
551 info[k].globalGroupMembership[globalID] = globalGroupIndex[groupOffset];
552 #else
553 globalID = info[k].atoms[tempI]->getIndex();
554 info[k].globalGroupMembership[globalID] = groupOffset;
555 #endif
556 myCutoffGroup->addAtom(info[k].atoms[tempI]);
557 cutoffAtomSet.insert(tempI);
558 }
559
560 molInfo.myCutoffGroups.push_back(myCutoffGroup);
561 groupOffset++;
562
563 }//end for (j = 0; j < molInfo.nCutoffGroups; j++)
564
565
566 // create a cutoff group for every atom in current molecule which
567 // does not belong to cutoffgroup defined at mdl file
568
569 for(j = 0; j < molInfo.nAtoms; j++){
570
571 if(cutoffAtomSet.find(molInfo.myAtoms[j]->getIndex()) == cutoffAtomSet.end()){
572 myCutoffGroup = new CutoffGroup();
573 myCutoffGroup->addAtom(molInfo.myAtoms[j]);
574
575 #ifdef IS_MPI
576 myCutoffGroup->setGlobalIndex(globalGroupIndex[groupOffset]);
577 globalID = info[k].atoms[atomOffset + j]->getGlobalIndex();
578 info[k].globalGroupMembership[globalID] = globalGroupIndex[groupOffset];
579 #else
580 myCutoffGroup->setGlobalIndex(groupOffset);
581 globalID = info[k].atoms[atomOffset + j]->getIndex();
582 info[k].globalGroupMembership[globalID] = groupOffset;
583 #endif
584 molInfo.myCutoffGroups.push_back(myCutoffGroup);
585 groupOffset++;
586 }
587 }
588
589 // After this is all set up, scan through the atoms to
590 // see if they can be added to the integrableObjects:
591
592 molInfo.myIntegrableObjects.clear();
593
594
595 for (j = 0; j < molInfo.nAtoms; j++){
596
597 #ifdef IS_MPI
598 slJ = molInfo.myAtoms[j]->getGlobalIndex();
599 #else
600 slJ = j+atomOffset;
601 #endif
602
603 // if they aren't on the skip list, then they can be integrated
604
605 if (skipList.find(slJ) == skipList.end()) {
606 mySD = (StuntDouble *) molInfo.myAtoms[j];
607 info[k].integrableObjects.push_back(mySD);
608 molInfo.myIntegrableObjects.push_back(mySD);
609 }
610 }
611
612 // all rigid bodies are integrated:
613
614 for (j = 0; j < molInfo.nRigidBodies; j++) {
615 mySD = (StuntDouble *) molInfo.myRigidBodies[j];
616 info[k].integrableObjects.push_back(mySD);
617 molInfo.myIntegrableObjects.push_back(mySD);
618 }
619
620 // send the arrays off to the forceField for init.
621
622 the_ff->initializeAtoms(molInfo.nAtoms, molInfo.myAtoms);
623 the_ff->initializeBonds(molInfo.nBonds, molInfo.myBonds, theBonds);
624 the_ff->initializeBends(molInfo.nBends, molInfo.myBends, theBends);
625 the_ff->initializeTorsions(molInfo.nTorsions, molInfo.myTorsions,
626 theTorsions);
627
628 info[k].molecules[i].initialize(molInfo);
629
630
631 atomOffset += molInfo.nAtoms;
632 delete[] theBonds;
633 delete[] theBends;
634 delete[] theTorsions;
635 }
636
637
638
639 #ifdef IS_MPI
640 // Since the globalGroupMembership has been zero filled and we've only
641 // poked values into the atoms we know, we can do an Allreduce
642 // to get the full globalGroupMembership array (We think).
643 // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
644 // docs said we could.
645
646 int* ggMjunk = new int[mpiSim->getNAtomsGlobal()];
647
648 MPI_Allreduce(info[k].globalGroupMembership,
649 ggMjunk,
650 mpiSim->getNAtomsGlobal(),
651 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
652
653 for (i = 0; i < mpiSim->getNAtomsGlobal(); i++)
654 info[k].globalGroupMembership[i] = ggMjunk[i];
655
656 delete[] ggMjunk;
657
658 #endif
659
660
661
662 }
663
664 #ifdef IS_MPI
665 sprintf(checkPointMsg, "all molecules initialized succesfully");
666 MPIcheckPoint();
667 #endif // is_mpi
668
669 }
670
671 void SimSetup::gatherInfo(void){
672 int i;
673
674 ensembleCase = -1;
675 ffCase = -1;
676
677 // set the easy ones first
678
679 for (i = 0; i < nInfo; i++){
680 if (globals->haveTargetTemp()) {
681 info[i].target_temp = globals->getTargetTemp();
682 info[i].have_target_temp = 1;
683 } else {
684 info[i].have_target_temp = 0;
685 }
686 if (globals->haveDt()) {
687 info[i].dt = globals->getDt();
688 }
689 if (globals->haveRunTime()) {
690 info[i].run_time = globals->getRunTime();
691 }
692 }
693 n_components = globals->getNComponents();
694
695
696 // get the forceField
697
698 strcpy(force_field, globals->getForceField());
699
700 if (!strcasecmp(force_field, "DUFF")){
701 ffCase = FF_DUFF;
702 }
703 else if (!strcasecmp(force_field, "LJ")){
704 ffCase = FF_LJ;
705 }
706 else if (!strcasecmp(force_field, "EAM")){
707 ffCase = FF_EAM;
708 }
709 else if (!strcasecmp(force_field, "WATER")){
710 ffCase = FF_H2O;
711 }
712 else{
713 sprintf(painCave.errMsg, "SimSetup Error. Unrecognized force field -> %s\n",
714 force_field);
715 painCave.isFatal = 1;
716 simError();
717 }
718 if (globals->haveForceFieldVariant()) {
719 strcpy(forcefield_variant, globals->getForceFieldVariant());
720 has_forcefield_variant = 1;
721 }
722
723 // get the ensemble
724
725
726 if (globals->haveEnsemble()) {
727
728 strcpy(ensemble, globals->getEnsemble());
729
730 if (!strcasecmp(ensemble, "NVE")){
731 ensembleCase = NVE_ENS;
732 }
733 else if (!strcasecmp(ensemble, "NVT")){
734 ensembleCase = NVT_ENS;
735 }
736 else if (!strcasecmp(ensemble, "NPTi") || !strcasecmp(ensemble, "NPT")){
737 ensembleCase = NPTi_ENS;
738 }
739 else if (!strcasecmp(ensemble, "NPTf")){
740 ensembleCase = NPTf_ENS;
741 }
742 else if (!strcasecmp(ensemble, "NPTxyz")){
743 ensembleCase = NPTxyz_ENS;
744 }
745 else{
746 sprintf(painCave.errMsg,
747 "SimSetup Warning. Unrecognized Ensemble -> %s \n"
748 "\treverting to NVE for this simulation.\n",
749 ensemble);
750 painCave.isFatal = 0;
751 simError();
752 strcpy(ensemble, "NVE");
753 ensembleCase = NVE_ENS;
754 }
755
756 for (i = 0; i < nInfo; i++)
757 strcpy(info[i].ensemble, ensemble);
758
759
760 //check whether sample time, status time, thermal time and reset time are divisble by dt
761 if (globals->haveSampleTime() && !isDivisible(globals->getSampleTime(), globals->getDt())){
762 sprintf(painCave.errMsg,
763 "Sample time is not divisible by dt.\n"
764 "\tThis will result in samples that are not uniformly\n"
765 "\tdistributed in time. If this is a problem, change\n"
766 "\tyour sampleTime variable.\n");
767 painCave.isFatal = 0;
768 simError();
769 }
770
771 if (globals->haveStatusTime() && !isDivisible(globals->getStatusTime(), globals->getDt())){
772 sprintf(painCave.errMsg,
773 "Status time is not divisible by dt.\n"
774 "\tThis will result in status reports that are not uniformly\n"
775 "\tdistributed in time. If this is a problem, change \n"
776 "\tyour statusTime variable.\n");
777 painCave.isFatal = 0;
778 simError();
779 }
780
781 if (globals->haveThermalTime() && !isDivisible(globals->getThermalTime(), globals->getDt())){
782 sprintf(painCave.errMsg,
783 "Thermal time is not divisible by dt.\n"
784 "\tThis will result in thermalizations that are not uniformly\n"
785 "\tdistributed in time. If this is a problem, change \n"
786 "\tyour thermalTime variable.\n");
787 painCave.isFatal = 0;
788 simError();
789 }
790
791 if (globals->haveResetTime() && !isDivisible(globals->getResetTime(), globals->getDt())){
792 sprintf(painCave.errMsg,
793 "Reset time is not divisible by dt.\n"
794 "\tThis will result in integrator resets that are not uniformly\n"
795 "\tdistributed in time. If this is a problem, change\n"
796 "\tyour resetTime variable.\n");
797 painCave.isFatal = 0;
798 simError();
799 }
800
801 // set the status, sample, and thermal kick times
802
803 for (i = 0; i < nInfo; i++){
804 if (globals->haveSampleTime()){
805 info[i].sampleTime = globals->getSampleTime();
806 info[i].statusTime = info[i].sampleTime;
807 }
808 else{
809 info[i].sampleTime = globals->getRunTime();
810 info[i].statusTime = info[i].sampleTime;
811 }
812
813 if (globals->haveStatusTime()){
814 info[i].statusTime = globals->getStatusTime();
815 }
816
817 if (globals->haveThermalTime()){
818 info[i].thermalTime = globals->getThermalTime();
819 } else {
820 info[i].thermalTime = globals->getRunTime();
821 }
822
823 info[i].resetIntegrator = 0;
824 if( globals->haveResetTime() ){
825 info[i].resetTime = globals->getResetTime();
826 info[i].resetIntegrator = 1;
827 }
828 }
829
830
831 for (i=0; i < nInfo; i++) {
832
833 // check for the temperature set flag
834
835 if (globals->haveTempSet())
836 info[i].setTemp = globals->getTempSet();
837
838 // check for the extended State init
839
840 info[i].useInitXSstate = globals->getUseInitXSstate();
841 info[i].orthoTolerance = globals->getOrthoBoxTolerance();
842
843 // check for thermodynamic integration
844 if (globals->getUseSolidThermInt() && !globals->getUseLiquidThermInt()) {
845 if (globals->haveThermIntLambda() && globals->haveThermIntK()) {
846 info[i].useSolidThermInt = globals->getUseSolidThermInt();
847 info[i].thermIntLambda = globals->getThermIntLambda();
848 info[i].thermIntK = globals->getThermIntK();
849
850 Restraints *myRestraint = new Restraints(tot_nmol, info[i].thermIntLambda, info[i].thermIntK);
851 info[i].restraint = myRestraint;
852 }
853 else {
854 sprintf(painCave.errMsg,
855 "SimSetup Error:\n"
856 "\tKeyword useSolidThermInt was set to 'true' but\n"
857 "\tthermodynamicIntegrationLambda (and/or\n"
858 "\tthermodynamicIntegrationK) was not specified.\n"
859 "\tPlease provide a lambda value and k value in your meta-data file.\n");
860 painCave.isFatal = 1;
861 simError();
862 }
863 }
864 else if(globals->getUseLiquidThermInt()) {
865 if (globals->getUseSolidThermInt()) {
866 sprintf( painCave.errMsg,
867 "SimSetup Warning: It appears that you have both solid and\n"
868 "\tliquid thermodynamic integration activated in your meta-data\n"
869 "\tfile. To avoid confusion, specify only one technique in\n"
870 "\tyour meta-data file. Liquid-state thermodynamic integration\n"
871 "\twill be assumed for the current simulation. If this is not\n"
872 "\twhat you desire, set useSolidThermInt to 'true' and\n"
873 "\tuseLiquidThermInt to 'false' in your meta-data file.\n");
874 painCave.isFatal = 0;
875 simError();
876 }
877 if (globals->haveThermIntLambda() && globals->haveThermIntK()) {
878 info[i].useLiquidThermInt = globals->getUseLiquidThermInt();
879 info[i].thermIntLambda = globals->getThermIntLambda();
880 info[i].thermIntK = globals->getThermIntK();
881 }
882 else {
883 sprintf(painCave.errMsg,
884 "SimSetup Error:\n"
885 "\tKeyword useLiquidThermInt was set to 'true' but\n"
886 "\tthermodynamicIntegrationLambda (and/or\n"
887 "\tthermodynamicIntegrationK) was not specified.\n"
888 "\tPlease provide a lambda value and k value in your meta-data file.\n");
889 painCave.isFatal = 1;
890 simError();
891 }
892 }
893 else if(globals->haveThermIntLambda() || globals->haveThermIntK()){
894 sprintf(painCave.errMsg,
895 "SimSetup Warning: If you want to use Thermodynamic\n"
896 "\tIntegration, set useSolidThermInt or useLiquidThermInt to\n"
897 "\t'true' in your meta-data file. These keywords are set to\n"
898 "\t'false' by default, so your lambda and/or k values are\n"
899 "\tbeing ignored.\n");
900 painCave.isFatal = 0;
901 simError();
902 }
903 }
904 }
905
906 for (i = 0; i < nInfo; i++) {
907 // get the mixing rule
908
909 strcpy(info[i].mixingRule, globals->getMixingRule());
910 info[i].usePBC = globals->getPBC();
911 }
912
913 // get the components and calculate the tot_nMol and indvidual n_mol
914
915 the_components = globals->getComponents();
916 components_nmol = new int[n_components];
917
918 if (!globals->haveNMol()){
919 // we don't have the total number of molecules, so we assume it is
920 // given in each component
921
922 tot_nmol = 0;
923 for (i = 0; i < n_components; i++){
924 if (!the_components[i]->haveNMol()){
925 // we have a problem
926 sprintf(painCave.errMsg,
927 "SimSetup Error. No global NMol or component NMol given.\n"
928 "\tCannot calculate the number of atoms.\n");
929 painCave.isFatal = 1;
930 simError();
931 }
932
933 tot_nmol += the_components[i]->getNMol();
934 components_nmol[i] = the_components[i]->getNMol();
935 }
936 }
937 else{
938 sprintf(painCave.errMsg,
939 "SimSetup error.\n"
940 "\tSorry, the ability to specify total"
941 " nMols and then give molfractions in the components\n"
942 "\tis not currently supported."
943 " Please give nMol in the components.\n");
944 painCave.isFatal = 1;
945 simError();
946 }
947
948
949
950
951 //setup seed for random number generator
952 int seedValue;
953
954 if (globals->haveSeed()){
955 seedValue = globals->getSeed();
956
957 if(seedValue / 1E9 == 0){
958 sprintf(painCave.errMsg,
959 "Seed for sprng library should contain at least 9 digits\n"
960 "OOPSE will generate a seed for user\n");
961 painCave.isFatal = 0;
962 simError();
963
964 //using seed generated by system instead of invalid seed set by user
965 #ifndef IS_MPI
966 seedValue = make_sprng_seed();
967 #else
968 if (worldRank == 0){
969 seedValue = make_sprng_seed();
970 }
971 MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD);
972 #endif
973 }
974 }//end of if branch of globals->haveSeed()
975 else{
976
977 #ifndef IS_MPI
978 seedValue = make_sprng_seed();
979 #else
980 if (worldRank == 0){
981 seedValue = make_sprng_seed();
982 }
983 MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD);
984 #endif
985 }//end of globals->haveSeed()
986
987 for (int i = 0; i < nInfo; i++){
988 info[i].setSeed(seedValue);
989 }
990
991 #ifdef IS_MPI
992 strcpy(checkPointMsg, "Successfully gathered all information from meta-data file\n");
993 MPIcheckPoint();
994 #endif // is_mpi
995 }
996
997
998 void SimSetup::finalInfoCheck(void){
999 int index;
1000 int usesDipoles;
1001 int usesCharges;
1002 int i;
1003
1004 for (i = 0; i < nInfo; i++){
1005 // check electrostatic parameters
1006
1007 index = 0;
1008 usesDipoles = 0;
1009 while ((index < info[i].n_atoms) && !usesDipoles){
1010 usesDipoles = (info[i].atoms[index])->hasDipole();
1011 index++;
1012 }
1013 index = 0;
1014 usesCharges = 0;
1015 while ((index < info[i].n_atoms) && !usesCharges){
1016 usesCharges= (info[i].atoms[index])->hasCharge();
1017 index++;
1018 }
1019 #ifdef IS_MPI
1020 int myUse = usesDipoles;
1021 MPI_Allreduce(&myUse, &usesDipoles, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
1022 #endif //is_mpi
1023
1024 double theRcut, theRsw;
1025
1026 if (globals->haveRcut()) {
1027 theRcut = globals->getRcut();
1028
1029 if (globals->haveRsw())
1030 theRsw = globals->getRsw();
1031 else
1032 theRsw = theRcut;
1033
1034 info[i].setDefaultRcut(theRcut, theRsw);
1035
1036 } else {
1037
1038 the_ff->calcRcut();
1039 theRcut = info[i].getRcut();
1040
1041 if (globals->haveRsw())
1042 theRsw = globals->getRsw();
1043 else
1044 theRsw = theRcut;
1045
1046 info[i].setDefaultRcut(theRcut, theRsw);
1047 }
1048
1049 if (globals->getUseRF()){
1050 info[i].useReactionField = 1;
1051
1052 if (!globals->haveRcut()){
1053 sprintf(painCave.errMsg,
1054 "SimSetup Warning: No value was set for the cutoffRadius.\n"
1055 "\tOOPSE will use a default value of 15.0 angstroms"
1056 "\tfor the cutoffRadius.\n");
1057 painCave.isFatal = 0;
1058 simError();
1059 theRcut = 15.0;
1060 }
1061 else{
1062 theRcut = globals->getRcut();
1063 }
1064
1065 if (!globals->haveRsw()){
1066 sprintf(painCave.errMsg,
1067 "SimSetup Warning: No value was set for switchingRadius.\n"
1068 "\tOOPSE will use a default value of\n"
1069 "\t0.95 * cutoffRadius for the switchingRadius\n");
1070 painCave.isFatal = 0;
1071 simError();
1072 theRsw = 0.95 * theRcut;
1073 }
1074 else{
1075 theRsw = globals->getRsw();
1076 }
1077
1078 info[i].setDefaultRcut(theRcut, theRsw);
1079
1080 if (!globals->haveDielectric()){
1081 sprintf(painCave.errMsg,
1082 "SimSetup Error: No Dielectric constant was set.\n"
1083 "\tYou are trying to use Reaction Field without"
1084 "\tsetting a dielectric constant!\n");
1085 painCave.isFatal = 1;
1086 simError();
1087 }
1088 info[i].dielectric = globals->getDielectric();
1089 }
1090 else{
1091 if (usesDipoles || usesCharges){
1092
1093 if (!globals->haveRcut()){
1094 sprintf(painCave.errMsg,
1095 "SimSetup Warning: No value was set for the cutoffRadius.\n"
1096 "\tOOPSE will use a default value of 15.0 angstroms"
1097 "\tfor the cutoffRadius.\n");
1098 painCave.isFatal = 0;
1099 simError();
1100 theRcut = 15.0;
1101 }
1102 else{
1103 theRcut = globals->getRcut();
1104 }
1105
1106 if (!globals->haveRsw()){
1107 sprintf(painCave.errMsg,
1108 "SimSetup Warning: No value was set for switchingRadius.\n"
1109 "\tOOPSE will use a default value of\n"
1110 "\t0.95 * cutoffRadius for the switchingRadius\n");
1111 painCave.isFatal = 0;
1112 simError();
1113 theRsw = 0.95 * theRcut;
1114 }
1115 else{
1116 theRsw = globals->getRsw();
1117 }
1118
1119 info[i].setDefaultRcut(theRcut, theRsw);
1120
1121 }
1122 }
1123 }
1124 #ifdef IS_MPI
1125 strcpy(checkPointMsg, "post processing checks out");
1126 MPIcheckPoint();
1127 #endif // is_mpi
1128
1129 // clean up the forcefield
1130 the_ff->cleanMe();
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 if (!strcmp(info[0].mixingRule, "standard")){
1766 the_ff->initForceField(LB_MIXING_RULE);
1767 }
1768 else if (!strcmp(info[0].mixingRule, "explicit")){
1769 the_ff->initForceField(EXPLICIT_MIXING_RULE);
1770 }
1771 else{
1772 sprintf(painCave.errMsg, "SimSetup Error: unknown mixing rule -> \"%s\"\n",
1773 info[0].mixingRule);
1774 painCave.isFatal = 1;
1775 simError();
1776 }
1777
1778
1779 #ifdef IS_MPI
1780 strcpy(checkPointMsg, "Successfully intialized the mixingRule for Fortran.");
1781 MPIcheckPoint();
1782 #endif // is_mpi
1783 }
1784
1785 void SimSetup::setupZConstraint(SimInfo& theInfo){
1786 int nZConstraints;
1787 ZconStamp** zconStamp;
1788
1789 if (globals->haveZconstraintTime()){
1790 //add sample time of z-constraint into SimInfo's property list
1791 DoubleData* zconsTimeProp = new DoubleData();
1792 zconsTimeProp->setID(ZCONSTIME_ID);
1793 zconsTimeProp->setData(globals->getZconsTime());
1794 theInfo.addProperty(zconsTimeProp);
1795 }
1796 else{
1797 sprintf(painCave.errMsg,
1798 "ZConstraint error: If you use a ZConstraint,\n"
1799 "\tyou must set zconsTime.\n");
1800 painCave.isFatal = 1;
1801 simError();
1802 }
1803
1804 //push zconsTol into siminfo, if user does not specify
1805 //value for zconsTol, a default value will be used
1806 DoubleData* zconsTol = new DoubleData();
1807 zconsTol->setID(ZCONSTOL_ID);
1808 if (globals->haveZconsTol()){
1809 zconsTol->setData(globals->getZconsTol());
1810 }
1811 else{
1812 double defaultZConsTol = 0.01;
1813 sprintf(painCave.errMsg,
1814 "ZConstraint Warning: Tolerance for z-constraint method is not specified.\n"
1815 "\tOOPSE will use a default value of %f.\n"
1816 "\tTo set the tolerance, use the zconsTol variable.\n",
1817 defaultZConsTol);
1818 painCave.isFatal = 0;
1819 simError();
1820
1821 zconsTol->setData(defaultZConsTol);
1822 }
1823 theInfo.addProperty(zconsTol);
1824
1825 //set Force Subtraction Policy
1826 StringData* zconsForcePolicy = new StringData();
1827 zconsForcePolicy->setID(ZCONSFORCEPOLICY_ID);
1828
1829 if (globals->haveZconsForcePolicy()){
1830 zconsForcePolicy->setData(globals->getZconsForcePolicy());
1831 }
1832 else{
1833 sprintf(painCave.errMsg,
1834 "ZConstraint Warning: No force subtraction policy was set.\n"
1835 "\tOOPSE will use PolicyByMass.\n"
1836 "\tTo set the policy, use the zconsForcePolicy variable.\n");
1837 painCave.isFatal = 0;
1838 simError();
1839 zconsForcePolicy->setData("BYMASS");
1840 }
1841
1842 theInfo.addProperty(zconsForcePolicy);
1843
1844 //set zcons gap
1845 DoubleData* zconsGap = new DoubleData();
1846 zconsGap->setID(ZCONSGAP_ID);
1847
1848 if (globals->haveZConsGap()){
1849 zconsGap->setData(globals->getZconsGap());
1850 theInfo.addProperty(zconsGap);
1851 }
1852
1853 //set zcons fixtime
1854 DoubleData* zconsFixtime = new DoubleData();
1855 zconsFixtime->setID(ZCONSFIXTIME_ID);
1856
1857 if (globals->haveZConsFixTime()){
1858 zconsFixtime->setData(globals->getZconsFixtime());
1859 theInfo.addProperty(zconsFixtime);
1860 }
1861
1862 //set zconsUsingSMD
1863 IntData* zconsUsingSMD = new IntData();
1864 zconsUsingSMD->setID(ZCONSUSINGSMD_ID);
1865
1866 if (globals->haveZConsUsingSMD()){
1867 zconsUsingSMD->setData(globals->getZconsUsingSMD());
1868 theInfo.addProperty(zconsUsingSMD);
1869 }
1870
1871 //Determine the name of ouput file and add it into SimInfo's property list
1872 //Be careful, do not use inFileName, since it is a pointer which
1873 //point to a string at master node, and slave nodes do not contain that string
1874
1875 string zconsOutput(theInfo.finalName);
1876
1877 zconsOutput = zconsOutput.substr(0, zconsOutput.rfind(".")) + ".fz";
1878
1879 StringData* zconsFilename = new StringData();
1880 zconsFilename->setID(ZCONSFILENAME_ID);
1881 zconsFilename->setData(zconsOutput);
1882
1883 theInfo.addProperty(zconsFilename);
1884
1885 //setup index, pos and other parameters of z-constraint molecules
1886 nZConstraints = globals->getNzConstraints();
1887 theInfo.nZconstraints = nZConstraints;
1888
1889 zconStamp = globals->getZconStamp();
1890 ZConsParaItem tempParaItem;
1891
1892 ZConsParaData* zconsParaData = new ZConsParaData();
1893 zconsParaData->setID(ZCONSPARADATA_ID);
1894
1895 for (int i = 0; i < nZConstraints; i++){
1896 tempParaItem.havingZPos = zconStamp[i]->haveZpos();
1897 tempParaItem.zPos = zconStamp[i]->getZpos();
1898 tempParaItem.zconsIndex = zconStamp[i]->getMolIndex();
1899 tempParaItem.kRatio = zconStamp[i]->getKratio();
1900 tempParaItem.havingCantVel = zconStamp[i]->haveCantVel();
1901 tempParaItem.cantVel = zconStamp[i]->getCantVel();
1902 zconsParaData->addItem(tempParaItem);
1903 }
1904
1905 //check the uniqueness of index
1906 if(!zconsParaData->isIndexUnique()){
1907 sprintf(painCave.errMsg,
1908 "ZConstraint Error: molIndex is not unique!\n");
1909 painCave.isFatal = 1;
1910 simError();
1911 }
1912
1913 //sort the parameters by index of molecules
1914 zconsParaData->sortByIndex();
1915
1916 //push data into siminfo, therefore, we can retrieve later
1917 theInfo.addProperty(zconsParaData);
1918 }
1919
1920 void SimSetup::makeMinimizer(){
1921
1922 OOPSEMinimizer* myOOPSEMinimizer;
1923 MinimizerParameterSet* param;
1924 char minimizerName[100];
1925
1926 for (int i = 0; i < nInfo; i++){
1927
1928 //prepare parameter set for minimizer
1929 param = new MinimizerParameterSet();
1930 param->setDefaultParameter();
1931
1932 if (globals->haveMinimizer()){
1933 param->setFTol(globals->getMinFTol());
1934 }
1935
1936 if (globals->haveMinGTol()){
1937 param->setGTol(globals->getMinGTol());
1938 }
1939
1940 if (globals->haveMinMaxIter()){
1941 param->setMaxIteration(globals->getMinMaxIter());
1942 }
1943
1944 if (globals->haveMinWriteFrq()){
1945 param->setMaxIteration(globals->getMinMaxIter());
1946 }
1947
1948 if (globals->haveMinWriteFrq()){
1949 param->setWriteFrq(globals->getMinWriteFrq());
1950 }
1951
1952 if (globals->haveMinStepSize()){
1953 param->setStepSize(globals->getMinStepSize());
1954 }
1955
1956 if (globals->haveMinLSMaxIter()){
1957 param->setLineSearchMaxIteration(globals->getMinLSMaxIter());
1958 }
1959
1960 if (globals->haveMinLSTol()){
1961 param->setLineSearchTol(globals->getMinLSTol());
1962 }
1963
1964 strcpy(minimizerName, globals->getMinimizer());
1965
1966 if (!strcasecmp(minimizerName, "CG")){
1967 myOOPSEMinimizer = new PRCGMinimizer(&(info[i]), the_ff, param);
1968 }
1969 else if (!strcasecmp(minimizerName, "SD")){
1970 //myOOPSEMinimizer = MinimizerFactory.creatMinimizer("", &(info[i]), the_ff, param);
1971 myOOPSEMinimizer = new SDMinimizer(&(info[i]), the_ff, param);
1972 }
1973 else{
1974 sprintf(painCave.errMsg,
1975 "SimSetup error: Unrecognized Minimizer, use Conjugate Gradient \n");
1976 painCave.isFatal = 0;
1977 simError();
1978
1979 myOOPSEMinimizer = new PRCGMinimizer(&(info[i]), the_ff, param);
1980 }
1981 info[i].the_integrator = myOOPSEMinimizer;
1982
1983 //store the minimizer into simInfo
1984 info[i].the_minimizer = myOOPSEMinimizer;
1985 info[i].has_minimizer = true;
1986 }
1987
1988 }