ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/SimSetup.cpp
Revision: 1452
Committed: Mon Aug 23 15:11:36 2004 UTC (19 years, 10 months ago) by tim
File size: 63570 byte(s)
Log Message:
*** empty log message ***

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