ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 1108
Committed: Wed Apr 14 15:37:41 2004 UTC (20 years, 2 months ago) by tim
File size: 51259 byte(s)
Log Message:
Change DumpWriter and InitFromFile

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