ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.hpp
Revision: 1587
Committed: Fri Jul 8 20:25:32 2011 UTC (13 years, 9 months ago) by gezelter
File size: 23355 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.hpp
44     * @author tlin
45     * @date 11/02/2004
46     * @version 1.0
47     */
48 gezelter 2
49 gezelter 246 #ifndef BRAINS_SIMMODEL_HPP
50     #define BRAINS_SIMMODEL_HPP
51    
52     #include <iostream>
53     #include <set>
54     #include <utility>
55 gezelter 2 #include <vector>
56    
57 gezelter 1287 #include "brains/PairList.hpp"
58 gezelter 246 #include "io/Globals.hpp"
59     #include "math/Vector3.hpp"
60 chuckv 555 #include "math/SquareMatrix3.hpp"
61 gezelter 246 #include "types/MoleculeStamp.hpp"
62     #include "UseTheForce/ForceField.hpp"
63     #include "utils/PropertyMap.hpp"
64     #include "utils/LocalIndexManager.hpp"
65 gezelter 1530 #include "nonbonded/SwitchingFunction.hpp"
66 tim 316
67 gezelter 1528 using namespace std;
68 gezelter 1390 namespace OpenMD{
69 gezelter 1553 //forward declaration
70 gezelter 507 class SnapshotManager;
71     class Molecule;
72     class SelectionManager;
73 tim 1024 class StuntDouble;
74 gezelter 1528
75 gezelter 507 /**
76 gezelter 1528 * @class SimInfo SimInfo.hpp "brains/SimInfo.hpp"
77     *
78     * @brief One of the heavy-weight classes of OpenMD, SimInfo
79     * maintains objects and variables relating to the current
80     * simulation. This includes the master list of Molecules. The
81     * Molecule class maintains all of the concrete objects (Atoms,
82     * Bond, Bend, Torsions, Inversions, RigidBodies, CutoffGroups,
83     * Constraints). In both the single and parallel versions, Atoms and
84     * RigidBodies have both global and local indices.
85     */
86 gezelter 507 class SimInfo {
87     public:
88 gezelter 1528 typedef map<int, Molecule*>::iterator MoleculeIterator;
89    
90 gezelter 507 /**
91     * Constructor of SimInfo
92 gezelter 1528 *
93     * @param molStampPairs MoleculeStamp Array. The first element of
94     * the pair is molecule stamp, the second element is the total
95     * number of molecules with the same molecule stamp in the system
96     *
97 gezelter 507 * @param ff pointer of a concrete ForceField instance
98 gezelter 1528 *
99 gezelter 507 * @param simParams
100     */
101 tim 770 SimInfo(ForceField* ff, Globals* simParams);
102 gezelter 507 virtual ~SimInfo();
103 gezelter 2
104 gezelter 507 /**
105     * Adds a molecule
106 gezelter 1528 *
107     * @return return true if adding successfully, return false if the
108     * molecule is already in SimInfo
109     *
110 gezelter 507 * @param mol molecule to be added
111     */
112     bool addMolecule(Molecule* mol);
113 gezelter 2
114 gezelter 507 /**
115     * Removes a molecule from SimInfo
116 gezelter 1528 *
117     * @return true if removing successfully, return false if molecule
118     * is not in this SimInfo
119 gezelter 507 */
120     bool removeMolecule(Molecule* mol);
121 gezelter 2
122 gezelter 507 /** Returns the total number of molecules in the system. */
123     int getNGlobalMolecules() {
124     return nGlobalMols_;
125     }
126 gezelter 2
127 gezelter 507 /** Returns the total number of atoms in the system. */
128     int getNGlobalAtoms() {
129     return nGlobalAtoms_;
130     }
131 gezelter 2
132 gezelter 507 /** Returns the total number of cutoff groups in the system. */
133     int getNGlobalCutoffGroups() {
134     return nGlobalCutoffGroups_;
135     }
136 gezelter 2
137 gezelter 507 /**
138 gezelter 1528 * Returns the total number of integrable objects (total number of
139     * rigid bodies plus the total number of atoms which do not belong
140     * to the rigid bodies) in the system
141 gezelter 507 */
142     int getNGlobalIntegrableObjects() {
143     return nGlobalIntegrableObjects_;
144     }
145 gezelter 2
146 gezelter 507 /**
147 gezelter 1528 * Returns the total number of integrable objects (total number of
148     * rigid bodies plus the total number of atoms which do not belong
149     * to the rigid bodies) in the system
150 gezelter 507 */
151     int getNGlobalRigidBodies() {
152     return nGlobalRigidBodies_;
153     }
154 gezelter 2
155 gezelter 507 int getNGlobalConstraints();
156     /**
157     * Returns the number of local molecules.
158     * @return the number of local molecules
159     */
160     int getNMolecules() {
161     return molecules_.size();
162     }
163 gezelter 2
164 gezelter 507 /** Returns the number of local atoms */
165     unsigned int getNAtoms() {
166     return nAtoms_;
167     }
168 gezelter 2
169 gezelter 1577 /** Returns the number of effective cutoff groups on local processor */
170     unsigned int getNLocalCutoffGroups();
171    
172 gezelter 507 /** Returns the number of local bonds */
173     unsigned int getNBonds(){
174     return nBonds_;
175     }
176 gezelter 2
177 gezelter 507 /** Returns the number of local bends */
178     unsigned int getNBends() {
179     return nBends_;
180     }
181 gezelter 2
182 gezelter 507 /** Returns the number of local torsions */
183     unsigned int getNTorsions() {
184     return nTorsions_;
185     }
186 gezelter 2
187 gezelter 1277 /** Returns the number of local torsions */
188     unsigned int getNInversions() {
189     return nInversions_;
190     }
191 gezelter 507 /** Returns the number of local rigid bodies */
192     unsigned int getNRigidBodies() {
193     return nRigidBodies_;
194     }
195 gezelter 2
196 gezelter 507 /** Returns the number of local integrable objects */
197     unsigned int getNIntegrableObjects() {
198     return nIntegrableObjects_;
199     }
200 gezelter 2
201 gezelter 507 /** Returns the number of local cutoff groups */
202     unsigned int getNCutoffGroups() {
203     return nCutoffGroups_;
204     }
205 gezelter 2
206 gezelter 507 /** Returns the total number of constraints in this SimInfo */
207     unsigned int getNConstraints() {
208     return nConstraints_;
209     }
210 gezelter 246
211 gezelter 507 /**
212     * Returns the first molecule in this SimInfo and intialize the iterator.
213     * @return the first molecule, return NULL if there is not molecule in this SimInfo
214     * @param i the iterator of molecule array (user shouldn't change it)
215     */
216     Molecule* beginMolecule(MoleculeIterator& i);
217 gezelter 2
218 gezelter 507 /**
219     * Returns the next avaliable Molecule based on the iterator.
220     * @return the next avaliable molecule, return NULL if reaching the end of the array
221     * @param i the iterator of molecule array
222     */
223     Molecule* nextMolecule(MoleculeIterator& i);
224 gezelter 2
225 gezelter 507 /** Returns the number of degrees of freedom */
226     int getNdf() {
227 gezelter 945 return ndf_ - getFdf();
228 gezelter 507 }
229 gezelter 2
230 gezelter 507 /** Returns the number of raw degrees of freedom */
231     int getNdfRaw() {
232     return ndfRaw_;
233     }
234 gezelter 2
235 gezelter 507 /** Returns the number of translational degrees of freedom */
236     int getNdfTrans() {
237     return ndfTrans_;
238     }
239 gezelter 2
240 gezelter 945 /** sets the current number of frozen degrees of freedom */
241     void setFdf(int fdf) {
242     fdf_local = fdf;
243     }
244    
245     int getFdf();
246    
247 gezelter 1528 //getNZconstraint and setNZconstraint ruin the coherence of
248     //SimInfo class, need refactoring
249 gezelter 246
250 gezelter 507 /** Returns the total number of z-constraint molecules in the system */
251     int getNZconstraint() {
252     return nZconstraint_;
253     }
254 gezelter 2
255 gezelter 507 /**
256     * Sets the number of z-constraint molecules in the system.
257     */
258     void setNZconstraint(int nZconstraint) {
259     nZconstraint_ = nZconstraint;
260     }
261 gezelter 246
262 gezelter 507 /** Returns the snapshot manager. */
263     SnapshotManager* getSnapshotManager() {
264     return sman_;
265     }
266 gezelter 2
267 gezelter 507 /** Sets the snapshot manager. */
268     void setSnapshotManager(SnapshotManager* sman);
269 gezelter 246
270 gezelter 507 /** Returns the force field */
271     ForceField* getForceField() {
272     return forceField_;
273     }
274 gezelter 2
275 gezelter 507 Globals* getSimParams() {
276     return simParams_;
277     }
278 gezelter 2
279 gezelter 507 /** Returns the velocity of center of mass of the whole system.*/
280     Vector3d getComVel();
281 gezelter 2
282 gezelter 507 /** Returns the center of the mass of the whole system.*/
283     Vector3d getCom();
284 gezelter 1528 /** Returns the center of the mass and Center of Mass velocity of
285     the whole system.*/
286 chuckv 555 void getComAll(Vector3d& com,Vector3d& comVel);
287 gezelter 2
288 gezelter 1528 /** Returns intertia tensor for the entire system and system
289     Angular Momentum.*/
290 chuckv 555 void getInertiaTensor(Mat3x3d &intertiaTensor,Vector3d &angularMomentum);
291    
292     /** Returns system angular momentum */
293     Vector3d getAngularMomentum();
294    
295 gezelter 1528 /** Returns volume of system as estimated by an ellipsoid defined
296     by the radii of gyration*/
297 chuckv 1103 void getGyrationalVolume(RealType &vol);
298 gezelter 1528 /** Overloaded version of gyrational volume that also returns
299     det(I) so dV/dr can be calculated*/
300 chuckv 1103 void getGyrationalVolume(RealType &vol, RealType &detI);
301 gezelter 1535
302 gezelter 507 void update();
303 gezelter 1535 /**
304 gezelter 1569 * Do final bookkeeping before Force managers need their data.
305 gezelter 1535 */
306 gezelter 1569 void prepareTopology();
307 gezelter 2
308 gezelter 1535
309 gezelter 507 /** Returns the local index manager */
310     LocalIndexManager* getLocalIndexManager() {
311     return &localIndexMan_;
312     }
313 gezelter 2
314 gezelter 507 int getMoleculeStampId(int globalIndex) {
315     //assert(globalIndex < molStampIds_.size())
316     return molStampIds_[globalIndex];
317     }
318 gezelter 2
319 gezelter 507 /** Returns the molecule stamp */
320     MoleculeStamp* getMoleculeStamp(int id) {
321     return moleculeStamps_[id];
322     }
323 gezelter 2
324 gezelter 507 /** Return the total number of the molecule stamps */
325     int getNMoleculeStamp() {
326     return moleculeStamps_.size();
327     }
328     /**
329     * Finds a molecule with a specified global index
330     * @return a pointer point to found molecule
331     * @param index
332     */
333     Molecule* getMoleculeByGlobalIndex(int index) {
334     MoleculeIterator i;
335     i = molecules_.find(index);
336 gezelter 2
337 gezelter 507 return i != molecules_.end() ? i->second : NULL;
338     }
339 gezelter 2
340 chuckv 1292 int getGlobalMolMembership(int id){
341     return globalMolMembership_[id];
342     }
343 gezelter 1549
344     /**
345     * returns a vector which maps the local atom index on this
346     * processor to the global atom index. With only one processor,
347     * these should be identical.
348     */
349     vector<int> getGlobalAtomIndices();
350    
351     /**
352     * returns a vector which maps the local cutoff group index on
353     * this processor to the global cutoff group index. With only one
354     * processor, these should be identical.
355     */
356     vector<int> getGlobalGroupIndices();
357 gezelter 1569
358 gezelter 246
359 gezelter 1528 string getFinalConfigFileName() {
360 gezelter 507 return finalConfigFileName_;
361     }
362 tim 1024
363 gezelter 1528 void setFinalConfigFileName(const string& fileName) {
364 gezelter 507 finalConfigFileName_ = fileName;
365     }
366 gezelter 2
367 gezelter 1528 string getRawMetaData() {
368 tim 1024 return rawMetaData_;
369     }
370 gezelter 1528 void setRawMetaData(const string& rawMetaData) {
371 tim 1024 rawMetaData_ = rawMetaData;
372     }
373    
374 gezelter 1528 string getDumpFileName() {
375 gezelter 507 return dumpFileName_;
376     }
377 gezelter 246
378 gezelter 1528 void setDumpFileName(const string& fileName) {
379 gezelter 507 dumpFileName_ = fileName;
380     }
381 gezelter 2
382 gezelter 1528 string getStatFileName() {
383 gezelter 507 return statFileName_;
384     }
385 gezelter 246
386 gezelter 1528 void setStatFileName(const string& fileName) {
387 gezelter 507 statFileName_ = fileName;
388     }
389 chrisfen 417
390 gezelter 1528 string getRestFileName() {
391 gezelter 507 return restFileName_;
392     }
393 chrisfen 417
394 gezelter 1528 void setRestFileName(const string& fileName) {
395 gezelter 507 restFileName_ = fileName;
396     }
397 gezelter 2
398 gezelter 507 /**
399     * Sets GlobalGroupMembership
400     * @see #SimCreator::setGlobalIndex
401     */
402 gezelter 1528 void setGlobalGroupMembership(const vector<int>& globalGroupMembership) {
403 gezelter 1287 assert(globalGroupMembership.size() == static_cast<size_t>(nGlobalAtoms_));
404 gezelter 507 globalGroupMembership_ = globalGroupMembership;
405     }
406 gezelter 2
407 gezelter 507 /**
408     * Sets GlobalMolMembership
409     * @see #SimCreator::setGlobalIndex
410     */
411 gezelter 1528 void setGlobalMolMembership(const vector<int>& globalMolMembership) {
412 gezelter 1287 assert(globalMolMembership.size() == static_cast<size_t>(nGlobalAtoms_));
413 gezelter 507 globalMolMembership_ = globalMolMembership;
414     }
415 gezelter 246
416    
417 gezelter 1569 bool isTopologyDone() {
418     return topologyDone_;
419 gezelter 507 }
420 gezelter 246
421 chrisfen 998 bool getCalcBoxDipole() {
422     return calcBoxDipole_;
423     }
424    
425 gezelter 1126 bool getUseAtomicVirial() {
426     return useAtomicVirial_;
427     }
428    
429 gezelter 507 /**
430     * Adds property into property map
431     * @param genData GenericData to be added into PropertyMap
432     */
433     void addProperty(GenericData* genData);
434 gezelter 246
435 gezelter 507 /**
436     * Removes property from PropertyMap by name
437     * @param propName the name of property to be removed
438     */
439 gezelter 1528 void removeProperty(const string& propName);
440 gezelter 246
441 gezelter 507 /**
442     * clear all of the properties
443     */
444     void clearProperties();
445 gezelter 246
446 gezelter 507 /**
447     * Returns all names of properties
448     * @return all names of properties
449     */
450 gezelter 1528 vector<string> getPropertyNames();
451 gezelter 246
452 gezelter 507 /**
453     * Returns all of the properties in PropertyMap
454     * @return all of the properties in PropertyMap
455     */
456 gezelter 1528 vector<GenericData*> getProperties();
457 gezelter 246
458 gezelter 507 /**
459     * Returns property
460     * @param propName name of property
461     * @return a pointer point to property with propName. If no property named propName
462     * exists, return NULL
463     */
464 gezelter 1528 GenericData* getPropertyByName(const string& propName);
465 gezelter 246
466 gezelter 507 /**
467 gezelter 1287 * add all special interaction pairs (including excluded
468     * interactions) in a molecule into the appropriate lists.
469 gezelter 507 */
470 gezelter 1287 void addInteractionPairs(Molecule* mol);
471 gezelter 246
472 gezelter 507 /**
473 gezelter 1287 * remove all special interaction pairs which belong to a molecule
474     * from the appropriate lists.
475 gezelter 507 */
476 gezelter 1287 void removeInteractionPairs(Molecule* mol);
477 gezelter 246
478 gezelter 1528 /** Returns the set of atom types present in this simulation */
479     set<AtomType*> getSimulatedAtomTypes();
480 tim 292
481 gezelter 1528 friend ostream& operator <<(ostream& o, SimInfo& info);
482 tim 326
483 tim 963 void getCutoff(RealType& rcut, RealType& rsw);
484 gezelter 246
485 gezelter 507 private:
486 gezelter 246
487 gezelter 1530 /** fill up the simtype struct and other simulation-related variables */
488     void setupSimVariables();
489 gezelter 246
490    
491 chrisfen 998 /** Determine if we need to accumulate the simulation box dipole */
492     void setupAccumulateBoxDipole();
493    
494 gezelter 507 /** Calculates the number of degress of freedom in the whole system */
495     void calcNdf();
496     void calcNdfRaw();
497     void calcNdfTrans();
498 gezelter 246
499 gezelter 507 /**
500 gezelter 1528 * Adds molecule stamp and the total number of the molecule with
501     * same molecule stamp in the whole system.
502 gezelter 507 */
503     void addMoleculeStamp(MoleculeStamp* molStamp, int nmol);
504 gezelter 246
505 gezelter 1528 // Other classes holdingn important information
506     ForceField* forceField_; /**< provides access to defined atom types, bond types, etc. */
507     Globals* simParams_; /**< provides access to simulation parameters set by user */
508 gezelter 246
509 gezelter 1528 /// Counts of local objects
510 gezelter 1277 int nAtoms_; /**< number of atoms in local processor */
511     int nBonds_; /**< number of bonds in local processor */
512     int nBends_; /**< number of bends in local processor */
513     int nTorsions_; /**< number of torsions in local processor */
514     int nInversions_; /**< number of inversions in local processor */
515     int nRigidBodies_; /**< number of rigid bodies in local processor */
516     int nIntegrableObjects_; /**< number of integrable objects in local processor */
517     int nCutoffGroups_; /**< number of cutoff groups in local processor */
518     int nConstraints_; /**< number of constraints in local processors */
519 gezelter 1528
520     /// Counts of global objects
521     int nGlobalMols_; /**< number of molecules in the system (GLOBAL) */
522     int nGlobalAtoms_; /**< number of atoms in the system (GLOBAL) */
523     int nGlobalCutoffGroups_; /**< number of cutoff groups in this system (GLOBAL) */
524     int nGlobalIntegrableObjects_; /**< number of integrable objects in this system */
525     int nGlobalRigidBodies_; /**< number of rigid bodies in this system (GLOBAL) */
526    
527     /// Degress of freedom
528     int ndf_; /**< number of degress of freedom (excludes constraints) (LOCAL) */
529     int fdf_local; /**< number of frozen degrees of freedom (LOCAL) */
530     int fdf_; /**< number of frozen degrees of freedom (GLOBAL) */
531     int ndfRaw_; /**< number of degress of freedom (includes constraints), (LOCAL) */
532     int ndfTrans_; /**< number of translation degress of freedom, (LOCAL) */
533     int nZconstraint_; /**< number of z-constraint molecules (GLOBAL) */
534 gezelter 246
535 gezelter 1528 /// logicals
536     bool usesPeriodicBoundaries_; /**< use periodic boundary conditions? */
537     bool usesDirectionalAtoms_; /**< are there atoms with position AND orientation? */
538     bool usesMetallicAtoms_; /**< are there transition metal atoms? */
539     bool usesElectrostaticAtoms_; /**< are there electrostatic atoms? */
540     bool usesAtomicVirial_; /**< are we computing atomic virials? */
541     bool requiresPrepair_; /**< does this simulation require a pre-pair loop? */
542     bool requiresSkipCorrection_; /**< does this simulation require a skip-correction? */
543     bool requiresSelfCorrection_; /**< does this simulation require a self-correction? */
544 gezelter 246
545 gezelter 1535 public:
546     bool usesElectrostaticAtoms() { return usesElectrostaticAtoms_; }
547     bool usesDirectionalAtoms() { return usesDirectionalAtoms_; }
548     bool usesMetallicAtoms() { return usesMetallicAtoms_; }
549 gezelter 1546 bool usesAtomicVirial() { return usesAtomicVirial_; }
550     bool requiresPrepair() { return requiresPrepair_; }
551     bool requiresSkipCorrection() { return requiresSkipCorrection_;}
552     bool requiresSelfCorrection() { return requiresSelfCorrection_;}
553 gezelter 1535
554     private:
555 gezelter 1528 /// Data structures holding primary simulation objects
556     map<int, Molecule*> molecules_; /**< map holding pointers to LOCAL molecules */
557 gezelter 1535
558 gezelter 1528 /// Stamps are templates for objects that are then used to create
559     /// groups of objects. For example, a molecule stamp contains
560     /// information on how to build that molecule (i.e. the topology,
561     /// the atoms, the bonds, etc.) Once the system is built, the
562     /// stamps are no longer useful.
563     vector<int> molStampIds_; /**< stamp id for molecules in the system */
564     vector<MoleculeStamp*> moleculeStamps_; /**< molecule stamps array */
565    
566     /**
567     * A vector that maps between the global index of an atom, and the
568     * global index of cutoff group the atom belong to. It is filled
569     * by SimCreator once and only once, since it never changed during
570     * the simulation. It should be nGlobalAtoms_ in size.
571     */
572     vector<int> globalGroupMembership_;
573 gezelter 1547 public:
574     vector<int> getGlobalGroupMembership() { return globalGroupMembership_; }
575     private:
576 gezelter 1528
577     /**
578     * A vector that maps between the global index of an atom and the
579     * global index of the molecule the atom belongs to. It is filled
580     * by SimCreator once and only once, since it is never changed
581     * during the simulation. It shoudl be nGlobalAtoms_ in size.
582     */
583 gezelter 1544 vector<int> globalMolMembership_;
584    
585     /**
586     * A vector that maps between the local index of an atom and the
587     * index of the AtomType.
588     */
589     vector<int> identArray_;
590 gezelter 1545 public:
591 gezelter 1544 vector<int> getIdentArray() { return identArray_; }
592 gezelter 1545 private:
593 gezelter 1569
594     /**
595     * A vector which contains the fractional contribution of an
596     * atom's mass to the total mass of the cutoffGroup that atom
597     * belongs to. In the case of single atom cutoff groups, the mass
598     * factor for that atom is 1. For massless atoms, the factor is
599     * also 1.
600     */
601     vector<RealType> massFactors_;
602     public:
603     vector<RealType> getMassFactors() { return massFactors_; }
604 gezelter 1570
605 gezelter 1587 PairList* getExcludedInteractions() { return &excludedInteractions_; }
606     PairList* getOneTwoInteractions() { return &oneTwoInteractions_; }
607     PairList* getOneThreeInteractions() { return &oneThreeInteractions_; }
608     PairList* getOneFourInteractions() { return &oneFourInteractions_; }
609 gezelter 1570
610 gezelter 1569 private:
611 gezelter 1528
612     /// lists to handle atoms needing special treatment in the non-bonded interactions
613     PairList excludedInteractions_; /**< atoms excluded from interacting with each other */
614     PairList oneTwoInteractions_; /**< atoms that are directly Bonded */
615     PairList oneThreeInteractions_; /**< atoms sharing a Bend */
616     PairList oneFourInteractions_; /**< atoms sharing a Torsion */
617    
618     PropertyMap properties_; /**< Generic Properties can be added */
619     SnapshotManager* sman_; /**< SnapshotManager (handles particle positions, etc.) */
620    
621 gezelter 507 /**
622 gezelter 1528 * The reason to have a local index manager is that when molecule
623     * is migrating to other processors, the atoms and the
624     * rigid-bodies will release their local indices to
625     * LocalIndexManager. Combining the information of molecule
626     * migrating to current processor, Migrator class can query the
627     * LocalIndexManager to make a efficient data moving plan.
628 gezelter 507 */
629     LocalIndexManager localIndexMan_;
630 gezelter 246
631 tim 1024 // unparsed MetaData block for storing in Dump and EOR files:
632 gezelter 1528 string rawMetaData_;
633 tim 1024
634 gezelter 1528 // file names
635     string finalConfigFileName_;
636     string dumpFileName_;
637     string statFileName_;
638     string restFileName_;
639 chrisfen 417
640 gezelter 246
641 gezelter 1569 bool topologyDone_; /** flag to indicate whether the topology has
642     been scanned and all the relevant
643     bookkeeping has been done*/
644 gezelter 1126
645     bool calcBoxDipole_; /**< flag to indicate whether or not we calculate
646     the simulation box dipole moment */
647    
648     bool useAtomicVirial_; /**< flag to indicate whether or not we use
649     Atomic Virials to calculate the pressure */
650 gezelter 1528
651     public:
652     /**
653     * return an integral objects by its global index. In MPI
654     * version, if the StuntDouble with specified global index does
655     * not belong to local processor, a NULL will be return.
656 tim 1024 */
657 gezelter 1528 StuntDouble* getIOIndexToIntegrableObject(int index);
658     void setIOIndexToIntegrableObject(const vector<StuntDouble*>& v);
659 tim 1024
660 gezelter 1528 private:
661     vector<StuntDouble*> IOIndexToIntegrableObject;
662    
663 gezelter 507 public:
664 gezelter 246
665 gezelter 507 /**
666     * Finds the processor where a molecule resides
667     * @return the id of the processor which contains the molecule
668     * @param globalIndex global Index of the molecule
669     */
670     int getMolToProc(int globalIndex) {
671     //assert(globalIndex < molToProcMap_.size());
672     return molToProcMap_[globalIndex];
673     }
674 gezelter 1528
675 gezelter 507 /**
676     * Set MolToProcMap array
677     * @see #SimCreator::divideMolecules
678     */
679 gezelter 1528 void setMolToProcMap(const vector<int>& molToProcMap) {
680 gezelter 507 molToProcMap_ = molToProcMap;
681     }
682 gezelter 246
683 gezelter 507 private:
684 gezelter 246
685 gezelter 507 /**
686 gezelter 1241 * The size of molToProcMap_ is equal to total number of molecules
687     * in the system. It maps a molecule to the processor on which it
688     * resides. it is filled by SimCreator once and only once.
689 gezelter 507 */
690 gezelter 1528 vector<int> molToProcMap_;
691 tim 292
692 gezelter 507 };
693 gezelter 2
694 gezelter 1390 } //namespace OpenMD
695 gezelter 246 #endif //BRAINS_SIMMODEL_HPP
696 gezelter 2

Properties

Name Value
svn:keywords Author Id Revision Date