ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 1116
Committed: Thu Apr 15 22:15:21 2004 UTC (21 years ago) by tim
File size: 51625 byte(s)
Log Message:
fix a bug in setting exclude list

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