ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/brains/SimInfo.cpp
Revision: 2448
Committed: Wed Nov 16 23:10:02 2005 UTC (18 years, 7 months ago) by tim
File size: 43157 byte(s)
Log Message:
OptionSectionParser get compiled

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 #include <map>
52
53 #include "brains/SimInfo.hpp"
54 #include "math/Vector3.hpp"
55 #include "primitives/Molecule.hpp"
56 #include "UseTheForce/fCutoffPolicy.h"
57 #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
58 #include "UseTheForce/DarkSide/fElectrostaticScreeningMethod.h"
59 #include "UseTheForce/DarkSide/fSwitchingFunctionType.h"
60 #include "UseTheForce/doForces_interface.h"
61 #include "UseTheForce/DarkSide/electrostatic_interface.h"
62 #include "UseTheForce/notifyCutoffs_interface.h"
63 #include "UseTheForce/DarkSide/switcheroo_interface.h"
64 #include "utils/MemoryUtils.hpp"
65 #include "utils/simError.h"
66 #include "selection/SelectionManager.hpp"
67
68 #ifdef IS_MPI
69 #include "UseTheForce/mpiComponentPlan.h"
70 #include "UseTheForce/DarkSide/simParallel_interface.h"
71 #endif
72
73 namespace oopse {
74 std::set<int> getRigidSet(int index, std::map<int, std::set<int> >& container) {
75 std::map<int, std::set<int> >::iterator i = container.find(index);
76 std::set<int> result;
77 if (i != container.end()) {
78 result = i->second;
79 }
80
81 return result;
82 }
83
84 SimInfo::SimInfo(MakeStamps* stamps, std::vector<std::pair<MoleculeStamp*, int> >& molStampPairs,
85 ForceField* ff, Globals* simParams) :
86 stamps_(stamps), forceField_(ff), simParams_(simParams),
87 ndf_(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0),
88 nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0),
89 nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
90 nAtoms_(0), nBonds_(0), nBends_(0), nTorsions_(0), nRigidBodies_(0),
91 nIntegrableObjects_(0), nCutoffGroups_(0), nConstraints_(0),
92 sman_(NULL), fortranInitialized_(false) {
93
94
95 std::vector<std::pair<MoleculeStamp*, int> >::iterator i;
96 MoleculeStamp* molStamp;
97 int nMolWithSameStamp;
98 int nCutoffAtoms = 0; // number of atoms belong to cutoff groups
99 int nGroups = 0; //total cutoff groups defined in meta-data file
100 CutoffGroupStamp* cgStamp;
101 RigidBodyStamp* rbStamp;
102 int nRigidAtoms = 0;
103
104 for (i = molStampPairs.begin(); i !=molStampPairs.end(); ++i) {
105 molStamp = i->first;
106 nMolWithSameStamp = i->second;
107
108 addMoleculeStamp(molStamp, nMolWithSameStamp);
109
110 //calculate atoms in molecules
111 nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp;
112
113
114 //calculate atoms in cutoff groups
115 int nAtomsInGroups = 0;
116 int nCutoffGroupsInStamp = molStamp->getNCutoffGroups();
117
118 for (int j=0; j < nCutoffGroupsInStamp; j++) {
119 cgStamp = molStamp->getCutoffGroup(j);
120 nAtomsInGroups += cgStamp->getNMembers();
121 }
122
123 nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
124
125 nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;
126
127 //calculate atoms in rigid bodies
128 int nAtomsInRigidBodies = 0;
129 int nRigidBodiesInStamp = molStamp->getNRigidBodies();
130
131 for (int j=0; j < nRigidBodiesInStamp; j++) {
132 rbStamp = molStamp->getRigidBody(j);
133 nAtomsInRigidBodies += rbStamp->getNMembers();
134 }
135
136 nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
137 nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;
138
139 }
140
141 //every free atom (atom does not belong to cutoff groups) is a cutoff
142 //group therefore the total number of cutoff groups in the system is
143 //equal to the total number of atoms minus number of atoms belong to
144 //cutoff group defined in meta-data file plus the number of cutoff
145 //groups defined in meta-data file
146 nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
147
148 //every free atom (atom does not belong to rigid bodies) is an
149 //integrable object therefore the total number of integrable objects
150 //in the system is equal to the total number of atoms minus number of
151 //atoms belong to rigid body defined in meta-data file plus the number
152 //of rigid bodies defined in meta-data file
153 nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms
154 + nGlobalRigidBodies_;
155
156 nGlobalMols_ = molStampIds_.size();
157
158 #ifdef IS_MPI
159 molToProcMap_.resize(nGlobalMols_);
160 #endif
161
162 }
163
164 SimInfo::~SimInfo() {
165 std::map<int, Molecule*>::iterator i;
166 for (i = molecules_.begin(); i != molecules_.end(); ++i) {
167 delete i->second;
168 }
169 molecules_.clear();
170
171 delete stamps_;
172 delete sman_;
173 delete simParams_;
174 delete forceField_;
175 }
176
177 int SimInfo::getNGlobalConstraints() {
178 int nGlobalConstraints;
179 #ifdef IS_MPI
180 MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM,
181 MPI_COMM_WORLD);
182 #else
183 nGlobalConstraints = nConstraints_;
184 #endif
185 return nGlobalConstraints;
186 }
187
188 bool SimInfo::addMolecule(Molecule* mol) {
189 MoleculeIterator i;
190
191 i = molecules_.find(mol->getGlobalIndex());
192 if (i == molecules_.end() ) {
193
194 molecules_.insert(std::make_pair(mol->getGlobalIndex(), mol));
195
196 nAtoms_ += mol->getNAtoms();
197 nBonds_ += mol->getNBonds();
198 nBends_ += mol->getNBends();
199 nTorsions_ += mol->getNTorsions();
200 nRigidBodies_ += mol->getNRigidBodies();
201 nIntegrableObjects_ += mol->getNIntegrableObjects();
202 nCutoffGroups_ += mol->getNCutoffGroups();
203 nConstraints_ += mol->getNConstraintPairs();
204
205 addExcludePairs(mol);
206
207 return true;
208 } else {
209 return false;
210 }
211 }
212
213 bool SimInfo::removeMolecule(Molecule* mol) {
214 MoleculeIterator i;
215 i = molecules_.find(mol->getGlobalIndex());
216
217 if (i != molecules_.end() ) {
218
219 assert(mol == i->second);
220
221 nAtoms_ -= mol->getNAtoms();
222 nBonds_ -= mol->getNBonds();
223 nBends_ -= mol->getNBends();
224 nTorsions_ -= mol->getNTorsions();
225 nRigidBodies_ -= mol->getNRigidBodies();
226 nIntegrableObjects_ -= mol->getNIntegrableObjects();
227 nCutoffGroups_ -= mol->getNCutoffGroups();
228 nConstraints_ -= mol->getNConstraintPairs();
229
230 removeExcludePairs(mol);
231 molecules_.erase(mol->getGlobalIndex());
232
233 delete mol;
234
235 return true;
236 } else {
237 return false;
238 }
239
240
241 }
242
243
244 Molecule* SimInfo::beginMolecule(MoleculeIterator& i) {
245 i = molecules_.begin();
246 return i == molecules_.end() ? NULL : i->second;
247 }
248
249 Molecule* SimInfo::nextMolecule(MoleculeIterator& i) {
250 ++i;
251 return i == molecules_.end() ? NULL : i->second;
252 }
253
254
255 void SimInfo::calcNdf() {
256 int ndf_local;
257 MoleculeIterator i;
258 std::vector<StuntDouble*>::iterator j;
259 Molecule* mol;
260 StuntDouble* integrableObject;
261
262 ndf_local = 0;
263
264 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
265 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
266 integrableObject = mol->nextIntegrableObject(j)) {
267
268 ndf_local += 3;
269
270 if (integrableObject->isDirectional()) {
271 if (integrableObject->isLinear()) {
272 ndf_local += 2;
273 } else {
274 ndf_local += 3;
275 }
276 }
277
278 }//end for (integrableObject)
279 }// end for (mol)
280
281 // n_constraints is local, so subtract them on each processor
282 ndf_local -= nConstraints_;
283
284 #ifdef IS_MPI
285 MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
286 #else
287 ndf_ = ndf_local;
288 #endif
289
290 // nZconstraints_ is global, as are the 3 COM translations for the
291 // entire system:
292 ndf_ = ndf_ - 3 - nZconstraint_;
293
294 }
295
296 void SimInfo::calcNdfRaw() {
297 int ndfRaw_local;
298
299 MoleculeIterator i;
300 std::vector<StuntDouble*>::iterator j;
301 Molecule* mol;
302 StuntDouble* integrableObject;
303
304 // Raw degrees of freedom that we have to set
305 ndfRaw_local = 0;
306
307 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
308 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
309 integrableObject = mol->nextIntegrableObject(j)) {
310
311 ndfRaw_local += 3;
312
313 if (integrableObject->isDirectional()) {
314 if (integrableObject->isLinear()) {
315 ndfRaw_local += 2;
316 } else {
317 ndfRaw_local += 3;
318 }
319 }
320
321 }
322 }
323
324 #ifdef IS_MPI
325 MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
326 #else
327 ndfRaw_ = ndfRaw_local;
328 #endif
329 }
330
331 void SimInfo::calcNdfTrans() {
332 int ndfTrans_local;
333
334 ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
335
336
337 #ifdef IS_MPI
338 MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
339 #else
340 ndfTrans_ = ndfTrans_local;
341 #endif
342
343 ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
344
345 }
346
347 void SimInfo::addExcludePairs(Molecule* mol) {
348 std::vector<Bond*>::iterator bondIter;
349 std::vector<Bend*>::iterator bendIter;
350 std::vector<Torsion*>::iterator torsionIter;
351 Bond* bond;
352 Bend* bend;
353 Torsion* torsion;
354 int a;
355 int b;
356 int c;
357 int d;
358
359 std::map<int, std::set<int> > atomGroups;
360
361 Molecule::RigidBodyIterator rbIter;
362 RigidBody* rb;
363 Molecule::IntegrableObjectIterator ii;
364 StuntDouble* integrableObject;
365
366 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
367 integrableObject = mol->nextIntegrableObject(ii)) {
368
369 if (integrableObject->isRigidBody()) {
370 rb = static_cast<RigidBody*>(integrableObject);
371 std::vector<Atom*> atoms = rb->getAtoms();
372 std::set<int> rigidAtoms;
373 for (int i = 0; i < atoms.size(); ++i) {
374 rigidAtoms.insert(atoms[i]->getGlobalIndex());
375 }
376 for (int i = 0; i < atoms.size(); ++i) {
377 atomGroups.insert(std::map<int, std::set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
378 }
379 } else {
380 std::set<int> oneAtomSet;
381 oneAtomSet.insert(integrableObject->getGlobalIndex());
382 atomGroups.insert(std::map<int, std::set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));
383 }
384 }
385
386
387
388 for (bond= mol->beginBond(bondIter); bond != NULL; bond = mol->nextBond(bondIter)) {
389 a = bond->getAtomA()->getGlobalIndex();
390 b = bond->getAtomB()->getGlobalIndex();
391 exclude_.addPair(a, b);
392 }
393
394 for (bend= mol->beginBend(bendIter); bend != NULL; bend = mol->nextBend(bendIter)) {
395 a = bend->getAtomA()->getGlobalIndex();
396 b = bend->getAtomB()->getGlobalIndex();
397 c = bend->getAtomC()->getGlobalIndex();
398 std::set<int> rigidSetA = getRigidSet(a, atomGroups);
399 std::set<int> rigidSetB = getRigidSet(b, atomGroups);
400 std::set<int> rigidSetC = getRigidSet(c, atomGroups);
401
402 exclude_.addPairs(rigidSetA, rigidSetB);
403 exclude_.addPairs(rigidSetA, rigidSetC);
404 exclude_.addPairs(rigidSetB, rigidSetC);
405
406 //exclude_.addPair(a, b);
407 //exclude_.addPair(a, c);
408 //exclude_.addPair(b, c);
409 }
410
411 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; torsion = mol->nextTorsion(torsionIter)) {
412 a = torsion->getAtomA()->getGlobalIndex();
413 b = torsion->getAtomB()->getGlobalIndex();
414 c = torsion->getAtomC()->getGlobalIndex();
415 d = torsion->getAtomD()->getGlobalIndex();
416 std::set<int> rigidSetA = getRigidSet(a, atomGroups);
417 std::set<int> rigidSetB = getRigidSet(b, atomGroups);
418 std::set<int> rigidSetC = getRigidSet(c, atomGroups);
419 std::set<int> rigidSetD = getRigidSet(d, atomGroups);
420
421 exclude_.addPairs(rigidSetA, rigidSetB);
422 exclude_.addPairs(rigidSetA, rigidSetC);
423 exclude_.addPairs(rigidSetA, rigidSetD);
424 exclude_.addPairs(rigidSetB, rigidSetC);
425 exclude_.addPairs(rigidSetB, rigidSetD);
426 exclude_.addPairs(rigidSetC, rigidSetD);
427
428 /*
429 exclude_.addPairs(rigidSetA.begin(), rigidSetA.end(), rigidSetB.begin(), rigidSetB.end());
430 exclude_.addPairs(rigidSetA.begin(), rigidSetA.end(), rigidSetC.begin(), rigidSetC.end());
431 exclude_.addPairs(rigidSetA.begin(), rigidSetA.end(), rigidSetD.begin(), rigidSetD.end());
432 exclude_.addPairs(rigidSetB.begin(), rigidSetB.end(), rigidSetC.begin(), rigidSetC.end());
433 exclude_.addPairs(rigidSetB.begin(), rigidSetB.end(), rigidSetD.begin(), rigidSetD.end());
434 exclude_.addPairs(rigidSetC.begin(), rigidSetC.end(), rigidSetD.begin(), rigidSetD.end());
435
436
437 exclude_.addPair(a, b);
438 exclude_.addPair(a, c);
439 exclude_.addPair(a, d);
440 exclude_.addPair(b, c);
441 exclude_.addPair(b, d);
442 exclude_.addPair(c, d);
443 */
444 }
445
446 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
447 std::vector<Atom*> atoms = rb->getAtoms();
448 for (int i = 0; i < atoms.size() -1 ; ++i) {
449 for (int j = i + 1; j < atoms.size(); ++j) {
450 a = atoms[i]->getGlobalIndex();
451 b = atoms[j]->getGlobalIndex();
452 exclude_.addPair(a, b);
453 }
454 }
455 }
456
457 }
458
459 void SimInfo::removeExcludePairs(Molecule* mol) {
460 std::vector<Bond*>::iterator bondIter;
461 std::vector<Bend*>::iterator bendIter;
462 std::vector<Torsion*>::iterator torsionIter;
463 Bond* bond;
464 Bend* bend;
465 Torsion* torsion;
466 int a;
467 int b;
468 int c;
469 int d;
470
471 std::map<int, std::set<int> > atomGroups;
472
473 Molecule::RigidBodyIterator rbIter;
474 RigidBody* rb;
475 Molecule::IntegrableObjectIterator ii;
476 StuntDouble* integrableObject;
477
478 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
479 integrableObject = mol->nextIntegrableObject(ii)) {
480
481 if (integrableObject->isRigidBody()) {
482 rb = static_cast<RigidBody*>(integrableObject);
483 std::vector<Atom*> atoms = rb->getAtoms();
484 std::set<int> rigidAtoms;
485 for (int i = 0; i < atoms.size(); ++i) {
486 rigidAtoms.insert(atoms[i]->getGlobalIndex());
487 }
488 for (int i = 0; i < atoms.size(); ++i) {
489 atomGroups.insert(std::map<int, std::set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
490 }
491 } else {
492 std::set<int> oneAtomSet;
493 oneAtomSet.insert(integrableObject->getGlobalIndex());
494 atomGroups.insert(std::map<int, std::set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));
495 }
496 }
497
498
499 for (bond= mol->beginBond(bondIter); bond != NULL; bond = mol->nextBond(bondIter)) {
500 a = bond->getAtomA()->getGlobalIndex();
501 b = bond->getAtomB()->getGlobalIndex();
502 exclude_.removePair(a, b);
503 }
504
505 for (bend= mol->beginBend(bendIter); bend != NULL; bend = mol->nextBend(bendIter)) {
506 a = bend->getAtomA()->getGlobalIndex();
507 b = bend->getAtomB()->getGlobalIndex();
508 c = bend->getAtomC()->getGlobalIndex();
509
510 std::set<int> rigidSetA = getRigidSet(a, atomGroups);
511 std::set<int> rigidSetB = getRigidSet(b, atomGroups);
512 std::set<int> rigidSetC = getRigidSet(c, atomGroups);
513
514 exclude_.removePairs(rigidSetA, rigidSetB);
515 exclude_.removePairs(rigidSetA, rigidSetC);
516 exclude_.removePairs(rigidSetB, rigidSetC);
517
518 //exclude_.removePair(a, b);
519 //exclude_.removePair(a, c);
520 //exclude_.removePair(b, c);
521 }
522
523 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; torsion = mol->nextTorsion(torsionIter)) {
524 a = torsion->getAtomA()->getGlobalIndex();
525 b = torsion->getAtomB()->getGlobalIndex();
526 c = torsion->getAtomC()->getGlobalIndex();
527 d = torsion->getAtomD()->getGlobalIndex();
528
529 std::set<int> rigidSetA = getRigidSet(a, atomGroups);
530 std::set<int> rigidSetB = getRigidSet(b, atomGroups);
531 std::set<int> rigidSetC = getRigidSet(c, atomGroups);
532 std::set<int> rigidSetD = getRigidSet(d, atomGroups);
533
534 exclude_.removePairs(rigidSetA, rigidSetB);
535 exclude_.removePairs(rigidSetA, rigidSetC);
536 exclude_.removePairs(rigidSetA, rigidSetD);
537 exclude_.removePairs(rigidSetB, rigidSetC);
538 exclude_.removePairs(rigidSetB, rigidSetD);
539 exclude_.removePairs(rigidSetC, rigidSetD);
540
541 /*
542 exclude_.removePairs(rigidSetA.begin(), rigidSetA.end(), rigidSetB.begin(), rigidSetB.end());
543 exclude_.removePairs(rigidSetA.begin(), rigidSetA.end(), rigidSetC.begin(), rigidSetC.end());
544 exclude_.removePairs(rigidSetA.begin(), rigidSetA.end(), rigidSetD.begin(), rigidSetD.end());
545 exclude_.removePairs(rigidSetB.begin(), rigidSetB.end(), rigidSetC.begin(), rigidSetC.end());
546 exclude_.removePairs(rigidSetB.begin(), rigidSetB.end(), rigidSetD.begin(), rigidSetD.end());
547 exclude_.removePairs(rigidSetC.begin(), rigidSetC.end(), rigidSetD.begin(), rigidSetD.end());
548
549
550 exclude_.removePair(a, b);
551 exclude_.removePair(a, c);
552 exclude_.removePair(a, d);
553 exclude_.removePair(b, c);
554 exclude_.removePair(b, d);
555 exclude_.removePair(c, d);
556 */
557 }
558
559 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
560 std::vector<Atom*> atoms = rb->getAtoms();
561 for (int i = 0; i < atoms.size() -1 ; ++i) {
562 for (int j = i + 1; j < atoms.size(); ++j) {
563 a = atoms[i]->getGlobalIndex();
564 b = atoms[j]->getGlobalIndex();
565 exclude_.removePair(a, b);
566 }
567 }
568 }
569
570 }
571
572
573 void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) {
574 int curStampId;
575
576 //index from 0
577 curStampId = moleculeStamps_.size();
578
579 moleculeStamps_.push_back(molStamp);
580 molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
581 }
582
583 void SimInfo::update() {
584
585 setupSimType();
586
587 #ifdef IS_MPI
588 setupFortranParallel();
589 #endif
590
591 setupFortranSim();
592
593 //setup fortran force field
594 /** @deprecate */
595 int isError = 0;
596
597 setupElectrostaticSummationMethod( isError );
598 setupSwitchingFunction();
599
600 if(isError){
601 sprintf( painCave.errMsg,
602 "ForceField error: There was an error initializing the forceField in fortran.\n" );
603 painCave.isFatal = 1;
604 simError();
605 }
606
607
608 setupCutoff();
609
610 calcNdf();
611 calcNdfRaw();
612 calcNdfTrans();
613
614 fortranInitialized_ = true;
615 }
616
617 std::set<AtomType*> SimInfo::getUniqueAtomTypes() {
618 SimInfo::MoleculeIterator mi;
619 Molecule* mol;
620 Molecule::AtomIterator ai;
621 Atom* atom;
622 std::set<AtomType*> atomTypes;
623
624 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
625
626 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
627 atomTypes.insert(atom->getAtomType());
628 }
629
630 }
631
632 return atomTypes;
633 }
634
635 void SimInfo::setupSimType() {
636 std::set<AtomType*>::iterator i;
637 std::set<AtomType*> atomTypes;
638 atomTypes = getUniqueAtomTypes();
639
640 int useLennardJones = 0;
641 int useElectrostatic = 0;
642 int useEAM = 0;
643 int useSC = 0;
644 int useCharge = 0;
645 int useDirectional = 0;
646 int useDipole = 0;
647 int useGayBerne = 0;
648 int useSticky = 0;
649 int useStickyPower = 0;
650 int useShape = 0;
651 int useFLARB = 0; //it is not in AtomType yet
652 int useDirectionalAtom = 0;
653 int useElectrostatics = 0;
654 //usePBC and useRF are from simParams
655 int usePBC = simParams_->getUsePeriodicBoundaryConditions();
656 int useRF;
657 int useSF;
658 std::string myMethod;
659
660 // set the useRF logical
661 useRF = 0;
662 useSF = 0;
663
664
665 if (simParams_->haveElectrostaticSummationMethod()) {
666 std::string myMethod = simParams_->getElectrostaticSummationMethod();
667 toUpper(myMethod);
668 if (myMethod == "REACTION_FIELD") {
669 useRF=1;
670 } else {
671 if (myMethod == "SHIFTED_FORCE") {
672 useSF = 1;
673 }
674 }
675 }
676
677 //loop over all of the atom types
678 for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
679 useLennardJones |= (*i)->isLennardJones();
680 useElectrostatic |= (*i)->isElectrostatic();
681 useEAM |= (*i)->isEAM();
682 useSC |= (*i)->isSC();
683 useCharge |= (*i)->isCharge();
684 useDirectional |= (*i)->isDirectional();
685 useDipole |= (*i)->isDipole();
686 useGayBerne |= (*i)->isGayBerne();
687 useSticky |= (*i)->isSticky();
688 useStickyPower |= (*i)->isStickyPower();
689 useShape |= (*i)->isShape();
690 }
691
692 if (useSticky || useStickyPower || useDipole || useGayBerne || useShape) {
693 useDirectionalAtom = 1;
694 }
695
696 if (useCharge || useDipole) {
697 useElectrostatics = 1;
698 }
699
700 #ifdef IS_MPI
701 int temp;
702
703 temp = usePBC;
704 MPI_Allreduce(&temp, &usePBC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
705
706 temp = useDirectionalAtom;
707 MPI_Allreduce(&temp, &useDirectionalAtom, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
708
709 temp = useLennardJones;
710 MPI_Allreduce(&temp, &useLennardJones, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
711
712 temp = useElectrostatics;
713 MPI_Allreduce(&temp, &useElectrostatics, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
714
715 temp = useCharge;
716 MPI_Allreduce(&temp, &useCharge, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
717
718 temp = useDipole;
719 MPI_Allreduce(&temp, &useDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
720
721 temp = useSticky;
722 MPI_Allreduce(&temp, &useSticky, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
723
724 temp = useStickyPower;
725 MPI_Allreduce(&temp, &useStickyPower, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
726
727 temp = useGayBerne;
728 MPI_Allreduce(&temp, &useGayBerne, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
729
730 temp = useEAM;
731 MPI_Allreduce(&temp, &useEAM, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
732
733 temp = useSC;
734 MPI_Allreduce(&temp, &useSC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
735
736 temp = useShape;
737 MPI_Allreduce(&temp, &useShape, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
738
739 temp = useFLARB;
740 MPI_Allreduce(&temp, &useFLARB, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
741
742 temp = useRF;
743 MPI_Allreduce(&temp, &useRF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
744
745 temp = useSF;
746 MPI_Allreduce(&temp, &useSF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
747
748 #endif
749
750 fInfo_.SIM_uses_PBC = usePBC;
751 fInfo_.SIM_uses_DirectionalAtoms = useDirectionalAtom;
752 fInfo_.SIM_uses_LennardJones = useLennardJones;
753 fInfo_.SIM_uses_Electrostatics = useElectrostatics;
754 fInfo_.SIM_uses_Charges = useCharge;
755 fInfo_.SIM_uses_Dipoles = useDipole;
756 fInfo_.SIM_uses_Sticky = useSticky;
757 fInfo_.SIM_uses_StickyPower = useStickyPower;
758 fInfo_.SIM_uses_GayBerne = useGayBerne;
759 fInfo_.SIM_uses_EAM = useEAM;
760 fInfo_.SIM_uses_SC = useSC;
761 fInfo_.SIM_uses_Shapes = useShape;
762 fInfo_.SIM_uses_FLARB = useFLARB;
763 fInfo_.SIM_uses_RF = useRF;
764 fInfo_.SIM_uses_SF = useSF;
765
766 if( myMethod == "REACTION_FIELD") {
767
768 if (simParams_->haveDielectric()) {
769 fInfo_.dielect = simParams_->getDielectric();
770 } else {
771 sprintf(painCave.errMsg,
772 "SimSetup Error: No Dielectric constant was set.\n"
773 "\tYou are trying to use Reaction Field without"
774 "\tsetting a dielectric constant!\n");
775 painCave.isFatal = 1;
776 simError();
777 }
778 }
779
780 }
781
782 void SimInfo::setupFortranSim() {
783 int isError;
784 int nExclude;
785 std::vector<int> fortranGlobalGroupMembership;
786
787 nExclude = exclude_.getSize();
788 isError = 0;
789
790 //globalGroupMembership_ is filled by SimCreator
791 for (int i = 0; i < nGlobalAtoms_; i++) {
792 fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
793 }
794
795 //calculate mass ratio of cutoff group
796 std::vector<double> mfact;
797 SimInfo::MoleculeIterator mi;
798 Molecule* mol;
799 Molecule::CutoffGroupIterator ci;
800 CutoffGroup* cg;
801 Molecule::AtomIterator ai;
802 Atom* atom;
803 double totalMass;
804
805 //to avoid memory reallocation, reserve enough space for mfact
806 mfact.reserve(getNCutoffGroups());
807
808 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
809 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
810
811 totalMass = cg->getMass();
812 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
813 // Check for massless groups - set mfact to 1 if true
814 if (totalMass != 0)
815 mfact.push_back(atom->getMass()/totalMass);
816 else
817 mfact.push_back( 1.0 );
818 }
819
820 }
821 }
822
823 //fill ident array of local atoms (it is actually ident of AtomType, it is so confusing !!!)
824 std::vector<int> identArray;
825
826 //to avoid memory reallocation, reserve enough space identArray
827 identArray.reserve(getNAtoms());
828
829 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
830 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
831 identArray.push_back(atom->getIdent());
832 }
833 }
834
835 //fill molMembershipArray
836 //molMembershipArray is filled by SimCreator
837 std::vector<int> molMembershipArray(nGlobalAtoms_);
838 for (int i = 0; i < nGlobalAtoms_; i++) {
839 molMembershipArray[i] = globalMolMembership_[i] + 1;
840 }
841
842 //setup fortran simulation
843 int nGlobalExcludes = 0;
844 int* globalExcludes = NULL;
845 int* excludeList = exclude_.getExcludeList();
846 setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray[0], &nExclude, excludeList ,
847 &nGlobalExcludes, globalExcludes, &molMembershipArray[0],
848 &mfact[0], &nCutoffGroups_, &fortranGlobalGroupMembership[0], &isError);
849
850 if( isError ){
851
852 sprintf( painCave.errMsg,
853 "There was an error setting the simulation information in fortran.\n" );
854 painCave.isFatal = 1;
855 painCave.severity = OOPSE_ERROR;
856 simError();
857 }
858
859 #ifdef IS_MPI
860 sprintf( checkPointMsg,
861 "succesfully sent the simulation information to fortran.\n");
862 MPIcheckPoint();
863 #endif // is_mpi
864 }
865
866
867 #ifdef IS_MPI
868 void SimInfo::setupFortranParallel() {
869
870 //SimInfo is responsible for creating localToGlobalAtomIndex and localToGlobalGroupIndex
871 std::vector<int> localToGlobalAtomIndex(getNAtoms(), 0);
872 std::vector<int> localToGlobalCutoffGroupIndex;
873 SimInfo::MoleculeIterator mi;
874 Molecule::AtomIterator ai;
875 Molecule::CutoffGroupIterator ci;
876 Molecule* mol;
877 Atom* atom;
878 CutoffGroup* cg;
879 mpiSimData parallelData;
880 int isError;
881
882 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
883
884 //local index(index in DataStorge) of atom is important
885 for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
886 localToGlobalAtomIndex[atom->getLocalIndex()] = atom->getGlobalIndex() + 1;
887 }
888
889 //local index of cutoff group is trivial, it only depends on the order of travesing
890 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
891 localToGlobalCutoffGroupIndex.push_back(cg->getGlobalIndex() + 1);
892 }
893
894 }
895
896 //fill up mpiSimData struct
897 parallelData.nMolGlobal = getNGlobalMolecules();
898 parallelData.nMolLocal = getNMolecules();
899 parallelData.nAtomsGlobal = getNGlobalAtoms();
900 parallelData.nAtomsLocal = getNAtoms();
901 parallelData.nGroupsGlobal = getNGlobalCutoffGroups();
902 parallelData.nGroupsLocal = getNCutoffGroups();
903 parallelData.myNode = worldRank;
904 MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors));
905
906 //pass mpiSimData struct and index arrays to fortran
907 setFsimParallel(&parallelData, &(parallelData.nAtomsLocal),
908 &localToGlobalAtomIndex[0], &(parallelData.nGroupsLocal),
909 &localToGlobalCutoffGroupIndex[0], &isError);
910
911 if (isError) {
912 sprintf(painCave.errMsg,
913 "mpiRefresh errror: fortran didn't like something we gave it.\n");
914 painCave.isFatal = 1;
915 simError();
916 }
917
918 sprintf(checkPointMsg, " mpiRefresh successful.\n");
919 MPIcheckPoint();
920
921
922 }
923
924 #endif
925
926 double SimInfo::calcMaxCutoffRadius() {
927
928
929 std::set<AtomType*> atomTypes;
930 std::set<AtomType*>::iterator i;
931 std::vector<double> cutoffRadius;
932
933 //get the unique atom types
934 atomTypes = getUniqueAtomTypes();
935
936 //query the max cutoff radius among these atom types
937 for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
938 cutoffRadius.push_back(forceField_->getRcutFromAtomType(*i));
939 }
940
941 double maxCutoffRadius = *(std::max_element(cutoffRadius.begin(), cutoffRadius.end()));
942 #ifdef IS_MPI
943 //pick the max cutoff radius among the processors
944 #endif
945
946 return maxCutoffRadius;
947 }
948
949 void SimInfo::getCutoff(double& rcut, double& rsw) {
950
951 if (fInfo_.SIM_uses_Charges | fInfo_.SIM_uses_Dipoles | fInfo_.SIM_uses_RF) {
952
953 if (!simParams_->haveCutoffRadius()){
954 sprintf(painCave.errMsg,
955 "SimCreator Warning: No value was set for the cutoffRadius.\n"
956 "\tOOPSE will use a default value of 15.0 angstroms"
957 "\tfor the cutoffRadius.\n");
958 painCave.isFatal = 0;
959 simError();
960 rcut = 15.0;
961 } else{
962 rcut = simParams_->getCutoffRadius();
963 }
964
965 if (!simParams_->haveSwitchingRadius()){
966 sprintf(painCave.errMsg,
967 "SimCreator Warning: No value was set for switchingRadius.\n"
968 "\tOOPSE will use a default value of\n"
969 "\t0.85 * cutoffRadius for the switchingRadius\n");
970 painCave.isFatal = 0;
971 simError();
972 rsw = 0.85 * rcut;
973 } else{
974 rsw = simParams_->getSwitchingRadius();
975 }
976
977 } else {
978 // if charge, dipole or reaction field is not used and the cutofff radius is not specified in
979 //meta-data file, the maximum cutoff radius calculated from forcefiled will be used
980
981 if (simParams_->haveCutoffRadius()) {
982 rcut = simParams_->getCutoffRadius();
983 } else {
984 //set cutoff radius to the maximum cutoff radius based on atom types in the whole system
985 rcut = calcMaxCutoffRadius();
986 }
987
988 if (simParams_->haveSwitchingRadius()) {
989 rsw = simParams_->getSwitchingRadius();
990 } else {
991 rsw = rcut;
992 }
993
994 }
995 }
996
997 void SimInfo::setupCutoff() {
998 getCutoff(rcut_, rsw_);
999 double rnblist = rcut_ + 1; // skin of neighbor list
1000
1001 //Pass these cutoff radius etc. to fortran. This function should be called once and only once
1002
1003 int cp = TRADITIONAL_CUTOFF_POLICY;
1004 if (simParams_->haveCutoffPolicy()) {
1005 std::string myPolicy = simParams_->getCutoffPolicy();
1006 toUpper(myPolicy);
1007 if (myPolicy == "MIX") {
1008 cp = MIX_CUTOFF_POLICY;
1009 } else {
1010 if (myPolicy == "MAX") {
1011 cp = MAX_CUTOFF_POLICY;
1012 } else {
1013 if (myPolicy == "TRADITIONAL") {
1014 cp = TRADITIONAL_CUTOFF_POLICY;
1015 } else {
1016 // throw error
1017 sprintf( painCave.errMsg,
1018 "SimInfo error: Unknown cutoffPolicy. (Input file specified %s .)\n\tcutoffPolicy must be one of: \"Mix\", \"Max\", or \"Traditional\".", myPolicy.c_str() );
1019 painCave.isFatal = 1;
1020 simError();
1021 }
1022 }
1023 }
1024 }
1025
1026
1027 if (simParams_->haveSkinThickness()) {
1028 double skinThickness = simParams_->getSkinThickness();
1029 }
1030
1031 notifyFortranCutoffs(&rcut_, &rsw_, &rnblist, &cp);
1032 // also send cutoff notification to electrostatics
1033 setElectrostaticCutoffRadius(&rcut_, &rsw_);
1034 }
1035
1036 void SimInfo::setupElectrostaticSummationMethod( int isError ) {
1037
1038 int errorOut;
1039 int esm = NONE;
1040 int sm = UNDAMPED;
1041 double alphaVal;
1042 double dielectric;
1043
1044 errorOut = isError;
1045 alphaVal = simParams_->getDampingAlpha();
1046 dielectric = simParams_->getDielectric();
1047
1048 if (simParams_->haveElectrostaticSummationMethod()) {
1049 std::string myMethod = simParams_->getElectrostaticSummationMethod();
1050 toUpper(myMethod);
1051 if (myMethod == "NONE") {
1052 esm = NONE;
1053 } else {
1054 if (myMethod == "SWITCHING_FUNCTION") {
1055 esm = SWITCHING_FUNCTION;
1056 } else {
1057 if (myMethod == "SHIFTED_POTENTIAL") {
1058 esm = SHIFTED_POTENTIAL;
1059 } else {
1060 if (myMethod == "SHIFTED_FORCE") {
1061 esm = SHIFTED_FORCE;
1062 } else {
1063 if (myMethod == "REACTION_FIELD") {
1064 esm = REACTION_FIELD;
1065 } else {
1066 // throw error
1067 sprintf( painCave.errMsg,
1068 "SimInfo error: Unknown electrostaticSummationMethod. (Input file specified %s .)\n\telectrostaticSummationMethod must be one of: \"none\", \"shifted_potential\", \"shifted_force\", or \"reaction_field\".", myMethod.c_str() );
1069 painCave.isFatal = 1;
1070 simError();
1071 }
1072 }
1073 }
1074 }
1075 }
1076 }
1077
1078 if (simParams_->haveElectrostaticScreeningMethod()) {
1079 std::string myScreen = simParams_->getElectrostaticScreeningMethod();
1080 toUpper(myScreen);
1081 if (myScreen == "UNDAMPED") {
1082 sm = UNDAMPED;
1083 } else {
1084 if (myScreen == "DAMPED") {
1085 sm = DAMPED;
1086 if (!simParams_->haveDampingAlpha()) {
1087 //throw error
1088 sprintf( painCave.errMsg,
1089 "SimInfo warning: dampingAlpha was not specified in the input file. A default value of %f (1/ang) will be used.", alphaVal);
1090 painCave.isFatal = 0;
1091 simError();
1092 }
1093 } else {
1094 // throw error
1095 sprintf( painCave.errMsg,
1096 "SimInfo error: Unknown electrostaticScreeningMethod. (Input file specified %s .)\n\telectrostaticScreeningMethod must be one of: \"undamped\" or \"damped\".", myScreen.c_str() );
1097 painCave.isFatal = 1;
1098 simError();
1099 }
1100 }
1101 }
1102
1103 // let's pass some summation method variables to fortran
1104 setElectrostaticSummationMethod( &esm );
1105 setScreeningMethod( &sm );
1106 setDampingAlpha( &alphaVal );
1107 setReactionFieldDielectric( &dielectric );
1108 initFortranFF( &esm, &errorOut );
1109 }
1110
1111 void SimInfo::setupSwitchingFunction() {
1112 int ft = CUBIC;
1113
1114 if (simParams_->haveSwitchingFunctionType()) {
1115 std::string funcType = simParams_->getSwitchingFunctionType();
1116 toUpper(funcType);
1117 if (funcType == "CUBIC") {
1118 ft = CUBIC;
1119 } else {
1120 if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
1121 ft = FIFTH_ORDER_POLY;
1122 } else {
1123 // throw error
1124 sprintf( painCave.errMsg,
1125 "SimInfo error: Unknown switchingFunctionType. (Input file specified %s .)\n\tswitchingFunctionType must be one of: \"cubic\" or \"fifth_order_polynomial\".", funcType.c_str() );
1126 painCave.isFatal = 1;
1127 simError();
1128 }
1129 }
1130 }
1131
1132 // send switching function notification to switcheroo
1133 setFunctionType(&ft);
1134
1135 }
1136
1137 void SimInfo::addProperty(GenericData* genData) {
1138 properties_.addProperty(genData);
1139 }
1140
1141 void SimInfo::removeProperty(const std::string& propName) {
1142 properties_.removeProperty(propName);
1143 }
1144
1145 void SimInfo::clearProperties() {
1146 properties_.clearProperties();
1147 }
1148
1149 std::vector<std::string> SimInfo::getPropertyNames() {
1150 return properties_.getPropertyNames();
1151 }
1152
1153 std::vector<GenericData*> SimInfo::getProperties() {
1154 return properties_.getProperties();
1155 }
1156
1157 GenericData* SimInfo::getPropertyByName(const std::string& propName) {
1158 return properties_.getPropertyByName(propName);
1159 }
1160
1161 void SimInfo::setSnapshotManager(SnapshotManager* sman) {
1162 if (sman_ == sman) {
1163 return;
1164 }
1165 delete sman_;
1166 sman_ = sman;
1167
1168 Molecule* mol;
1169 RigidBody* rb;
1170 Atom* atom;
1171 SimInfo::MoleculeIterator mi;
1172 Molecule::RigidBodyIterator rbIter;
1173 Molecule::AtomIterator atomIter;;
1174
1175 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
1176
1177 for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) {
1178 atom->setSnapshotManager(sman_);
1179 }
1180
1181 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
1182 rb->setSnapshotManager(sman_);
1183 }
1184 }
1185
1186 }
1187
1188 Vector3d SimInfo::getComVel(){
1189 SimInfo::MoleculeIterator i;
1190 Molecule* mol;
1191
1192 Vector3d comVel(0.0);
1193 double totalMass = 0.0;
1194
1195
1196 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1197 double mass = mol->getMass();
1198 totalMass += mass;
1199 comVel += mass * mol->getComVel();
1200 }
1201
1202 #ifdef IS_MPI
1203 double tmpMass = totalMass;
1204 Vector3d tmpComVel(comVel);
1205 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1206 MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1207 #endif
1208
1209 comVel /= totalMass;
1210
1211 return comVel;
1212 }
1213
1214 Vector3d SimInfo::getCom(){
1215 SimInfo::MoleculeIterator i;
1216 Molecule* mol;
1217
1218 Vector3d com(0.0);
1219 double totalMass = 0.0;
1220
1221 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1222 double mass = mol->getMass();
1223 totalMass += mass;
1224 com += mass * mol->getCom();
1225 }
1226
1227 #ifdef IS_MPI
1228 double tmpMass = totalMass;
1229 Vector3d tmpCom(com);
1230 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1231 MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1232 #endif
1233
1234 com /= totalMass;
1235
1236 return com;
1237
1238 }
1239
1240 std::ostream& operator <<(std::ostream& o, SimInfo& info) {
1241
1242 return o;
1243 }
1244
1245
1246 /*
1247 Returns center of mass and center of mass velocity in one function call.
1248 */
1249
1250 void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){
1251 SimInfo::MoleculeIterator i;
1252 Molecule* mol;
1253
1254
1255 double totalMass = 0.0;
1256
1257
1258 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1259 double mass = mol->getMass();
1260 totalMass += mass;
1261 com += mass * mol->getCom();
1262 comVel += mass * mol->getComVel();
1263 }
1264
1265 #ifdef IS_MPI
1266 double tmpMass = totalMass;
1267 Vector3d tmpCom(com);
1268 Vector3d tmpComVel(comVel);
1269 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1270 MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1271 MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1272 #endif
1273
1274 com /= totalMass;
1275 comVel /= totalMass;
1276 }
1277
1278 /*
1279 Return intertia tensor for entire system and angular momentum Vector.
1280
1281
1282 [ Ixx -Ixy -Ixz ]
1283 J =| -Iyx Iyy -Iyz |
1284 [ -Izx -Iyz Izz ]
1285 */
1286
1287 void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){
1288
1289
1290 double xx = 0.0;
1291 double yy = 0.0;
1292 double zz = 0.0;
1293 double xy = 0.0;
1294 double xz = 0.0;
1295 double yz = 0.0;
1296 Vector3d com(0.0);
1297 Vector3d comVel(0.0);
1298
1299 getComAll(com, comVel);
1300
1301 SimInfo::MoleculeIterator i;
1302 Molecule* mol;
1303
1304 Vector3d thisq(0.0);
1305 Vector3d thisv(0.0);
1306
1307 double thisMass = 0.0;
1308
1309
1310
1311
1312 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1313
1314 thisq = mol->getCom()-com;
1315 thisv = mol->getComVel()-comVel;
1316 thisMass = mol->getMass();
1317 // Compute moment of intertia coefficients.
1318 xx += thisq[0]*thisq[0]*thisMass;
1319 yy += thisq[1]*thisq[1]*thisMass;
1320 zz += thisq[2]*thisq[2]*thisMass;
1321
1322 // compute products of intertia
1323 xy += thisq[0]*thisq[1]*thisMass;
1324 xz += thisq[0]*thisq[2]*thisMass;
1325 yz += thisq[1]*thisq[2]*thisMass;
1326
1327 angularMomentum += cross( thisq, thisv ) * thisMass;
1328
1329 }
1330
1331
1332 inertiaTensor(0,0) = yy + zz;
1333 inertiaTensor(0,1) = -xy;
1334 inertiaTensor(0,2) = -xz;
1335 inertiaTensor(1,0) = -xy;
1336 inertiaTensor(1,1) = xx + zz;
1337 inertiaTensor(1,2) = -yz;
1338 inertiaTensor(2,0) = -xz;
1339 inertiaTensor(2,1) = -yz;
1340 inertiaTensor(2,2) = xx + yy;
1341
1342 #ifdef IS_MPI
1343 Mat3x3d tmpI(inertiaTensor);
1344 Vector3d tmpAngMom;
1345 MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1346 MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1347 #endif
1348
1349 return;
1350 }
1351
1352 //Returns the angular momentum of the system
1353 Vector3d SimInfo::getAngularMomentum(){
1354
1355 Vector3d com(0.0);
1356 Vector3d comVel(0.0);
1357 Vector3d angularMomentum(0.0);
1358
1359 getComAll(com,comVel);
1360
1361 SimInfo::MoleculeIterator i;
1362 Molecule* mol;
1363
1364 Vector3d thisr(0.0);
1365 Vector3d thisp(0.0);
1366
1367 double thisMass;
1368
1369 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1370 thisMass = mol->getMass();
1371 thisr = mol->getCom()-com;
1372 thisp = (mol->getComVel()-comVel)*thisMass;
1373
1374 angularMomentum += cross( thisr, thisp );
1375
1376 }
1377
1378 #ifdef IS_MPI
1379 Vector3d tmpAngMom;
1380 MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1381 #endif
1382
1383 return angularMomentum;
1384 }
1385
1386
1387 }//end namespace oopse
1388