ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1549
Committed: Wed Apr 27 18:38:15 2011 UTC (14 years ago) by gezelter
File size: 40769 byte(s)
Log Message:
a few more tweaks   We're getting somewhat closer to deleting fortran.

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

Properties

Name Value
svn:keywords Author Id Revision Date