ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/brains/SimSetup.cpp
Revision: 1614
Committed: Wed Oct 20 04:55:21 2004 UTC (19 years, 8 months ago) by gezelter
File size: 56545 byte(s)
Log Message:
namespace problem prevented linking

File Contents

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