ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/brains/SimInfo.cpp
Revision: 2310
Committed: Mon Sep 19 23:21:46 2005 UTC (18 years, 10 months ago) by chrisfen
File size: 35255 byte(s)
Log Message:
Fixed bugs in reaction field, now it appears as though it really is working...

File Contents

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