ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1529
Committed: Mon Dec 27 18:35:59 2010 UTC (14 years, 4 months ago) by gezelter
File size: 43121 byte(s)
Log Message:
fixes

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date