ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-1.0/libmdtools/SimSetup.cpp
Revision: 1334
Committed: Fri Jul 16 18:58:03 2004 UTC (19 years, 11 months ago) by gezelter
File size: 60534 byte(s)
Log Message:
Initial import of OOPSE-1.0 source tree

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 info[i].target_temp = globals->getTargetTemp();
828 info[i].dt = globals->getDt();
829 info[i].run_time = globals->getRunTime();
830 }
831 n_components = globals->getNComponents();
832
833
834 // get the forceField
835
836 strcpy(force_field, globals->getForceField());
837
838 if (!strcasecmp(force_field, "DUFF")){
839 ffCase = FF_DUFF;
840 }
841 else if (!strcasecmp(force_field, "LJ")){
842 ffCase = FF_LJ;
843 }
844 else if (!strcasecmp(force_field, "EAM")){
845 ffCase = FF_EAM;
846 }
847 else if (!strcasecmp(force_field, "WATER")){
848 ffCase = FF_H2O;
849 }
850 else{
851 sprintf(painCave.errMsg, "SimSetup Error. Unrecognized force field -> %s\n",
852 force_field);
853 painCave.isFatal = 1;
854 simError();
855 }
856 if (globals->haveForceFieldVariant()) {
857 strcpy(forcefield_variant, globals->getForceFieldVariant());
858 has_forcefield_variant = 1;
859 }
860
861 // get the ensemble
862
863 strcpy(ensemble, globals->getEnsemble());
864
865 if (!strcasecmp(ensemble, "NVE")){
866 ensembleCase = NVE_ENS;
867 }
868 else if (!strcasecmp(ensemble, "NVT")){
869 ensembleCase = NVT_ENS;
870 }
871 else if (!strcasecmp(ensemble, "NPTi") || !strcasecmp(ensemble, "NPT")){
872 ensembleCase = NPTi_ENS;
873 }
874 else if (!strcasecmp(ensemble, "NPTf")){
875 ensembleCase = NPTf_ENS;
876 }
877 else if (!strcasecmp(ensemble, "NPTxyz")){
878 ensembleCase = NPTxyz_ENS;
879 }
880 else{
881 sprintf(painCave.errMsg,
882 "SimSetup Warning. Unrecognized Ensemble -> %s \n"
883 "\treverting to NVE for this simulation.\n",
884 ensemble);
885 painCave.isFatal = 0;
886 simError();
887 strcpy(ensemble, "NVE");
888 ensembleCase = NVE_ENS;
889 }
890
891 for (i = 0; i < nInfo; i++){
892 strcpy(info[i].ensemble, ensemble);
893
894 // get the mixing rule
895
896 strcpy(info[i].mixingRule, globals->getMixingRule());
897 info[i].usePBC = globals->getPBC();
898 }
899
900 // get the components and calculate the tot_nMol and indvidual n_mol
901
902 the_components = globals->getComponents();
903 components_nmol = new int[n_components];
904
905
906 if (!globals->haveNMol()){
907 // we don't have the total number of molecules, so we assume it is
908 // given in each component
909
910 tot_nmol = 0;
911 for (i = 0; i < n_components; i++){
912 if (!the_components[i]->haveNMol()){
913 // we have a problem
914 sprintf(painCave.errMsg,
915 "SimSetup Error. No global NMol or component NMol given.\n"
916 "\tCannot calculate the number of atoms.\n");
917 painCave.isFatal = 1;
918 simError();
919 }
920
921 tot_nmol += the_components[i]->getNMol();
922 components_nmol[i] = the_components[i]->getNMol();
923 }
924 }
925 else{
926 sprintf(painCave.errMsg,
927 "SimSetup error.\n"
928 "\tSorry, the ability to specify total"
929 " nMols and then give molfractions in the components\n"
930 "\tis not currently supported."
931 " Please give nMol in the components.\n");
932 painCave.isFatal = 1;
933 simError();
934 }
935
936 //check whether sample time, status time, thermal time and reset time are divisble by dt
937 if (globals->haveSampleTime() && !isDivisible(globals->getSampleTime(), globals->getDt())){
938 sprintf(painCave.errMsg,
939 "Sample time is not divisible by dt.\n"
940 "\tThis will result in samples that are not uniformly\n"
941 "\tdistributed in time. If this is a problem, change\n"
942 "\tyour sampleTime variable.\n");
943 painCave.isFatal = 0;
944 simError();
945 }
946
947 if (globals->haveStatusTime() && !isDivisible(globals->getStatusTime(), globals->getDt())){
948 sprintf(painCave.errMsg,
949 "Status time is not divisible by dt.\n"
950 "\tThis will result in status reports that are not uniformly\n"
951 "\tdistributed in time. If this is a problem, change \n"
952 "\tyour statusTime variable.\n");
953 painCave.isFatal = 0;
954 simError();
955 }
956
957 if (globals->haveThermalTime() && !isDivisible(globals->getThermalTime(), globals->getDt())){
958 sprintf(painCave.errMsg,
959 "Thermal time is not divisible by dt.\n"
960 "\tThis will result in thermalizations that are not uniformly\n"
961 "\tdistributed in time. If this is a problem, change \n"
962 "\tyour thermalTime variable.\n");
963 painCave.isFatal = 0;
964 simError();
965 }
966
967 if (globals->haveResetTime() && !isDivisible(globals->getResetTime(), globals->getDt())){
968 sprintf(painCave.errMsg,
969 "Reset time is not divisible by dt.\n"
970 "\tThis will result in integrator resets that are not uniformly\n"
971 "\tdistributed in time. If this is a problem, change\n"
972 "\tyour resetTime variable.\n");
973 painCave.isFatal = 0;
974 simError();
975 }
976
977 // set the status, sample, and thermal kick times
978
979 for (i = 0; i < nInfo; i++){
980 if (globals->haveSampleTime()){
981 info[i].sampleTime = globals->getSampleTime();
982 info[i].statusTime = info[i].sampleTime;
983 }
984 else{
985 info[i].sampleTime = globals->getRunTime();
986 info[i].statusTime = info[i].sampleTime;
987 }
988
989 if (globals->haveStatusTime()){
990 info[i].statusTime = globals->getStatusTime();
991 }
992
993 if (globals->haveThermalTime()){
994 info[i].thermalTime = globals->getThermalTime();
995 } else {
996 info[i].thermalTime = globals->getRunTime();
997 }
998
999 info[i].resetIntegrator = 0;
1000 if( globals->haveResetTime() ){
1001 info[i].resetTime = globals->getResetTime();
1002 info[i].resetIntegrator = 1;
1003 }
1004
1005 // check for the temperature set flag
1006
1007 if (globals->haveTempSet())
1008 info[i].setTemp = globals->getTempSet();
1009
1010 // check for the extended State init
1011
1012 info[i].useInitXSstate = globals->getUseInitXSstate();
1013 info[i].orthoTolerance = globals->getOrthoBoxTolerance();
1014
1015 // check for thermodynamic integration
1016 if (globals->getUseSolidThermInt() && !globals->getUseLiquidThermInt()) {
1017 if (globals->haveThermIntLambda() && globals->haveThermIntK()) {
1018 info[i].useSolidThermInt = globals->getUseSolidThermInt();
1019 info[i].thermIntLambda = globals->getThermIntLambda();
1020 info[i].thermIntK = globals->getThermIntK();
1021
1022 Restraints *myRestraint = new Restraints(tot_nmol, info[i].thermIntLambda, info[i].thermIntK);
1023 info[i].restraint = myRestraint;
1024 }
1025 else {
1026 sprintf(painCave.errMsg,
1027 "SimSetup Error:\n"
1028 "\tKeyword useSolidThermInt was set to 'true' but\n"
1029 "\tthermodynamicIntegrationLambda (and/or\n"
1030 "\tthermodynamicIntegrationK) was not specified.\n"
1031 "\tPlease provide a lambda value and k value in your .bass file.\n");
1032 painCave.isFatal = 1;
1033 simError();
1034 }
1035 }
1036 else if(globals->getUseLiquidThermInt()) {
1037 if (globals->getUseSolidThermInt()) {
1038 sprintf( painCave.errMsg,
1039 "SimSetup Warning: It appears that you have both solid and\n"
1040 "\tliquid thermodynamic integration activated in your .bass\n"
1041 "\tfile. To avoid confusion, specify only one technique in\n"
1042 "\tyour .bass file. Liquid-state thermodynamic integration\n"
1043 "\twill be assumed for the current simulation. If this is not\n"
1044 "\twhat you desire, set useSolidThermInt to 'true' and\n"
1045 "\tuseLiquidThermInt to 'false' in your .bass file.\n");
1046 painCave.isFatal = 0;
1047 simError();
1048 }
1049 if (globals->haveThermIntLambda() && globals->haveThermIntK()) {
1050 info[i].useLiquidThermInt = globals->getUseLiquidThermInt();
1051 info[i].thermIntLambda = globals->getThermIntLambda();
1052 info[i].thermIntK = globals->getThermIntK();
1053 }
1054 else {
1055 sprintf(painCave.errMsg,
1056 "SimSetup Error:\n"
1057 "\tKeyword useLiquidThermInt was set to 'true' but\n"
1058 "\tthermodynamicIntegrationLambda (and/or\n"
1059 "\tthermodynamicIntegrationK) was not specified.\n"
1060 "\tPlease provide a lambda value and k value in your .bass file.\n");
1061 painCave.isFatal = 1;
1062 simError();
1063 }
1064 }
1065 else if(globals->haveThermIntLambda() || globals->haveThermIntK()){
1066 sprintf(painCave.errMsg,
1067 "SimSetup Warning: If you want to use Thermodynamic\n"
1068 "\tIntegration, set useSolidThermInt or useLiquidThermInt to\n"
1069 "\t'true' in your .bass file. These keywords are set to\n"
1070 "\t'false' by default, so your lambda and/or k values are\n"
1071 "\tbeing ignored.\n");
1072 painCave.isFatal = 0;
1073 simError();
1074 }
1075 }
1076
1077 //setup seed for random number generator
1078 int seedValue;
1079
1080 if (globals->haveSeed()){
1081 seedValue = globals->getSeed();
1082
1083 if(seedValue / 1E9 == 0){
1084 sprintf(painCave.errMsg,
1085 "Seed for sprng library should contain at least 9 digits\n"
1086 "OOPSE will generate a seed for user\n");
1087 painCave.isFatal = 0;
1088 simError();
1089
1090 //using seed generated by system instead of invalid seed set by user
1091 #ifndef IS_MPI
1092 seedValue = make_sprng_seed();
1093 #else
1094 if (worldRank == 0){
1095 seedValue = make_sprng_seed();
1096 }
1097 MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD);
1098 #endif
1099 }
1100 }//end of if branch of globals->haveSeed()
1101 else{
1102
1103 #ifndef IS_MPI
1104 seedValue = make_sprng_seed();
1105 #else
1106 if (worldRank == 0){
1107 seedValue = make_sprng_seed();
1108 }
1109 MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD);
1110 #endif
1111 }//end of globals->haveSeed()
1112
1113 for (int i = 0; i < nInfo; i++){
1114 info[i].setSeed(seedValue);
1115 }
1116
1117 #ifdef IS_MPI
1118 strcpy(checkPointMsg, "Successfully gathered all information from Bass\n");
1119 MPIcheckPoint();
1120 #endif // is_mpi
1121 }
1122
1123
1124 void SimSetup::finalInfoCheck(void){
1125 int index;
1126 int usesDipoles;
1127 int usesCharges;
1128 int i;
1129
1130 for (i = 0; i < nInfo; i++){
1131 // check electrostatic parameters
1132
1133 index = 0;
1134 usesDipoles = 0;
1135 while ((index < info[i].n_atoms) && !usesDipoles){
1136 usesDipoles = (info[i].atoms[index])->hasDipole();
1137 index++;
1138 }
1139 index = 0;
1140 usesCharges = 0;
1141 while ((index < info[i].n_atoms) && !usesCharges){
1142 usesCharges= (info[i].atoms[index])->hasCharge();
1143 index++;
1144 }
1145 #ifdef IS_MPI
1146 int myUse = usesDipoles;
1147 MPI_Allreduce(&myUse, &usesDipoles, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
1148 #endif //is_mpi
1149
1150 double theRcut, theRsw;
1151
1152 if (globals->haveRcut()) {
1153 theRcut = globals->getRcut();
1154
1155 if (globals->haveRsw())
1156 theRsw = globals->getRsw();
1157 else
1158 theRsw = theRcut;
1159
1160 info[i].setDefaultRcut(theRcut, theRsw);
1161
1162 } else {
1163
1164 the_ff->calcRcut();
1165 theRcut = info[i].getRcut();
1166
1167 if (globals->haveRsw())
1168 theRsw = globals->getRsw();
1169 else
1170 theRsw = theRcut;
1171
1172 info[i].setDefaultRcut(theRcut, theRsw);
1173 }
1174
1175 if (globals->getUseRF()){
1176 info[i].useReactionField = 1;
1177
1178 if (!globals->haveRcut()){
1179 sprintf(painCave.errMsg,
1180 "SimSetup Warning: No value was set for the cutoffRadius.\n"
1181 "\tOOPSE will use a default value of 15.0 angstroms"
1182 "\tfor the cutoffRadius.\n");
1183 painCave.isFatal = 0;
1184 simError();
1185 theRcut = 15.0;
1186 }
1187 else{
1188 theRcut = globals->getRcut();
1189 }
1190
1191 if (!globals->haveRsw()){
1192 sprintf(painCave.errMsg,
1193 "SimSetup Warning: No value was set for switchingRadius.\n"
1194 "\tOOPSE will use a default value of\n"
1195 "\t0.95 * cutoffRadius for the switchingRadius\n");
1196 painCave.isFatal = 0;
1197 simError();
1198 theRsw = 0.95 * theRcut;
1199 }
1200 else{
1201 theRsw = globals->getRsw();
1202 }
1203
1204 info[i].setDefaultRcut(theRcut, theRsw);
1205
1206 if (!globals->haveDielectric()){
1207 sprintf(painCave.errMsg,
1208 "SimSetup Error: No Dielectric constant was set.\n"
1209 "\tYou are trying to use Reaction Field without"
1210 "\tsetting a dielectric constant!\n");
1211 painCave.isFatal = 1;
1212 simError();
1213 }
1214 info[i].dielectric = globals->getDielectric();
1215 }
1216 else{
1217 if (usesDipoles || usesCharges){
1218
1219 if (!globals->haveRcut()){
1220 sprintf(painCave.errMsg,
1221 "SimSetup Warning: No value was set for the cutoffRadius.\n"
1222 "\tOOPSE will use a default value of 15.0 angstroms"
1223 "\tfor the cutoffRadius.\n");
1224 painCave.isFatal = 0;
1225 simError();
1226 theRcut = 15.0;
1227 }
1228 else{
1229 theRcut = globals->getRcut();
1230 }
1231
1232 if (!globals->haveRsw()){
1233 sprintf(painCave.errMsg,
1234 "SimSetup Warning: No value was set for switchingRadius.\n"
1235 "\tOOPSE will use a default value of\n"
1236 "\t0.95 * cutoffRadius for the switchingRadius\n");
1237 painCave.isFatal = 0;
1238 simError();
1239 theRsw = 0.95 * theRcut;
1240 }
1241 else{
1242 theRsw = globals->getRsw();
1243 }
1244
1245 info[i].setDefaultRcut(theRcut, theRsw);
1246
1247 }
1248 }
1249 }
1250 #ifdef IS_MPI
1251 strcpy(checkPointMsg, "post processing checks out");
1252 MPIcheckPoint();
1253 #endif // is_mpi
1254
1255 // clean up the forcefield
1256 the_ff->cleanMe();
1257 }
1258
1259 void SimSetup::initSystemCoords(void){
1260 int i;
1261
1262 char* inName;
1263
1264 (info[0].getConfiguration())->createArrays(info[0].n_atoms);
1265
1266 for (i = 0; i < info[0].n_atoms; i++)
1267 info[0].atoms[i]->setCoords();
1268
1269 if (globals->haveInitialConfig()){
1270 InitializeFromFile* fileInit;
1271 #ifdef IS_MPI // is_mpi
1272 if (worldRank == 0){
1273 #endif //is_mpi
1274 inName = globals->getInitialConfig();
1275 fileInit = new InitializeFromFile(inName);
1276 #ifdef IS_MPI
1277 }
1278 else
1279 fileInit = new InitializeFromFile(NULL);
1280 #endif
1281 fileInit->readInit(info); // default velocities on
1282
1283 delete fileInit;
1284 }
1285 else{
1286
1287 // no init from bass
1288
1289 sprintf(painCave.errMsg,
1290 "Cannot intialize a simulation without an initial configuration file.\n");
1291 painCave.isFatal = 1;;
1292 simError();
1293
1294 }
1295
1296 #ifdef IS_MPI
1297 strcpy(checkPointMsg, "Successfully read in the initial configuration");
1298 MPIcheckPoint();
1299 #endif // is_mpi
1300 }
1301
1302
1303 void SimSetup::makeOutNames(void){
1304 int k;
1305
1306
1307 for (k = 0; k < nInfo; k++){
1308 #ifdef IS_MPI
1309 if (worldRank == 0){
1310 #endif // is_mpi
1311
1312 if (globals->haveFinalConfig()){
1313 strcpy(info[k].finalName, globals->getFinalConfig());
1314 }
1315 else{
1316 strcpy(info[k].finalName, inFileName);
1317 char* endTest;
1318 int nameLength = strlen(info[k].finalName);
1319 endTest = &(info[k].finalName[nameLength - 5]);
1320 if (!strcmp(endTest, ".bass")){
1321 strcpy(endTest, ".eor");
1322 }
1323 else if (!strcmp(endTest, ".BASS")){
1324 strcpy(endTest, ".eor");
1325 }
1326 else{
1327 endTest = &(info[k].finalName[nameLength - 4]);
1328 if (!strcmp(endTest, ".bss")){
1329 strcpy(endTest, ".eor");
1330 }
1331 else if (!strcmp(endTest, ".mdl")){
1332 strcpy(endTest, ".eor");
1333 }
1334 else{
1335 strcat(info[k].finalName, ".eor");
1336 }
1337 }
1338 }
1339
1340 // make the sample and status out names
1341
1342 strcpy(info[k].sampleName, inFileName);
1343 char* endTest;
1344 int nameLength = strlen(info[k].sampleName);
1345 endTest = &(info[k].sampleName[nameLength - 5]);
1346 if (!strcmp(endTest, ".bass")){
1347 strcpy(endTest, ".dump");
1348 }
1349 else if (!strcmp(endTest, ".BASS")){
1350 strcpy(endTest, ".dump");
1351 }
1352 else{
1353 endTest = &(info[k].sampleName[nameLength - 4]);
1354 if (!strcmp(endTest, ".bss")){
1355 strcpy(endTest, ".dump");
1356 }
1357 else if (!strcmp(endTest, ".mdl")){
1358 strcpy(endTest, ".dump");
1359 }
1360 else{
1361 strcat(info[k].sampleName, ".dump");
1362 }
1363 }
1364
1365 strcpy(info[k].statusName, inFileName);
1366 nameLength = strlen(info[k].statusName);
1367 endTest = &(info[k].statusName[nameLength - 5]);
1368 if (!strcmp(endTest, ".bass")){
1369 strcpy(endTest, ".stat");
1370 }
1371 else if (!strcmp(endTest, ".BASS")){
1372 strcpy(endTest, ".stat");
1373 }
1374 else{
1375 endTest = &(info[k].statusName[nameLength - 4]);
1376 if (!strcmp(endTest, ".bss")){
1377 strcpy(endTest, ".stat");
1378 }
1379 else if (!strcmp(endTest, ".mdl")){
1380 strcpy(endTest, ".stat");
1381 }
1382 else{
1383 strcat(info[k].statusName, ".stat");
1384 }
1385 }
1386
1387 strcpy(info[k].rawPotName, inFileName);
1388 nameLength = strlen(info[k].rawPotName);
1389 endTest = &(info[k].rawPotName[nameLength - 5]);
1390 if (!strcmp(endTest, ".bass")){
1391 strcpy(endTest, ".raw");
1392 }
1393 else if (!strcmp(endTest, ".BASS")){
1394 strcpy(endTest, ".raw");
1395 }
1396 else{
1397 endTest = &(info[k].rawPotName[nameLength - 4]);
1398 if (!strcmp(endTest, ".bss")){
1399 strcpy(endTest, ".raw");
1400 }
1401 else if (!strcmp(endTest, ".mdl")){
1402 strcpy(endTest, ".raw");
1403 }
1404 else{
1405 strcat(info[k].rawPotName, ".raw");
1406 }
1407 }
1408
1409 #ifdef IS_MPI
1410
1411 }
1412 #endif // is_mpi
1413 }
1414 }
1415
1416
1417 void SimSetup::sysObjectsCreation(void){
1418 int i, k;
1419
1420 // create the forceField
1421
1422 createFF();
1423
1424 // extract componentList
1425
1426 compList();
1427
1428 // calc the number of atoms, bond, bends, and torsions
1429
1430 calcSysValues();
1431
1432 #ifdef IS_MPI
1433 // divide the molecules among the processors
1434
1435 mpiMolDivide();
1436 #endif //is_mpi
1437
1438 // create the atom and SRI arrays. Also initialize Molecule Stamp ID's
1439
1440 makeSysArrays();
1441
1442 // make and initialize the molecules (all but atomic coordinates)
1443
1444 makeMolecules();
1445
1446 for (k = 0; k < nInfo; k++){
1447 info[k].identArray = new int[info[k].n_atoms];
1448 for (i = 0; i < info[k].n_atoms; i++){
1449 info[k].identArray[i] = info[k].atoms[i]->getIdent();
1450 }
1451 }
1452 }
1453
1454
1455 void SimSetup::createFF(void){
1456 switch (ffCase){
1457 case FF_DUFF:
1458 the_ff = new DUFF();
1459 break;
1460
1461 case FF_LJ:
1462 the_ff = new LJFF();
1463 break;
1464
1465 case FF_EAM:
1466 if (has_forcefield_variant)
1467 the_ff = new EAM_FF(forcefield_variant);
1468 else
1469 the_ff = new EAM_FF();
1470 break;
1471
1472 case FF_H2O:
1473 the_ff = new WATER();
1474 break;
1475
1476 default:
1477 sprintf(painCave.errMsg,
1478 "SimSetup Error. Unrecognized force field in case statement.\n");
1479 painCave.isFatal = 1;
1480 simError();
1481 }
1482
1483
1484 #ifdef IS_MPI
1485 strcpy(checkPointMsg, "ForceField creation successful");
1486 MPIcheckPoint();
1487 #endif // is_mpi
1488 }
1489
1490
1491 void SimSetup::compList(void){
1492 int i;
1493 char* id;
1494 LinkedMolStamp* headStamp = new LinkedMolStamp();
1495 LinkedMolStamp* currentStamp = NULL;
1496 comp_stamps = new MoleculeStamp * [n_components];
1497 bool haveCutoffGroups;
1498
1499 haveCutoffGroups = false;
1500
1501 // make an array of molecule stamps that match the components used.
1502 // also extract the used stamps out into a separate linked list
1503
1504 for (i = 0; i < nInfo; i++){
1505 info[i].nComponents = n_components;
1506 info[i].componentsNmol = components_nmol;
1507 info[i].compStamps = comp_stamps;
1508 info[i].headStamp = headStamp;
1509 }
1510
1511
1512 for (i = 0; i < n_components; i++){
1513 id = the_components[i]->getType();
1514 comp_stamps[i] = NULL;
1515
1516 // check to make sure the component isn't already in the list
1517
1518 comp_stamps[i] = headStamp->match(id);
1519 if (comp_stamps[i] == NULL){
1520 // extract the component from the list;
1521
1522 currentStamp = stamps->extractMolStamp(id);
1523 if (currentStamp == NULL){
1524 sprintf(painCave.errMsg,
1525 "SimSetup error: Component \"%s\" was not found in the "
1526 "list of declared molecules\n",
1527 id);
1528 painCave.isFatal = 1;
1529 simError();
1530 }
1531
1532 headStamp->add(currentStamp);
1533 comp_stamps[i] = headStamp->match(id);
1534 }
1535
1536 if(comp_stamps[i]->getNCutoffGroups() > 0)
1537 haveCutoffGroups = true;
1538 }
1539
1540 for (i = 0; i < nInfo; i++)
1541 info[i].haveCutoffGroups = haveCutoffGroups;
1542
1543 #ifdef IS_MPI
1544 strcpy(checkPointMsg, "Component stamps successfully extracted\n");
1545 MPIcheckPoint();
1546 #endif // is_mpi
1547 }
1548
1549 void SimSetup::calcSysValues(void){
1550 int i, j;
1551 int ncutgroups, atomsingroups, ngroupsinstamp;
1552
1553 int* molMembershipArray;
1554 CutoffGroupStamp* cg;
1555
1556 tot_atoms = 0;
1557 tot_bonds = 0;
1558 tot_bends = 0;
1559 tot_torsions = 0;
1560 tot_rigid = 0;
1561 tot_groups = 0;
1562 for (i = 0; i < n_components; i++){
1563 tot_atoms += components_nmol[i] * comp_stamps[i]->getNAtoms();
1564 tot_bonds += components_nmol[i] * comp_stamps[i]->getNBonds();
1565 tot_bends += components_nmol[i] * comp_stamps[i]->getNBends();
1566 tot_torsions += components_nmol[i] * comp_stamps[i]->getNTorsions();
1567 tot_rigid += components_nmol[i] * comp_stamps[i]->getNRigidBodies();
1568
1569 ncutgroups = comp_stamps[i]->getNCutoffGroups();
1570 atomsingroups = 0;
1571 for (j=0; j < ncutgroups; j++) {
1572 cg = comp_stamps[i]->getCutoffGroup(j);
1573 atomsingroups += cg->getNMembers();
1574 }
1575 ngroupsinstamp = comp_stamps[i]->getNAtoms() - atomsingroups + ncutgroups;
1576 tot_groups += components_nmol[i] * ngroupsinstamp;
1577 }
1578
1579 tot_SRI = tot_bonds + tot_bends + tot_torsions;
1580 molMembershipArray = new int[tot_atoms];
1581
1582 for (i = 0; i < nInfo; i++){
1583 info[i].n_atoms = tot_atoms;
1584 info[i].n_bonds = tot_bonds;
1585 info[i].n_bends = tot_bends;
1586 info[i].n_torsions = tot_torsions;
1587 info[i].n_SRI = tot_SRI;
1588 info[i].n_mol = tot_nmol;
1589 info[i].ngroup = tot_groups;
1590 info[i].molMembershipArray = molMembershipArray;
1591 }
1592 }
1593
1594 #ifdef IS_MPI
1595
1596 void SimSetup::mpiMolDivide(void){
1597 int i, j, k;
1598 int localMol, allMol;
1599 int local_atoms, local_bonds, local_bends, local_torsions, local_SRI;
1600 int local_rigid, local_groups;
1601 vector<int> globalMolIndex;
1602 int ncutgroups, atomsingroups, ngroupsinstamp;
1603 CutoffGroupStamp* cg;
1604
1605 mpiSim = new mpiSimulation(info);
1606
1607 mpiSim->divideLabor();
1608 globalAtomIndex = mpiSim->getGlobalAtomIndex();
1609 globalGroupIndex = mpiSim->getGlobalGroupIndex();
1610 //globalMolIndex = mpiSim->getGlobalMolIndex();
1611
1612 // set up the local variables
1613
1614 mol2proc = mpiSim->getMolToProcMap();
1615 molCompType = mpiSim->getMolComponentType();
1616
1617 allMol = 0;
1618 localMol = 0;
1619 local_atoms = 0;
1620 local_bonds = 0;
1621 local_bends = 0;
1622 local_torsions = 0;
1623 local_rigid = 0;
1624 local_groups = 0;
1625 globalAtomCounter = 0;
1626
1627 for (i = 0; i < n_components; i++){
1628 for (j = 0; j < components_nmol[i]; j++){
1629 if (mol2proc[allMol] == worldRank){
1630 local_atoms += comp_stamps[i]->getNAtoms();
1631 local_bonds += comp_stamps[i]->getNBonds();
1632 local_bends += comp_stamps[i]->getNBends();
1633 local_torsions += comp_stamps[i]->getNTorsions();
1634 local_rigid += comp_stamps[i]->getNRigidBodies();
1635
1636 ncutgroups = comp_stamps[i]->getNCutoffGroups();
1637 atomsingroups = 0;
1638 for (k=0; k < ncutgroups; k++) {
1639 cg = comp_stamps[i]->getCutoffGroup(k);
1640 atomsingroups += cg->getNMembers();
1641 }
1642 ngroupsinstamp = comp_stamps[i]->getNAtoms() - atomsingroups +
1643 ncutgroups;
1644 local_groups += ngroupsinstamp;
1645
1646 localMol++;
1647 }
1648 for (k = 0; k < comp_stamps[i]->getNAtoms(); k++){
1649 info[0].molMembershipArray[globalAtomCounter] = allMol;
1650 globalAtomCounter++;
1651 }
1652
1653 allMol++;
1654 }
1655 }
1656 local_SRI = local_bonds + local_bends + local_torsions;
1657
1658 info[0].n_atoms = mpiSim->getNAtomsLocal();
1659
1660 if (local_atoms != info[0].n_atoms){
1661 sprintf(painCave.errMsg,
1662 "SimSetup error: mpiSim's localAtom (%d) and SimSetup's\n"
1663 "\tlocalAtom (%d) are not equal.\n",
1664 info[0].n_atoms, local_atoms);
1665 painCave.isFatal = 1;
1666 simError();
1667 }
1668
1669 info[0].ngroup = mpiSim->getNGroupsLocal();
1670 if (local_groups != info[0].ngroup){
1671 sprintf(painCave.errMsg,
1672 "SimSetup error: mpiSim's localGroups (%d) and SimSetup's\n"
1673 "\tlocalGroups (%d) are not equal.\n",
1674 info[0].ngroup, local_groups);
1675 painCave.isFatal = 1;
1676 simError();
1677 }
1678
1679 info[0].n_bonds = local_bonds;
1680 info[0].n_bends = local_bends;
1681 info[0].n_torsions = local_torsions;
1682 info[0].n_SRI = local_SRI;
1683 info[0].n_mol = localMol;
1684
1685 strcpy(checkPointMsg, "Passed nlocal consistency check.");
1686 MPIcheckPoint();
1687 }
1688
1689 #endif // is_mpi
1690
1691
1692 void SimSetup::makeSysArrays(void){
1693
1694 #ifndef IS_MPI
1695 int k, j;
1696 #endif // is_mpi
1697 int i, l;
1698
1699 Atom** the_atoms;
1700 Molecule* the_molecules;
1701
1702 for (l = 0; l < nInfo; l++){
1703 // create the atom and short range interaction arrays
1704
1705 the_atoms = new Atom * [info[l].n_atoms];
1706 the_molecules = new Molecule[info[l].n_mol];
1707 int molIndex;
1708
1709 // initialize the molecule's stampID's
1710
1711 #ifdef IS_MPI
1712
1713
1714 molIndex = 0;
1715 for (i = 0; i < mpiSim->getNMolGlobal(); i++){
1716 if (mol2proc[i] == worldRank){
1717 the_molecules[molIndex].setStampID(molCompType[i]);
1718 the_molecules[molIndex].setMyIndex(molIndex);
1719 the_molecules[molIndex].setGlobalIndex(i);
1720 molIndex++;
1721 }
1722 }
1723
1724 #else // is_mpi
1725
1726 molIndex = 0;
1727 globalAtomCounter = 0;
1728 for (i = 0; i < n_components; i++){
1729 for (j = 0; j < components_nmol[i]; j++){
1730 the_molecules[molIndex].setStampID(i);
1731 the_molecules[molIndex].setMyIndex(molIndex);
1732 the_molecules[molIndex].setGlobalIndex(molIndex);
1733 for (k = 0; k < comp_stamps[i]->getNAtoms(); k++){
1734 info[l].molMembershipArray[globalAtomCounter] = molIndex;
1735 globalAtomCounter++;
1736 }
1737 molIndex++;
1738 }
1739 }
1740
1741
1742 #endif // is_mpi
1743
1744 info[l].globalExcludes = new int;
1745 info[l].globalExcludes[0] = 0;
1746
1747 // set the arrays into the SimInfo object
1748
1749 info[l].atoms = the_atoms;
1750 info[l].molecules = the_molecules;
1751 info[l].nGlobalExcludes = 0;
1752
1753 the_ff->setSimInfo(info);
1754 }
1755 }
1756
1757 void SimSetup::makeIntegrator(void){
1758 int k;
1759
1760 NVE<Integrator<BaseIntegrator> >* myNVE = NULL;
1761 NVT<Integrator<BaseIntegrator> >* myNVT = NULL;
1762 NPTi<NPT<Integrator<BaseIntegrator> > >* myNPTi = NULL;
1763 NPTf<NPT<Integrator<BaseIntegrator> > >* myNPTf = NULL;
1764 NPTxyz<NPT<Integrator<BaseIntegrator> > >* myNPTxyz = NULL;
1765
1766 for (k = 0; k < nInfo; k++){
1767 switch (ensembleCase){
1768 case NVE_ENS:
1769 if (globals->haveZconstraints()){
1770 setupZConstraint(info[k]);
1771 myNVE = new ZConstraint<NVE<RealIntegrator> >(&(info[k]), the_ff);
1772 }
1773 else{
1774 myNVE = new NVE<RealIntegrator>(&(info[k]), the_ff);
1775 }
1776
1777 info->the_integrator = myNVE;
1778 break;
1779
1780 case NVT_ENS:
1781 if (globals->haveZconstraints()){
1782 setupZConstraint(info[k]);
1783 myNVT = new ZConstraint<NVT<RealIntegrator> >(&(info[k]), the_ff);
1784 }
1785 else
1786 myNVT = new NVT<RealIntegrator>(&(info[k]), the_ff);
1787
1788 myNVT->setTargetTemp(globals->getTargetTemp());
1789
1790 if (globals->haveTauThermostat())
1791 myNVT->setTauThermostat(globals->getTauThermostat());
1792 else{
1793 sprintf(painCave.errMsg,
1794 "SimSetup error: If you use the NVT\n"
1795 "\tensemble, you must set tauThermostat.\n");
1796 painCave.isFatal = 1;
1797 simError();
1798 }
1799
1800 info->the_integrator = myNVT;
1801 break;
1802
1803 case NPTi_ENS:
1804 if (globals->haveZconstraints()){
1805 setupZConstraint(info[k]);
1806 myNPTi = new ZConstraint<NPTi<NPT <RealIntegrator> > >(&(info[k]), the_ff);
1807 }
1808 else
1809 myNPTi = new NPTi<NPT<RealIntegrator> >(&(info[k]), the_ff);
1810
1811 myNPTi->setTargetTemp(globals->getTargetTemp());
1812
1813 if (globals->haveTargetPressure())
1814 myNPTi->setTargetPressure(globals->getTargetPressure());
1815 else{
1816 sprintf(painCave.errMsg,
1817 "SimSetup error: If you use a constant pressure\n"
1818 "\tensemble, you must set targetPressure in the BASS file.\n");
1819 painCave.isFatal = 1;
1820 simError();
1821 }
1822
1823 if (globals->haveTauThermostat())
1824 myNPTi->setTauThermostat(globals->getTauThermostat());
1825 else{
1826 sprintf(painCave.errMsg,
1827 "SimSetup error: If you use an NPT\n"
1828 "\tensemble, you must set tauThermostat.\n");
1829 painCave.isFatal = 1;
1830 simError();
1831 }
1832
1833 if (globals->haveTauBarostat())
1834 myNPTi->setTauBarostat(globals->getTauBarostat());
1835 else{
1836 sprintf(painCave.errMsg,
1837 "SimSetup error: If you use an NPT\n"
1838 "\tensemble, you must set tauBarostat.\n");
1839 painCave.isFatal = 1;
1840 simError();
1841 }
1842
1843 info->the_integrator = myNPTi;
1844 break;
1845
1846 case NPTf_ENS:
1847 if (globals->haveZconstraints()){
1848 setupZConstraint(info[k]);
1849 myNPTf = new ZConstraint<NPTf<NPT <RealIntegrator> > >(&(info[k]), the_ff);
1850 }
1851 else
1852 myNPTf = new NPTf<NPT <RealIntegrator> >(&(info[k]), the_ff);
1853
1854 myNPTf->setTargetTemp(globals->getTargetTemp());
1855
1856 if (globals->haveTargetPressure())
1857 myNPTf->setTargetPressure(globals->getTargetPressure());
1858 else{
1859 sprintf(painCave.errMsg,
1860 "SimSetup error: If you use a constant pressure\n"
1861 "\tensemble, you must set targetPressure in the BASS file.\n");
1862 painCave.isFatal = 1;
1863 simError();
1864 }
1865
1866 if (globals->haveTauThermostat())
1867 myNPTf->setTauThermostat(globals->getTauThermostat());
1868
1869 else{
1870 sprintf(painCave.errMsg,
1871 "SimSetup error: If you use an NPT\n"
1872 "\tensemble, you must set tauThermostat.\n");
1873 painCave.isFatal = 1;
1874 simError();
1875 }
1876
1877 if (globals->haveTauBarostat())
1878 myNPTf->setTauBarostat(globals->getTauBarostat());
1879
1880 else{
1881 sprintf(painCave.errMsg,
1882 "SimSetup error: If you use an NPT\n"
1883 "\tensemble, you must set tauBarostat.\n");
1884 painCave.isFatal = 1;
1885 simError();
1886 }
1887
1888 info->the_integrator = myNPTf;
1889 break;
1890
1891 case NPTxyz_ENS:
1892 if (globals->haveZconstraints()){
1893 setupZConstraint(info[k]);
1894 myNPTxyz = new ZConstraint<NPTxyz<NPT <RealIntegrator> > >(&(info[k]), the_ff);
1895 }
1896 else
1897 myNPTxyz = new NPTxyz<NPT <RealIntegrator> >(&(info[k]), the_ff);
1898
1899 myNPTxyz->setTargetTemp(globals->getTargetTemp());
1900
1901 if (globals->haveTargetPressure())
1902 myNPTxyz->setTargetPressure(globals->getTargetPressure());
1903 else{
1904 sprintf(painCave.errMsg,
1905 "SimSetup error: If you use a constant pressure\n"
1906 "\tensemble, you must set targetPressure in the BASS file.\n");
1907 painCave.isFatal = 1;
1908 simError();
1909 }
1910
1911 if (globals->haveTauThermostat())
1912 myNPTxyz->setTauThermostat(globals->getTauThermostat());
1913 else{
1914 sprintf(painCave.errMsg,
1915 "SimSetup error: If you use an NPT\n"
1916 "\tensemble, you must set tauThermostat.\n");
1917 painCave.isFatal = 1;
1918 simError();
1919 }
1920
1921 if (globals->haveTauBarostat())
1922 myNPTxyz->setTauBarostat(globals->getTauBarostat());
1923 else{
1924 sprintf(painCave.errMsg,
1925 "SimSetup error: If you use an NPT\n"
1926 "\tensemble, you must set tauBarostat.\n");
1927 painCave.isFatal = 1;
1928 simError();
1929 }
1930
1931 info->the_integrator = myNPTxyz;
1932 break;
1933
1934 default:
1935 sprintf(painCave.errMsg,
1936 "SimSetup Error. Unrecognized ensemble in case statement.\n");
1937 painCave.isFatal = 1;
1938 simError();
1939 }
1940 }
1941 }
1942
1943 void SimSetup::initFortran(void){
1944 info[0].refreshSim();
1945
1946 if (!strcmp(info[0].mixingRule, "standard")){
1947 the_ff->initForceField(LB_MIXING_RULE);
1948 }
1949 else if (!strcmp(info[0].mixingRule, "explicit")){
1950 the_ff->initForceField(EXPLICIT_MIXING_RULE);
1951 }
1952 else{
1953 sprintf(painCave.errMsg, "SimSetup Error: unknown mixing rule -> \"%s\"\n",
1954 info[0].mixingRule);
1955 painCave.isFatal = 1;
1956 simError();
1957 }
1958
1959
1960 #ifdef IS_MPI
1961 strcpy(checkPointMsg, "Successfully intialized the mixingRule for Fortran.");
1962 MPIcheckPoint();
1963 #endif // is_mpi
1964 }
1965
1966 void SimSetup::setupZConstraint(SimInfo& theInfo){
1967 int nZConstraints;
1968 ZconStamp** zconStamp;
1969
1970 if (globals->haveZconstraintTime()){
1971 //add sample time of z-constraint into SimInfo's property list
1972 DoubleData* zconsTimeProp = new DoubleData();
1973 zconsTimeProp->setID(ZCONSTIME_ID);
1974 zconsTimeProp->setData(globals->getZconsTime());
1975 theInfo.addProperty(zconsTimeProp);
1976 }
1977 else{
1978 sprintf(painCave.errMsg,
1979 "ZConstraint error: If you use a ZConstraint,\n"
1980 "\tyou must set zconsTime.\n");
1981 painCave.isFatal = 1;
1982 simError();
1983 }
1984
1985 //push zconsTol into siminfo, if user does not specify
1986 //value for zconsTol, a default value will be used
1987 DoubleData* zconsTol = new DoubleData();
1988 zconsTol->setID(ZCONSTOL_ID);
1989 if (globals->haveZconsTol()){
1990 zconsTol->setData(globals->getZconsTol());
1991 }
1992 else{
1993 double defaultZConsTol = 0.01;
1994 sprintf(painCave.errMsg,
1995 "ZConstraint Warning: Tolerance for z-constraint method is not specified.\n"
1996 "\tOOPSE will use a default value of %f.\n"
1997 "\tTo set the tolerance, use the zconsTol variable.\n",
1998 defaultZConsTol);
1999 painCave.isFatal = 0;
2000 simError();
2001
2002 zconsTol->setData(defaultZConsTol);
2003 }
2004 theInfo.addProperty(zconsTol);
2005
2006 //set Force Subtraction Policy
2007 StringData* zconsForcePolicy = new StringData();
2008 zconsForcePolicy->setID(ZCONSFORCEPOLICY_ID);
2009
2010 if (globals->haveZconsForcePolicy()){
2011 zconsForcePolicy->setData(globals->getZconsForcePolicy());
2012 }
2013 else{
2014 sprintf(painCave.errMsg,
2015 "ZConstraint Warning: No force subtraction policy was set.\n"
2016 "\tOOPSE will use PolicyByMass.\n"
2017 "\tTo set the policy, use the zconsForcePolicy variable.\n");
2018 painCave.isFatal = 0;
2019 simError();
2020 zconsForcePolicy->setData("BYMASS");
2021 }
2022
2023 theInfo.addProperty(zconsForcePolicy);
2024
2025 //set zcons gap
2026 DoubleData* zconsGap = new DoubleData();
2027 zconsGap->setID(ZCONSGAP_ID);
2028
2029 if (globals->haveZConsGap()){
2030 zconsGap->setData(globals->getZconsGap());
2031 theInfo.addProperty(zconsGap);
2032 }
2033
2034 //set zcons fixtime
2035 DoubleData* zconsFixtime = new DoubleData();
2036 zconsFixtime->setID(ZCONSFIXTIME_ID);
2037
2038 if (globals->haveZConsFixTime()){
2039 zconsFixtime->setData(globals->getZconsFixtime());
2040 theInfo.addProperty(zconsFixtime);
2041 }
2042
2043 //set zconsUsingSMD
2044 IntData* zconsUsingSMD = new IntData();
2045 zconsUsingSMD->setID(ZCONSUSINGSMD_ID);
2046
2047 if (globals->haveZConsUsingSMD()){
2048 zconsUsingSMD->setData(globals->getZconsUsingSMD());
2049 theInfo.addProperty(zconsUsingSMD);
2050 }
2051
2052 //Determine the name of ouput file and add it into SimInfo's property list
2053 //Be careful, do not use inFileName, since it is a pointer which
2054 //point to a string at master node, and slave nodes do not contain that string
2055
2056 string zconsOutput(theInfo.finalName);
2057
2058 zconsOutput = zconsOutput.substr(0, zconsOutput.rfind(".")) + ".fz";
2059
2060 StringData* zconsFilename = new StringData();
2061 zconsFilename->setID(ZCONSFILENAME_ID);
2062 zconsFilename->setData(zconsOutput);
2063
2064 theInfo.addProperty(zconsFilename);
2065
2066 //setup index, pos and other parameters of z-constraint molecules
2067 nZConstraints = globals->getNzConstraints();
2068 theInfo.nZconstraints = nZConstraints;
2069
2070 zconStamp = globals->getZconStamp();
2071 ZConsParaItem tempParaItem;
2072
2073 ZConsParaData* zconsParaData = new ZConsParaData();
2074 zconsParaData->setID(ZCONSPARADATA_ID);
2075
2076 for (int i = 0; i < nZConstraints; i++){
2077 tempParaItem.havingZPos = zconStamp[i]->haveZpos();
2078 tempParaItem.zPos = zconStamp[i]->getZpos();
2079 tempParaItem.zconsIndex = zconStamp[i]->getMolIndex();
2080 tempParaItem.kRatio = zconStamp[i]->getKratio();
2081 tempParaItem.havingCantVel = zconStamp[i]->haveCantVel();
2082 tempParaItem.cantVel = zconStamp[i]->getCantVel();
2083 zconsParaData->addItem(tempParaItem);
2084 }
2085
2086 //check the uniqueness of index
2087 if(!zconsParaData->isIndexUnique()){
2088 sprintf(painCave.errMsg,
2089 "ZConstraint Error: molIndex is not unique!\n");
2090 painCave.isFatal = 1;
2091 simError();
2092 }
2093
2094 //sort the parameters by index of molecules
2095 zconsParaData->sortByIndex();
2096
2097 //push data into siminfo, therefore, we can retrieve later
2098 theInfo.addProperty(zconsParaData);
2099 }
2100
2101 void SimSetup::makeMinimizer(){
2102
2103 OOPSEMinimizer* myOOPSEMinimizer;
2104 MinimizerParameterSet* param;
2105 char minimizerName[100];
2106
2107 for (int i = 0; i < nInfo; i++){
2108
2109 //prepare parameter set for minimizer
2110 param = new MinimizerParameterSet();
2111 param->setDefaultParameter();
2112
2113 if (globals->haveMinimizer()){
2114 param->setFTol(globals->getMinFTol());
2115 }
2116
2117 if (globals->haveMinGTol()){
2118 param->setGTol(globals->getMinGTol());
2119 }
2120
2121 if (globals->haveMinMaxIter()){
2122 param->setMaxIteration(globals->getMinMaxIter());
2123 }
2124
2125 if (globals->haveMinWriteFrq()){
2126 param->setMaxIteration(globals->getMinMaxIter());
2127 }
2128
2129 if (globals->haveMinWriteFrq()){
2130 param->setWriteFrq(globals->getMinWriteFrq());
2131 }
2132
2133 if (globals->haveMinStepSize()){
2134 param->setStepSize(globals->getMinStepSize());
2135 }
2136
2137 if (globals->haveMinLSMaxIter()){
2138 param->setLineSearchMaxIteration(globals->getMinLSMaxIter());
2139 }
2140
2141 if (globals->haveMinLSTol()){
2142 param->setLineSearchTol(globals->getMinLSTol());
2143 }
2144
2145 strcpy(minimizerName, globals->getMinimizer());
2146
2147 if (!strcasecmp(minimizerName, "CG")){
2148 myOOPSEMinimizer = new PRCGMinimizer(&(info[i]), the_ff, param);
2149 }
2150 else if (!strcasecmp(minimizerName, "SD")){
2151 //myOOPSEMinimizer = MinimizerFactory.creatMinimizer("", &(info[i]), the_ff, param);
2152 myOOPSEMinimizer = new SDMinimizer(&(info[i]), the_ff, param);
2153 }
2154 else{
2155 sprintf(painCave.errMsg,
2156 "SimSetup error: Unrecognized Minimizer, use Conjugate Gradient \n");
2157 painCave.isFatal = 0;
2158 simError();
2159
2160 myOOPSEMinimizer = new PRCGMinimizer(&(info[i]), the_ff, param);
2161 }
2162 info[i].the_integrator = myOOPSEMinimizer;
2163
2164 //store the minimizer into simInfo
2165 info[i].the_minimizer = myOOPSEMinimizer;
2166 info[i].has_minimizer = true;
2167 }
2168
2169 }