ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1583
Committed: Thu Jun 16 22:00:08 2011 UTC (13 years, 10 months ago) by gezelter
File size: 38472 byte(s)
Log Message:
Bug squashing

File Contents

# User Rev Content
1 gezelter 507 /*
2 gezelter 246 * 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 gezelter 1390 * 1. Redistributions of source code must retain the above copyright
10 gezelter 246 * notice, this list of conditions and the following disclaimer.
11     *
12 gezelter 1390 * 2. Redistributions in binary form must reproduce the above copyright
13 gezelter 246 * notice, this list of conditions and the following disclaimer in the
14     * documentation and/or other materials provided with the
15     * distribution.
16     *
17     * This software is provided "AS IS," without a warranty of any
18     * kind. All express or implied conditions, representations and
19     * warranties, including any implied warranty of merchantability,
20     * fitness for a particular purpose or non-infringement, are hereby
21     * excluded. The University of Notre Dame and its licensors shall not
22     * be liable for any damages suffered by licensee as a result of
23     * using, modifying or distributing the software or its
24     * derivatives. In no event will the University of Notre Dame or its
25     * licensors be liable for any lost revenue, profit or data, or for
26     * direct, indirect, special, consequential, incidental or punitive
27     * damages, however caused and regardless of the theory of liability,
28     * arising out of the use of or inability to use software, even if the
29     * University of Notre Dame has been advised of the possibility of
30     * such damages.
31 gezelter 1390 *
32     * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33     * research, please cite the appropriate papers when you publish your
34     * work. Good starting points are:
35     *
36     * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37     * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38     * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39     * [4] Vardeman & Gezelter, in progress (2009).
40 gezelter 246 */
41    
42     /**
43     * @file SimInfo.cpp
44     * @author tlin
45     * @date 11/02/2004
46     * @version 1.0
47     */
48 gezelter 2
49 gezelter 246 #include <algorithm>
50     #include <set>
51 tim 749 #include <map>
52 gezelter 2
53 tim 3 #include "brains/SimInfo.hpp"
54 gezelter 246 #include "math/Vector3.hpp"
55     #include "primitives/Molecule.hpp"
56 tim 1024 #include "primitives/StuntDouble.hpp"
57 gezelter 246 #include "utils/MemoryUtils.hpp"
58 tim 3 #include "utils/simError.h"
59 tim 316 #include "selection/SelectionManager.hpp"
60 chuckv 834 #include "io/ForceFieldOptions.hpp"
61     #include "UseTheForce/ForceField.hpp"
62 gezelter 1530 #include "nonbonded/SwitchingFunction.hpp"
63 gezelter 2
64 gezelter 1528 using namespace std;
65 gezelter 1390 namespace OpenMD {
66 tim 749
67 tim 770 SimInfo::SimInfo(ForceField* ff, Globals* simParams) :
68     forceField_(ff), simParams_(simParams),
69 gezelter 945 ndf_(0), fdf_local(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0),
70 gezelter 507 nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0),
71     nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
72 gezelter 1277 nAtoms_(0), nBonds_(0), nBends_(0), nTorsions_(0), nInversions_(0),
73     nRigidBodies_(0), nIntegrableObjects_(0), nCutoffGroups_(0),
74 gezelter 1569 nConstraints_(0), sman_(NULL), topologyDone_(false),
75 gezelter 1528 calcBoxDipole_(false), useAtomicVirial_(true) {
76    
77     MoleculeStamp* molStamp;
78     int nMolWithSameStamp;
79     int nCutoffAtoms = 0; // number of atoms belong to cutoff groups
80     int nGroups = 0; //total cutoff groups defined in meta-data file
81     CutoffGroupStamp* cgStamp;
82     RigidBodyStamp* rbStamp;
83     int nRigidAtoms = 0;
84    
85     vector<Component*> components = simParams->getComponents();
86    
87     for (vector<Component*>::iterator i = components.begin(); i !=components.end(); ++i) {
88     molStamp = (*i)->getMoleculeStamp();
89     nMolWithSameStamp = (*i)->getNMol();
90 tim 770
91 gezelter 1528 addMoleculeStamp(molStamp, nMolWithSameStamp);
92    
93     //calculate atoms in molecules
94     nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp;
95    
96     //calculate atoms in cutoff groups
97     int nAtomsInGroups = 0;
98     int nCutoffGroupsInStamp = molStamp->getNCutoffGroups();
99    
100     for (int j=0; j < nCutoffGroupsInStamp; j++) {
101     cgStamp = molStamp->getCutoffGroupStamp(j);
102     nAtomsInGroups += cgStamp->getNMembers();
103 gezelter 507 }
104 gezelter 1528
105     nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
106    
107     nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;
108    
109     //calculate atoms in rigid bodies
110     int nAtomsInRigidBodies = 0;
111     int nRigidBodiesInStamp = molStamp->getNRigidBodies();
112    
113     for (int j=0; j < nRigidBodiesInStamp; j++) {
114     rbStamp = molStamp->getRigidBodyStamp(j);
115     nAtomsInRigidBodies += rbStamp->getNMembers();
116     }
117    
118     nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
119     nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;
120    
121     }
122    
123     //every free atom (atom does not belong to cutoff groups) is a cutoff
124     //group therefore the total number of cutoff groups in the system is
125     //equal to the total number of atoms minus number of atoms belong to
126     //cutoff group defined in meta-data file plus the number of cutoff
127     //groups defined in meta-data file
128 gezelter 1540
129 gezelter 1528 nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
130    
131     //every free atom (atom does not belong to rigid bodies) is an
132     //integrable object therefore the total number of integrable objects
133     //in the system is equal to the total number of atoms minus number of
134     //atoms belong to rigid body defined in meta-data file plus the number
135     //of rigid bodies defined in meta-data file
136     nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms
137     + nGlobalRigidBodies_;
138    
139     nGlobalMols_ = molStampIds_.size();
140     molToProcMap_.resize(nGlobalMols_);
141     }
142 chrisfen 645
143 gezelter 507 SimInfo::~SimInfo() {
144 gezelter 1528 map<int, Molecule*>::iterator i;
145 tim 398 for (i = molecules_.begin(); i != molecules_.end(); ++i) {
146 gezelter 507 delete i->second;
147 tim 398 }
148     molecules_.clear();
149 tim 490
150 gezelter 246 delete sman_;
151     delete simParams_;
152     delete forceField_;
153 gezelter 507 }
154 gezelter 2
155    
156 gezelter 507 bool SimInfo::addMolecule(Molecule* mol) {
157 gezelter 246 MoleculeIterator i;
158 gezelter 1528
159 gezelter 246 i = molecules_.find(mol->getGlobalIndex());
160     if (i == molecules_.end() ) {
161 gezelter 1528
162     molecules_.insert(make_pair(mol->getGlobalIndex(), mol));
163    
164 gezelter 507 nAtoms_ += mol->getNAtoms();
165     nBonds_ += mol->getNBonds();
166     nBends_ += mol->getNBends();
167     nTorsions_ += mol->getNTorsions();
168 gezelter 1277 nInversions_ += mol->getNInversions();
169 gezelter 507 nRigidBodies_ += mol->getNRigidBodies();
170     nIntegrableObjects_ += mol->getNIntegrableObjects();
171     nCutoffGroups_ += mol->getNCutoffGroups();
172     nConstraints_ += mol->getNConstraintPairs();
173 gezelter 1528
174 gezelter 1287 addInteractionPairs(mol);
175 gezelter 1528
176 gezelter 507 return true;
177 gezelter 246 } else {
178 gezelter 507 return false;
179 gezelter 246 }
180 gezelter 507 }
181 gezelter 1528
182 gezelter 507 bool SimInfo::removeMolecule(Molecule* mol) {
183 gezelter 246 MoleculeIterator i;
184     i = molecules_.find(mol->getGlobalIndex());
185 gezelter 2
186 gezelter 246 if (i != molecules_.end() ) {
187 gezelter 2
188 gezelter 507 assert(mol == i->second);
189 gezelter 246
190 gezelter 507 nAtoms_ -= mol->getNAtoms();
191     nBonds_ -= mol->getNBonds();
192     nBends_ -= mol->getNBends();
193     nTorsions_ -= mol->getNTorsions();
194 gezelter 1277 nInversions_ -= mol->getNInversions();
195 gezelter 507 nRigidBodies_ -= mol->getNRigidBodies();
196     nIntegrableObjects_ -= mol->getNIntegrableObjects();
197     nCutoffGroups_ -= mol->getNCutoffGroups();
198     nConstraints_ -= mol->getNConstraintPairs();
199 gezelter 2
200 gezelter 1287 removeInteractionPairs(mol);
201 gezelter 507 molecules_.erase(mol->getGlobalIndex());
202 gezelter 2
203 gezelter 507 delete mol;
204 gezelter 246
205 gezelter 507 return true;
206 gezelter 246 } else {
207 gezelter 507 return false;
208 gezelter 246 }
209 gezelter 507 }
210 gezelter 246
211    
212 gezelter 507 Molecule* SimInfo::beginMolecule(MoleculeIterator& i) {
213 gezelter 246 i = molecules_.begin();
214     return i == molecules_.end() ? NULL : i->second;
215 gezelter 507 }
216 gezelter 246
217 gezelter 507 Molecule* SimInfo::nextMolecule(MoleculeIterator& i) {
218 gezelter 246 ++i;
219     return i == molecules_.end() ? NULL : i->second;
220 gezelter 507 }
221 gezelter 2
222    
223 gezelter 507 void SimInfo::calcNdf() {
224 gezelter 246 int ndf_local;
225     MoleculeIterator i;
226 gezelter 1528 vector<StuntDouble*>::iterator j;
227 gezelter 246 Molecule* mol;
228     StuntDouble* integrableObject;
229 gezelter 2
230 gezelter 246 ndf_local = 0;
231    
232     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
233 gezelter 507 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
234     integrableObject = mol->nextIntegrableObject(j)) {
235 gezelter 2
236 gezelter 507 ndf_local += 3;
237 gezelter 2
238 gezelter 507 if (integrableObject->isDirectional()) {
239     if (integrableObject->isLinear()) {
240     ndf_local += 2;
241     } else {
242     ndf_local += 3;
243     }
244     }
245 gezelter 246
246 tim 770 }
247     }
248 gezelter 246
249     // n_constraints is local, so subtract them on each processor
250     ndf_local -= nConstraints_;
251    
252     #ifdef IS_MPI
253     MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
254     #else
255     ndf_ = ndf_local;
256     #endif
257    
258     // nZconstraints_ is global, as are the 3 COM translations for the
259     // entire system:
260     ndf_ = ndf_ - 3 - nZconstraint_;
261    
262 gezelter 507 }
263 gezelter 2
264 gezelter 945 int SimInfo::getFdf() {
265     #ifdef IS_MPI
266     MPI_Allreduce(&fdf_local,&fdf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
267     #else
268     fdf_ = fdf_local;
269     #endif
270     return fdf_;
271     }
272 gezelter 1577
273     unsigned int SimInfo::getNLocalCutoffGroups(){
274     int nLocalCutoffAtoms = 0;
275     Molecule* mol;
276     MoleculeIterator mi;
277     CutoffGroup* cg;
278     Molecule::CutoffGroupIterator ci;
279 gezelter 945
280 gezelter 1577 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
281    
282     for (cg = mol->beginCutoffGroup(ci); cg != NULL;
283     cg = mol->nextCutoffGroup(ci)) {
284     nLocalCutoffAtoms += cg->getNumAtom();
285    
286     }
287     }
288    
289     return nAtoms_ - nLocalCutoffAtoms + nCutoffGroups_;
290     }
291    
292 gezelter 507 void SimInfo::calcNdfRaw() {
293 gezelter 246 int ndfRaw_local;
294 gezelter 2
295 gezelter 246 MoleculeIterator i;
296 gezelter 1528 vector<StuntDouble*>::iterator j;
297 gezelter 246 Molecule* mol;
298     StuntDouble* integrableObject;
299    
300     // Raw degrees of freedom that we have to set
301     ndfRaw_local = 0;
302    
303     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
304 gezelter 507 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
305     integrableObject = mol->nextIntegrableObject(j)) {
306 gezelter 246
307 gezelter 507 ndfRaw_local += 3;
308 gezelter 246
309 gezelter 507 if (integrableObject->isDirectional()) {
310     if (integrableObject->isLinear()) {
311     ndfRaw_local += 2;
312     } else {
313     ndfRaw_local += 3;
314     }
315     }
316 gezelter 246
317 gezelter 507 }
318 gezelter 246 }
319    
320     #ifdef IS_MPI
321     MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
322     #else
323     ndfRaw_ = ndfRaw_local;
324     #endif
325 gezelter 507 }
326 gezelter 2
327 gezelter 507 void SimInfo::calcNdfTrans() {
328 gezelter 246 int ndfTrans_local;
329 gezelter 2
330 gezelter 246 ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
331 gezelter 2
332    
333 gezelter 246 #ifdef IS_MPI
334     MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
335     #else
336     ndfTrans_ = ndfTrans_local;
337     #endif
338 gezelter 2
339 gezelter 246 ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
340    
341 gezelter 507 }
342 gezelter 2
343 gezelter 1287 void SimInfo::addInteractionPairs(Molecule* mol) {
344     ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
345 gezelter 1528 vector<Bond*>::iterator bondIter;
346     vector<Bend*>::iterator bendIter;
347     vector<Torsion*>::iterator torsionIter;
348     vector<Inversion*>::iterator inversionIter;
349 gezelter 246 Bond* bond;
350     Bend* bend;
351     Torsion* torsion;
352 gezelter 1277 Inversion* inversion;
353 gezelter 246 int a;
354     int b;
355     int c;
356     int d;
357 tim 749
358 gezelter 1287 // atomGroups can be used to add special interaction maps between
359     // groups of atoms that are in two separate rigid bodies.
360     // However, most site-site interactions between two rigid bodies
361     // are probably not special, just the ones between the physically
362     // bonded atoms. Interactions *within* a single rigid body should
363     // always be excluded. These are done at the bottom of this
364     // function.
365    
366 gezelter 1528 map<int, set<int> > atomGroups;
367 tim 749 Molecule::RigidBodyIterator rbIter;
368     RigidBody* rb;
369     Molecule::IntegrableObjectIterator ii;
370     StuntDouble* integrableObject;
371 gezelter 246
372 gezelter 1287 for (integrableObject = mol->beginIntegrableObject(ii);
373     integrableObject != NULL;
374     integrableObject = mol->nextIntegrableObject(ii)) {
375    
376 tim 749 if (integrableObject->isRigidBody()) {
377 gezelter 1287 rb = static_cast<RigidBody*>(integrableObject);
378 gezelter 1528 vector<Atom*> atoms = rb->getAtoms();
379     set<int> rigidAtoms;
380 gezelter 1287 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
381     rigidAtoms.insert(atoms[i]->getGlobalIndex());
382     }
383     for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
384 gezelter 1528 atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
385 gezelter 1287 }
386 tim 749 } else {
387 gezelter 1528 set<int> oneAtomSet;
388 tim 749 oneAtomSet.insert(integrableObject->getGlobalIndex());
389 gezelter 1528 atomGroups.insert(map<int, set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));
390 tim 749 }
391     }
392 gezelter 1287
393     for (bond= mol->beginBond(bondIter); bond != NULL;
394     bond = mol->nextBond(bondIter)) {
395 tim 749
396 gezelter 1287 a = bond->getAtomA()->getGlobalIndex();
397     b = bond->getAtomB()->getGlobalIndex();
398 tim 749
399 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
400     oneTwoInteractions_.addPair(a, b);
401     } else {
402     excludedInteractions_.addPair(a, b);
403     }
404 gezelter 246 }
405 gezelter 2
406 gezelter 1287 for (bend= mol->beginBend(bendIter); bend != NULL;
407     bend = mol->nextBend(bendIter)) {
408    
409 gezelter 507 a = bend->getAtomA()->getGlobalIndex();
410     b = bend->getAtomB()->getGlobalIndex();
411     c = bend->getAtomC()->getGlobalIndex();
412 gezelter 1287
413     if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
414     oneTwoInteractions_.addPair(a, b);
415     oneTwoInteractions_.addPair(b, c);
416     } else {
417     excludedInteractions_.addPair(a, b);
418     excludedInteractions_.addPair(b, c);
419     }
420 gezelter 2
421 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
422     oneThreeInteractions_.addPair(a, c);
423     } else {
424     excludedInteractions_.addPair(a, c);
425     }
426 gezelter 246 }
427 gezelter 2
428 gezelter 1287 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
429     torsion = mol->nextTorsion(torsionIter)) {
430    
431 gezelter 507 a = torsion->getAtomA()->getGlobalIndex();
432     b = torsion->getAtomB()->getGlobalIndex();
433     c = torsion->getAtomC()->getGlobalIndex();
434 gezelter 1287 d = torsion->getAtomD()->getGlobalIndex();
435 cli2 1290
436 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
437     oneTwoInteractions_.addPair(a, b);
438     oneTwoInteractions_.addPair(b, c);
439     oneTwoInteractions_.addPair(c, d);
440     } else {
441     excludedInteractions_.addPair(a, b);
442     excludedInteractions_.addPair(b, c);
443     excludedInteractions_.addPair(c, d);
444     }
445 gezelter 2
446 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
447     oneThreeInteractions_.addPair(a, c);
448     oneThreeInteractions_.addPair(b, d);
449     } else {
450     excludedInteractions_.addPair(a, c);
451     excludedInteractions_.addPair(b, d);
452     }
453 tim 749
454 gezelter 1287 if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
455     oneFourInteractions_.addPair(a, d);
456     } else {
457     excludedInteractions_.addPair(a, d);
458     }
459 gezelter 2 }
460    
461 gezelter 1277 for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
462     inversion = mol->nextInversion(inversionIter)) {
463 gezelter 1287
464 gezelter 1277 a = inversion->getAtomA()->getGlobalIndex();
465     b = inversion->getAtomB()->getGlobalIndex();
466     c = inversion->getAtomC()->getGlobalIndex();
467     d = inversion->getAtomD()->getGlobalIndex();
468    
469 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
470     oneTwoInteractions_.addPair(a, b);
471     oneTwoInteractions_.addPair(a, c);
472     oneTwoInteractions_.addPair(a, d);
473     } else {
474     excludedInteractions_.addPair(a, b);
475     excludedInteractions_.addPair(a, c);
476     excludedInteractions_.addPair(a, d);
477     }
478 gezelter 1277
479 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
480     oneThreeInteractions_.addPair(b, c);
481     oneThreeInteractions_.addPair(b, d);
482     oneThreeInteractions_.addPair(c, d);
483     } else {
484     excludedInteractions_.addPair(b, c);
485     excludedInteractions_.addPair(b, d);
486     excludedInteractions_.addPair(c, d);
487     }
488 gezelter 1277 }
489    
490 gezelter 1287 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
491     rb = mol->nextRigidBody(rbIter)) {
492 gezelter 1528 vector<Atom*> atoms = rb->getAtoms();
493 gezelter 1287 for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
494     for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
495 gezelter 507 a = atoms[i]->getGlobalIndex();
496     b = atoms[j]->getGlobalIndex();
497 gezelter 1287 excludedInteractions_.addPair(a, b);
498 gezelter 507 }
499     }
500 tim 430 }
501    
502 gezelter 507 }
503 gezelter 246
504 gezelter 1287 void SimInfo::removeInteractionPairs(Molecule* mol) {
505     ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
506 gezelter 1528 vector<Bond*>::iterator bondIter;
507     vector<Bend*>::iterator bendIter;
508     vector<Torsion*>::iterator torsionIter;
509     vector<Inversion*>::iterator inversionIter;
510 gezelter 246 Bond* bond;
511     Bend* bend;
512     Torsion* torsion;
513 gezelter 1277 Inversion* inversion;
514 gezelter 246 int a;
515     int b;
516     int c;
517     int d;
518 tim 749
519 gezelter 1528 map<int, set<int> > atomGroups;
520 tim 749 Molecule::RigidBodyIterator rbIter;
521     RigidBody* rb;
522     Molecule::IntegrableObjectIterator ii;
523     StuntDouble* integrableObject;
524 gezelter 246
525 gezelter 1287 for (integrableObject = mol->beginIntegrableObject(ii);
526     integrableObject != NULL;
527     integrableObject = mol->nextIntegrableObject(ii)) {
528    
529 tim 749 if (integrableObject->isRigidBody()) {
530 gezelter 1287 rb = static_cast<RigidBody*>(integrableObject);
531 gezelter 1528 vector<Atom*> atoms = rb->getAtoms();
532     set<int> rigidAtoms;
533 gezelter 1287 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
534     rigidAtoms.insert(atoms[i]->getGlobalIndex());
535     }
536     for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
537 gezelter 1528 atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
538 gezelter 1287 }
539 tim 749 } else {
540 gezelter 1528 set<int> oneAtomSet;
541 tim 749 oneAtomSet.insert(integrableObject->getGlobalIndex());
542 gezelter 1528 atomGroups.insert(map<int, set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));
543 tim 749 }
544     }
545    
546 gezelter 1287 for (bond= mol->beginBond(bondIter); bond != NULL;
547     bond = mol->nextBond(bondIter)) {
548    
549     a = bond->getAtomA()->getGlobalIndex();
550     b = bond->getAtomB()->getGlobalIndex();
551 tim 749
552 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
553     oneTwoInteractions_.removePair(a, b);
554     } else {
555     excludedInteractions_.removePair(a, b);
556     }
557 gezelter 2 }
558 gezelter 246
559 gezelter 1287 for (bend= mol->beginBend(bendIter); bend != NULL;
560     bend = mol->nextBend(bendIter)) {
561    
562 gezelter 507 a = bend->getAtomA()->getGlobalIndex();
563     b = bend->getAtomB()->getGlobalIndex();
564     c = bend->getAtomC()->getGlobalIndex();
565 gezelter 1287
566     if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
567     oneTwoInteractions_.removePair(a, b);
568     oneTwoInteractions_.removePair(b, c);
569     } else {
570     excludedInteractions_.removePair(a, b);
571     excludedInteractions_.removePair(b, c);
572     }
573 gezelter 246
574 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
575     oneThreeInteractions_.removePair(a, c);
576     } else {
577     excludedInteractions_.removePair(a, c);
578     }
579 gezelter 2 }
580 gezelter 246
581 gezelter 1287 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
582     torsion = mol->nextTorsion(torsionIter)) {
583    
584 gezelter 507 a = torsion->getAtomA()->getGlobalIndex();
585     b = torsion->getAtomB()->getGlobalIndex();
586     c = torsion->getAtomC()->getGlobalIndex();
587 gezelter 1287 d = torsion->getAtomD()->getGlobalIndex();
588    
589     if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
590     oneTwoInteractions_.removePair(a, b);
591     oneTwoInteractions_.removePair(b, c);
592     oneTwoInteractions_.removePair(c, d);
593     } else {
594     excludedInteractions_.removePair(a, b);
595     excludedInteractions_.removePair(b, c);
596     excludedInteractions_.removePair(c, d);
597     }
598 gezelter 246
599 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
600     oneThreeInteractions_.removePair(a, c);
601     oneThreeInteractions_.removePair(b, d);
602     } else {
603     excludedInteractions_.removePair(a, c);
604     excludedInteractions_.removePair(b, d);
605     }
606 tim 749
607 gezelter 1287 if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
608     oneFourInteractions_.removePair(a, d);
609     } else {
610     excludedInteractions_.removePair(a, d);
611     }
612     }
613 tim 749
614 gezelter 1287 for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
615     inversion = mol->nextInversion(inversionIter)) {
616 tim 749
617 gezelter 1277 a = inversion->getAtomA()->getGlobalIndex();
618     b = inversion->getAtomB()->getGlobalIndex();
619     c = inversion->getAtomC()->getGlobalIndex();
620     d = inversion->getAtomD()->getGlobalIndex();
621    
622 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
623     oneTwoInteractions_.removePair(a, b);
624     oneTwoInteractions_.removePair(a, c);
625     oneTwoInteractions_.removePair(a, d);
626     } else {
627     excludedInteractions_.removePair(a, b);
628     excludedInteractions_.removePair(a, c);
629     excludedInteractions_.removePair(a, d);
630     }
631 gezelter 1277
632 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
633     oneThreeInteractions_.removePair(b, c);
634     oneThreeInteractions_.removePair(b, d);
635     oneThreeInteractions_.removePair(c, d);
636     } else {
637     excludedInteractions_.removePair(b, c);
638     excludedInteractions_.removePair(b, d);
639     excludedInteractions_.removePair(c, d);
640     }
641 gezelter 1277 }
642    
643 gezelter 1287 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
644     rb = mol->nextRigidBody(rbIter)) {
645 gezelter 1528 vector<Atom*> atoms = rb->getAtoms();
646 gezelter 1287 for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
647     for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
648 gezelter 507 a = atoms[i]->getGlobalIndex();
649     b = atoms[j]->getGlobalIndex();
650 gezelter 1287 excludedInteractions_.removePair(a, b);
651 gezelter 507 }
652     }
653 tim 430 }
654 gezelter 1287
655 gezelter 507 }
656 gezelter 1287
657    
658 gezelter 507 void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) {
659 gezelter 246 int curStampId;
660 gezelter 1287
661 gezelter 246 //index from 0
662     curStampId = moleculeStamps_.size();
663 gezelter 2
664 gezelter 246 moleculeStamps_.push_back(molStamp);
665     molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
666 gezelter 507 }
667 gezelter 2
668 gezelter 1530
669     /**
670     * update
671     *
672 gezelter 1535 * Performs the global checks and variable settings after the
673     * objects have been created.
674 gezelter 1530 *
675     */
676 gezelter 1535 void SimInfo::update() {
677 gezelter 1530 setupSimVariables();
678 gezelter 246 calcNdf();
679     calcNdfRaw();
680     calcNdfTrans();
681 gezelter 507 }
682 gezelter 1528
683 gezelter 1535 /**
684     * getSimulatedAtomTypes
685     *
686     * Returns an STL set of AtomType* that are actually present in this
687     * simulation. Must query all processors to assemble this information.
688     *
689     */
690 gezelter 1528 set<AtomType*> SimInfo::getSimulatedAtomTypes() {
691 gezelter 246 SimInfo::MoleculeIterator mi;
692     Molecule* mol;
693     Molecule::AtomIterator ai;
694     Atom* atom;
695 gezelter 1528 set<AtomType*> atomTypes;
696    
697 gezelter 1529 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
698 gezelter 507 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
699     atomTypes.insert(atom->getAtomType());
700 gezelter 1529 }
701     }
702 gezelter 2
703 gezelter 1535 #ifdef IS_MPI
704    
705     // loop over the found atom types on this processor, and add their
706     // numerical idents to a vector:
707    
708     vector<int> foundTypes;
709     set<AtomType*>::iterator i;
710     for (i = atomTypes.begin(); i != atomTypes.end(); ++i)
711     foundTypes.push_back( (*i)->getIdent() );
712    
713     // count_local holds the number of found types on this processor
714     int count_local = foundTypes.size();
715    
716     // count holds the total number of found types on all processors
717     // (some will be redundant with the ones found locally):
718     int count;
719     MPI::COMM_WORLD.Allreduce(&count_local, &count, 1, MPI::INT, MPI::SUM);
720    
721     // create a vector to hold the globally found types, and resize it:
722     vector<int> ftGlobal;
723     ftGlobal.resize(count);
724     vector<int> counts;
725    
726     int nproc = MPI::COMM_WORLD.Get_size();
727     counts.resize(nproc);
728     vector<int> disps;
729     disps.resize(nproc);
730    
731     // now spray out the foundTypes to all the other processors:
732 gezelter 2
733 gezelter 1535 MPI::COMM_WORLD.Allgatherv(&foundTypes[0], count_local, MPI::INT,
734     &ftGlobal[0], &counts[0], &disps[0], MPI::INT);
735 gezelter 1126
736 gezelter 1535 // foundIdents is a stl set, so inserting an already found ident
737     // will have no effect.
738     set<int> foundIdents;
739     vector<int>::iterator j;
740     for (j = ftGlobal.begin(); j != ftGlobal.end(); ++j)
741     foundIdents.insert((*j));
742 gezelter 1528
743 gezelter 1535 // now iterate over the foundIdents and get the actual atom types
744     // that correspond to these:
745     set<int>::iterator it;
746     for (it = foundIdents.begin(); it != foundIdents.end(); ++it)
747     atomTypes.insert( forceField_->getAtomType((*it)) );
748    
749     #endif
750 gezelter 1530
751 gezelter 1535 return atomTypes;
752 gezelter 1528 }
753 chrisfen 611
754 gezelter 1534 void SimInfo::setupSimVariables() {
755     useAtomicVirial_ = simParams_->getUseAtomicVirial();
756     // we only call setAccumulateBoxDipole if the accumulateBoxDipole parameter is true
757     calcBoxDipole_ = false;
758     if ( simParams_->haveAccumulateBoxDipole() )
759     if ( simParams_->getAccumulateBoxDipole() ) {
760     calcBoxDipole_ = true;
761     }
762 gezelter 1583
763 gezelter 1528 set<AtomType*>::iterator i;
764     set<AtomType*> atomTypes;
765 gezelter 1534 atomTypes = getSimulatedAtomTypes();
766 gezelter 1528 int usesElectrostatic = 0;
767     int usesMetallic = 0;
768     int usesDirectional = 0;
769 gezelter 246 //loop over all of the atom types
770     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
771 gezelter 1528 usesElectrostatic |= (*i)->isElectrostatic();
772     usesMetallic |= (*i)->isMetal();
773     usesDirectional |= (*i)->isDirectional();
774 gezelter 246 }
775 gezelter 1583
776 gezelter 246 #ifdef IS_MPI
777     int temp;
778 gezelter 1528 temp = usesDirectional;
779     MPI_Allreduce(&temp, &usesDirectionalAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
780 gezelter 1583
781 gezelter 1528 temp = usesMetallic;
782     MPI_Allreduce(&temp, &usesMetallicAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
783 gezelter 1583
784 gezelter 1528 temp = usesElectrostatic;
785     MPI_Allreduce(&temp, &usesElectrostaticAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
786 gezelter 1583 #else
787    
788     usesDirectionalAtoms_ = usesDirectional;
789     usesMetallicAtoms_ = usesMetallic;
790     usesElectrostaticAtoms_ = usesElectrostatic;
791    
792 gezelter 2 #endif
793 gezelter 1583
794     requiresPrepair_ = usesMetallicAtoms_ ? true : false;
795     requiresSkipCorrection_ = usesElectrostaticAtoms_ ? true : false;
796     requiresSelfCorrection_ = usesElectrostaticAtoms_ ? true : false;
797 gezelter 507 }
798 gezelter 2
799 gezelter 1549
800     vector<int> SimInfo::getGlobalAtomIndices() {
801     SimInfo::MoleculeIterator mi;
802     Molecule* mol;
803     Molecule::AtomIterator ai;
804     Atom* atom;
805    
806     vector<int> GlobalAtomIndices(getNAtoms(), 0);
807    
808     for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
809    
810     for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
811     GlobalAtomIndices[atom->getLocalIndex()] = atom->getGlobalIndex();
812     }
813     }
814     return GlobalAtomIndices;
815     }
816    
817    
818     vector<int> SimInfo::getGlobalGroupIndices() {
819     SimInfo::MoleculeIterator mi;
820     Molecule* mol;
821     Molecule::CutoffGroupIterator ci;
822     CutoffGroup* cg;
823    
824     vector<int> GlobalGroupIndices;
825    
826     for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
827    
828     //local index of cutoff group is trivial, it only depends on the
829     //order of travesing
830     for (cg = mol->beginCutoffGroup(ci); cg != NULL;
831     cg = mol->nextCutoffGroup(ci)) {
832     GlobalGroupIndices.push_back(cg->getGlobalIndex());
833     }
834     }
835     return GlobalGroupIndices;
836     }
837    
838    
839 gezelter 1569 void SimInfo::prepareTopology() {
840 gezelter 1287 int nExclude, nOneTwo, nOneThree, nOneFour;
841 gezelter 2
842 gezelter 246 //calculate mass ratio of cutoff group
843     SimInfo::MoleculeIterator mi;
844     Molecule* mol;
845     Molecule::CutoffGroupIterator ci;
846     CutoffGroup* cg;
847     Molecule::AtomIterator ai;
848     Atom* atom;
849 tim 963 RealType totalMass;
850 gezelter 246
851 gezelter 1581 /**
852     * The mass factor is the relative mass of an atom to the total
853     * mass of the cutoff group it belongs to. By default, all atoms
854     * are their own cutoff groups, and therefore have mass factors of
855     * 1. We need some special handling for massless atoms, which
856     * will be treated as carrying the entire mass of the cutoff
857     * group.
858     */
859 gezelter 1569 massFactors_.clear();
860 gezelter 1581 massFactors_.resize(getNAtoms(), 1.0);
861 gezelter 2
862 gezelter 246 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
863 gezelter 1569 for (cg = mol->beginCutoffGroup(ci); cg != NULL;
864     cg = mol->nextCutoffGroup(ci)) {
865 gezelter 2
866 gezelter 507 totalMass = cg->getMass();
867     for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
868 chrisfen 645 // Check for massless groups - set mfact to 1 if true
869 gezelter 1581 if (totalMass != 0)
870     massFactors_[atom->getLocalIndex()] = atom->getMass()/totalMass;
871 chrisfen 645 else
872 gezelter 1581 massFactors_[atom->getLocalIndex()] = 1.0;
873 gezelter 507 }
874     }
875 gezelter 246 }
876 gezelter 2
877 gezelter 1544 // Build the identArray_
878 gezelter 2
879 gezelter 1544 identArray_.clear();
880     identArray_.reserve(getNAtoms());
881 gezelter 246 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
882 gezelter 507 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
883 gezelter 1544 identArray_.push_back(atom->getIdent());
884 gezelter 507 }
885 gezelter 246 }
886    
887 gezelter 1569 //scan topology
888 gezelter 1287
889     nExclude = excludedInteractions_.getSize();
890     nOneTwo = oneTwoInteractions_.getSize();
891     nOneThree = oneThreeInteractions_.getSize();
892     nOneFour = oneFourInteractions_.getSize();
893    
894     int* excludeList = excludedInteractions_.getPairList();
895     int* oneTwoList = oneTwoInteractions_.getPairList();
896     int* oneThreeList = oneThreeInteractions_.getPairList();
897     int* oneFourList = oneFourInteractions_.getPairList();
898    
899 gezelter 1569 topologyDone_ = true;
900 gezelter 507 }
901 chrisfen 143
902 gezelter 507 void SimInfo::addProperty(GenericData* genData) {
903 gezelter 246 properties_.addProperty(genData);
904 gezelter 507 }
905 gezelter 2
906 gezelter 1528 void SimInfo::removeProperty(const string& propName) {
907 gezelter 246 properties_.removeProperty(propName);
908 gezelter 507 }
909 gezelter 2
910 gezelter 507 void SimInfo::clearProperties() {
911 gezelter 246 properties_.clearProperties();
912 gezelter 507 }
913 gezelter 2
914 gezelter 1528 vector<string> SimInfo::getPropertyNames() {
915 gezelter 246 return properties_.getPropertyNames();
916 gezelter 507 }
917 gezelter 246
918 gezelter 1528 vector<GenericData*> SimInfo::getProperties() {
919 gezelter 246 return properties_.getProperties();
920 gezelter 507 }
921 gezelter 2
922 gezelter 1528 GenericData* SimInfo::getPropertyByName(const string& propName) {
923 gezelter 246 return properties_.getPropertyByName(propName);
924 gezelter 507 }
925 gezelter 2
926 gezelter 507 void SimInfo::setSnapshotManager(SnapshotManager* sman) {
927 tim 432 if (sman_ == sman) {
928 gezelter 507 return;
929 tim 432 }
930     delete sman_;
931 gezelter 246 sman_ = sman;
932 gezelter 2
933 gezelter 246 Molecule* mol;
934     RigidBody* rb;
935     Atom* atom;
936 gezelter 1540 CutoffGroup* cg;
937 gezelter 246 SimInfo::MoleculeIterator mi;
938     Molecule::RigidBodyIterator rbIter;
939 gezelter 1540 Molecule::AtomIterator atomIter;
940     Molecule::CutoffGroupIterator cgIter;
941 gezelter 246
942     for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
943    
944 gezelter 507 for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) {
945     atom->setSnapshotManager(sman_);
946     }
947 gezelter 246
948 gezelter 507 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
949     rb->setSnapshotManager(sman_);
950     }
951 gezelter 1540
952     for (cg = mol->beginCutoffGroup(cgIter); cg != NULL; cg = mol->nextCutoffGroup(cgIter)) {
953     cg->setSnapshotManager(sman_);
954     }
955 gezelter 246 }
956 gezelter 2
957 gezelter 507 }
958 gezelter 2
959 gezelter 507 Vector3d SimInfo::getComVel(){
960 gezelter 246 SimInfo::MoleculeIterator i;
961     Molecule* mol;
962 gezelter 2
963 gezelter 246 Vector3d comVel(0.0);
964 tim 963 RealType totalMass = 0.0;
965 gezelter 2
966 gezelter 246
967     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
968 tim 963 RealType mass = mol->getMass();
969 gezelter 507 totalMass += mass;
970     comVel += mass * mol->getComVel();
971 gezelter 246 }
972 gezelter 2
973 gezelter 246 #ifdef IS_MPI
974 tim 963 RealType tmpMass = totalMass;
975 gezelter 246 Vector3d tmpComVel(comVel);
976 tim 963 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
977     MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
978 gezelter 246 #endif
979    
980     comVel /= totalMass;
981    
982     return comVel;
983 gezelter 507 }
984 gezelter 2
985 gezelter 507 Vector3d SimInfo::getCom(){
986 gezelter 246 SimInfo::MoleculeIterator i;
987     Molecule* mol;
988 gezelter 2
989 gezelter 246 Vector3d com(0.0);
990 tim 963 RealType totalMass = 0.0;
991 gezelter 246
992     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
993 tim 963 RealType mass = mol->getMass();
994 gezelter 507 totalMass += mass;
995     com += mass * mol->getCom();
996 gezelter 246 }
997 gezelter 2
998     #ifdef IS_MPI
999 tim 963 RealType tmpMass = totalMass;
1000 gezelter 246 Vector3d tmpCom(com);
1001 tim 963 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1002     MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1003 gezelter 2 #endif
1004    
1005 gezelter 246 com /= totalMass;
1006 gezelter 2
1007 gezelter 246 return com;
1008 gezelter 2
1009 gezelter 507 }
1010 gezelter 246
1011 gezelter 1528 ostream& operator <<(ostream& o, SimInfo& info) {
1012 gezelter 246
1013     return o;
1014 gezelter 507 }
1015 chuckv 555
1016    
1017     /*
1018     Returns center of mass and center of mass velocity in one function call.
1019     */
1020    
1021     void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){
1022     SimInfo::MoleculeIterator i;
1023     Molecule* mol;
1024    
1025    
1026 tim 963 RealType totalMass = 0.0;
1027 chuckv 555
1028 gezelter 246
1029 chuckv 555 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1030 tim 963 RealType mass = mol->getMass();
1031 chuckv 555 totalMass += mass;
1032     com += mass * mol->getCom();
1033     comVel += mass * mol->getComVel();
1034     }
1035    
1036     #ifdef IS_MPI
1037 tim 963 RealType tmpMass = totalMass;
1038 chuckv 555 Vector3d tmpCom(com);
1039     Vector3d tmpComVel(comVel);
1040 tim 963 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1041     MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1042     MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1043 chuckv 555 #endif
1044    
1045     com /= totalMass;
1046     comVel /= totalMass;
1047     }
1048    
1049     /*
1050     Return intertia tensor for entire system and angular momentum Vector.
1051 chuckv 557
1052    
1053     [ Ixx -Ixy -Ixz ]
1054 gezelter 1505 J =| -Iyx Iyy -Iyz |
1055 chuckv 557 [ -Izx -Iyz Izz ]
1056 chuckv 555 */
1057    
1058     void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){
1059    
1060    
1061 tim 963 RealType xx = 0.0;
1062     RealType yy = 0.0;
1063     RealType zz = 0.0;
1064     RealType xy = 0.0;
1065     RealType xz = 0.0;
1066     RealType yz = 0.0;
1067 chuckv 555 Vector3d com(0.0);
1068     Vector3d comVel(0.0);
1069    
1070     getComAll(com, comVel);
1071    
1072     SimInfo::MoleculeIterator i;
1073     Molecule* mol;
1074    
1075     Vector3d thisq(0.0);
1076     Vector3d thisv(0.0);
1077    
1078 tim 963 RealType thisMass = 0.0;
1079 chuckv 555
1080    
1081    
1082    
1083     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1084    
1085     thisq = mol->getCom()-com;
1086     thisv = mol->getComVel()-comVel;
1087     thisMass = mol->getMass();
1088     // Compute moment of intertia coefficients.
1089     xx += thisq[0]*thisq[0]*thisMass;
1090     yy += thisq[1]*thisq[1]*thisMass;
1091     zz += thisq[2]*thisq[2]*thisMass;
1092    
1093     // compute products of intertia
1094     xy += thisq[0]*thisq[1]*thisMass;
1095     xz += thisq[0]*thisq[2]*thisMass;
1096     yz += thisq[1]*thisq[2]*thisMass;
1097    
1098     angularMomentum += cross( thisq, thisv ) * thisMass;
1099    
1100     }
1101    
1102    
1103     inertiaTensor(0,0) = yy + zz;
1104     inertiaTensor(0,1) = -xy;
1105     inertiaTensor(0,2) = -xz;
1106     inertiaTensor(1,0) = -xy;
1107 chuckv 557 inertiaTensor(1,1) = xx + zz;
1108 chuckv 555 inertiaTensor(1,2) = -yz;
1109     inertiaTensor(2,0) = -xz;
1110     inertiaTensor(2,1) = -yz;
1111     inertiaTensor(2,2) = xx + yy;
1112    
1113     #ifdef IS_MPI
1114     Mat3x3d tmpI(inertiaTensor);
1115     Vector3d tmpAngMom;
1116 tim 963 MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1117     MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1118 chuckv 555 #endif
1119    
1120     return;
1121     }
1122    
1123     //Returns the angular momentum of the system
1124     Vector3d SimInfo::getAngularMomentum(){
1125    
1126     Vector3d com(0.0);
1127     Vector3d comVel(0.0);
1128     Vector3d angularMomentum(0.0);
1129    
1130     getComAll(com,comVel);
1131    
1132     SimInfo::MoleculeIterator i;
1133     Molecule* mol;
1134    
1135 chuckv 557 Vector3d thisr(0.0);
1136     Vector3d thisp(0.0);
1137 chuckv 555
1138 tim 963 RealType thisMass;
1139 chuckv 555
1140     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1141 chuckv 557 thisMass = mol->getMass();
1142     thisr = mol->getCom()-com;
1143     thisp = (mol->getComVel()-comVel)*thisMass;
1144 chuckv 555
1145 chuckv 557 angularMomentum += cross( thisr, thisp );
1146    
1147 chuckv 555 }
1148    
1149     #ifdef IS_MPI
1150     Vector3d tmpAngMom;
1151 tim 963 MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1152 chuckv 555 #endif
1153    
1154     return angularMomentum;
1155     }
1156    
1157 tim 1024 StuntDouble* SimInfo::getIOIndexToIntegrableObject(int index) {
1158     return IOIndexToIntegrableObject.at(index);
1159     }
1160    
1161 gezelter 1528 void SimInfo::setIOIndexToIntegrableObject(const vector<StuntDouble*>& v) {
1162 tim 1024 IOIndexToIntegrableObject= v;
1163     }
1164    
1165 chuckv 1103 /* Returns the Volume of the simulation based on a ellipsoid with semi-axes
1166     based on the radius of gyration V=4/3*Pi*R_1*R_2*R_3
1167     where R_i are related to the principle inertia moments R_i = sqrt(C*I_i/N), this reduces to
1168     V = 4/3*Pi*(C/N)^3/2*sqrt(det(I)). See S.E. Baltazar et. al. Comp. Mat. Sci. 37 (2006) 526-536.
1169     */
1170     void SimInfo::getGyrationalVolume(RealType &volume){
1171     Mat3x3d intTensor;
1172     RealType det;
1173     Vector3d dummyAngMom;
1174     RealType sysconstants;
1175     RealType geomCnst;
1176    
1177     geomCnst = 3.0/2.0;
1178     /* Get the inertial tensor and angular momentum for free*/
1179     getInertiaTensor(intTensor,dummyAngMom);
1180    
1181     det = intTensor.determinant();
1182     sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1183     volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(det);
1184     return;
1185     }
1186    
1187     void SimInfo::getGyrationalVolume(RealType &volume, RealType &detI){
1188     Mat3x3d intTensor;
1189     Vector3d dummyAngMom;
1190     RealType sysconstants;
1191     RealType geomCnst;
1192    
1193     geomCnst = 3.0/2.0;
1194     /* Get the inertial tensor and angular momentum for free*/
1195     getInertiaTensor(intTensor,dummyAngMom);
1196    
1197     detI = intTensor.determinant();
1198     sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1199     volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(detI);
1200     return;
1201     }
1202 tim 1024 /*
1203 gezelter 1528 void SimInfo::setStuntDoubleFromGlobalIndex(vector<StuntDouble*> v) {
1204 tim 1024 assert( v.size() == nAtoms_ + nRigidBodies_);
1205     sdByGlobalIndex_ = v;
1206     }
1207    
1208     StuntDouble* SimInfo::getStuntDoubleFromGlobalIndex(int index) {
1209     //assert(index < nAtoms_ + nRigidBodies_);
1210     return sdByGlobalIndex_.at(index);
1211     }
1212     */
1213 gezelter 1528 int SimInfo::getNGlobalConstraints() {
1214     int nGlobalConstraints;
1215     #ifdef IS_MPI
1216     MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM,
1217     MPI_COMM_WORLD);
1218     #else
1219     nGlobalConstraints = nConstraints_;
1220     #endif
1221     return nGlobalConstraints;
1222     }
1223    
1224 gezelter 1390 }//end namespace OpenMD
1225 gezelter 246

Properties

Name Value
svn:keywords Author Id Revision Date