ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 682
Committed: Tue Aug 12 17:51:33 2003 UTC (20 years, 10 months ago) by tim
File size: 40284 byte(s)
Log Message:
added harmonical potential to z-constraint method

File Contents

# Content
1 #include <algorithm>
2 #include <cstdlib>
3 #include <iostream>
4 #include <cmath>
5 #include <string>
6
7 #include "SimSetup.hpp"
8 #include "parse_me.h"
9 #include "Integrator.hpp"
10 #include "simError.h"
11
12 #ifdef IS_MPI
13 #include "mpiBASS.h"
14 #include "mpiSimulation.hpp"
15 #endif
16
17 // some defines for ensemble and Forcefield cases
18
19 #define NVE_ENS 0
20 #define NVT_ENS 1
21 #define NPTi_ENS 2
22 #define NPTf_ENS 3
23 #define NPTim_ENS 4
24 #define NPTfm_ENS 5
25
26 #define FF_DUFF 0
27 #define FF_LJ 1
28 #define FF_EAM 2
29
30 using namespace std;
31
32 SimSetup::SimSetup(){
33
34 isInfoArray = 0;
35 nInfo = 1;
36
37 stamps = new MakeStamps();
38 globals = new Globals();
39
40
41 #ifdef IS_MPI
42 strcpy( checkPointMsg, "SimSetup creation successful" );
43 MPIcheckPoint();
44 #endif // IS_MPI
45 }
46
47 SimSetup::~SimSetup(){
48 delete stamps;
49 delete globals;
50 }
51
52 void SimSetup::setSimInfo( SimInfo* the_info, int theNinfo ) {
53 info = the_info;
54 nInfo = theNinfo;
55 isInfoArray = 1;
56 }
57
58
59 void SimSetup::parseFile( char* fileName ){
60
61 #ifdef IS_MPI
62 if( worldRank == 0 ){
63 #endif // is_mpi
64
65 inFileName = fileName;
66 set_interface_stamps( stamps, globals );
67
68 #ifdef IS_MPI
69 mpiEventInit();
70 #endif
71
72 yacc_BASS( fileName );
73
74 #ifdef IS_MPI
75 throwMPIEvent(NULL);
76 }
77 else receiveParse();
78 #endif
79
80 }
81
82 #ifdef IS_MPI
83 void SimSetup::receiveParse(void){
84
85 set_interface_stamps( stamps, globals );
86 mpiEventInit();
87 MPIcheckPoint();
88 mpiEventLoop();
89
90 }
91
92 #endif // is_mpi
93
94 void SimSetup::createSim(void){
95
96 int i, j, k, globalAtomIndex;
97
98 // gather all of the information from the Bass file
99
100 gatherInfo();
101
102 // creation of complex system objects
103
104 sysObjectsCreation();
105
106 // check on the post processing info
107
108 finalInfoCheck();
109
110 // initialize the system coordinates
111
112 if( !isInfoArray ) initSystemCoords();
113
114 // make the output filenames
115
116 makeOutNames();
117
118 // make the integrator
119
120 makeIntegrator();
121
122 #ifdef IS_MPI
123 mpiSim->mpiRefresh();
124 #endif
125
126 // initialize the Fortran
127
128 initFortran();
129
130
131
132 }
133
134
135 void SimSetup::makeMolecules( void ){
136
137 int k,l;
138 int i, j, exI, exJ, tempEx, stampID, atomOffset, excludeOffset;
139 molInit molInfo;
140 DirectionalAtom* dAtom;
141 LinkedAssign* extras;
142 LinkedAssign* current_extra;
143 AtomStamp* currentAtom;
144 BondStamp* currentBond;
145 BendStamp* currentBend;
146 TorsionStamp* currentTorsion;
147
148 bond_pair* theBonds;
149 bend_set* theBends;
150 torsion_set* theTorsions;
151
152
153 //init the forceField paramters
154
155 the_ff->readParams();
156
157
158 // init the atoms
159
160 double ux, uy, uz, u, uSqr;
161
162 for(k=0; k<nInfo; k++){
163
164 the_ff->setSimInfo( &(info[k]) );
165
166 atomOffset = 0;
167 excludeOffset = 0;
168 for(i=0; i<info[k].n_mol; i++){
169
170 stampID = info[k].molecules[i].getStampID();
171
172 molInfo.nAtoms = comp_stamps[stampID]->getNAtoms();
173 molInfo.nBonds = comp_stamps[stampID]->getNBonds();
174 molInfo.nBends = comp_stamps[stampID]->getNBends();
175 molInfo.nTorsions = comp_stamps[stampID]->getNTorsions();
176 molInfo.nExcludes = molInfo.nBonds + molInfo.nBends + molInfo.nTorsions;
177
178 molInfo.myAtoms = &(info[k].atoms[atomOffset]);
179 molInfo.myExcludes = &(info[k].excludes[excludeOffset]);
180 molInfo.myBonds = new Bond*[molInfo.nBonds];
181 molInfo.myBends = new Bend*[molInfo.nBends];
182 molInfo.myTorsions = new Torsion*[molInfo.nTorsions];
183
184 theBonds = new bond_pair[molInfo.nBonds];
185 theBends = new bend_set[molInfo.nBends];
186 theTorsions = new torsion_set[molInfo.nTorsions];
187
188 // make the Atoms
189
190 for(j=0; j<molInfo.nAtoms; j++){
191
192 currentAtom = comp_stamps[stampID]->getAtom( j );
193 if( currentAtom->haveOrientation() ){
194
195 dAtom = new DirectionalAtom( (j + atomOffset),
196 info[k].getConfiguration() );
197 info[k].n_oriented++;
198 molInfo.myAtoms[j] = dAtom;
199
200 ux = currentAtom->getOrntX();
201 uy = currentAtom->getOrntY();
202 uz = currentAtom->getOrntZ();
203
204 uSqr = (ux * ux) + (uy * uy) + (uz * uz);
205
206 u = sqrt( uSqr );
207 ux = ux / u;
208 uy = uy / u;
209 uz = uz / u;
210
211 dAtom->setSUx( ux );
212 dAtom->setSUy( uy );
213 dAtom->setSUz( uz );
214 }
215 else{
216 molInfo.myAtoms[j] = new GeneralAtom( (j + atomOffset),
217 info[k].getConfiguration() );
218 }
219 molInfo.myAtoms[j]->setType( currentAtom->getType() );
220
221 #ifdef IS_MPI
222
223 molInfo.myAtoms[j]->setGlobalIndex( globalIndex[j+atomOffset] );
224
225 #endif // is_mpi
226 }
227
228 // make the bonds
229 for(j=0; j<molInfo.nBonds; j++){
230
231 currentBond = comp_stamps[stampID]->getBond( j );
232 theBonds[j].a = currentBond->getA() + atomOffset;
233 theBonds[j].b = currentBond->getB() + atomOffset;
234
235 exI = theBonds[j].a;
236 exJ = theBonds[j].b;
237
238 // exclude_I must always be the smaller of the pair
239 if( exI > exJ ){
240 tempEx = exI;
241 exI = exJ;
242 exJ = tempEx;
243 }
244 #ifdef IS_MPI
245 tempEx = exI;
246 exI = info[k].atoms[tempEx]->getGlobalIndex() + 1;
247 tempEx = exJ;
248 exJ = info[k].atoms[tempEx]->getGlobalIndex() + 1;
249
250 info[k].excludes[j+excludeOffset]->setPair( exI, exJ );
251 #else // isn't MPI
252
253 info[k].excludes[j+excludeOffset]->setPair( (exI+1), (exJ+1) );
254 #endif //is_mpi
255 }
256 excludeOffset += molInfo.nBonds;
257
258 //make the bends
259 for(j=0; j<molInfo.nBends; j++){
260
261 currentBend = comp_stamps[stampID]->getBend( j );
262 theBends[j].a = currentBend->getA() + atomOffset;
263 theBends[j].b = currentBend->getB() + atomOffset;
264 theBends[j].c = currentBend->getC() + atomOffset;
265
266 if( currentBend->haveExtras() ){
267
268 extras = currentBend->getExtras();
269 current_extra = extras;
270
271 while( current_extra != NULL ){
272 if( !strcmp( current_extra->getlhs(), "ghostVectorSource" )){
273
274 switch( current_extra->getType() ){
275
276 case 0:
277 theBends[j].ghost =
278 current_extra->getInt() + atomOffset;
279 theBends[j].isGhost = 1;
280 break;
281
282 case 1:
283 theBends[j].ghost =
284 (int)current_extra->getDouble() + atomOffset;
285 theBends[j].isGhost = 1;
286 break;
287
288 default:
289 sprintf( painCave.errMsg,
290 "SimSetup Error: ghostVectorSource was neither a "
291 "double nor an int.\n"
292 "-->Bend[%d] in %s\n",
293 j, comp_stamps[stampID]->getID() );
294 painCave.isFatal = 1;
295 simError();
296 }
297 }
298
299 else{
300
301 sprintf( painCave.errMsg,
302 "SimSetup Error: unhandled bend assignment:\n"
303 " -->%s in Bend[%d] in %s\n",
304 current_extra->getlhs(),
305 j, comp_stamps[stampID]->getID() );
306 painCave.isFatal = 1;
307 simError();
308 }
309
310 current_extra = current_extra->getNext();
311 }
312 }
313
314 if( !theBends[j].isGhost ){
315
316 exI = theBends[j].a;
317 exJ = theBends[j].c;
318 }
319 else{
320
321 exI = theBends[j].a;
322 exJ = theBends[j].b;
323 }
324
325 // exclude_I must always be the smaller of the pair
326 if( exI > exJ ){
327 tempEx = exI;
328 exI = exJ;
329 exJ = tempEx;
330 }
331 #ifdef IS_MPI
332 tempEx = exI;
333 exI = info[k].atoms[tempEx]->getGlobalIndex() + 1;
334 tempEx = exJ;
335 exJ = info[k].atoms[tempEx]->getGlobalIndex() + 1;
336
337 info[k].excludes[j+excludeOffset]->setPair( exI, exJ );
338 #else // isn't MPI
339 info[k].excludes[j+excludeOffset]->setPair( (exI+1), (exJ+1) );
340 #endif //is_mpi
341 }
342 excludeOffset += molInfo.nBends;
343
344 for(j=0; j<molInfo.nTorsions; j++){
345
346 currentTorsion = comp_stamps[stampID]->getTorsion( j );
347 theTorsions[j].a = currentTorsion->getA() + atomOffset;
348 theTorsions[j].b = currentTorsion->getB() + atomOffset;
349 theTorsions[j].c = currentTorsion->getC() + atomOffset;
350 theTorsions[j].d = currentTorsion->getD() + atomOffset;
351
352 exI = theTorsions[j].a;
353 exJ = theTorsions[j].d;
354
355 // exclude_I must always be the smaller of the pair
356 if( exI > exJ ){
357 tempEx = exI;
358 exI = exJ;
359 exJ = tempEx;
360 }
361 #ifdef IS_MPI
362 tempEx = exI;
363 exI = info[k].atoms[tempEx]->getGlobalIndex() + 1;
364 tempEx = exJ;
365 exJ = info[k].atoms[tempEx]->getGlobalIndex() + 1;
366
367 info[k].excludes[j+excludeOffset]->setPair( exI, exJ );
368 #else // isn't MPI
369 info[k].excludes[j+excludeOffset]->setPair( (exI+1), (exJ+1) );
370 #endif //is_mpi
371 }
372 excludeOffset += molInfo.nTorsions;
373
374
375 // send the arrays off to the forceField for init.
376
377 the_ff->initializeAtoms( molInfo.nAtoms, molInfo.myAtoms );
378 the_ff->initializeBonds( molInfo.nBonds, molInfo.myBonds, theBonds );
379 the_ff->initializeBends( molInfo.nBends, molInfo.myBends, theBends );
380 the_ff->initializeTorsions( molInfo.nTorsions, molInfo.myTorsions, theTorsions );
381
382
383 info[k].molecules[i].initialize( molInfo );
384
385
386 atomOffset += molInfo.nAtoms;
387 delete[] theBonds;
388 delete[] theBends;
389 delete[] theTorsions;
390 }
391 }
392
393 #ifdef IS_MPI
394 sprintf( checkPointMsg, "all molecules initialized succesfully" );
395 MPIcheckPoint();
396 #endif // is_mpi
397
398 // clean up the forcefield
399
400 the_ff->calcRcut();
401 the_ff->cleanMe();
402
403 }
404
405 void SimSetup::initFromBass( void ){
406
407 int i, j, k;
408 int n_cells;
409 double cellx, celly, cellz;
410 double temp1, temp2, temp3;
411 int n_per_extra;
412 int n_extra;
413 int have_extra, done;
414
415 double vel[3];
416 vel[0] = 0.0;
417 vel[1] = 0.0;
418 vel[2] = 0.0;
419
420 temp1 = (double)tot_nmol / 4.0;
421 temp2 = pow( temp1, ( 1.0 / 3.0 ) );
422 temp3 = ceil( temp2 );
423
424 have_extra =0;
425 if( temp2 < temp3 ){ // we have a non-complete lattice
426 have_extra =1;
427
428 n_cells = (int)temp3 - 1;
429 cellx = info[0].boxL[0] / temp3;
430 celly = info[0].boxL[1] / temp3;
431 cellz = info[0].boxL[2] / temp3;
432 n_extra = tot_nmol - ( 4 * n_cells * n_cells * n_cells );
433 temp1 = ((double)n_extra) / ( pow( temp3, 3.0 ) - pow( n_cells, 3.0 ) );
434 n_per_extra = (int)ceil( temp1 );
435
436 if( n_per_extra > 4){
437 sprintf( painCave.errMsg,
438 "SimSetup error. There has been an error in constructing"
439 " the non-complete lattice.\n" );
440 painCave.isFatal = 1;
441 simError();
442 }
443 }
444 else{
445 n_cells = (int)temp3;
446 cellx = info[0].boxL[0] / temp3;
447 celly = info[0].boxL[1] / temp3;
448 cellz = info[0].boxL[2] / temp3;
449 }
450
451 current_mol = 0;
452 current_comp_mol = 0;
453 current_comp = 0;
454 current_atom_ndx = 0;
455
456 for( i=0; i < n_cells ; i++ ){
457 for( j=0; j < n_cells; j++ ){
458 for( k=0; k < n_cells; k++ ){
459
460 makeElement( i * cellx,
461 j * celly,
462 k * cellz );
463
464 makeElement( i * cellx + 0.5 * cellx,
465 j * celly + 0.5 * celly,
466 k * cellz );
467
468 makeElement( i * cellx,
469 j * celly + 0.5 * celly,
470 k * cellz + 0.5 * cellz );
471
472 makeElement( i * cellx + 0.5 * cellx,
473 j * celly,
474 k * cellz + 0.5 * cellz );
475 }
476 }
477 }
478
479 if( have_extra ){
480 done = 0;
481
482 int start_ndx;
483 for( i=0; i < (n_cells+1) && !done; i++ ){
484 for( j=0; j < (n_cells+1) && !done; j++ ){
485
486 if( i < n_cells ){
487
488 if( j < n_cells ){
489 start_ndx = n_cells;
490 }
491 else start_ndx = 0;
492 }
493 else start_ndx = 0;
494
495 for( k=start_ndx; k < (n_cells+1) && !done; k++ ){
496
497 makeElement( i * cellx,
498 j * celly,
499 k * cellz );
500 done = ( current_mol >= tot_nmol );
501
502 if( !done && n_per_extra > 1 ){
503 makeElement( i * cellx + 0.5 * cellx,
504 j * celly + 0.5 * celly,
505 k * cellz );
506 done = ( current_mol >= tot_nmol );
507 }
508
509 if( !done && n_per_extra > 2){
510 makeElement( i * cellx,
511 j * celly + 0.5 * celly,
512 k * cellz + 0.5 * cellz );
513 done = ( current_mol >= tot_nmol );
514 }
515
516 if( !done && n_per_extra > 3){
517 makeElement( i * cellx + 0.5 * cellx,
518 j * celly,
519 k * cellz + 0.5 * cellz );
520 done = ( current_mol >= tot_nmol );
521 }
522 }
523 }
524 }
525 }
526
527 for( i=0; i<info[0].n_atoms; i++ ){
528 info[0].atoms[i]->setVel( vel );
529 }
530 }
531
532 void SimSetup::makeElement( double x, double y, double z ){
533
534 int k;
535 AtomStamp* current_atom;
536 DirectionalAtom* dAtom;
537 double rotMat[3][3];
538 double pos[3];
539
540 for( k=0; k<comp_stamps[current_comp]->getNAtoms(); k++ ){
541
542 current_atom = comp_stamps[current_comp]->getAtom( k );
543 if( !current_atom->havePosition() ){
544 sprintf( painCave.errMsg,
545 "SimSetup:initFromBass error.\n"
546 "\tComponent %s, atom %s does not have a position specified.\n"
547 "\tThe initialization routine is unable to give a start"
548 " position.\n",
549 comp_stamps[current_comp]->getID(),
550 current_atom->getType() );
551 painCave.isFatal = 1;
552 simError();
553 }
554
555 pos[0] = x + current_atom->getPosX();
556 pos[1] = y + current_atom->getPosY();
557 pos[2] = z + current_atom->getPosZ();
558
559 info[0].atoms[current_atom_ndx]->setPos( pos );
560
561 if( info[0].atoms[current_atom_ndx]->isDirectional() ){
562
563 dAtom = (DirectionalAtom *)info[0].atoms[current_atom_ndx];
564
565 rotMat[0][0] = 1.0;
566 rotMat[0][1] = 0.0;
567 rotMat[0][2] = 0.0;
568
569 rotMat[1][0] = 0.0;
570 rotMat[1][1] = 1.0;
571 rotMat[1][2] = 0.0;
572
573 rotMat[2][0] = 0.0;
574 rotMat[2][1] = 0.0;
575 rotMat[2][2] = 1.0;
576
577 dAtom->setA( rotMat );
578 }
579
580 current_atom_ndx++;
581 }
582
583 current_mol++;
584 current_comp_mol++;
585
586 if( current_comp_mol >= components_nmol[current_comp] ){
587
588 current_comp_mol = 0;
589 current_comp++;
590 }
591 }
592
593
594 void SimSetup::gatherInfo( void ){
595 int i,j,k;
596
597 ensembleCase = -1;
598 ffCase = -1;
599
600 // set the easy ones first
601
602 for( i=0; i<nInfo; i++){
603 info[i].target_temp = globals->getTargetTemp();
604 info[i].dt = globals->getDt();
605 info[i].run_time = globals->getRunTime();
606 }
607 n_components = globals->getNComponents();
608
609
610 // get the forceField
611
612 strcpy( force_field, globals->getForceField() );
613
614 if( !strcasecmp( force_field, "DUFF" )) ffCase = FF_DUFF;
615 else if( !strcasecmp( force_field, "LJ" )) ffCase = FF_LJ;
616 else if( !strcasecmp( force_field, "EAM" )) ffCase = FF_EAM;
617 else{
618 sprintf( painCave.errMsg,
619 "SimSetup Error. Unrecognized force field -> %s\n",
620 force_field );
621 painCave.isFatal = 1;
622 simError();
623 }
624
625 // get the ensemble
626
627 strcpy( ensemble, globals->getEnsemble() );
628
629 if( !strcasecmp( ensemble, "NVE" )) ensembleCase = NVE_ENS;
630 else if( !strcasecmp( ensemble, "NVT" )) ensembleCase = NVT_ENS;
631 else if( !strcasecmp( ensemble, "NPTi" ) || !strcasecmp( ensemble, "NPT") )
632 ensembleCase = NPTi_ENS;
633 else if( !strcasecmp( ensemble, "NPTf" )) ensembleCase = NPTf_ENS;
634 else if( !strcasecmp( ensemble, "NPTim" )) ensembleCase = NPTim_ENS;
635 else if( !strcasecmp( ensemble, "NPTfm" )) ensembleCase = NPTfm_ENS;
636 else{
637 sprintf( painCave.errMsg,
638 "SimSetup Warning. Unrecognized Ensemble -> %s, "
639 "reverting to NVE for this simulation.\n",
640 ensemble );
641 painCave.isFatal = 0;
642 simError();
643 strcpy( ensemble, "NVE" );
644 ensembleCase = NVE_ENS;
645 }
646
647 for(i=0; i<nInfo; i++){
648
649 strcpy( info[i].ensemble, ensemble );
650
651 // get the mixing rule
652
653 strcpy( info[i].mixingRule, globals->getMixingRule() );
654 info[i].usePBC = globals->getPBC();
655 }
656
657 // get the components and calculate the tot_nMol and indvidual n_mol
658
659 the_components = globals->getComponents();
660 components_nmol = new int[n_components];
661
662
663 if( !globals->haveNMol() ){
664 // we don't have the total number of molecules, so we assume it is
665 // given in each component
666
667 tot_nmol = 0;
668 for( i=0; i<n_components; i++ ){
669
670 if( !the_components[i]->haveNMol() ){
671 // we have a problem
672 sprintf( painCave.errMsg,
673 "SimSetup Error. No global NMol or component NMol"
674 " given. Cannot calculate the number of atoms.\n" );
675 painCave.isFatal = 1;
676 simError();
677 }
678
679 tot_nmol += the_components[i]->getNMol();
680 components_nmol[i] = the_components[i]->getNMol();
681 }
682 }
683 else{
684 sprintf( painCave.errMsg,
685 "SimSetup error.\n"
686 "\tSorry, the ability to specify total"
687 " nMols and then give molfractions in the components\n"
688 "\tis not currently supported."
689 " Please give nMol in the components.\n" );
690 painCave.isFatal = 1;
691 simError();
692 }
693
694 // set the status, sample, and thermal kick times
695
696 for(i=0; i<nInfo; i++){
697
698 if( globals->haveSampleTime() ){
699 info[i].sampleTime = globals->getSampleTime();
700 info[i].statusTime = info[i].sampleTime;
701 info[i].thermalTime = info[i].sampleTime;
702 }
703 else{
704 info[i].sampleTime = globals->getRunTime();
705 info[i].statusTime = info[i].sampleTime;
706 info[i].thermalTime = info[i].sampleTime;
707 }
708
709 if( globals->haveStatusTime() ){
710 info[i].statusTime = globals->getStatusTime();
711 }
712
713 if( globals->haveThermalTime() ){
714 info[i].thermalTime = globals->getThermalTime();
715 }
716
717 // check for the temperature set flag
718
719 if( globals->haveTempSet() ) info[i].setTemp = globals->getTempSet();
720
721 // get some of the tricky things that may still be in the globals
722
723 double boxVector[3];
724 if( globals->haveBox() ){
725 boxVector[0] = globals->getBox();
726 boxVector[1] = globals->getBox();
727 boxVector[2] = globals->getBox();
728
729 info[i].setBox( boxVector );
730 }
731 else if( globals->haveDensity() ){
732
733 double vol;
734 vol = (double)tot_nmol / globals->getDensity();
735 boxVector[0] = pow( vol, ( 1.0 / 3.0 ) );
736 boxVector[1] = boxVector[0];
737 boxVector[2] = boxVector[0];
738
739 info[i].setBox( boxVector );
740 }
741 else{
742 if( !globals->haveBoxX() ){
743 sprintf( painCave.errMsg,
744 "SimSetup error, no periodic BoxX size given.\n" );
745 painCave.isFatal = 1;
746 simError();
747 }
748 boxVector[0] = globals->getBoxX();
749
750 if( !globals->haveBoxY() ){
751 sprintf( painCave.errMsg,
752 "SimSetup error, no periodic BoxY size given.\n" );
753 painCave.isFatal = 1;
754 simError();
755 }
756 boxVector[1] = globals->getBoxY();
757
758 if( !globals->haveBoxZ() ){
759 sprintf( painCave.errMsg,
760 "SimSetup error, no periodic BoxZ size given.\n" );
761 painCave.isFatal = 1;
762 simError();
763 }
764 boxVector[2] = globals->getBoxZ();
765
766 info[i].setBox( boxVector );
767 }
768
769 }
770
771 #ifdef IS_MPI
772 strcpy( checkPointMsg, "Succesfully gathered all information from Bass\n" );
773 MPIcheckPoint();
774 #endif // is_mpi
775
776 }
777
778
779 void SimSetup::finalInfoCheck( void ){
780 int index;
781 int usesDipoles;
782 int i;
783
784 for(i=0; i<nInfo; i++){
785 // check electrostatic parameters
786
787 index = 0;
788 usesDipoles = 0;
789 while( (index < info[i].n_atoms) && !usesDipoles ){
790 usesDipoles = (info[i].atoms[index])->hasDipole();
791 index++;
792 }
793
794 #ifdef IS_MPI
795 int myUse = usesDipoles;
796 MPI_Allreduce( &myUse, &usesDipoles, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD );
797 #endif //is_mpi
798
799 double theEcr, theEst;
800
801 if (globals->getUseRF() ) {
802 info[i].useReactionField = 1;
803
804 if( !globals->haveECR() ){
805 sprintf( painCave.errMsg,
806 "SimSetup Warning: using default value of 1/2 the smallest "
807 "box length for the electrostaticCutoffRadius.\n"
808 "I hope you have a very fast processor!\n");
809 painCave.isFatal = 0;
810 simError();
811 double smallest;
812 smallest = info[i].boxL[0];
813 if (info[i].boxL[1] <= smallest) smallest = info[i].boxL[1];
814 if (info[i].boxL[2] <= smallest) smallest = info[i].boxL[2];
815 theEcr = 0.5 * smallest;
816 } else {
817 theEcr = globals->getECR();
818 }
819
820 if( !globals->haveEST() ){
821 sprintf( painCave.errMsg,
822 "SimSetup Warning: using default value of 0.05 * the "
823 "electrostaticCutoffRadius for the electrostaticSkinThickness\n"
824 );
825 painCave.isFatal = 0;
826 simError();
827 theEst = 0.05 * theEcr;
828 } else {
829 theEst= globals->getEST();
830 }
831
832 info[i].setEcr( theEcr, theEst );
833
834 if(!globals->haveDielectric() ){
835 sprintf( painCave.errMsg,
836 "SimSetup Error: You are trying to use Reaction Field without"
837 "setting a dielectric constant!\n"
838 );
839 painCave.isFatal = 1;
840 simError();
841 }
842 info[i].dielectric = globals->getDielectric();
843 }
844 else {
845 if (usesDipoles) {
846
847 if( !globals->haveECR() ){
848 sprintf( painCave.errMsg,
849 "SimSetup Warning: using default value of 1/2 the smallest "
850 "box length for the electrostaticCutoffRadius.\n"
851 "I hope you have a very fast processor!\n");
852 painCave.isFatal = 0;
853 simError();
854 double smallest;
855 smallest = info[i].boxL[0];
856 if (info[i].boxL[1] <= smallest) smallest = info[i].boxL[1];
857 if (info[i].boxL[2] <= smallest) smallest = info[i].boxL[2];
858 theEcr = 0.5 * smallest;
859 } else {
860 theEcr = globals->getECR();
861 }
862
863 if( !globals->haveEST() ){
864 sprintf( painCave.errMsg,
865 "SimSetup Warning: using default value of 0.05 * the "
866 "electrostaticCutoffRadius for the "
867 "electrostaticSkinThickness\n"
868 );
869 painCave.isFatal = 0;
870 simError();
871 theEst = 0.05 * theEcr;
872 } else {
873 theEst= globals->getEST();
874 }
875
876 info[i].setEcr( theEcr, theEst );
877 }
878 }
879 }
880
881 #ifdef IS_MPI
882 strcpy( checkPointMsg, "post processing checks out" );
883 MPIcheckPoint();
884 #endif // is_mpi
885
886 }
887
888 void SimSetup::initSystemCoords( void ){
889 int i;
890
891 std::cerr << "Setting atom Coords\n";
892
893 (info[0].getConfiguration())->createArrays( info[0].n_atoms );
894
895 for(i=0; i<info[0].n_atoms; i++) info[0].atoms[i]->setCoords();
896
897 if( globals->haveInitialConfig() ){
898
899 InitializeFromFile* fileInit;
900 #ifdef IS_MPI // is_mpi
901 if( worldRank == 0 ){
902 #endif //is_mpi
903 fileInit = new InitializeFromFile( globals->getInitialConfig() );
904 #ifdef IS_MPI
905 }else fileInit = new InitializeFromFile( NULL );
906 #endif
907 fileInit->readInit( info ); // default velocities on
908
909 delete fileInit;
910 }
911 else{
912
913 #ifdef IS_MPI
914
915 // no init from bass
916
917 sprintf( painCave.errMsg,
918 "Cannot intialize a parallel simulation without an initial configuration file.\n" );
919 painCave.isFatal;
920 simError();
921
922 #else
923
924 initFromBass();
925
926
927 #endif
928 }
929
930 #ifdef IS_MPI
931 strcpy( checkPointMsg, "Successfully read in the initial configuration" );
932 MPIcheckPoint();
933 #endif // is_mpi
934
935 }
936
937
938 void SimSetup::makeOutNames( void ){
939
940 int k;
941
942
943 for(k=0; k<nInfo; k++){
944
945 #ifdef IS_MPI
946 if( worldRank == 0 ){
947 #endif // is_mpi
948
949 if( globals->haveFinalConfig() ){
950 strcpy( info[k].finalName, globals->getFinalConfig() );
951 }
952 else{
953 strcpy( info[k].finalName, inFileName );
954 char* endTest;
955 int nameLength = strlen( info[k].finalName );
956 endTest = &(info[k].finalName[nameLength - 5]);
957 if( !strcmp( endTest, ".bass" ) ){
958 strcpy( endTest, ".eor" );
959 }
960 else if( !strcmp( endTest, ".BASS" ) ){
961 strcpy( endTest, ".eor" );
962 }
963 else{
964 endTest = &(info[k].finalName[nameLength - 4]);
965 if( !strcmp( endTest, ".bss" ) ){
966 strcpy( endTest, ".eor" );
967 }
968 else if( !strcmp( endTest, ".mdl" ) ){
969 strcpy( endTest, ".eor" );
970 }
971 else{
972 strcat( info[k].finalName, ".eor" );
973 }
974 }
975 }
976
977 // make the sample and status out names
978
979 strcpy( info[k].sampleName, inFileName );
980 char* endTest;
981 int nameLength = strlen( info[k].sampleName );
982 endTest = &(info[k].sampleName[nameLength - 5]);
983 if( !strcmp( endTest, ".bass" ) ){
984 strcpy( endTest, ".dump" );
985 }
986 else if( !strcmp( endTest, ".BASS" ) ){
987 strcpy( endTest, ".dump" );
988 }
989 else{
990 endTest = &(info[k].sampleName[nameLength - 4]);
991 if( !strcmp( endTest, ".bss" ) ){
992 strcpy( endTest, ".dump" );
993 }
994 else if( !strcmp( endTest, ".mdl" ) ){
995 strcpy( endTest, ".dump" );
996 }
997 else{
998 strcat( info[k].sampleName, ".dump" );
999 }
1000 }
1001
1002 strcpy( info[k].statusName, inFileName );
1003 nameLength = strlen( info[k].statusName );
1004 endTest = &(info[k].statusName[nameLength - 5]);
1005 if( !strcmp( endTest, ".bass" ) ){
1006 strcpy( endTest, ".stat" );
1007 }
1008 else if( !strcmp( endTest, ".BASS" ) ){
1009 strcpy( endTest, ".stat" );
1010 }
1011 else{
1012 endTest = &(info[k].statusName[nameLength - 4]);
1013 if( !strcmp( endTest, ".bss" ) ){
1014 strcpy( endTest, ".stat" );
1015 }
1016 else if( !strcmp( endTest, ".mdl" ) ){
1017 strcpy( endTest, ".stat" );
1018 }
1019 else{
1020 strcat( info[k].statusName, ".stat" );
1021 }
1022 }
1023
1024 #ifdef IS_MPI
1025 }
1026 #endif // is_mpi
1027 }
1028 }
1029
1030
1031 void SimSetup::sysObjectsCreation( void ){
1032
1033 int i,k;
1034
1035 // create the forceField
1036
1037 createFF();
1038
1039 // extract componentList
1040
1041 compList();
1042
1043 // calc the number of atoms, bond, bends, and torsions
1044
1045 calcSysValues();
1046
1047 #ifdef IS_MPI
1048 // divide the molecules among the processors
1049
1050 mpiMolDivide();
1051 #endif //is_mpi
1052
1053 // create the atom and SRI arrays. Also initialize Molecule Stamp ID's
1054
1055 makeSysArrays();
1056
1057 // make and initialize the molecules (all but atomic coordinates)
1058
1059 makeMolecules();
1060
1061 for(k=0; k<nInfo; k++){
1062 info[k].identArray = new int[info[k].n_atoms];
1063 for(i=0; i<info[k].n_atoms; i++){
1064 info[k].identArray[i] = info[k].atoms[i]->getIdent();
1065 }
1066 }
1067 }
1068
1069
1070 void SimSetup::createFF( void ){
1071
1072 switch( ffCase ){
1073
1074 case FF_DUFF:
1075 the_ff = new DUFF();
1076 break;
1077
1078 case FF_LJ:
1079 the_ff = new LJFF();
1080 break;
1081
1082 case FF_EAM:
1083 the_ff = new EAM_FF();
1084 break;
1085
1086 default:
1087 sprintf( painCave.errMsg,
1088 "SimSetup Error. Unrecognized force field in case statement.\n");
1089 painCave.isFatal = 1;
1090 simError();
1091 }
1092
1093 #ifdef IS_MPI
1094 strcpy( checkPointMsg, "ForceField creation successful" );
1095 MPIcheckPoint();
1096 #endif // is_mpi
1097
1098 }
1099
1100
1101 void SimSetup::compList( void ){
1102
1103 int i;
1104 char* id;
1105 LinkedMolStamp* headStamp = new LinkedMolStamp();
1106 LinkedMolStamp* currentStamp = NULL;
1107 comp_stamps = new MoleculeStamp*[n_components];
1108
1109 // make an array of molecule stamps that match the components used.
1110 // also extract the used stamps out into a separate linked list
1111
1112 for(i=0; i<nInfo; i++){
1113 info[i].nComponents = n_components;
1114 info[i].componentsNmol = components_nmol;
1115 info[i].compStamps = comp_stamps;
1116 info[i].headStamp = headStamp;
1117 }
1118
1119
1120 for( i=0; i<n_components; i++ ){
1121
1122 id = the_components[i]->getType();
1123 comp_stamps[i] = NULL;
1124
1125 // check to make sure the component isn't already in the list
1126
1127 comp_stamps[i] = headStamp->match( id );
1128 if( comp_stamps[i] == NULL ){
1129
1130 // extract the component from the list;
1131
1132 currentStamp = stamps->extractMolStamp( id );
1133 if( currentStamp == NULL ){
1134 sprintf( painCave.errMsg,
1135 "SimSetup error: Component \"%s\" was not found in the "
1136 "list of declared molecules\n",
1137 id );
1138 painCave.isFatal = 1;
1139 simError();
1140 }
1141
1142 headStamp->add( currentStamp );
1143 comp_stamps[i] = headStamp->match( id );
1144 }
1145 }
1146
1147 #ifdef IS_MPI
1148 strcpy( checkPointMsg, "Component stamps successfully extracted\n" );
1149 MPIcheckPoint();
1150 #endif // is_mpi
1151
1152
1153 }
1154
1155 void SimSetup::calcSysValues( void ){
1156 int i, j, k;
1157
1158 int *molMembershipArray;
1159
1160 tot_atoms = 0;
1161 tot_bonds = 0;
1162 tot_bends = 0;
1163 tot_torsions = 0;
1164 for( i=0; i<n_components; i++ ){
1165
1166 tot_atoms += components_nmol[i] * comp_stamps[i]->getNAtoms();
1167 tot_bonds += components_nmol[i] * comp_stamps[i]->getNBonds();
1168 tot_bends += components_nmol[i] * comp_stamps[i]->getNBends();
1169 tot_torsions += components_nmol[i] * comp_stamps[i]->getNTorsions();
1170 }
1171
1172 tot_SRI = tot_bonds + tot_bends + tot_torsions;
1173 molMembershipArray = new int[tot_atoms];
1174
1175 for(i=0; i<nInfo; i++){
1176 info[i].n_atoms = tot_atoms;
1177 info[i].n_bonds = tot_bonds;
1178 info[i].n_bends = tot_bends;
1179 info[i].n_torsions = tot_torsions;
1180 info[i].n_SRI = tot_SRI;
1181 info[i].n_mol = tot_nmol;
1182
1183 info[i].molMembershipArray = molMembershipArray;
1184 }
1185 }
1186
1187 #ifdef IS_MPI
1188
1189 void SimSetup::mpiMolDivide( void ){
1190
1191 int i, j, k;
1192 int localMol, allMol;
1193 int local_atoms, local_bonds, local_bends, local_torsions, local_SRI;
1194
1195 mpiSim = new mpiSimulation( info );
1196
1197 globalIndex = mpiSim->divideLabor();
1198
1199 // set up the local variables
1200
1201 mol2proc = mpiSim->getMolToProcMap();
1202 molCompType = mpiSim->getMolComponentType();
1203
1204 allMol = 0;
1205 localMol = 0;
1206 local_atoms = 0;
1207 local_bonds = 0;
1208 local_bends = 0;
1209 local_torsions = 0;
1210 globalAtomIndex = 0;
1211
1212
1213 for( i=0; i<n_components; i++ ){
1214
1215 for( j=0; j<components_nmol[i]; j++ ){
1216
1217 if( mol2proc[allMol] == worldRank ){
1218
1219 local_atoms += comp_stamps[i]->getNAtoms();
1220 local_bonds += comp_stamps[i]->getNBonds();
1221 local_bends += comp_stamps[i]->getNBends();
1222 local_torsions += comp_stamps[i]->getNTorsions();
1223 localMol++;
1224 }
1225 for (k = 0; k < comp_stamps[i]->getNAtoms(); k++) {
1226 info[0].molMembershipArray[globalAtomIndex] = allMol;
1227 globalAtomIndex++;
1228 }
1229
1230 allMol++;
1231 }
1232 }
1233 local_SRI = local_bonds + local_bends + local_torsions;
1234
1235 info[0].n_atoms = mpiSim->getMyNlocal();
1236
1237 if( local_atoms != info[0].n_atoms ){
1238 sprintf( painCave.errMsg,
1239 "SimSetup error: mpiSim's localAtom (%d) and SimSetup's"
1240 " localAtom (%d) are not equal.\n",
1241 info[0].n_atoms,
1242 local_atoms );
1243 painCave.isFatal = 1;
1244 simError();
1245 }
1246
1247 info[0].n_bonds = local_bonds;
1248 info[0].n_bends = local_bends;
1249 info[0].n_torsions = local_torsions;
1250 info[0].n_SRI = local_SRI;
1251 info[0].n_mol = localMol;
1252
1253 strcpy( checkPointMsg, "Passed nlocal consistency check." );
1254 MPIcheckPoint();
1255 }
1256
1257 #endif // is_mpi
1258
1259
1260 void SimSetup::makeSysArrays( void ){
1261 int i, j, k, l;
1262
1263 Atom** the_atoms;
1264 Molecule* the_molecules;
1265 Exclude** the_excludes;
1266
1267
1268 for(l=0; l<nInfo; l++){
1269
1270 // create the atom and short range interaction arrays
1271
1272 the_atoms = new Atom*[info[l].n_atoms];
1273 the_molecules = new Molecule[info[l].n_mol];
1274 int molIndex;
1275
1276 // initialize the molecule's stampID's
1277
1278 #ifdef IS_MPI
1279
1280
1281 molIndex = 0;
1282 for(i=0; i<mpiSim->getTotNmol(); i++){
1283
1284 if(mol2proc[i] == worldRank ){
1285 the_molecules[molIndex].setStampID( molCompType[i] );
1286 the_molecules[molIndex].setMyIndex( molIndex );
1287 the_molecules[molIndex].setGlobalIndex( i );
1288 molIndex++;
1289 }
1290 }
1291
1292 #else // is_mpi
1293
1294 molIndex = 0;
1295 globalAtomIndex = 0;
1296 for(i=0; i<n_components; i++){
1297 for(j=0; j<components_nmol[i]; j++ ){
1298 the_molecules[molIndex].setStampID( i );
1299 the_molecules[molIndex].setMyIndex( molIndex );
1300 the_molecules[molIndex].setGlobalIndex( molIndex );
1301 for (k = 0; k < comp_stamps[i]->getNAtoms(); k++) {
1302 info[l].molMembershipArray[globalAtomIndex] = molIndex;
1303 globalAtomIndex++;
1304 }
1305 molIndex++;
1306 }
1307 }
1308
1309
1310 #endif // is_mpi
1311
1312
1313 if( info[l].n_SRI ){
1314
1315 Exclude::createArray(info[l].n_SRI);
1316 the_excludes = new Exclude*[info[l].n_SRI];
1317 for( int ex=0; ex<info[l].n_SRI; ex++){
1318 the_excludes[ex] = new Exclude(ex);
1319 }
1320 info[l].globalExcludes = new int;
1321 info[l].n_exclude = info[l].n_SRI;
1322 }
1323 else{
1324
1325 Exclude::createArray( 1 );
1326 the_excludes = new Exclude*;
1327 the_excludes[0] = new Exclude(0);
1328 the_excludes[0]->setPair( 0,0 );
1329 info[l].globalExcludes = new int;
1330 info[l].globalExcludes[0] = 0;
1331 info[l].n_exclude = 0;
1332 }
1333
1334 // set the arrays into the SimInfo object
1335
1336 info[l].atoms = the_atoms;
1337 info[l].molecules = the_molecules;
1338 info[l].nGlobalExcludes = 0;
1339 info[l].excludes = the_excludes;
1340
1341 the_ff->setSimInfo( info );
1342
1343 }
1344 }
1345
1346 void SimSetup::makeIntegrator( void ){
1347
1348 int k;
1349
1350 NVT<RealIntegrator>* myNVT = NULL;
1351 NPTi<RealIntegrator>* myNPTi = NULL;
1352 NPTf<RealIntegrator>* myNPTf = NULL;
1353 NPTim<RealIntegrator>* myNPTim = NULL;
1354 NPTfm<RealIntegrator>* myNPTfm = NULL;
1355
1356 for(k=0; k<nInfo; k++){
1357
1358 switch( ensembleCase ){
1359
1360 case NVE_ENS:
1361 if (globals->haveZconstraints()){
1362 setupZConstraint(info[k]);
1363 new ZConstraint<NVE<RealIntegrator> >( &(info[k]), the_ff );
1364 }
1365
1366 else
1367 new NVE<RealIntegrator>( &(info[k]), the_ff );
1368 break;
1369
1370 case NVT_ENS:
1371 if (globals->haveZconstraints()){
1372 setupZConstraint(info[k]);
1373 myNVT = new ZConstraint<NVT<RealIntegrator> >( &(info[k]), the_ff );
1374 }
1375 else
1376 myNVT = new NVT<RealIntegrator>( &(info[k]), the_ff );
1377
1378 myNVT->setTargetTemp(globals->getTargetTemp());
1379
1380 if (globals->haveTauThermostat())
1381 myNVT->setTauThermostat(globals->getTauThermostat());
1382
1383 else {
1384 sprintf( painCave.errMsg,
1385 "SimSetup error: If you use the NVT\n"
1386 " ensemble, you must set tauThermostat.\n");
1387 painCave.isFatal = 1;
1388 simError();
1389 }
1390 break;
1391
1392 case NPTi_ENS:
1393 if (globals->haveZconstraints()){
1394 setupZConstraint(info[k]);
1395 myNPTi = new ZConstraint<NPTi<RealIntegrator> >( &(info[k]), the_ff );
1396 }
1397 else
1398 myNPTi = new NPTi<RealIntegrator>( &(info[k]), the_ff );
1399
1400 myNPTi->setTargetTemp( globals->getTargetTemp() );
1401
1402 if (globals->haveTargetPressure())
1403 myNPTi->setTargetPressure(globals->getTargetPressure());
1404 else {
1405 sprintf( painCave.errMsg,
1406 "SimSetup error: If you use a constant pressure\n"
1407 " ensemble, you must set targetPressure in the BASS file.\n");
1408 painCave.isFatal = 1;
1409 simError();
1410 }
1411
1412 if( globals->haveTauThermostat() )
1413 myNPTi->setTauThermostat( globals->getTauThermostat() );
1414 else{
1415 sprintf( painCave.errMsg,
1416 "SimSetup error: If you use an NPT\n"
1417 " ensemble, you must set tauThermostat.\n");
1418 painCave.isFatal = 1;
1419 simError();
1420 }
1421
1422 if( globals->haveTauBarostat() )
1423 myNPTi->setTauBarostat( globals->getTauBarostat() );
1424 else{
1425 sprintf( painCave.errMsg,
1426 "SimSetup error: If you use an NPT\n"
1427 " ensemble, you must set tauBarostat.\n");
1428 painCave.isFatal = 1;
1429 simError();
1430 }
1431 break;
1432
1433 case NPTf_ENS:
1434 if (globals->haveZconstraints()){
1435 setupZConstraint(info[k]);
1436 myNPTf = new ZConstraint<NPTf<RealIntegrator> >( &(info[k]), the_ff );
1437 }
1438 else
1439 myNPTf = new NPTf<RealIntegrator>( &(info[k]), the_ff );
1440
1441 myNPTf->setTargetTemp( globals->getTargetTemp());
1442
1443 if (globals->haveTargetPressure())
1444 myNPTf->setTargetPressure(globals->getTargetPressure());
1445 else {
1446 sprintf( painCave.errMsg,
1447 "SimSetup error: If you use a constant pressure\n"
1448 " ensemble, you must set targetPressure in the BASS file.\n");
1449 painCave.isFatal = 1;
1450 simError();
1451 }
1452
1453 if( globals->haveTauThermostat() )
1454 myNPTf->setTauThermostat( globals->getTauThermostat() );
1455 else{
1456 sprintf( painCave.errMsg,
1457 "SimSetup error: If you use an NPT\n"
1458 " ensemble, you must set tauThermostat.\n");
1459 painCave.isFatal = 1;
1460 simError();
1461 }
1462
1463 if( globals->haveTauBarostat() )
1464 myNPTf->setTauBarostat( globals->getTauBarostat() );
1465 else{
1466 sprintf( painCave.errMsg,
1467 "SimSetup error: If you use an NPT\n"
1468 " ensemble, you must set tauBarostat.\n");
1469 painCave.isFatal = 1;
1470 simError();
1471 }
1472 break;
1473
1474 case NPTim_ENS:
1475 if (globals->haveZconstraints()){
1476 setupZConstraint(info[k]);
1477 myNPTim = new ZConstraint<NPTim<RealIntegrator> >( &(info[k]), the_ff );
1478 }
1479 else
1480 myNPTim = new NPTim<RealIntegrator>( &(info[k]), the_ff );
1481
1482 myNPTim->setTargetTemp( globals->getTargetTemp());
1483
1484 if (globals->haveTargetPressure())
1485 myNPTim->setTargetPressure(globals->getTargetPressure());
1486 else {
1487 sprintf( painCave.errMsg,
1488 "SimSetup error: If you use a constant pressure\n"
1489 " ensemble, you must set targetPressure in the BASS file.\n");
1490 painCave.isFatal = 1;
1491 simError();
1492 }
1493
1494 if( globals->haveTauThermostat() )
1495 myNPTim->setTauThermostat( globals->getTauThermostat() );
1496 else{
1497 sprintf( painCave.errMsg,
1498 "SimSetup error: If you use an NPT\n"
1499 " ensemble, you must set tauThermostat.\n");
1500 painCave.isFatal = 1;
1501 simError();
1502 }
1503
1504 if( globals->haveTauBarostat() )
1505 myNPTim->setTauBarostat( globals->getTauBarostat() );
1506 else{
1507 sprintf( painCave.errMsg,
1508 "SimSetup error: If you use an NPT\n"
1509 " ensemble, you must set tauBarostat.\n");
1510 painCave.isFatal = 1;
1511 simError();
1512 }
1513 break;
1514
1515 case NPTfm_ENS:
1516 if (globals->haveZconstraints()){
1517 setupZConstraint(info[k]);
1518 myNPTfm = new ZConstraint<NPTfm<RealIntegrator> >( &(info[k]), the_ff );
1519 }
1520 else
1521 myNPTfm = new NPTfm<RealIntegrator>( &(info[k]), the_ff );
1522
1523 myNPTfm->setTargetTemp( globals->getTargetTemp());
1524
1525 if (globals->haveTargetPressure())
1526 myNPTfm->setTargetPressure(globals->getTargetPressure());
1527 else {
1528 sprintf( painCave.errMsg,
1529 "SimSetup error: If you use a constant pressure\n"
1530 " ensemble, you must set targetPressure in the BASS file.\n");
1531 painCave.isFatal = 1;
1532 simError();
1533 }
1534
1535 if( globals->haveTauThermostat() )
1536 myNPTfm->setTauThermostat( globals->getTauThermostat() );
1537 else{
1538 sprintf( painCave.errMsg,
1539 "SimSetup error: If you use an NPT\n"
1540 " ensemble, you must set tauThermostat.\n");
1541 painCave.isFatal = 1;
1542 simError();
1543 }
1544
1545 if( globals->haveTauBarostat() )
1546 myNPTfm->setTauBarostat( globals->getTauBarostat() );
1547 else{
1548 sprintf( painCave.errMsg,
1549 "SimSetup error: If you use an NPT\n"
1550 " ensemble, you must set tauBarostat.\n");
1551 painCave.isFatal = 1;
1552 simError();
1553 }
1554 break;
1555
1556 default:
1557 sprintf( painCave.errMsg,
1558 "SimSetup Error. Unrecognized ensemble in case statement.\n");
1559 painCave.isFatal = 1;
1560 simError();
1561 }
1562 }
1563 }
1564
1565 void SimSetup::initFortran( void ){
1566
1567 info[0].refreshSim();
1568
1569 if( !strcmp( info[0].mixingRule, "standard") ){
1570 the_ff->initForceField( LB_MIXING_RULE );
1571 }
1572 else if( !strcmp( info[0].mixingRule, "explicit") ){
1573 the_ff->initForceField( EXPLICIT_MIXING_RULE );
1574 }
1575 else{
1576 sprintf( painCave.errMsg,
1577 "SimSetup Error: unknown mixing rule -> \"%s\"\n",
1578 info[0].mixingRule );
1579 painCave.isFatal = 1;
1580 simError();
1581 }
1582
1583
1584 #ifdef IS_MPI
1585 strcpy( checkPointMsg,
1586 "Successfully intialized the mixingRule for Fortran." );
1587 MPIcheckPoint();
1588 #endif // is_mpi
1589
1590 }
1591
1592 void SimSetup::setupZConstraint(SimInfo& theInfo)
1593 {
1594 int nZConstraints;
1595 ZconStamp** zconStamp;
1596
1597 if(globals->haveZconstraintTime()){
1598
1599 //add sample time of z-constraint into SimInfo's property list
1600 DoubleData* zconsTimeProp = new DoubleData();
1601 zconsTimeProp->setID(ZCONSTIME_ID);
1602 zconsTimeProp->setData(globals->getZconsTime());
1603 theInfo.addProperty(zconsTimeProp);
1604 }
1605 else{
1606 sprintf( painCave.errMsg,
1607 "ZConstraint error: If you use an ZConstraint\n"
1608 " , you must set sample time.\n");
1609 painCave.isFatal = 1;
1610 simError();
1611 }
1612
1613 //
1614 nZConstraints = globals->getNzConstraints();
1615 zconStamp = globals->getZconStamp();
1616 ZConsParaItem tempParaItem;
1617
1618 ZConsParaData* zconsParaData = new ZConsParaData();
1619 zconsParaData->setID(ZCONSPARADATA_ID);
1620
1621 for(int i = 0; i < nZConstraints; i++){
1622 tempParaItem.havingZPos = zconStamp[i]->haveZpos();
1623 tempParaItem.zPos = zconStamp[i]->getZpos();
1624 tempParaItem.zconsIndex = zconStamp[i]->getMolIndex();
1625 tempParaItem.kRatio = zconStamp[i]->getKratio();
1626
1627 zconsParaData->addItem(tempParaItem);
1628 }
1629
1630 //sort the parameters by index of molecules
1631 zconsParaData->sortByIndex();
1632
1633 //push data into siminfo, therefore, we can retrieve later
1634 theInfo.addProperty(zconsParaData);
1635
1636 //push zconsTol into siminfo, if user does not specify
1637 //value for zconsTol, a default value will be used
1638 DoubleData* zconsTol = new DoubleData();
1639 zconsTol->setID(ZCONSTOL_ID);
1640 if(globals->haveZconsTol()){
1641 zconsTol->setData(globals->getZconsTol());
1642 }
1643 else{
1644 double defaultZConsTol = 1E-6;
1645 sprintf( painCave.errMsg,
1646 "ZConstraint Waring: Tolerance for z-constraint methodl is not specified\n"
1647 " , default value %f is used.\n", defaultZConsTol);
1648 painCave.isFatal = 0;
1649 simError();
1650
1651 zconsTol->setData(defaultZConsTol);
1652 }
1653 theInfo.addProperty(zconsTol);
1654
1655 //Determine the name of ouput file and add it into SimInfo's property list
1656 //Be careful, do not use inFileName, since it is a pointer which
1657 //point to a string at master node, and slave nodes do not contain that string
1658
1659 string zconsOutput(theInfo.finalName);
1660
1661 zconsOutput = zconsOutput.substr(0, zconsOutput.rfind(".")) + ".fz";
1662
1663 StringData* zconsFilename = new StringData();
1664 zconsFilename->setID(ZCONSFILENAME_ID);
1665 zconsFilename->setData(zconsOutput);
1666
1667 theInfo.addProperty(zconsFilename);
1668 }