ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1535
Committed: Fri Dec 31 18:31:56 2010 UTC (14 years, 4 months ago) by gezelter
File size: 39948 byte(s)
Log Message:
Well, it compiles and builds, but still has a bus error at runtime.

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

Properties

Name Value
svn:keywords Author Id Revision Date