ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/brains/SimInfo.cpp
Revision: 2309
Committed: Sun Sep 18 20:45:38 2005 UTC (18 years, 9 months ago) by chrisfen
File size: 34919 byte(s)
Log Message:
reaction field seems to work now, still need to do some testing...

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 gezelter 1490
527 gezelter 1930 //loop over all of the atom types
528     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
529 gezelter 2204 useLennardJones |= (*i)->isLennardJones();
530     useElectrostatic |= (*i)->isElectrostatic();
531     useEAM |= (*i)->isEAM();
532     useCharge |= (*i)->isCharge();
533     useDirectional |= (*i)->isDirectional();
534     useDipole |= (*i)->isDipole();
535     useGayBerne |= (*i)->isGayBerne();
536     useSticky |= (*i)->isSticky();
537 chrisfen 2220 useStickyPower |= (*i)->isStickyPower();
538 gezelter 2204 useShape |= (*i)->isShape();
539 gezelter 1930 }
540 gezelter 1490
541 chrisfen 2220 if (useSticky || useStickyPower || useDipole || useGayBerne || useShape) {
542 gezelter 2204 useDirectionalAtom = 1;
543 gezelter 1930 }
544 gezelter 1490
545 gezelter 1930 if (useCharge || useDipole) {
546 gezelter 2204 useElectrostatics = 1;
547 gezelter 1930 }
548 gezelter 1490
549 gezelter 1930 #ifdef IS_MPI
550     int temp;
551 gezelter 1490
552 gezelter 1930 temp = usePBC;
553     MPI_Allreduce(&temp, &usePBC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
554 gezelter 1490
555 gezelter 1930 temp = useDirectionalAtom;
556     MPI_Allreduce(&temp, &useDirectionalAtom, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
557 gezelter 1490
558 gezelter 1930 temp = useLennardJones;
559     MPI_Allreduce(&temp, &useLennardJones, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
560 gezelter 1490
561 gezelter 1930 temp = useElectrostatics;
562     MPI_Allreduce(&temp, &useElectrostatics, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
563 gezelter 1490
564 gezelter 1930 temp = useCharge;
565     MPI_Allreduce(&temp, &useCharge, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
566 gezelter 1490
567 gezelter 1930 temp = useDipole;
568     MPI_Allreduce(&temp, &useDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
569 gezelter 1490
570 gezelter 1930 temp = useSticky;
571     MPI_Allreduce(&temp, &useSticky, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
572 gezelter 1490
573 chrisfen 2220 temp = useStickyPower;
574     MPI_Allreduce(&temp, &useStickyPower, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
575    
576 gezelter 1930 temp = useGayBerne;
577     MPI_Allreduce(&temp, &useGayBerne, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
578 gezelter 1490
579 gezelter 1930 temp = useEAM;
580     MPI_Allreduce(&temp, &useEAM, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
581 gezelter 1490
582 gezelter 1930 temp = useShape;
583     MPI_Allreduce(&temp, &useShape, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
584    
585     temp = useFLARB;
586     MPI_Allreduce(&temp, &useFLARB, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
587    
588 gezelter 1490 #endif
589    
590 gezelter 1930 fInfo_.SIM_uses_PBC = usePBC;
591     fInfo_.SIM_uses_DirectionalAtoms = useDirectionalAtom;
592     fInfo_.SIM_uses_LennardJones = useLennardJones;
593     fInfo_.SIM_uses_Electrostatics = useElectrostatics;
594     fInfo_.SIM_uses_Charges = useCharge;
595     fInfo_.SIM_uses_Dipoles = useDipole;
596     fInfo_.SIM_uses_Sticky = useSticky;
597 chrisfen 2220 fInfo_.SIM_uses_StickyPower = useStickyPower;
598 gezelter 1930 fInfo_.SIM_uses_GayBerne = useGayBerne;
599     fInfo_.SIM_uses_EAM = useEAM;
600     fInfo_.SIM_uses_Shapes = useShape;
601     fInfo_.SIM_uses_FLARB = useFLARB;
602 gezelter 1490
603 gezelter 1930 if( fInfo_.SIM_uses_Dipoles && fInfo_.SIM_uses_RF) {
604 gezelter 1490
605 gezelter 2204 if (simParams_->haveDielectric()) {
606     fInfo_.dielect = simParams_->getDielectric();
607     } else {
608     sprintf(painCave.errMsg,
609     "SimSetup Error: No Dielectric constant was set.\n"
610     "\tYou are trying to use Reaction Field without"
611     "\tsetting a dielectric constant!\n");
612     painCave.isFatal = 1;
613     simError();
614     }
615 gezelter 1930
616     } else {
617 gezelter 2204 fInfo_.dielect = 0.0;
618 gezelter 1930 }
619    
620 gezelter 2204 }
621 gezelter 1490
622 gezelter 2204 void SimInfo::setupFortranSim() {
623 gezelter 1930 int isError;
624     int nExclude;
625     std::vector<int> fortranGlobalGroupMembership;
626    
627     nExclude = exclude_.getSize();
628     isError = 0;
629 gezelter 1490
630 gezelter 1930 //globalGroupMembership_ is filled by SimCreator
631     for (int i = 0; i < nGlobalAtoms_; i++) {
632 gezelter 2204 fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
633 gezelter 1930 }
634 gezelter 1490
635 gezelter 1930 //calculate mass ratio of cutoff group
636     std::vector<double> mfact;
637     SimInfo::MoleculeIterator mi;
638     Molecule* mol;
639     Molecule::CutoffGroupIterator ci;
640     CutoffGroup* cg;
641     Molecule::AtomIterator ai;
642     Atom* atom;
643     double totalMass;
644    
645     //to avoid memory reallocation, reserve enough space for mfact
646     mfact.reserve(getNCutoffGroups());
647 gezelter 1490
648 gezelter 1930 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
649 gezelter 2204 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
650 gezelter 1490
651 gezelter 2204 totalMass = cg->getMass();
652     for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
653     mfact.push_back(atom->getMass()/totalMass);
654     }
655 gezelter 1490
656 gezelter 2204 }
657 gezelter 1930 }
658 gezelter 1490
659 gezelter 1930 //fill ident array of local atoms (it is actually ident of AtomType, it is so confusing !!!)
660     std::vector<int> identArray;
661 gezelter 1490
662 gezelter 1930 //to avoid memory reallocation, reserve enough space identArray
663     identArray.reserve(getNAtoms());
664    
665     for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
666 gezelter 2204 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
667     identArray.push_back(atom->getIdent());
668     }
669 gezelter 1930 }
670 gezelter 1490
671 gezelter 1930 //fill molMembershipArray
672     //molMembershipArray is filled by SimCreator
673     std::vector<int> molMembershipArray(nGlobalAtoms_);
674     for (int i = 0; i < nGlobalAtoms_; i++) {
675 gezelter 2204 molMembershipArray[i] = globalMolMembership_[i] + 1;
676 gezelter 1930 }
677    
678     //setup fortran simulation
679     int nGlobalExcludes = 0;
680     int* globalExcludes = NULL;
681     int* excludeList = exclude_.getExcludeList();
682     setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray[0], &nExclude, excludeList ,
683 gezelter 2204 &nGlobalExcludes, globalExcludes, &molMembershipArray[0],
684     &mfact[0], &nCutoffGroups_, &fortranGlobalGroupMembership[0], &isError);
685 gezelter 1490
686 gezelter 1930 if( isError ){
687 gezelter 1490
688 gezelter 2204 sprintf( painCave.errMsg,
689     "There was an error setting the simulation information in fortran.\n" );
690     painCave.isFatal = 1;
691     painCave.severity = OOPSE_ERROR;
692     simError();
693 gezelter 1930 }
694    
695     #ifdef IS_MPI
696     sprintf( checkPointMsg,
697 gezelter 2204 "succesfully sent the simulation information to fortran.\n");
698 gezelter 1930 MPIcheckPoint();
699     #endif // is_mpi
700 gezelter 2204 }
701 gezelter 1490
702    
703 gezelter 1930 #ifdef IS_MPI
704 gezelter 2204 void SimInfo::setupFortranParallel() {
705 gezelter 1930
706     //SimInfo is responsible for creating localToGlobalAtomIndex and localToGlobalGroupIndex
707     std::vector<int> localToGlobalAtomIndex(getNAtoms(), 0);
708     std::vector<int> localToGlobalCutoffGroupIndex;
709     SimInfo::MoleculeIterator mi;
710     Molecule::AtomIterator ai;
711     Molecule::CutoffGroupIterator ci;
712     Molecule* mol;
713     Atom* atom;
714     CutoffGroup* cg;
715     mpiSimData parallelData;
716     int isError;
717 gezelter 1490
718 gezelter 1930 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
719 gezelter 1490
720 gezelter 2204 //local index(index in DataStorge) of atom is important
721     for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
722     localToGlobalAtomIndex[atom->getLocalIndex()] = atom->getGlobalIndex() + 1;
723     }
724 gezelter 1490
725 gezelter 2204 //local index of cutoff group is trivial, it only depends on the order of travesing
726     for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
727     localToGlobalCutoffGroupIndex.push_back(cg->getGlobalIndex() + 1);
728     }
729 gezelter 1930
730     }
731 gezelter 1490
732 gezelter 1930 //fill up mpiSimData struct
733     parallelData.nMolGlobal = getNGlobalMolecules();
734     parallelData.nMolLocal = getNMolecules();
735     parallelData.nAtomsGlobal = getNGlobalAtoms();
736     parallelData.nAtomsLocal = getNAtoms();
737     parallelData.nGroupsGlobal = getNGlobalCutoffGroups();
738     parallelData.nGroupsLocal = getNCutoffGroups();
739     parallelData.myNode = worldRank;
740     MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors));
741 gezelter 1490
742 gezelter 1930 //pass mpiSimData struct and index arrays to fortran
743     setFsimParallel(&parallelData, &(parallelData.nAtomsLocal),
744     &localToGlobalAtomIndex[0], &(parallelData.nGroupsLocal),
745     &localToGlobalCutoffGroupIndex[0], &isError);
746 gezelter 1490
747 gezelter 1930 if (isError) {
748 gezelter 2204 sprintf(painCave.errMsg,
749     "mpiRefresh errror: fortran didn't like something we gave it.\n");
750     painCave.isFatal = 1;
751     simError();
752 gezelter 1930 }
753 gezelter 1490
754 gezelter 1930 sprintf(checkPointMsg, " mpiRefresh successful.\n");
755     MPIcheckPoint();
756 gezelter 1490
757    
758 gezelter 2204 }
759 chrisfen 1636
760 gezelter 1930 #endif
761 chrisfen 1636
762 gezelter 2204 double SimInfo::calcMaxCutoffRadius() {
763 chrisfen 1636
764    
765 gezelter 1930 std::set<AtomType*> atomTypes;
766     std::set<AtomType*>::iterator i;
767     std::vector<double> cutoffRadius;
768 gezelter 1490
769 gezelter 1930 //get the unique atom types
770     atomTypes = getUniqueAtomTypes();
771    
772     //query the max cutoff radius among these atom types
773     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
774 gezelter 2204 cutoffRadius.push_back(forceField_->getRcutFromAtomType(*i));
775 gezelter 1930 }
776    
777     double maxCutoffRadius = *(std::max_element(cutoffRadius.begin(), cutoffRadius.end()));
778 gezelter 1490 #ifdef IS_MPI
779 gezelter 1930 //pick the max cutoff radius among the processors
780 gezelter 1490 #endif
781    
782 gezelter 1930 return maxCutoffRadius;
783 gezelter 2204 }
784 gezelter 1930
785 gezelter 2204 void SimInfo::getCutoff(double& rcut, double& rsw) {
786 gezelter 1490
787 gezelter 1930 if (fInfo_.SIM_uses_Charges | fInfo_.SIM_uses_Dipoles | fInfo_.SIM_uses_RF) {
788    
789 gezelter 2204 if (!simParams_->haveRcut()){
790     sprintf(painCave.errMsg,
791 gezelter 1930 "SimCreator Warning: No value was set for the cutoffRadius.\n"
792     "\tOOPSE will use a default value of 15.0 angstroms"
793     "\tfor the cutoffRadius.\n");
794 gezelter 2204 painCave.isFatal = 0;
795     simError();
796     rcut = 15.0;
797     } else{
798     rcut = simParams_->getRcut();
799     }
800 gezelter 1930
801 gezelter 2204 if (!simParams_->haveRsw()){
802     sprintf(painCave.errMsg,
803 gezelter 1930 "SimCreator Warning: No value was set for switchingRadius.\n"
804     "\tOOPSE will use a default value of\n"
805     "\t0.95 * cutoffRadius for the switchingRadius\n");
806 gezelter 2204 painCave.isFatal = 0;
807     simError();
808     rsw = 0.95 * rcut;
809     } else{
810     rsw = simParams_->getRsw();
811     }
812 gezelter 1930
813     } else {
814 gezelter 2204 // if charge, dipole or reaction field is not used and the cutofff radius is not specified in
815     //meta-data file, the maximum cutoff radius calculated from forcefiled will be used
816 gezelter 1930
817 gezelter 2204 if (simParams_->haveRcut()) {
818     rcut = simParams_->getRcut();
819     } else {
820     //set cutoff radius to the maximum cutoff radius based on atom types in the whole system
821     rcut = calcMaxCutoffRadius();
822     }
823 gezelter 1930
824 gezelter 2204 if (simParams_->haveRsw()) {
825     rsw = simParams_->getRsw();
826     } else {
827     rsw = rcut;
828     }
829 gezelter 1930
830     }
831 gezelter 2204 }
832 tim 2010
833 gezelter 2285 void SimInfo::setupCutoff() {
834 tim 2010 getCutoff(rcut_, rsw_);
835 gezelter 1930 double rnblist = rcut_ + 1; // skin of neighbor list
836    
837     //Pass these cutoff radius etc. to fortran. This function should be called once and only once
838 gezelter 2285
839     int cp = TRADITIONAL_CUTOFF_POLICY;
840     if (simParams_->haveCutoffPolicy()) {
841     std::string myPolicy = simParams_->getCutoffPolicy();
842     if (myPolicy == "MIX") {
843     cp = MIX_CUTOFF_POLICY;
844     } else {
845     if (myPolicy == "MAX") {
846     cp = MAX_CUTOFF_POLICY;
847     } else {
848     if (myPolicy == "TRADITIONAL") {
849     cp = TRADITIONAL_CUTOFF_POLICY;
850     } else {
851     // throw error
852     sprintf( painCave.errMsg,
853     "SimInfo error: Unknown cutoffPolicy. (Input file specified %s .)\n\tcutoffPolicy must be one of: \"Mix\", \"Max\", or \"Traditional\".", myPolicy.c_str() );
854     painCave.isFatal = 1;
855     simError();
856     }
857     }
858     }
859     }
860     notifyFortranCutoffs(&rcut_, &rsw_, &rnblist, &cp);
861 chrisfen 2309 // also send cutoff notification to electrostatics
862     setElectrostaticCutoffRadius(&rcut_);
863 gezelter 2204 }
864 gezelter 1490
865 chrisfen 2302 void SimInfo::setupElectrostaticSummationMethod( int isError ) {
866 chrisfen 2297
867     int errorOut;
868 chrisfen 2302 int esm = NONE;
869 chrisfen 2297 double alphaVal;
870 chrisfen 2309 double dielectric;
871 chrisfen 2297
872     errorOut = isError;
873 chrisfen 2309 alphaVal = simParams_->getDampingAlpha();
874     dielectric = simParams_->getDielectric();
875 chrisfen 2297
876 chrisfen 2302 if (simParams_->haveElectrostaticSummationMethod()) {
877 chrisfen 2303 std::string myMethod = simParams_->getElectrostaticSummationMethod();
878 chrisfen 2302 if (myMethod == "NONE") {
879     esm = NONE;
880 chrisfen 2297 } else {
881 chrisfen 2302 if (myMethod == "UNDAMPED_WOLF") {
882     esm = UNDAMPED_WOLF;
883 chrisfen 2297 } else {
884 chrisfen 2302 if (myMethod == "DAMPED_WOLF") {
885 chrisfen 2303 esm = DAMPED_WOLF;
886 chrisfen 2297 if (!simParams_->haveDampingAlpha()) {
887     //throw error
888     sprintf( painCave.errMsg,
889 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);
890 chrisfen 2297 painCave.isFatal = 0;
891     simError();
892     }
893     } else {
894 chrisfen 2302 if (myMethod == "REACTION_FIELD") {
895     esm = REACTION_FIELD;
896 chrisfen 2297 } else {
897     // throw error
898     sprintf( painCave.errMsg,
899 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() );
900 chrisfen 2297 painCave.isFatal = 1;
901     simError();
902     }
903     }
904     }
905     }
906     }
907 chrisfen 2309 // let's pass some summation method variables to fortran
908     setElectrostaticSummationMethod( &esm );
909     setDampedWolfAlpha( &alphaVal );
910     setReactionFieldDielectric( &dielectric );
911     initFortranFF( &esm, &errorOut );
912 chrisfen 2297 }
913    
914 gezelter 2204 void SimInfo::addProperty(GenericData* genData) {
915 gezelter 1930 properties_.addProperty(genData);
916 gezelter 2204 }
917 gezelter 1490
918 gezelter 2204 void SimInfo::removeProperty(const std::string& propName) {
919 gezelter 1930 properties_.removeProperty(propName);
920 gezelter 2204 }
921 gezelter 1490
922 gezelter 2204 void SimInfo::clearProperties() {
923 gezelter 1930 properties_.clearProperties();
924 gezelter 2204 }
925 gezelter 1490
926 gezelter 2204 std::vector<std::string> SimInfo::getPropertyNames() {
927 gezelter 1930 return properties_.getPropertyNames();
928 gezelter 2204 }
929 gezelter 1930
930 gezelter 2204 std::vector<GenericData*> SimInfo::getProperties() {
931 gezelter 1930 return properties_.getProperties();
932 gezelter 2204 }
933 gezelter 1490
934 gezelter 2204 GenericData* SimInfo::getPropertyByName(const std::string& propName) {
935 gezelter 1930 return properties_.getPropertyByName(propName);
936 gezelter 2204 }
937 gezelter 1490
938 gezelter 2204 void SimInfo::setSnapshotManager(SnapshotManager* sman) {
939 tim 2116 if (sman_ == sman) {
940 gezelter 2204 return;
941 tim 2116 }
942     delete sman_;
943 gezelter 1930 sman_ = sman;
944 gezelter 1490
945 gezelter 1930 Molecule* mol;
946     RigidBody* rb;
947     Atom* atom;
948     SimInfo::MoleculeIterator mi;
949     Molecule::RigidBodyIterator rbIter;
950     Molecule::AtomIterator atomIter;;
951    
952     for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
953    
954 gezelter 2204 for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) {
955     atom->setSnapshotManager(sman_);
956     }
957 gezelter 1930
958 gezelter 2204 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
959     rb->setSnapshotManager(sman_);
960     }
961 gezelter 1930 }
962 gezelter 1490
963 gezelter 2204 }
964 gezelter 1490
965 gezelter 2204 Vector3d SimInfo::getComVel(){
966 gezelter 1930 SimInfo::MoleculeIterator i;
967     Molecule* mol;
968 gezelter 1490
969 gezelter 1930 Vector3d comVel(0.0);
970     double totalMass = 0.0;
971 gezelter 1490
972 gezelter 1930
973     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
974 gezelter 2204 double mass = mol->getMass();
975     totalMass += mass;
976     comVel += mass * mol->getComVel();
977 gezelter 1930 }
978 gezelter 1490
979 gezelter 1930 #ifdef IS_MPI
980     double tmpMass = totalMass;
981     Vector3d tmpComVel(comVel);
982     MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
983     MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
984     #endif
985    
986     comVel /= totalMass;
987    
988     return comVel;
989 gezelter 2204 }
990 gezelter 1490
991 gezelter 2204 Vector3d SimInfo::getCom(){
992 gezelter 1930 SimInfo::MoleculeIterator i;
993     Molecule* mol;
994 gezelter 1490
995 gezelter 1930 Vector3d com(0.0);
996     double totalMass = 0.0;
997    
998     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
999 gezelter 2204 double mass = mol->getMass();
1000     totalMass += mass;
1001     com += mass * mol->getCom();
1002 gezelter 1930 }
1003 gezelter 1490
1004     #ifdef IS_MPI
1005 gezelter 1930 double tmpMass = totalMass;
1006     Vector3d tmpCom(com);
1007     MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1008     MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1009 gezelter 1490 #endif
1010    
1011 gezelter 1930 com /= totalMass;
1012 gezelter 1490
1013 gezelter 1930 return com;
1014 gezelter 1490
1015 gezelter 2204 }
1016 gezelter 1930
1017 gezelter 2204 std::ostream& operator <<(std::ostream& o, SimInfo& info) {
1018 gezelter 1930
1019     return o;
1020 gezelter 2204 }
1021 chuckv 2252
1022    
1023     /*
1024     Returns center of mass and center of mass velocity in one function call.
1025     */
1026    
1027     void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){
1028     SimInfo::MoleculeIterator i;
1029     Molecule* mol;
1030    
1031    
1032     double totalMass = 0.0;
1033    
1034 gezelter 1930
1035 chuckv 2252 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1036     double mass = mol->getMass();
1037     totalMass += mass;
1038     com += mass * mol->getCom();
1039     comVel += mass * mol->getComVel();
1040     }
1041    
1042     #ifdef IS_MPI
1043     double tmpMass = totalMass;
1044     Vector3d tmpCom(com);
1045     Vector3d tmpComVel(comVel);
1046     MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1047     MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1048     MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1049     #endif
1050    
1051     com /= totalMass;
1052     comVel /= totalMass;
1053     }
1054    
1055     /*
1056     Return intertia tensor for entire system and angular momentum Vector.
1057 chuckv 2256
1058    
1059     [ Ixx -Ixy -Ixz ]
1060     J =| -Iyx Iyy -Iyz |
1061     [ -Izx -Iyz Izz ]
1062 chuckv 2252 */
1063    
1064     void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){
1065    
1066    
1067     double xx = 0.0;
1068     double yy = 0.0;
1069     double zz = 0.0;
1070     double xy = 0.0;
1071     double xz = 0.0;
1072     double yz = 0.0;
1073     Vector3d com(0.0);
1074     Vector3d comVel(0.0);
1075    
1076     getComAll(com, comVel);
1077    
1078     SimInfo::MoleculeIterator i;
1079     Molecule* mol;
1080    
1081     Vector3d thisq(0.0);
1082     Vector3d thisv(0.0);
1083    
1084     double thisMass = 0.0;
1085    
1086    
1087    
1088    
1089     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1090    
1091     thisq = mol->getCom()-com;
1092     thisv = mol->getComVel()-comVel;
1093     thisMass = mol->getMass();
1094     // Compute moment of intertia coefficients.
1095     xx += thisq[0]*thisq[0]*thisMass;
1096     yy += thisq[1]*thisq[1]*thisMass;
1097     zz += thisq[2]*thisq[2]*thisMass;
1098    
1099     // compute products of intertia
1100     xy += thisq[0]*thisq[1]*thisMass;
1101     xz += thisq[0]*thisq[2]*thisMass;
1102     yz += thisq[1]*thisq[2]*thisMass;
1103    
1104     angularMomentum += cross( thisq, thisv ) * thisMass;
1105    
1106     }
1107    
1108    
1109     inertiaTensor(0,0) = yy + zz;
1110     inertiaTensor(0,1) = -xy;
1111     inertiaTensor(0,2) = -xz;
1112     inertiaTensor(1,0) = -xy;
1113 chuckv 2256 inertiaTensor(1,1) = xx + zz;
1114 chuckv 2252 inertiaTensor(1,2) = -yz;
1115     inertiaTensor(2,0) = -xz;
1116     inertiaTensor(2,1) = -yz;
1117     inertiaTensor(2,2) = xx + yy;
1118    
1119     #ifdef IS_MPI
1120     Mat3x3d tmpI(inertiaTensor);
1121     Vector3d tmpAngMom;
1122     MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1123     MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1124     #endif
1125    
1126     return;
1127     }
1128    
1129     //Returns the angular momentum of the system
1130     Vector3d SimInfo::getAngularMomentum(){
1131    
1132     Vector3d com(0.0);
1133     Vector3d comVel(0.0);
1134     Vector3d angularMomentum(0.0);
1135    
1136     getComAll(com,comVel);
1137    
1138     SimInfo::MoleculeIterator i;
1139     Molecule* mol;
1140    
1141 chuckv 2256 Vector3d thisr(0.0);
1142     Vector3d thisp(0.0);
1143 chuckv 2252
1144 chuckv 2256 double thisMass;
1145 chuckv 2252
1146     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1147 chuckv 2256 thisMass = mol->getMass();
1148     thisr = mol->getCom()-com;
1149     thisp = (mol->getComVel()-comVel)*thisMass;
1150 chuckv 2252
1151 chuckv 2256 angularMomentum += cross( thisr, thisp );
1152    
1153 chuckv 2252 }
1154    
1155     #ifdef IS_MPI
1156     Vector3d tmpAngMom;
1157     MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1158     #endif
1159    
1160     return angularMomentum;
1161     }
1162    
1163    
1164 gezelter 1930 }//end namespace oopse
1165