ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/brains/SimInfo.cpp
Revision: 2344
Committed: Tue Oct 4 19:34:03 2005 UTC (18 years, 9 months ago) by chrisfen
File size: 35547 byte(s)
Log Message:
fixed an annoying mass ratio bug that results in simulation failure with massless particles

File Contents

# Content
1 /*
2 * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3 *
4 * The University of Notre Dame grants you ("Licensee") a
5 * non-exclusive, royalty free, license to use, modify and
6 * redistribute this software in source and binary code form, provided
7 * that the following conditions are met:
8 *
9 * 1. Acknowledgement of the program authors must be made in any
10 * publication of scientific results based in part on use of the
11 * program. An acceptable form of acknowledgement is citation of
12 * the article in which the program was described (Matthew
13 * A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14 * J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15 * Parallel Simulation Engine for Molecular Dynamics,"
16 * J. Comput. Chem. 26, pp. 252-271 (2005))
17 *
18 * 2. Redistributions of source code must retain the above copyright
19 * notice, this list of conditions and the following disclaimer.
20 *
21 * 3. Redistributions in binary form must reproduce the above copyright
22 * notice, this list of conditions and the following disclaimer in the
23 * documentation and/or other materials provided with the
24 * distribution.
25 *
26 * This software is provided "AS IS," without a warranty of any
27 * kind. All express or implied conditions, representations and
28 * warranties, including any implied warranty of merchantability,
29 * fitness for a particular purpose or non-infringement, are hereby
30 * excluded. The University of Notre Dame and its licensors shall not
31 * be liable for any damages suffered by licensee as a result of
32 * using, modifying or distributing the software or its
33 * derivatives. In no event will the University of Notre Dame or its
34 * licensors be liable for any lost revenue, profit or data, or for
35 * direct, indirect, special, consequential, incidental or punitive
36 * damages, however caused and regardless of the theory of liability,
37 * arising out of the use of or inability to use software, even if the
38 * University of Notre Dame has been advised of the possibility of
39 * such damages.
40 */
41
42 /**
43 * @file SimInfo.cpp
44 * @author tlin
45 * @date 11/02/2004
46 * @version 1.0
47 */
48
49 #include <algorithm>
50 #include <set>
51
52 #include "brains/SimInfo.hpp"
53 #include "math/Vector3.hpp"
54 #include "primitives/Molecule.hpp"
55 #include "UseTheForce/fCutoffPolicy.h"
56 #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
57 #include "UseTheForce/doForces_interface.h"
58 #include "UseTheForce/DarkSide/electrostatic_interface.h"
59 #include "UseTheForce/notifyCutoffs_interface.h"
60 #include "utils/MemoryUtils.hpp"
61 #include "utils/simError.h"
62 #include "selection/SelectionManager.hpp"
63
64 #ifdef IS_MPI
65 #include "UseTheForce/mpiComponentPlan.h"
66 #include "UseTheForce/DarkSide/simParallel_interface.h"
67 #endif
68
69 namespace oopse {
70
71 SimInfo::SimInfo(MakeStamps* stamps, std::vector<std::pair<MoleculeStamp*, int> >& molStampPairs,
72 ForceField* ff, Globals* simParams) :
73 stamps_(stamps), forceField_(ff), simParams_(simParams),
74 ndf_(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0),
75 nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0),
76 nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
77 nAtoms_(0), nBonds_(0), nBends_(0), nTorsions_(0), nRigidBodies_(0),
78 nIntegrableObjects_(0), nCutoffGroups_(0), nConstraints_(0),
79 sman_(NULL), fortranInitialized_(false) {
80
81
82 std::vector<std::pair<MoleculeStamp*, int> >::iterator i;
83 MoleculeStamp* molStamp;
84 int nMolWithSameStamp;
85 int nCutoffAtoms = 0; // number of atoms belong to cutoff groups
86 int nGroups = 0; //total cutoff groups defined in meta-data file
87 CutoffGroupStamp* cgStamp;
88 RigidBodyStamp* rbStamp;
89 int nRigidAtoms = 0;
90
91 for (i = molStampPairs.begin(); i !=molStampPairs.end(); ++i) {
92 molStamp = i->first;
93 nMolWithSameStamp = i->second;
94
95 addMoleculeStamp(molStamp, nMolWithSameStamp);
96
97 //calculate atoms in molecules
98 nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp;
99
100
101 //calculate atoms in cutoff groups
102 int nAtomsInGroups = 0;
103 int nCutoffGroupsInStamp = molStamp->getNCutoffGroups();
104
105 for (int j=0; j < nCutoffGroupsInStamp; j++) {
106 cgStamp = molStamp->getCutoffGroup(j);
107 nAtomsInGroups += cgStamp->getNMembers();
108 }
109
110 nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
111
112 nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;
113
114 //calculate atoms in rigid bodies
115 int nAtomsInRigidBodies = 0;
116 int nRigidBodiesInStamp = molStamp->getNRigidBodies();
117
118 for (int j=0; j < nRigidBodiesInStamp; j++) {
119 rbStamp = molStamp->getRigidBody(j);
120 nAtomsInRigidBodies += rbStamp->getNMembers();
121 }
122
123 nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
124 nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;
125
126 }
127
128 //every free atom (atom does not belong to cutoff groups) is a cutoff
129 //group therefore the total number of cutoff groups in the system is
130 //equal to the total number of atoms minus number of atoms belong to
131 //cutoff group defined in meta-data file plus the number of cutoff
132 //groups defined in meta-data file
133 nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
134
135 //every free atom (atom does not belong to rigid bodies) is an
136 //integrable object therefore the total number of integrable objects
137 //in the system is equal to the total number of atoms minus number of
138 //atoms belong to rigid body defined in meta-data file plus the number
139 //of rigid bodies defined in meta-data file
140 nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms
141 + nGlobalRigidBodies_;
142
143 nGlobalMols_ = molStampIds_.size();
144
145 #ifdef IS_MPI
146 molToProcMap_.resize(nGlobalMols_);
147 #endif
148
149 }
150
151 SimInfo::~SimInfo() {
152 std::map<int, Molecule*>::iterator i;
153 for (i = molecules_.begin(); i != molecules_.end(); ++i) {
154 delete i->second;
155 }
156 molecules_.clear();
157
158 delete stamps_;
159 delete sman_;
160 delete simParams_;
161 delete forceField_;
162 }
163
164 int SimInfo::getNGlobalConstraints() {
165 int nGlobalConstraints;
166 #ifdef IS_MPI
167 MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM,
168 MPI_COMM_WORLD);
169 #else
170 nGlobalConstraints = nConstraints_;
171 #endif
172 return nGlobalConstraints;
173 }
174
175 bool SimInfo::addMolecule(Molecule* mol) {
176 MoleculeIterator i;
177
178 i = molecules_.find(mol->getGlobalIndex());
179 if (i == molecules_.end() ) {
180
181 molecules_.insert(std::make_pair(mol->getGlobalIndex(), mol));
182
183 nAtoms_ += mol->getNAtoms();
184 nBonds_ += mol->getNBonds();
185 nBends_ += mol->getNBends();
186 nTorsions_ += mol->getNTorsions();
187 nRigidBodies_ += mol->getNRigidBodies();
188 nIntegrableObjects_ += mol->getNIntegrableObjects();
189 nCutoffGroups_ += mol->getNCutoffGroups();
190 nConstraints_ += mol->getNConstraintPairs();
191
192 addExcludePairs(mol);
193
194 return true;
195 } else {
196 return false;
197 }
198 }
199
200 bool SimInfo::removeMolecule(Molecule* mol) {
201 MoleculeIterator i;
202 i = molecules_.find(mol->getGlobalIndex());
203
204 if (i != molecules_.end() ) {
205
206 assert(mol == i->second);
207
208 nAtoms_ -= mol->getNAtoms();
209 nBonds_ -= mol->getNBonds();
210 nBends_ -= mol->getNBends();
211 nTorsions_ -= mol->getNTorsions();
212 nRigidBodies_ -= mol->getNRigidBodies();
213 nIntegrableObjects_ -= mol->getNIntegrableObjects();
214 nCutoffGroups_ -= mol->getNCutoffGroups();
215 nConstraints_ -= mol->getNConstraintPairs();
216
217 removeExcludePairs(mol);
218 molecules_.erase(mol->getGlobalIndex());
219
220 delete mol;
221
222 return true;
223 } else {
224 return false;
225 }
226
227
228 }
229
230
231 Molecule* SimInfo::beginMolecule(MoleculeIterator& i) {
232 i = molecules_.begin();
233 return i == molecules_.end() ? NULL : i->second;
234 }
235
236 Molecule* SimInfo::nextMolecule(MoleculeIterator& i) {
237 ++i;
238 return i == molecules_.end() ? NULL : i->second;
239 }
240
241
242 void SimInfo::calcNdf() {
243 int ndf_local;
244 MoleculeIterator i;
245 std::vector<StuntDouble*>::iterator j;
246 Molecule* mol;
247 StuntDouble* integrableObject;
248
249 ndf_local = 0;
250
251 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
252 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
253 integrableObject = mol->nextIntegrableObject(j)) {
254
255 ndf_local += 3;
256
257 if (integrableObject->isDirectional()) {
258 if (integrableObject->isLinear()) {
259 ndf_local += 2;
260 } else {
261 ndf_local += 3;
262 }
263 }
264
265 }//end for (integrableObject)
266 }// end for (mol)
267
268 // n_constraints is local, so subtract them on each processor
269 ndf_local -= nConstraints_;
270
271 #ifdef IS_MPI
272 MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
273 #else
274 ndf_ = ndf_local;
275 #endif
276
277 // nZconstraints_ is global, as are the 3 COM translations for the
278 // entire system:
279 ndf_ = ndf_ - 3 - nZconstraint_;
280
281 }
282
283 void SimInfo::calcNdfRaw() {
284 int ndfRaw_local;
285
286 MoleculeIterator i;
287 std::vector<StuntDouble*>::iterator j;
288 Molecule* mol;
289 StuntDouble* integrableObject;
290
291 // Raw degrees of freedom that we have to set
292 ndfRaw_local = 0;
293
294 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
295 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
296 integrableObject = mol->nextIntegrableObject(j)) {
297
298 ndfRaw_local += 3;
299
300 if (integrableObject->isDirectional()) {
301 if (integrableObject->isLinear()) {
302 ndfRaw_local += 2;
303 } else {
304 ndfRaw_local += 3;
305 }
306 }
307
308 }
309 }
310
311 #ifdef IS_MPI
312 MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
313 #else
314 ndfRaw_ = ndfRaw_local;
315 #endif
316 }
317
318 void SimInfo::calcNdfTrans() {
319 int ndfTrans_local;
320
321 ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
322
323
324 #ifdef IS_MPI
325 MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
326 #else
327 ndfTrans_ = ndfTrans_local;
328 #endif
329
330 ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
331
332 }
333
334 void SimInfo::addExcludePairs(Molecule* mol) {
335 std::vector<Bond*>::iterator bondIter;
336 std::vector<Bend*>::iterator bendIter;
337 std::vector<Torsion*>::iterator torsionIter;
338 Bond* bond;
339 Bend* bend;
340 Torsion* torsion;
341 int a;
342 int b;
343 int c;
344 int d;
345
346 for (bond= mol->beginBond(bondIter); bond != NULL; bond = mol->nextBond(bondIter)) {
347 a = bond->getAtomA()->getGlobalIndex();
348 b = bond->getAtomB()->getGlobalIndex();
349 exclude_.addPair(a, b);
350 }
351
352 for (bend= mol->beginBend(bendIter); bend != NULL; bend = mol->nextBend(bendIter)) {
353 a = bend->getAtomA()->getGlobalIndex();
354 b = bend->getAtomB()->getGlobalIndex();
355 c = bend->getAtomC()->getGlobalIndex();
356
357 exclude_.addPair(a, b);
358 exclude_.addPair(a, c);
359 exclude_.addPair(b, c);
360 }
361
362 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; torsion = mol->nextTorsion(torsionIter)) {
363 a = torsion->getAtomA()->getGlobalIndex();
364 b = torsion->getAtomB()->getGlobalIndex();
365 c = torsion->getAtomC()->getGlobalIndex();
366 d = torsion->getAtomD()->getGlobalIndex();
367
368 exclude_.addPair(a, b);
369 exclude_.addPair(a, c);
370 exclude_.addPair(a, d);
371 exclude_.addPair(b, c);
372 exclude_.addPair(b, d);
373 exclude_.addPair(c, d);
374 }
375
376 Molecule::RigidBodyIterator rbIter;
377 RigidBody* rb;
378 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
379 std::vector<Atom*> atoms = rb->getAtoms();
380 for (int i = 0; i < atoms.size() -1 ; ++i) {
381 for (int j = i + 1; j < atoms.size(); ++j) {
382 a = atoms[i]->getGlobalIndex();
383 b = atoms[j]->getGlobalIndex();
384 exclude_.addPair(a, b);
385 }
386 }
387 }
388
389 }
390
391 void SimInfo::removeExcludePairs(Molecule* mol) {
392 std::vector<Bond*>::iterator bondIter;
393 std::vector<Bend*>::iterator bendIter;
394 std::vector<Torsion*>::iterator torsionIter;
395 Bond* bond;
396 Bend* bend;
397 Torsion* torsion;
398 int a;
399 int b;
400 int c;
401 int d;
402
403 for (bond= mol->beginBond(bondIter); bond != NULL; bond = mol->nextBond(bondIter)) {
404 a = bond->getAtomA()->getGlobalIndex();
405 b = bond->getAtomB()->getGlobalIndex();
406 exclude_.removePair(a, b);
407 }
408
409 for (bend= mol->beginBend(bendIter); bend != NULL; bend = mol->nextBend(bendIter)) {
410 a = bend->getAtomA()->getGlobalIndex();
411 b = bend->getAtomB()->getGlobalIndex();
412 c = bend->getAtomC()->getGlobalIndex();
413
414 exclude_.removePair(a, b);
415 exclude_.removePair(a, c);
416 exclude_.removePair(b, c);
417 }
418
419 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; torsion = mol->nextTorsion(torsionIter)) {
420 a = torsion->getAtomA()->getGlobalIndex();
421 b = torsion->getAtomB()->getGlobalIndex();
422 c = torsion->getAtomC()->getGlobalIndex();
423 d = torsion->getAtomD()->getGlobalIndex();
424
425 exclude_.removePair(a, b);
426 exclude_.removePair(a, c);
427 exclude_.removePair(a, d);
428 exclude_.removePair(b, c);
429 exclude_.removePair(b, d);
430 exclude_.removePair(c, d);
431 }
432
433 Molecule::RigidBodyIterator rbIter;
434 RigidBody* rb;
435 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
436 std::vector<Atom*> atoms = rb->getAtoms();
437 for (int i = 0; i < atoms.size() -1 ; ++i) {
438 for (int j = i + 1; j < atoms.size(); ++j) {
439 a = atoms[i]->getGlobalIndex();
440 b = atoms[j]->getGlobalIndex();
441 exclude_.removePair(a, b);
442 }
443 }
444 }
445
446 }
447
448
449 void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) {
450 int curStampId;
451
452 //index from 0
453 curStampId = moleculeStamps_.size();
454
455 moleculeStamps_.push_back(molStamp);
456 molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
457 }
458
459 void SimInfo::update() {
460
461 setupSimType();
462
463 #ifdef IS_MPI
464 setupFortranParallel();
465 #endif
466
467 setupFortranSim();
468
469 //setup fortran force field
470 /** @deprecate */
471 int isError = 0;
472
473 setupElectrostaticSummationMethod( isError );
474
475 if(isError){
476 sprintf( painCave.errMsg,
477 "ForceField error: There was an error initializing the forceField in fortran.\n" );
478 painCave.isFatal = 1;
479 simError();
480 }
481
482
483 setupCutoff();
484
485 calcNdf();
486 calcNdfRaw();
487 calcNdfTrans();
488
489 fortranInitialized_ = true;
490 }
491
492 std::set<AtomType*> SimInfo::getUniqueAtomTypes() {
493 SimInfo::MoleculeIterator mi;
494 Molecule* mol;
495 Molecule::AtomIterator ai;
496 Atom* atom;
497 std::set<AtomType*> atomTypes;
498
499 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
500
501 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
502 atomTypes.insert(atom->getAtomType());
503 }
504
505 }
506
507 return atomTypes;
508 }
509
510 void SimInfo::setupSimType() {
511 std::set<AtomType*>::iterator i;
512 std::set<AtomType*> atomTypes;
513 atomTypes = getUniqueAtomTypes();
514
515 int useLennardJones = 0;
516 int useElectrostatic = 0;
517 int useEAM = 0;
518 int useCharge = 0;
519 int useDirectional = 0;
520 int useDipole = 0;
521 int useGayBerne = 0;
522 int useSticky = 0;
523 int useStickyPower = 0;
524 int useShape = 0;
525 int useFLARB = 0; //it is not in AtomType yet
526 int useDirectionalAtom = 0;
527 int useElectrostatics = 0;
528 //usePBC and useRF are from simParams
529 int usePBC = simParams_->getPBC();
530 int useRF;
531
532 // set the useRF logical
533 std::string myMethod = simParams_->getElectrostaticSummationMethod();
534 if (myMethod == "REACTION_FIELD")
535 useRF = 1;
536 else
537 useRF = 0;
538
539 //loop over all of the atom types
540 for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
541 useLennardJones |= (*i)->isLennardJones();
542 useElectrostatic |= (*i)->isElectrostatic();
543 useEAM |= (*i)->isEAM();
544 useCharge |= (*i)->isCharge();
545 useDirectional |= (*i)->isDirectional();
546 useDipole |= (*i)->isDipole();
547 useGayBerne |= (*i)->isGayBerne();
548 useSticky |= (*i)->isSticky();
549 useStickyPower |= (*i)->isStickyPower();
550 useShape |= (*i)->isShape();
551 }
552
553 if (useSticky || useStickyPower || useDipole || useGayBerne || useShape) {
554 useDirectionalAtom = 1;
555 }
556
557 if (useCharge || useDipole) {
558 useElectrostatics = 1;
559 }
560
561 #ifdef IS_MPI
562 int temp;
563
564 temp = usePBC;
565 MPI_Allreduce(&temp, &usePBC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
566
567 temp = useDirectionalAtom;
568 MPI_Allreduce(&temp, &useDirectionalAtom, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
569
570 temp = useLennardJones;
571 MPI_Allreduce(&temp, &useLennardJones, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
572
573 temp = useElectrostatics;
574 MPI_Allreduce(&temp, &useElectrostatics, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
575
576 temp = useCharge;
577 MPI_Allreduce(&temp, &useCharge, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
578
579 temp = useDipole;
580 MPI_Allreduce(&temp, &useDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
581
582 temp = useSticky;
583 MPI_Allreduce(&temp, &useSticky, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
584
585 temp = useStickyPower;
586 MPI_Allreduce(&temp, &useStickyPower, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
587
588 temp = useGayBerne;
589 MPI_Allreduce(&temp, &useGayBerne, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
590
591 temp = useEAM;
592 MPI_Allreduce(&temp, &useEAM, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
593
594 temp = useShape;
595 MPI_Allreduce(&temp, &useShape, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
596
597 temp = useFLARB;
598 MPI_Allreduce(&temp, &useFLARB, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
599
600 temp = useRF;
601 MPI_Allreduce(&temp, &useRF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
602
603 #endif
604
605 fInfo_.SIM_uses_PBC = usePBC;
606 fInfo_.SIM_uses_DirectionalAtoms = useDirectionalAtom;
607 fInfo_.SIM_uses_LennardJones = useLennardJones;
608 fInfo_.SIM_uses_Electrostatics = useElectrostatics;
609 fInfo_.SIM_uses_Charges = useCharge;
610 fInfo_.SIM_uses_Dipoles = useDipole;
611 fInfo_.SIM_uses_Sticky = useSticky;
612 fInfo_.SIM_uses_StickyPower = useStickyPower;
613 fInfo_.SIM_uses_GayBerne = useGayBerne;
614 fInfo_.SIM_uses_EAM = useEAM;
615 fInfo_.SIM_uses_Shapes = useShape;
616 fInfo_.SIM_uses_FLARB = useFLARB;
617 fInfo_.SIM_uses_RF = useRF;
618
619 if( fInfo_.SIM_uses_Dipoles && myMethod == "REACTION_FIELD") {
620
621 if (simParams_->haveDielectric()) {
622 fInfo_.dielect = simParams_->getDielectric();
623 } else {
624 sprintf(painCave.errMsg,
625 "SimSetup Error: No Dielectric constant was set.\n"
626 "\tYou are trying to use Reaction Field without"
627 "\tsetting a dielectric constant!\n");
628 painCave.isFatal = 1;
629 simError();
630 }
631
632 } else {
633 fInfo_.dielect = 0.0;
634 }
635
636 }
637
638 void SimInfo::setupFortranSim() {
639 int isError;
640 int nExclude;
641 std::vector<int> fortranGlobalGroupMembership;
642
643 nExclude = exclude_.getSize();
644 isError = 0;
645
646 //globalGroupMembership_ is filled by SimCreator
647 for (int i = 0; i < nGlobalAtoms_; i++) {
648 fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
649 }
650
651 //calculate mass ratio of cutoff group
652 std::vector<double> mfact;
653 SimInfo::MoleculeIterator mi;
654 Molecule* mol;
655 Molecule::CutoffGroupIterator ci;
656 CutoffGroup* cg;
657 Molecule::AtomIterator ai;
658 Atom* atom;
659 double totalMass;
660
661 //to avoid memory reallocation, reserve enough space for mfact
662 mfact.reserve(getNCutoffGroups());
663
664 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
665 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
666
667 totalMass = cg->getMass();
668 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
669 // Check for massless groups - set mfact to 1 if true
670 if (totalMass != 0)
671 mfact.push_back(atom->getMass()/totalMass);
672 else
673 mfact.push_back( 1.0 );
674 }
675
676 }
677 }
678
679 //fill ident array of local atoms (it is actually ident of AtomType, it is so confusing !!!)
680 std::vector<int> identArray;
681
682 //to avoid memory reallocation, reserve enough space identArray
683 identArray.reserve(getNAtoms());
684
685 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
686 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
687 identArray.push_back(atom->getIdent());
688 }
689 }
690
691 //fill molMembershipArray
692 //molMembershipArray is filled by SimCreator
693 std::vector<int> molMembershipArray(nGlobalAtoms_);
694 for (int i = 0; i < nGlobalAtoms_; i++) {
695 molMembershipArray[i] = globalMolMembership_[i] + 1;
696 }
697
698 //setup fortran simulation
699 int nGlobalExcludes = 0;
700 int* globalExcludes = NULL;
701 int* excludeList = exclude_.getExcludeList();
702 setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray[0], &nExclude, excludeList ,
703 &nGlobalExcludes, globalExcludes, &molMembershipArray[0],
704 &mfact[0], &nCutoffGroups_, &fortranGlobalGroupMembership[0], &isError);
705
706 if( isError ){
707
708 sprintf( painCave.errMsg,
709 "There was an error setting the simulation information in fortran.\n" );
710 painCave.isFatal = 1;
711 painCave.severity = OOPSE_ERROR;
712 simError();
713 }
714
715 #ifdef IS_MPI
716 sprintf( checkPointMsg,
717 "succesfully sent the simulation information to fortran.\n");
718 MPIcheckPoint();
719 #endif // is_mpi
720 }
721
722
723 #ifdef IS_MPI
724 void SimInfo::setupFortranParallel() {
725
726 //SimInfo is responsible for creating localToGlobalAtomIndex and localToGlobalGroupIndex
727 std::vector<int> localToGlobalAtomIndex(getNAtoms(), 0);
728 std::vector<int> localToGlobalCutoffGroupIndex;
729 SimInfo::MoleculeIterator mi;
730 Molecule::AtomIterator ai;
731 Molecule::CutoffGroupIterator ci;
732 Molecule* mol;
733 Atom* atom;
734 CutoffGroup* cg;
735 mpiSimData parallelData;
736 int isError;
737
738 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
739
740 //local index(index in DataStorge) of atom is important
741 for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
742 localToGlobalAtomIndex[atom->getLocalIndex()] = atom->getGlobalIndex() + 1;
743 }
744
745 //local index of cutoff group is trivial, it only depends on the order of travesing
746 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
747 localToGlobalCutoffGroupIndex.push_back(cg->getGlobalIndex() + 1);
748 }
749
750 }
751
752 //fill up mpiSimData struct
753 parallelData.nMolGlobal = getNGlobalMolecules();
754 parallelData.nMolLocal = getNMolecules();
755 parallelData.nAtomsGlobal = getNGlobalAtoms();
756 parallelData.nAtomsLocal = getNAtoms();
757 parallelData.nGroupsGlobal = getNGlobalCutoffGroups();
758 parallelData.nGroupsLocal = getNCutoffGroups();
759 parallelData.myNode = worldRank;
760 MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors));
761
762 //pass mpiSimData struct and index arrays to fortran
763 setFsimParallel(&parallelData, &(parallelData.nAtomsLocal),
764 &localToGlobalAtomIndex[0], &(parallelData.nGroupsLocal),
765 &localToGlobalCutoffGroupIndex[0], &isError);
766
767 if (isError) {
768 sprintf(painCave.errMsg,
769 "mpiRefresh errror: fortran didn't like something we gave it.\n");
770 painCave.isFatal = 1;
771 simError();
772 }
773
774 sprintf(checkPointMsg, " mpiRefresh successful.\n");
775 MPIcheckPoint();
776
777
778 }
779
780 #endif
781
782 double SimInfo::calcMaxCutoffRadius() {
783
784
785 std::set<AtomType*> atomTypes;
786 std::set<AtomType*>::iterator i;
787 std::vector<double> cutoffRadius;
788
789 //get the unique atom types
790 atomTypes = getUniqueAtomTypes();
791
792 //query the max cutoff radius among these atom types
793 for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
794 cutoffRadius.push_back(forceField_->getRcutFromAtomType(*i));
795 }
796
797 double maxCutoffRadius = *(std::max_element(cutoffRadius.begin(), cutoffRadius.end()));
798 #ifdef IS_MPI
799 //pick the max cutoff radius among the processors
800 #endif
801
802 return maxCutoffRadius;
803 }
804
805 void SimInfo::getCutoff(double& rcut, double& rsw) {
806
807 if (fInfo_.SIM_uses_Charges | fInfo_.SIM_uses_Dipoles | fInfo_.SIM_uses_RF) {
808
809 if (!simParams_->haveRcut()){
810 sprintf(painCave.errMsg,
811 "SimCreator Warning: No value was set for the cutoffRadius.\n"
812 "\tOOPSE will use a default value of 15.0 angstroms"
813 "\tfor the cutoffRadius.\n");
814 painCave.isFatal = 0;
815 simError();
816 rcut = 15.0;
817 } else{
818 rcut = simParams_->getRcut();
819 }
820
821 if (!simParams_->haveRsw()){
822 sprintf(painCave.errMsg,
823 "SimCreator Warning: No value was set for switchingRadius.\n"
824 "\tOOPSE will use a default value of\n"
825 "\t0.95 * cutoffRadius for the switchingRadius\n");
826 painCave.isFatal = 0;
827 simError();
828 rsw = 0.95 * rcut;
829 } else{
830 rsw = simParams_->getRsw();
831 }
832
833 } else {
834 // if charge, dipole or reaction field is not used and the cutofff radius is not specified in
835 //meta-data file, the maximum cutoff radius calculated from forcefiled will be used
836
837 if (simParams_->haveRcut()) {
838 rcut = simParams_->getRcut();
839 } else {
840 //set cutoff radius to the maximum cutoff radius based on atom types in the whole system
841 rcut = calcMaxCutoffRadius();
842 }
843
844 if (simParams_->haveRsw()) {
845 rsw = simParams_->getRsw();
846 } else {
847 rsw = rcut;
848 }
849
850 }
851 }
852
853 void SimInfo::setupCutoff() {
854 getCutoff(rcut_, rsw_);
855 double rnblist = rcut_ + 1; // skin of neighbor list
856
857 //Pass these cutoff radius etc. to fortran. This function should be called once and only once
858
859 int cp = TRADITIONAL_CUTOFF_POLICY;
860 if (simParams_->haveCutoffPolicy()) {
861 std::string myPolicy = simParams_->getCutoffPolicy();
862 if (myPolicy == "MIX") {
863 cp = MIX_CUTOFF_POLICY;
864 } else {
865 if (myPolicy == "MAX") {
866 cp = MAX_CUTOFF_POLICY;
867 } else {
868 if (myPolicy == "TRADITIONAL") {
869 cp = TRADITIONAL_CUTOFF_POLICY;
870 } else {
871 // throw error
872 sprintf( painCave.errMsg,
873 "SimInfo error: Unknown cutoffPolicy. (Input file specified %s .)\n\tcutoffPolicy must be one of: \"Mix\", \"Max\", or \"Traditional\".", myPolicy.c_str() );
874 painCave.isFatal = 1;
875 simError();
876 }
877 }
878 }
879 }
880
881
882 if (simParams_->haveSkinThickness()) {
883 double skinThickness = simParams_->getSkinThickness();
884 }
885
886 notifyFortranCutoffs(&rcut_, &rsw_, &rnblist, &cp);
887 // also send cutoff notification to electrostatics
888 setElectrostaticCutoffRadius(&rcut_);
889 }
890
891 void SimInfo::setupElectrostaticSummationMethod( int isError ) {
892
893 int errorOut;
894 int esm = NONE;
895 double alphaVal;
896 double dielectric;
897
898 errorOut = isError;
899 alphaVal = simParams_->getDampingAlpha();
900 dielectric = simParams_->getDielectric();
901
902 if (simParams_->haveElectrostaticSummationMethod()) {
903 std::string myMethod = simParams_->getElectrostaticSummationMethod();
904 if (myMethod == "NONE") {
905 esm = NONE;
906 } else {
907 if (myMethod == "UNDAMPED_WOLF") {
908 esm = UNDAMPED_WOLF;
909 } else {
910 if (myMethod == "DAMPED_WOLF") {
911 esm = DAMPED_WOLF;
912 if (!simParams_->haveDampingAlpha()) {
913 //throw error
914 sprintf( painCave.errMsg,
915 "SimInfo warning: dampingAlpha was not specified in the input file. A default value of %f (1/ang) will be used for the Damped Wolf Method.", alphaVal);
916 painCave.isFatal = 0;
917 simError();
918 }
919 } else {
920 if (myMethod == "REACTION_FIELD") {
921 esm = REACTION_FIELD;
922 } else {
923 // throw error
924 sprintf( painCave.errMsg,
925 "SimInfo error: Unknown electrostaticSummationMethod. (Input file specified %s .)\n\telectrostaticSummationMethod must be one of: \"none\", \"undamped_wolf\", \"damped_wolf\", or \"reaction_field\".", myMethod.c_str() );
926 painCave.isFatal = 1;
927 simError();
928 }
929 }
930 }
931 }
932 }
933 // let's pass some summation method variables to fortran
934 setElectrostaticSummationMethod( &esm );
935 setDampedWolfAlpha( &alphaVal );
936 setReactionFieldDielectric( &dielectric );
937 initFortranFF( &esm, &errorOut );
938 }
939
940 void SimInfo::addProperty(GenericData* genData) {
941 properties_.addProperty(genData);
942 }
943
944 void SimInfo::removeProperty(const std::string& propName) {
945 properties_.removeProperty(propName);
946 }
947
948 void SimInfo::clearProperties() {
949 properties_.clearProperties();
950 }
951
952 std::vector<std::string> SimInfo::getPropertyNames() {
953 return properties_.getPropertyNames();
954 }
955
956 std::vector<GenericData*> SimInfo::getProperties() {
957 return properties_.getProperties();
958 }
959
960 GenericData* SimInfo::getPropertyByName(const std::string& propName) {
961 return properties_.getPropertyByName(propName);
962 }
963
964 void SimInfo::setSnapshotManager(SnapshotManager* sman) {
965 if (sman_ == sman) {
966 return;
967 }
968 delete sman_;
969 sman_ = sman;
970
971 Molecule* mol;
972 RigidBody* rb;
973 Atom* atom;
974 SimInfo::MoleculeIterator mi;
975 Molecule::RigidBodyIterator rbIter;
976 Molecule::AtomIterator atomIter;;
977
978 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
979
980 for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) {
981 atom->setSnapshotManager(sman_);
982 }
983
984 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
985 rb->setSnapshotManager(sman_);
986 }
987 }
988
989 }
990
991 Vector3d SimInfo::getComVel(){
992 SimInfo::MoleculeIterator i;
993 Molecule* mol;
994
995 Vector3d comVel(0.0);
996 double totalMass = 0.0;
997
998
999 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1000 double mass = mol->getMass();
1001 totalMass += mass;
1002 comVel += mass * mol->getComVel();
1003 }
1004
1005 #ifdef IS_MPI
1006 double tmpMass = totalMass;
1007 Vector3d tmpComVel(comVel);
1008 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1009 MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1010 #endif
1011
1012 comVel /= totalMass;
1013
1014 return comVel;
1015 }
1016
1017 Vector3d SimInfo::getCom(){
1018 SimInfo::MoleculeIterator i;
1019 Molecule* mol;
1020
1021 Vector3d com(0.0);
1022 double totalMass = 0.0;
1023
1024 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1025 double mass = mol->getMass();
1026 totalMass += mass;
1027 com += mass * mol->getCom();
1028 }
1029
1030 #ifdef IS_MPI
1031 double tmpMass = totalMass;
1032 Vector3d tmpCom(com);
1033 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1034 MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1035 #endif
1036
1037 com /= totalMass;
1038
1039 return com;
1040
1041 }
1042
1043 std::ostream& operator <<(std::ostream& o, SimInfo& info) {
1044
1045 return o;
1046 }
1047
1048
1049 /*
1050 Returns center of mass and center of mass velocity in one function call.
1051 */
1052
1053 void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){
1054 SimInfo::MoleculeIterator i;
1055 Molecule* mol;
1056
1057
1058 double totalMass = 0.0;
1059
1060
1061 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1062 double mass = mol->getMass();
1063 totalMass += mass;
1064 com += mass * mol->getCom();
1065 comVel += mass * mol->getComVel();
1066 }
1067
1068 #ifdef IS_MPI
1069 double tmpMass = totalMass;
1070 Vector3d tmpCom(com);
1071 Vector3d tmpComVel(comVel);
1072 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1073 MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1074 MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1075 #endif
1076
1077 com /= totalMass;
1078 comVel /= totalMass;
1079 }
1080
1081 /*
1082 Return intertia tensor for entire system and angular momentum Vector.
1083
1084
1085 [ Ixx -Ixy -Ixz ]
1086 J =| -Iyx Iyy -Iyz |
1087 [ -Izx -Iyz Izz ]
1088 */
1089
1090 void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){
1091
1092
1093 double xx = 0.0;
1094 double yy = 0.0;
1095 double zz = 0.0;
1096 double xy = 0.0;
1097 double xz = 0.0;
1098 double yz = 0.0;
1099 Vector3d com(0.0);
1100 Vector3d comVel(0.0);
1101
1102 getComAll(com, comVel);
1103
1104 SimInfo::MoleculeIterator i;
1105 Molecule* mol;
1106
1107 Vector3d thisq(0.0);
1108 Vector3d thisv(0.0);
1109
1110 double thisMass = 0.0;
1111
1112
1113
1114
1115 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1116
1117 thisq = mol->getCom()-com;
1118 thisv = mol->getComVel()-comVel;
1119 thisMass = mol->getMass();
1120 // Compute moment of intertia coefficients.
1121 xx += thisq[0]*thisq[0]*thisMass;
1122 yy += thisq[1]*thisq[1]*thisMass;
1123 zz += thisq[2]*thisq[2]*thisMass;
1124
1125 // compute products of intertia
1126 xy += thisq[0]*thisq[1]*thisMass;
1127 xz += thisq[0]*thisq[2]*thisMass;
1128 yz += thisq[1]*thisq[2]*thisMass;
1129
1130 angularMomentum += cross( thisq, thisv ) * thisMass;
1131
1132 }
1133
1134
1135 inertiaTensor(0,0) = yy + zz;
1136 inertiaTensor(0,1) = -xy;
1137 inertiaTensor(0,2) = -xz;
1138 inertiaTensor(1,0) = -xy;
1139 inertiaTensor(1,1) = xx + zz;
1140 inertiaTensor(1,2) = -yz;
1141 inertiaTensor(2,0) = -xz;
1142 inertiaTensor(2,1) = -yz;
1143 inertiaTensor(2,2) = xx + yy;
1144
1145 #ifdef IS_MPI
1146 Mat3x3d tmpI(inertiaTensor);
1147 Vector3d tmpAngMom;
1148 MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1149 MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1150 #endif
1151
1152 return;
1153 }
1154
1155 //Returns the angular momentum of the system
1156 Vector3d SimInfo::getAngularMomentum(){
1157
1158 Vector3d com(0.0);
1159 Vector3d comVel(0.0);
1160 Vector3d angularMomentum(0.0);
1161
1162 getComAll(com,comVel);
1163
1164 SimInfo::MoleculeIterator i;
1165 Molecule* mol;
1166
1167 Vector3d thisr(0.0);
1168 Vector3d thisp(0.0);
1169
1170 double thisMass;
1171
1172 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1173 thisMass = mol->getMass();
1174 thisr = mol->getCom()-com;
1175 thisp = (mol->getComVel()-comVel)*thisMass;
1176
1177 angularMomentum += cross( thisr, thisp );
1178
1179 }
1180
1181 #ifdef IS_MPI
1182 Vector3d tmpAngMom;
1183 MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1184 #endif
1185
1186 return angularMomentum;
1187 }
1188
1189
1190 }//end namespace oopse
1191