ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-1.0/libmdtools/SimSetup.cpp
Revision: 1415
Committed: Mon Jul 26 17:50:57 2004 UTC (20 years, 9 months ago) by gezelter
File size: 62841 byte(s)
Log Message:
Fixing required keyword annoyances

File Contents

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