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

File Contents

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