ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.hpp
Revision: 1549
Committed: Wed Apr 27 18:38:15 2011 UTC (14 years ago) by gezelter
File size: 22793 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.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 246 //another nonsense macro declaration
68 gezelter 1390 #define __OPENMD_C
69 tim 3 #include "brains/fSimulation.h"
70 gezelter 2
71 gezelter 1528 using namespace std;
72 gezelter 1390 namespace OpenMD{
73 gezelter 507 //forward decalration
74     class SnapshotManager;
75     class Molecule;
76     class SelectionManager;
77 tim 1024 class StuntDouble;
78 gezelter 1528
79 gezelter 507 /**
80 gezelter 1528 * @class SimInfo SimInfo.hpp "brains/SimInfo.hpp"
81     *
82     * @brief One of the heavy-weight classes of OpenMD, SimInfo
83     * maintains objects and variables relating to the current
84     * simulation. This includes the master list of Molecules. The
85     * Molecule class maintains all of the concrete objects (Atoms,
86     * Bond, Bend, Torsions, Inversions, RigidBodies, CutoffGroups,
87     * Constraints). In both the single and parallel versions, Atoms and
88     * RigidBodies have both global and local indices.
89     */
90 gezelter 507 class SimInfo {
91     public:
92 gezelter 1528 typedef map<int, Molecule*>::iterator MoleculeIterator;
93    
94 gezelter 507 /**
95     * Constructor of SimInfo
96 gezelter 1528 *
97     * @param molStampPairs MoleculeStamp Array. The first element of
98     * the pair is molecule stamp, the second element is the total
99     * number of molecules with the same molecule stamp in the system
100     *
101 gezelter 507 * @param ff pointer of a concrete ForceField instance
102 gezelter 1528 *
103 gezelter 507 * @param simParams
104     */
105 tim 770 SimInfo(ForceField* ff, Globals* simParams);
106 gezelter 507 virtual ~SimInfo();
107 gezelter 2
108 gezelter 507 /**
109     * Adds a molecule
110 gezelter 1528 *
111     * @return return true if adding successfully, return false if the
112     * molecule is already in SimInfo
113     *
114 gezelter 507 * @param mol molecule to be added
115     */
116     bool addMolecule(Molecule* mol);
117 gezelter 2
118 gezelter 507 /**
119     * Removes a molecule from SimInfo
120 gezelter 1528 *
121     * @return true if removing successfully, return false if molecule
122     * is not in this SimInfo
123 gezelter 507 */
124     bool removeMolecule(Molecule* mol);
125 gezelter 2
126 gezelter 507 /** Returns the total number of molecules in the system. */
127     int getNGlobalMolecules() {
128     return nGlobalMols_;
129     }
130 gezelter 2
131 gezelter 507 /** Returns the total number of atoms in the system. */
132     int getNGlobalAtoms() {
133     return nGlobalAtoms_;
134     }
135 gezelter 2
136 gezelter 507 /** Returns the total number of cutoff groups in the system. */
137     int getNGlobalCutoffGroups() {
138     return nGlobalCutoffGroups_;
139     }
140 gezelter 2
141 gezelter 507 /**
142 gezelter 1528 * Returns the total number of integrable objects (total number of
143     * rigid bodies plus the total number of atoms which do not belong
144     * to the rigid bodies) in the system
145 gezelter 507 */
146     int getNGlobalIntegrableObjects() {
147     return nGlobalIntegrableObjects_;
148     }
149 gezelter 2
150 gezelter 507 /**
151 gezelter 1528 * Returns the total number of integrable objects (total number of
152     * rigid bodies plus the total number of atoms which do not belong
153     * to the rigid bodies) in the system
154 gezelter 507 */
155     int getNGlobalRigidBodies() {
156     return nGlobalRigidBodies_;
157     }
158 gezelter 2
159 gezelter 507 int getNGlobalConstraints();
160     /**
161     * Returns the number of local molecules.
162     * @return the number of local molecules
163     */
164     int getNMolecules() {
165     return molecules_.size();
166     }
167 gezelter 2
168 gezelter 507 /** Returns the number of local atoms */
169     unsigned int getNAtoms() {
170     return nAtoms_;
171     }
172 gezelter 2
173 gezelter 507 /** Returns the number of local bonds */
174     unsigned int getNBonds(){
175     return nBonds_;
176     }
177 gezelter 2
178 gezelter 507 /** Returns the number of local bends */
179     unsigned int getNBends() {
180     return nBends_;
181     }
182 gezelter 2
183 gezelter 507 /** Returns the number of local torsions */
184     unsigned int getNTorsions() {
185     return nTorsions_;
186     }
187 gezelter 2
188 gezelter 1277 /** Returns the number of local torsions */
189     unsigned int getNInversions() {
190     return nInversions_;
191     }
192 gezelter 507 /** Returns the number of local rigid bodies */
193     unsigned int getNRigidBodies() {
194     return nRigidBodies_;
195     }
196 gezelter 2
197 gezelter 507 /** Returns the number of local integrable objects */
198     unsigned int getNIntegrableObjects() {
199     return nIntegrableObjects_;
200     }
201 gezelter 2
202 gezelter 507 /** Returns the number of local cutoff groups */
203     unsigned int getNCutoffGroups() {
204     return nCutoffGroups_;
205     }
206 gezelter 2
207 gezelter 507 /** Returns the total number of constraints in this SimInfo */
208     unsigned int getNConstraints() {
209     return nConstraints_;
210     }
211 gezelter 246
212 gezelter 507 /**
213     * Returns the first molecule in this SimInfo and intialize the iterator.
214     * @return the first molecule, return NULL if there is not molecule in this SimInfo
215     * @param i the iterator of molecule array (user shouldn't change it)
216     */
217     Molecule* beginMolecule(MoleculeIterator& i);
218 gezelter 2
219 gezelter 507 /**
220     * Returns the next avaliable Molecule based on the iterator.
221     * @return the next avaliable molecule, return NULL if reaching the end of the array
222     * @param i the iterator of molecule array
223     */
224     Molecule* nextMolecule(MoleculeIterator& i);
225 gezelter 2
226 gezelter 507 /** Returns the number of degrees of freedom */
227     int getNdf() {
228 gezelter 945 return ndf_ - getFdf();
229 gezelter 507 }
230 gezelter 2
231 gezelter 507 /** Returns the number of raw degrees of freedom */
232     int getNdfRaw() {
233     return ndfRaw_;
234     }
235 gezelter 2
236 gezelter 507 /** Returns the number of translational degrees of freedom */
237     int getNdfTrans() {
238     return ndfTrans_;
239     }
240 gezelter 2
241 gezelter 945 /** sets the current number of frozen degrees of freedom */
242     void setFdf(int fdf) {
243     fdf_local = fdf;
244     }
245    
246     int getFdf();
247    
248 gezelter 1528 //getNZconstraint and setNZconstraint ruin the coherence of
249     //SimInfo class, need refactoring
250 gezelter 246
251 gezelter 507 /** Returns the total number of z-constraint molecules in the system */
252     int getNZconstraint() {
253     return nZconstraint_;
254     }
255 gezelter 2
256 gezelter 507 /**
257     * Sets the number of z-constraint molecules in the system.
258     */
259     void setNZconstraint(int nZconstraint) {
260     nZconstraint_ = nZconstraint;
261     }
262 gezelter 246
263 gezelter 507 /** Returns the snapshot manager. */
264     SnapshotManager* getSnapshotManager() {
265     return sman_;
266     }
267 gezelter 2
268 gezelter 507 /** Sets the snapshot manager. */
269     void setSnapshotManager(SnapshotManager* sman);
270 gezelter 246
271 gezelter 507 /** Returns the force field */
272     ForceField* getForceField() {
273     return forceField_;
274     }
275 gezelter 2
276 gezelter 507 Globals* getSimParams() {
277     return simParams_;
278     }
279 gezelter 2
280 gezelter 507 /** Returns the velocity of center of mass of the whole system.*/
281     Vector3d getComVel();
282 gezelter 2
283 gezelter 507 /** Returns the center of the mass of the whole system.*/
284     Vector3d getCom();
285 gezelter 1528 /** Returns the center of the mass and Center of Mass velocity of
286     the whole system.*/
287 chuckv 555 void getComAll(Vector3d& com,Vector3d& comVel);
288 gezelter 2
289 gezelter 1528 /** Returns intertia tensor for the entire system and system
290     Angular Momentum.*/
291 chuckv 555 void getInertiaTensor(Mat3x3d &intertiaTensor,Vector3d &angularMomentum);
292    
293     /** Returns system angular momentum */
294     Vector3d getAngularMomentum();
295    
296 gezelter 1528 /** Returns volume of system as estimated by an ellipsoid defined
297     by the radii of gyration*/
298 chuckv 1103 void getGyrationalVolume(RealType &vol);
299 gezelter 1528 /** Overloaded version of gyrational volume that also returns
300     det(I) so dV/dr can be calculated*/
301 chuckv 1103 void getGyrationalVolume(RealType &vol, RealType &detI);
302 gezelter 1535
303 gezelter 507 void update();
304 gezelter 1535 /**
305     * Setup Fortran Simulation
306     */
307     void setupFortran();
308 gezelter 2
309 gezelter 1535
310 gezelter 507 /** Returns the local index manager */
311     LocalIndexManager* getLocalIndexManager() {
312     return &localIndexMan_;
313     }
314 gezelter 2
315 gezelter 507 int getMoleculeStampId(int globalIndex) {
316     //assert(globalIndex < molStampIds_.size())
317     return molStampIds_[globalIndex];
318     }
319 gezelter 2
320 gezelter 507 /** Returns the molecule stamp */
321     MoleculeStamp* getMoleculeStamp(int id) {
322     return moleculeStamps_[id];
323     }
324 gezelter 2
325 gezelter 507 /** Return the total number of the molecule stamps */
326     int getNMoleculeStamp() {
327     return moleculeStamps_.size();
328     }
329     /**
330     * Finds a molecule with a specified global index
331     * @return a pointer point to found molecule
332     * @param index
333     */
334     Molecule* getMoleculeByGlobalIndex(int index) {
335     MoleculeIterator i;
336     i = molecules_.find(index);
337 gezelter 2
338 gezelter 507 return i != molecules_.end() ? i->second : NULL;
339     }
340 gezelter 2
341 chuckv 1292 int getGlobalMolMembership(int id){
342     return globalMolMembership_[id];
343     }
344 gezelter 1549
345     /**
346     * returns a vector which maps the local atom index on this
347     * processor to the global atom index. With only one processor,
348     * these should be identical.
349     */
350     vector<int> getGlobalAtomIndices();
351    
352     /**
353     * returns a vector which maps the local cutoff group index on
354     * this processor to the global cutoff group index. With only one
355     * processor, these should be identical.
356     */
357     vector<int> getGlobalGroupIndices();
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 507 bool isFortranInitialized() {
418     return fortranInitialized_;
419     }
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     simtype fInfo_; /**< A dual struct shared by C++
558     and Fortran to pass
559     information about what types
560     of calculation are
561     required */
562 gezelter 1535
563 gezelter 1528 /// Stamps are templates for objects that are then used to create
564     /// groups of objects. For example, a molecule stamp contains
565     /// information on how to build that molecule (i.e. the topology,
566     /// the atoms, the bonds, etc.) Once the system is built, the
567     /// stamps are no longer useful.
568     vector<int> molStampIds_; /**< stamp id for molecules in the system */
569     vector<MoleculeStamp*> moleculeStamps_; /**< molecule stamps array */
570    
571     /**
572     * A vector that maps between the global index of an atom, and the
573     * global index of cutoff group the atom belong to. It is filled
574     * by SimCreator once and only once, since it never changed during
575     * the simulation. It should be nGlobalAtoms_ in size.
576     */
577     vector<int> globalGroupMembership_;
578 gezelter 1547 public:
579     vector<int> getGlobalGroupMembership() { return globalGroupMembership_; }
580     private:
581 gezelter 1528
582     /**
583     * A vector that maps between the global index of an atom and the
584     * global index of the molecule the atom belongs to. It is filled
585     * by SimCreator once and only once, since it is never changed
586     * during the simulation. It shoudl be nGlobalAtoms_ in size.
587     */
588 gezelter 1544 vector<int> globalMolMembership_;
589    
590     /**
591     * A vector that maps between the local index of an atom and the
592     * index of the AtomType.
593     */
594     vector<int> identArray_;
595 gezelter 1545 public:
596 gezelter 1544 vector<int> getIdentArray() { return identArray_; }
597 gezelter 1545 private:
598 gezelter 1528
599     /// lists to handle atoms needing special treatment in the non-bonded interactions
600     PairList excludedInteractions_; /**< atoms excluded from interacting with each other */
601     PairList oneTwoInteractions_; /**< atoms that are directly Bonded */
602     PairList oneThreeInteractions_; /**< atoms sharing a Bend */
603     PairList oneFourInteractions_; /**< atoms sharing a Torsion */
604    
605     PropertyMap properties_; /**< Generic Properties can be added */
606     SnapshotManager* sman_; /**< SnapshotManager (handles particle positions, etc.) */
607    
608 gezelter 507 /**
609 gezelter 1528 * The reason to have a local index manager is that when molecule
610     * is migrating to other processors, the atoms and the
611     * rigid-bodies will release their local indices to
612     * LocalIndexManager. Combining the information of molecule
613     * migrating to current processor, Migrator class can query the
614     * LocalIndexManager to make a efficient data moving plan.
615 gezelter 507 */
616     LocalIndexManager localIndexMan_;
617 gezelter 246
618 tim 1024 // unparsed MetaData block for storing in Dump and EOR files:
619 gezelter 1528 string rawMetaData_;
620 tim 1024
621 gezelter 1528 // file names
622     string finalConfigFileName_;
623     string dumpFileName_;
624     string statFileName_;
625     string restFileName_;
626 chrisfen 417
627 gezelter 246
628 gezelter 1467 bool fortranInitialized_; /** flag to indicate whether the fortran side is initialized */
629 gezelter 1126
630     bool calcBoxDipole_; /**< flag to indicate whether or not we calculate
631     the simulation box dipole moment */
632    
633     bool useAtomicVirial_; /**< flag to indicate whether or not we use
634     Atomic Virials to calculate the pressure */
635 gezelter 1528
636     public:
637     /**
638     * return an integral objects by its global index. In MPI
639     * version, if the StuntDouble with specified global index does
640     * not belong to local processor, a NULL will be return.
641 tim 1024 */
642 gezelter 1528 StuntDouble* getIOIndexToIntegrableObject(int index);
643     void setIOIndexToIntegrableObject(const vector<StuntDouble*>& v);
644 tim 1024
645 gezelter 1528 private:
646     vector<StuntDouble*> IOIndexToIntegrableObject;
647    
648 gezelter 507 public:
649 gezelter 246
650 gezelter 507 /**
651     * Finds the processor where a molecule resides
652     * @return the id of the processor which contains the molecule
653     * @param globalIndex global Index of the molecule
654     */
655     int getMolToProc(int globalIndex) {
656     //assert(globalIndex < molToProcMap_.size());
657     return molToProcMap_[globalIndex];
658     }
659 gezelter 1528
660 gezelter 507 /**
661     * Set MolToProcMap array
662     * @see #SimCreator::divideMolecules
663     */
664 gezelter 1528 void setMolToProcMap(const vector<int>& molToProcMap) {
665 gezelter 507 molToProcMap_ = molToProcMap;
666     }
667 gezelter 246
668 gezelter 507 private:
669 gezelter 246
670 gezelter 507 /**
671 gezelter 1241 * The size of molToProcMap_ is equal to total number of molecules
672     * in the system. It maps a molecule to the processor on which it
673     * resides. it is filled by SimCreator once and only once.
674 gezelter 507 */
675 gezelter 1528 vector<int> molToProcMap_;
676 tim 292
677 gezelter 507 };
678 gezelter 2
679 gezelter 1390 } //namespace OpenMD
680 gezelter 246 #endif //BRAINS_SIMMODEL_HPP
681 gezelter 2

Properties

Name Value
svn:keywords Author Id Revision Date