ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-1.0/libmdtools/SimSetup.cpp
Revision: 1435
Committed: Thu Jul 29 18:16:16 2004 UTC (20 years, 9 months ago) by tim
File size: 56474 byte(s)
Log Message:
working version of simpleBuilder

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 // clean up the forcefield
100 the_ff->cleanMe();
101
102 delete stamps;
103 delete globals;
104 }
105
106 void SimSetup::setSimInfo(SimInfo* the_info, int theNinfo){
107 info = the_info;
108 nInfo = theNinfo;
109 isInfoArray = 1;
110 initSuspend = true;
111 }
112
113
114 void SimSetup::parseFile(char* fileName){
115 #ifdef IS_MPI
116 if (worldRank == 0){
117 #endif // is_mpi
118
119 inFileName = fileName;
120
121 globals->initalize();
122 set_interface_stamps(stamps, globals);
123
124 #ifdef IS_MPI
125 mpiEventInit();
126 #endif
127
128 yacc_BASS(fileName);
129
130 #ifdef IS_MPI
131 throwMPIEvent(NULL);
132 }
133 else{
134 receiveParse();
135 }
136 #endif
137
138 }
139
140 #ifdef IS_MPI
141 void SimSetup::receiveParse(void){
142 set_interface_stamps(stamps, globals);
143 mpiEventInit();
144 MPIcheckPoint();
145 mpiEventLoop();
146 }
147
148 #endif // is_mpi
149
150 void SimSetup::createSim(void){
151
152 // gather all of the information from the meta-data file
153
154 gatherInfo();
155
156 // creation of complex system objects
157
158 sysObjectsCreation();
159
160 // check on the post processing info
161
162 finalInfoCheck();
163
164 // initialize the system coordinates
165
166 if ( !initSuspend ){
167 initSystemCoords();
168
169 if( !(globals->getUseInitTime()) )
170 info[0].currentTime = 0.0;
171 }
172
173 // make the output filenames
174
175 makeOutNames();
176
177 #ifdef IS_MPI
178 mpiSim->mpiRefresh();
179 #endif
180
181 // initialize the Fortran
182
183 initFortran();
184
185 if (globals->haveMinimizer())
186 // make minimizer
187 makeMinimizer();
188 else
189 // make the integrator
190 makeIntegrator();
191
192 }
193
194
195 void SimSetup::makeMolecules(void){
196 int i, j, k;
197 int exI, exJ, exK, exL, slI, slJ;
198 int tempI, tempJ, tempK, tempL;
199 int molI, globalID;
200 int stampID, atomOffset, rbOffset, groupOffset;
201 molInit molInfo;
202 DirectionalAtom* dAtom;
203 RigidBody* myRB;
204 StuntDouble* mySD;
205 LinkedAssign* extras;
206 LinkedAssign* current_extra;
207 AtomStamp* currentAtom;
208 BondStamp* currentBond;
209 BendStamp* currentBend;
210 TorsionStamp* currentTorsion;
211 RigidBodyStamp* currentRigidBody;
212 CutoffGroupStamp* currentCutoffGroup;
213 CutoffGroup* myCutoffGroup;
214 int nCutoffGroups;// number of cutoff group of a molecule defined in mdl file
215 set<int> cutoffAtomSet; //atoms belong to cutoffgroup defined at mdl file
216
217 bond_pair* theBonds;
218 bend_set* theBends;
219 torsion_set* theTorsions;
220
221 set<int> skipList;
222
223 double phi, theta, psi;
224 char* molName;
225 char rbName[100];
226
227 int whichRigidBody;
228 int consAtomIndex; //index of constraint atom in rigid body's atom array
229 double bondLength2;
230 //init the forceField paramters
231
232 the_ff->readParams();
233
234 // init the atoms
235
236 int nMembers, nNew, rb1, rb2;
237
238 for (k = 0; k < nInfo; k++){
239 the_ff->setSimInfo(&(info[k]));
240
241 #ifdef IS_MPI
242 info[k].globalGroupMembership = new int[mpiSim->getNAtomsGlobal()];
243 for (i = 0; i < mpiSim->getNAtomsGlobal(); i++)
244 info[k].globalGroupMembership[i] = 0;
245 #else
246 info[k].globalGroupMembership = new int[info[k].n_atoms];
247 for (i = 0; i < info[k].n_atoms; i++)
248 info[k].globalGroupMembership[i] = 0;
249 #endif
250
251 atomOffset = 0;
252 groupOffset = 0;
253
254 for (i = 0; i < info[k].n_mol; i++){
255 stampID = info[k].molecules[i].getStampID();
256 molName = comp_stamps[stampID]->getID();
257
258 molInfo.nAtoms = comp_stamps[stampID]->getNAtoms();
259 molInfo.nBonds = comp_stamps[stampID]->getNBonds();
260 molInfo.nBends = comp_stamps[stampID]->getNBends();
261 molInfo.nTorsions = comp_stamps[stampID]->getNTorsions();
262 molInfo.nRigidBodies = comp_stamps[stampID]->getNRigidBodies();
263
264 nCutoffGroups = comp_stamps[stampID]->getNCutoffGroups();
265
266 molInfo.myAtoms = &(info[k].atoms[atomOffset]);
267
268 if (molInfo.nBonds > 0)
269 molInfo.myBonds = new Bond*[molInfo.nBonds];
270 else
271 molInfo.myBonds = NULL;
272
273 if (molInfo.nBends > 0)
274 molInfo.myBends = new Bend*[molInfo.nBends];
275 else
276 molInfo.myBends = NULL;
277
278 if (molInfo.nTorsions > 0)
279 molInfo.myTorsions = new Torsion *[molInfo.nTorsions];
280 else
281 molInfo.myTorsions = NULL;
282
283 theBonds = new bond_pair[molInfo.nBonds];
284 theBends = new bend_set[molInfo.nBends];
285 theTorsions = new torsion_set[molInfo.nTorsions];
286
287 // make the Atoms
288
289 for (j = 0; j < molInfo.nAtoms; j++){
290 currentAtom = comp_stamps[stampID]->getAtom(j);
291
292 if (currentAtom->haveOrientation()){
293 dAtom = new DirectionalAtom((j + atomOffset),
294 info[k].getConfiguration());
295 info[k].n_oriented++;
296 molInfo.myAtoms[j] = dAtom;
297
298 // Directional Atoms have standard unit vectors which are oriented
299 // in space using the three Euler angles. We assume the standard
300 // unit vector was originally along the z axis below.
301
302 phi = currentAtom->getEulerPhi() * M_PI / 180.0;
303 theta = currentAtom->getEulerTheta() * M_PI / 180.0;
304 psi = currentAtom->getEulerPsi()* M_PI / 180.0;
305
306 dAtom->setUnitFrameFromEuler(phi, theta, psi);
307
308 }
309 else{
310
311 molInfo.myAtoms[j] = new Atom((j + atomOffset), info[k].getConfiguration());
312
313 }
314
315 molInfo.myAtoms[j]->setType(currentAtom->getType());
316 #ifdef IS_MPI
317 molInfo.myAtoms[j]->setGlobalIndex(globalAtomIndex[j + atomOffset]);
318 #endif // is_mpi
319 }
320
321 // make the bonds
322 for (j = 0; j < molInfo.nBonds; j++){
323 currentBond = comp_stamps[stampID]->getBond(j);
324 theBonds[j].a = currentBond->getA() + atomOffset;
325 theBonds[j].b = currentBond->getB() + atomOffset;
326
327 tempI = theBonds[j].a;
328 tempJ = theBonds[j].b;
329
330 #ifdef IS_MPI
331 exI = info[k].atoms[tempI]->getGlobalIndex() + 1;
332 exJ = info[k].atoms[tempJ]->getGlobalIndex() + 1;
333 #else
334 exI = tempI + 1;
335 exJ = tempJ + 1;
336 #endif
337
338 info[k].excludes->addPair(exI, exJ);
339 }
340
341 //make the bends
342 for (j = 0; j < molInfo.nBends; j++){
343 currentBend = comp_stamps[stampID]->getBend(j);
344 theBends[j].a = currentBend->getA() + atomOffset;
345 theBends[j].b = currentBend->getB() + atomOffset;
346 theBends[j].c = currentBend->getC() + atomOffset;
347
348 if (currentBend->haveExtras()){
349 extras = currentBend->getExtras();
350 current_extra = extras;
351
352 while (current_extra != NULL){
353 if (!strcmp(current_extra->getlhs(), "ghostVectorSource")){
354 switch (current_extra->getType()){
355 case 0:
356 theBends[j].ghost = current_extra->getInt() + atomOffset;
357 theBends[j].isGhost = 1;
358 break;
359
360 case 1:
361 theBends[j].ghost = (int) current_extra->getDouble() +
362 atomOffset;
363 theBends[j].isGhost = 1;
364 break;
365
366 default:
367 sprintf(painCave.errMsg,
368 "SimSetup Error: ghostVectorSource was neither a "
369 "double nor an int.\n"
370 "-->Bend[%d] in %s\n",
371 j, comp_stamps[stampID]->getID());
372 painCave.isFatal = 1;
373 simError();
374 }
375 }
376 else{
377 sprintf(painCave.errMsg,
378 "SimSetup Error: unhandled bend assignment:\n"
379 " -->%s in Bend[%d] in %s\n",
380 current_extra->getlhs(), j, comp_stamps[stampID]->getID());
381 painCave.isFatal = 1;
382 simError();
383 }
384
385 current_extra = current_extra->getNext();
386 }
387 }
388
389 if (theBends[j].isGhost) {
390
391 tempI = theBends[j].a;
392 tempJ = theBends[j].b;
393
394 #ifdef IS_MPI
395 exI = info[k].atoms[tempI]->getGlobalIndex() + 1;
396 exJ = info[k].atoms[tempJ]->getGlobalIndex() + 1;
397 #else
398 exI = tempI + 1;
399 exJ = tempJ + 1;
400 #endif
401 info[k].excludes->addPair(exI, exJ);
402
403 } else {
404
405 tempI = theBends[j].a;
406 tempJ = theBends[j].b;
407 tempK = theBends[j].c;
408
409 #ifdef IS_MPI
410 exI = info[k].atoms[tempI]->getGlobalIndex() + 1;
411 exJ = info[k].atoms[tempJ]->getGlobalIndex() + 1;
412 exK = info[k].atoms[tempK]->getGlobalIndex() + 1;
413 #else
414 exI = tempI + 1;
415 exJ = tempJ + 1;
416 exK = tempK + 1;
417 #endif
418
419 info[k].excludes->addPair(exI, exK);
420 info[k].excludes->addPair(exI, exJ);
421 info[k].excludes->addPair(exJ, exK);
422 }
423 }
424
425 for (j = 0; j < molInfo.nTorsions; j++){
426 currentTorsion = comp_stamps[stampID]->getTorsion(j);
427 theTorsions[j].a = currentTorsion->getA() + atomOffset;
428 theTorsions[j].b = currentTorsion->getB() + atomOffset;
429 theTorsions[j].c = currentTorsion->getC() + atomOffset;
430 theTorsions[j].d = currentTorsion->getD() + atomOffset;
431
432 tempI = theTorsions[j].a;
433 tempJ = theTorsions[j].b;
434 tempK = theTorsions[j].c;
435 tempL = theTorsions[j].d;
436
437 #ifdef IS_MPI
438 exI = info[k].atoms[tempI]->getGlobalIndex() + 1;
439 exJ = info[k].atoms[tempJ]->getGlobalIndex() + 1;
440 exK = info[k].atoms[tempK]->getGlobalIndex() + 1;
441 exL = info[k].atoms[tempL]->getGlobalIndex() + 1;
442 #else
443 exI = tempI + 1;
444 exJ = tempJ + 1;
445 exK = tempK + 1;
446 exL = tempL + 1;
447 #endif
448
449 info[k].excludes->addPair(exI, exJ);
450 info[k].excludes->addPair(exI, exK);
451 info[k].excludes->addPair(exI, exL);
452 info[k].excludes->addPair(exJ, exK);
453 info[k].excludes->addPair(exJ, exL);
454 info[k].excludes->addPair(exK, exL);
455 }
456
457
458 molInfo.myRigidBodies.clear();
459
460 for (j = 0; j < molInfo.nRigidBodies; j++){
461
462 currentRigidBody = comp_stamps[stampID]->getRigidBody(j);
463 nMembers = currentRigidBody->getNMembers();
464
465 // Create the Rigid Body:
466
467 myRB = new RigidBody();
468
469 sprintf(rbName,"%s_RB_%d", molName, j);
470 myRB->setType(rbName);
471
472 for (rb1 = 0; rb1 < nMembers; rb1++) {
473
474 // molI is atom numbering inside this molecule
475 molI = currentRigidBody->getMember(rb1);
476
477 // tempI is atom numbering on local processor
478 tempI = molI + atomOffset;
479
480 // currentAtom is the AtomStamp (which we need for
481 // rigid body reference positions)
482 currentAtom = comp_stamps[stampID]->getAtom(molI);
483
484 // When we add to the rigid body, add the atom itself and
485 // the stamp info:
486
487 myRB->addAtom(info[k].atoms[tempI], currentAtom);
488
489 // Add this atom to the Skip List for the integrators
490 #ifdef IS_MPI
491 slI = info[k].atoms[tempI]->getGlobalIndex();
492 #else
493 slI = tempI;
494 #endif
495 skipList.insert(slI);
496
497 }
498
499 for(rb1 = 0; rb1 < nMembers - 1; rb1++) {
500 for(rb2 = rb1+1; rb2 < nMembers; rb2++) {
501
502 tempI = currentRigidBody->getMember(rb1);
503 tempJ = currentRigidBody->getMember(rb2);
504
505 // Some explanation is required here.
506 // Fortran indexing starts at 1, while c indexing starts at 0
507 // Also, in parallel computations, the GlobalIndex is
508 // used for the exclude list:
509
510 #ifdef IS_MPI
511 exI = molInfo.myAtoms[tempI]->getGlobalIndex() + 1;
512 exJ = molInfo.myAtoms[tempJ]->getGlobalIndex() + 1;
513 #else
514 exI = molInfo.myAtoms[tempI]->getIndex() + 1;
515 exJ = molInfo.myAtoms[tempJ]->getIndex() + 1;
516 #endif
517
518 info[k].excludes->addPair(exI, exJ);
519
520 }
521 }
522
523 molInfo.myRigidBodies.push_back(myRB);
524 info[k].rigidBodies.push_back(myRB);
525 }
526
527
528 //create cutoff group for molecule
529
530 cutoffAtomSet.clear();
531 molInfo.myCutoffGroups.clear();
532
533 for (j = 0; j < nCutoffGroups; j++){
534
535 currentCutoffGroup = comp_stamps[stampID]->getCutoffGroup(j);
536 nMembers = currentCutoffGroup->getNMembers();
537
538 myCutoffGroup = new CutoffGroup();
539
540 #ifdef IS_MPI
541 myCutoffGroup->setGlobalIndex(globalGroupIndex[groupOffset]);
542 #else
543 myCutoffGroup->setGlobalIndex(groupOffset);
544 #endif
545
546 for (int cg = 0; cg < nMembers; cg++) {
547
548 // molI is atom numbering inside this molecule
549 molI = currentCutoffGroup->getMember(cg);
550
551 // tempI is atom numbering on local processor
552 tempI = molI + atomOffset;
553
554 #ifdef IS_MPI
555 globalID = info[k].atoms[tempI]->getGlobalIndex();
556 info[k].globalGroupMembership[globalID] = globalGroupIndex[groupOffset];
557 #else
558 globalID = info[k].atoms[tempI]->getIndex();
559 info[k].globalGroupMembership[globalID] = groupOffset;
560 #endif
561 myCutoffGroup->addAtom(info[k].atoms[tempI]);
562 cutoffAtomSet.insert(tempI);
563 }
564
565 molInfo.myCutoffGroups.push_back(myCutoffGroup);
566 groupOffset++;
567
568 }//end for (j = 0; j < molInfo.nCutoffGroups; j++)
569
570
571 // create a cutoff group for every atom in current molecule which
572 // does not belong to cutoffgroup defined at mdl file
573
574 for(j = 0; j < molInfo.nAtoms; j++){
575
576 if(cutoffAtomSet.find(molInfo.myAtoms[j]->getIndex()) == cutoffAtomSet.end()){
577 myCutoffGroup = new CutoffGroup();
578 myCutoffGroup->addAtom(molInfo.myAtoms[j]);
579
580 #ifdef IS_MPI
581 myCutoffGroup->setGlobalIndex(globalGroupIndex[groupOffset]);
582 globalID = info[k].atoms[atomOffset + j]->getGlobalIndex();
583 info[k].globalGroupMembership[globalID] = globalGroupIndex[groupOffset];
584 #else
585 myCutoffGroup->setGlobalIndex(groupOffset);
586 globalID = info[k].atoms[atomOffset + j]->getIndex();
587 info[k].globalGroupMembership[globalID] = groupOffset;
588 #endif
589 molInfo.myCutoffGroups.push_back(myCutoffGroup);
590 groupOffset++;
591 }
592 }
593
594 // After this is all set up, scan through the atoms to
595 // see if they can be added to the integrableObjects:
596
597 molInfo.myIntegrableObjects.clear();
598
599
600 for (j = 0; j < molInfo.nAtoms; j++){
601
602 #ifdef IS_MPI
603 slJ = molInfo.myAtoms[j]->getGlobalIndex();
604 #else
605 slJ = j+atomOffset;
606 #endif
607
608 // if they aren't on the skip list, then they can be integrated
609
610 if (skipList.find(slJ) == skipList.end()) {
611 mySD = (StuntDouble *) molInfo.myAtoms[j];
612 info[k].integrableObjects.push_back(mySD);
613 molInfo.myIntegrableObjects.push_back(mySD);
614 }
615 }
616
617 // all rigid bodies are integrated:
618
619 for (j = 0; j < molInfo.nRigidBodies; j++) {
620 mySD = (StuntDouble *) molInfo.myRigidBodies[j];
621 info[k].integrableObjects.push_back(mySD);
622 molInfo.myIntegrableObjects.push_back(mySD);
623 }
624
625 // send the arrays off to the forceField for init.
626
627 the_ff->initializeAtoms(molInfo.nAtoms, molInfo.myAtoms);
628 the_ff->initializeBonds(molInfo.nBonds, molInfo.myBonds, theBonds);
629 the_ff->initializeBends(molInfo.nBends, molInfo.myBends, theBends);
630 the_ff->initializeTorsions(molInfo.nTorsions, molInfo.myTorsions,
631 theTorsions);
632
633 info[k].molecules[i].initialize(molInfo);
634
635
636 atomOffset += molInfo.nAtoms;
637 delete[] theBonds;
638 delete[] theBends;
639 delete[] theTorsions;
640 }
641
642
643
644 #ifdef IS_MPI
645 // Since the globalGroupMembership has been zero filled and we've only
646 // poked values into the atoms we know, we can do an Allreduce
647 // to get the full globalGroupMembership array (We think).
648 // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
649 // docs said we could.
650
651 int* ggMjunk = new int[mpiSim->getNAtomsGlobal()];
652
653 MPI_Allreduce(info[k].globalGroupMembership,
654 ggMjunk,
655 mpiSim->getNAtomsGlobal(),
656 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
657
658 for (i = 0; i < mpiSim->getNAtomsGlobal(); i++)
659 info[k].globalGroupMembership[i] = ggMjunk[i];
660
661 delete[] ggMjunk;
662
663 #endif
664
665
666
667 }
668
669 #ifdef IS_MPI
670 sprintf(checkPointMsg, "all molecules initialized succesfully");
671 MPIcheckPoint();
672 #endif // is_mpi
673
674 }
675
676 void SimSetup::gatherInfo(void){
677 int i;
678
679 ensembleCase = -1;
680 ffCase = -1;
681
682 // set the easy ones first
683
684 for (i = 0; i < nInfo; i++){
685 if (globals->haveTargetTemp()) {
686 info[i].target_temp = globals->getTargetTemp();
687 info[i].have_target_temp = 1;
688 } else {
689 info[i].have_target_temp = 0;
690 }
691 if (globals->haveDt()) {
692 info[i].dt = globals->getDt();
693 }
694 if (globals->haveRunTime()) {
695 info[i].run_time = globals->getRunTime();
696 }
697 }
698 n_components = globals->getNComponents();
699
700
701 // get the forceField
702
703 strcpy(force_field, globals->getForceField());
704
705 if (!strcasecmp(force_field, "DUFF")){
706 ffCase = FF_DUFF;
707 }
708 else if (!strcasecmp(force_field, "LJ")){
709 ffCase = FF_LJ;
710 }
711 else if (!strcasecmp(force_field, "EAM")){
712 ffCase = FF_EAM;
713 }
714 else if (!strcasecmp(force_field, "WATER")){
715 ffCase = FF_H2O;
716 }
717 else{
718 sprintf(painCave.errMsg, "SimSetup Error. Unrecognized force field -> %s\n",
719 force_field);
720 painCave.isFatal = 1;
721 simError();
722 }
723 if (globals->haveForceFieldVariant()) {
724 strcpy(forcefield_variant, globals->getForceFieldVariant());
725 has_forcefield_variant = 1;
726 }
727
728 // get the ensemble
729
730
731 if (globals->haveEnsemble()) {
732
733 strcpy(ensemble, globals->getEnsemble());
734
735 if (!strcasecmp(ensemble, "NVE")){
736 ensembleCase = NVE_ENS;
737 }
738 else if (!strcasecmp(ensemble, "NVT")){
739 ensembleCase = NVT_ENS;
740 }
741 else if (!strcasecmp(ensemble, "NPTi") || !strcasecmp(ensemble, "NPT")){
742 ensembleCase = NPTi_ENS;
743 }
744 else if (!strcasecmp(ensemble, "NPTf")){
745 ensembleCase = NPTf_ENS;
746 }
747 else if (!strcasecmp(ensemble, "NPTxyz")){
748 ensembleCase = NPTxyz_ENS;
749 }
750 else{
751 sprintf(painCave.errMsg,
752 "SimSetup Warning. Unrecognized Ensemble -> %s \n"
753 "\treverting to NVE for this simulation.\n",
754 ensemble);
755 painCave.isFatal = 0;
756 simError();
757 strcpy(ensemble, "NVE");
758 ensembleCase = NVE_ENS;
759 }
760
761 for (i = 0; i < nInfo; i++)
762 strcpy(info[i].ensemble, ensemble);
763
764
765 //check whether sample time, status time, thermal time and reset time are divisble by dt
766 if (globals->haveSampleTime() && !isDivisible(globals->getSampleTime(), globals->getDt())){
767 sprintf(painCave.errMsg,
768 "Sample time is not divisible by dt.\n"
769 "\tThis will result in samples that are not uniformly\n"
770 "\tdistributed in time. If this is a problem, change\n"
771 "\tyour sampleTime variable.\n");
772 painCave.isFatal = 0;
773 simError();
774 }
775
776 if (globals->haveStatusTime() && !isDivisible(globals->getStatusTime(), globals->getDt())){
777 sprintf(painCave.errMsg,
778 "Status time is not divisible by dt.\n"
779 "\tThis will result in status reports that are not uniformly\n"
780 "\tdistributed in time. If this is a problem, change \n"
781 "\tyour statusTime variable.\n");
782 painCave.isFatal = 0;
783 simError();
784 }
785
786 if (globals->haveThermalTime() && !isDivisible(globals->getThermalTime(), globals->getDt())){
787 sprintf(painCave.errMsg,
788 "Thermal time is not divisible by dt.\n"
789 "\tThis will result in thermalizations that are not uniformly\n"
790 "\tdistributed in time. If this is a problem, change \n"
791 "\tyour thermalTime variable.\n");
792 painCave.isFatal = 0;
793 simError();
794 }
795
796 if (globals->haveResetTime() && !isDivisible(globals->getResetTime(), globals->getDt())){
797 sprintf(painCave.errMsg,
798 "Reset time is not divisible by dt.\n"
799 "\tThis will result in integrator resets that are not uniformly\n"
800 "\tdistributed in time. If this is a problem, change\n"
801 "\tyour resetTime variable.\n");
802 painCave.isFatal = 0;
803 simError();
804 }
805
806 // set the status, sample, and thermal kick times
807
808 for (i = 0; i < nInfo; i++){
809 if (globals->haveSampleTime()){
810 info[i].sampleTime = globals->getSampleTime();
811 info[i].statusTime = info[i].sampleTime;
812 }
813 else{
814 info[i].sampleTime = globals->getRunTime();
815 info[i].statusTime = info[i].sampleTime;
816 }
817
818 if (globals->haveStatusTime()){
819 info[i].statusTime = globals->getStatusTime();
820 }
821
822 if (globals->haveThermalTime()){
823 info[i].thermalTime = globals->getThermalTime();
824 } else {
825 info[i].thermalTime = globals->getRunTime();
826 }
827
828 info[i].resetIntegrator = 0;
829 if( globals->haveResetTime() ){
830 info[i].resetTime = globals->getResetTime();
831 info[i].resetIntegrator = 1;
832 }
833 }
834
835
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(tot_nmol, 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 }