ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.hpp
Revision: 1553
Committed: Fri Apr 29 17:25:12 2011 UTC (14 years ago) by gezelter
File size: 22389 byte(s)
Log Message:
More fortran removal, and trimming out unused code

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 507 /** Returns the number of local bonds */
170     unsigned int getNBonds(){
171     return nBonds_;
172     }
173 gezelter 2
174 gezelter 507 /** Returns the number of local bends */
175     unsigned int getNBends() {
176     return nBends_;
177     }
178 gezelter 2
179 gezelter 507 /** Returns the number of local torsions */
180     unsigned int getNTorsions() {
181     return nTorsions_;
182     }
183 gezelter 2
184 gezelter 1277 /** Returns the number of local torsions */
185     unsigned int getNInversions() {
186     return nInversions_;
187     }
188 gezelter 507 /** Returns the number of local rigid bodies */
189     unsigned int getNRigidBodies() {
190     return nRigidBodies_;
191     }
192 gezelter 2
193 gezelter 507 /** Returns the number of local integrable objects */
194     unsigned int getNIntegrableObjects() {
195     return nIntegrableObjects_;
196     }
197 gezelter 2
198 gezelter 507 /** Returns the number of local cutoff groups */
199     unsigned int getNCutoffGroups() {
200     return nCutoffGroups_;
201     }
202 gezelter 2
203 gezelter 507 /** Returns the total number of constraints in this SimInfo */
204     unsigned int getNConstraints() {
205     return nConstraints_;
206     }
207 gezelter 246
208 gezelter 507 /**
209     * Returns the first molecule in this SimInfo and intialize the iterator.
210     * @return the first molecule, return NULL if there is not molecule in this SimInfo
211     * @param i the iterator of molecule array (user shouldn't change it)
212     */
213     Molecule* beginMolecule(MoleculeIterator& i);
214 gezelter 2
215 gezelter 507 /**
216     * Returns the next avaliable Molecule based on the iterator.
217     * @return the next avaliable molecule, return NULL if reaching the end of the array
218     * @param i the iterator of molecule array
219     */
220     Molecule* nextMolecule(MoleculeIterator& i);
221 gezelter 2
222 gezelter 507 /** Returns the number of degrees of freedom */
223     int getNdf() {
224 gezelter 945 return ndf_ - getFdf();
225 gezelter 507 }
226 gezelter 2
227 gezelter 507 /** Returns the number of raw degrees of freedom */
228     int getNdfRaw() {
229     return ndfRaw_;
230     }
231 gezelter 2
232 gezelter 507 /** Returns the number of translational degrees of freedom */
233     int getNdfTrans() {
234     return ndfTrans_;
235     }
236 gezelter 2
237 gezelter 945 /** sets the current number of frozen degrees of freedom */
238     void setFdf(int fdf) {
239     fdf_local = fdf;
240     }
241    
242     int getFdf();
243    
244 gezelter 1528 //getNZconstraint and setNZconstraint ruin the coherence of
245     //SimInfo class, need refactoring
246 gezelter 246
247 gezelter 507 /** Returns the total number of z-constraint molecules in the system */
248     int getNZconstraint() {
249     return nZconstraint_;
250     }
251 gezelter 2
252 gezelter 507 /**
253     * Sets the number of z-constraint molecules in the system.
254     */
255     void setNZconstraint(int nZconstraint) {
256     nZconstraint_ = nZconstraint;
257     }
258 gezelter 246
259 gezelter 507 /** Returns the snapshot manager. */
260     SnapshotManager* getSnapshotManager() {
261     return sman_;
262     }
263 gezelter 2
264 gezelter 507 /** Sets the snapshot manager. */
265     void setSnapshotManager(SnapshotManager* sman);
266 gezelter 246
267 gezelter 507 /** Returns the force field */
268     ForceField* getForceField() {
269     return forceField_;
270     }
271 gezelter 2
272 gezelter 507 Globals* getSimParams() {
273     return simParams_;
274     }
275 gezelter 2
276 gezelter 507 /** Returns the velocity of center of mass of the whole system.*/
277     Vector3d getComVel();
278 gezelter 2
279 gezelter 507 /** Returns the center of the mass of the whole system.*/
280     Vector3d getCom();
281 gezelter 1528 /** Returns the center of the mass and Center of Mass velocity of
282     the whole system.*/
283 chuckv 555 void getComAll(Vector3d& com,Vector3d& comVel);
284 gezelter 2
285 gezelter 1528 /** Returns intertia tensor for the entire system and system
286     Angular Momentum.*/
287 chuckv 555 void getInertiaTensor(Mat3x3d &intertiaTensor,Vector3d &angularMomentum);
288    
289     /** Returns system angular momentum */
290     Vector3d getAngularMomentum();
291    
292 gezelter 1528 /** Returns volume of system as estimated by an ellipsoid defined
293     by the radii of gyration*/
294 chuckv 1103 void getGyrationalVolume(RealType &vol);
295 gezelter 1528 /** Overloaded version of gyrational volume that also returns
296     det(I) so dV/dr can be calculated*/
297 chuckv 1103 void getGyrationalVolume(RealType &vol, RealType &detI);
298 gezelter 1535
299 gezelter 507 void update();
300 gezelter 1535 /**
301     * Setup Fortran Simulation
302     */
303     void setupFortran();
304 gezelter 2
305 gezelter 1535
306 gezelter 507 /** Returns the local index manager */
307     LocalIndexManager* getLocalIndexManager() {
308     return &localIndexMan_;
309     }
310 gezelter 2
311 gezelter 507 int getMoleculeStampId(int globalIndex) {
312     //assert(globalIndex < molStampIds_.size())
313     return molStampIds_[globalIndex];
314     }
315 gezelter 2
316 gezelter 507 /** Returns the molecule stamp */
317     MoleculeStamp* getMoleculeStamp(int id) {
318     return moleculeStamps_[id];
319     }
320 gezelter 2
321 gezelter 507 /** Return the total number of the molecule stamps */
322     int getNMoleculeStamp() {
323     return moleculeStamps_.size();
324     }
325     /**
326     * Finds a molecule with a specified global index
327     * @return a pointer point to found molecule
328     * @param index
329     */
330     Molecule* getMoleculeByGlobalIndex(int index) {
331     MoleculeIterator i;
332     i = molecules_.find(index);
333 gezelter 2
334 gezelter 507 return i != molecules_.end() ? i->second : NULL;
335     }
336 gezelter 2
337 chuckv 1292 int getGlobalMolMembership(int id){
338     return globalMolMembership_[id];
339     }
340 gezelter 1549
341     /**
342     * returns a vector which maps the local atom index on this
343     * processor to the global atom index. With only one processor,
344     * these should be identical.
345     */
346     vector<int> getGlobalAtomIndices();
347    
348     /**
349     * returns a vector which maps the local cutoff group index on
350     * this processor to the global cutoff group index. With only one
351     * processor, these should be identical.
352     */
353     vector<int> getGlobalGroupIndices();
354 gezelter 246
355 gezelter 1528 string getFinalConfigFileName() {
356 gezelter 507 return finalConfigFileName_;
357     }
358 tim 1024
359 gezelter 1528 void setFinalConfigFileName(const string& fileName) {
360 gezelter 507 finalConfigFileName_ = fileName;
361     }
362 gezelter 2
363 gezelter 1528 string getRawMetaData() {
364 tim 1024 return rawMetaData_;
365     }
366 gezelter 1528 void setRawMetaData(const string& rawMetaData) {
367 tim 1024 rawMetaData_ = rawMetaData;
368     }
369    
370 gezelter 1528 string getDumpFileName() {
371 gezelter 507 return dumpFileName_;
372     }
373 gezelter 246
374 gezelter 1528 void setDumpFileName(const string& fileName) {
375 gezelter 507 dumpFileName_ = fileName;
376     }
377 gezelter 2
378 gezelter 1528 string getStatFileName() {
379 gezelter 507 return statFileName_;
380     }
381 gezelter 246
382 gezelter 1528 void setStatFileName(const string& fileName) {
383 gezelter 507 statFileName_ = fileName;
384     }
385 chrisfen 417
386 gezelter 1528 string getRestFileName() {
387 gezelter 507 return restFileName_;
388     }
389 chrisfen 417
390 gezelter 1528 void setRestFileName(const string& fileName) {
391 gezelter 507 restFileName_ = fileName;
392     }
393 gezelter 2
394 gezelter 507 /**
395     * Sets GlobalGroupMembership
396     * @see #SimCreator::setGlobalIndex
397     */
398 gezelter 1528 void setGlobalGroupMembership(const vector<int>& globalGroupMembership) {
399 gezelter 1287 assert(globalGroupMembership.size() == static_cast<size_t>(nGlobalAtoms_));
400 gezelter 507 globalGroupMembership_ = globalGroupMembership;
401     }
402 gezelter 2
403 gezelter 507 /**
404     * Sets GlobalMolMembership
405     * @see #SimCreator::setGlobalIndex
406     */
407 gezelter 1528 void setGlobalMolMembership(const vector<int>& globalMolMembership) {
408 gezelter 1287 assert(globalMolMembership.size() == static_cast<size_t>(nGlobalAtoms_));
409 gezelter 507 globalMolMembership_ = globalMolMembership;
410     }
411 gezelter 246
412    
413 gezelter 507 bool isFortranInitialized() {
414     return fortranInitialized_;
415     }
416 gezelter 246
417 chrisfen 998 bool getCalcBoxDipole() {
418     return calcBoxDipole_;
419     }
420    
421 gezelter 1126 bool getUseAtomicVirial() {
422     return useAtomicVirial_;
423     }
424    
425 gezelter 507 /**
426     * Adds property into property map
427     * @param genData GenericData to be added into PropertyMap
428     */
429     void addProperty(GenericData* genData);
430 gezelter 246
431 gezelter 507 /**
432     * Removes property from PropertyMap by name
433     * @param propName the name of property to be removed
434     */
435 gezelter 1528 void removeProperty(const string& propName);
436 gezelter 246
437 gezelter 507 /**
438     * clear all of the properties
439     */
440     void clearProperties();
441 gezelter 246
442 gezelter 507 /**
443     * Returns all names of properties
444     * @return all names of properties
445     */
446 gezelter 1528 vector<string> getPropertyNames();
447 gezelter 246
448 gezelter 507 /**
449     * Returns all of the properties in PropertyMap
450     * @return all of the properties in PropertyMap
451     */
452 gezelter 1528 vector<GenericData*> getProperties();
453 gezelter 246
454 gezelter 507 /**
455     * Returns property
456     * @param propName name of property
457     * @return a pointer point to property with propName. If no property named propName
458     * exists, return NULL
459     */
460 gezelter 1528 GenericData* getPropertyByName(const string& propName);
461 gezelter 246
462 gezelter 507 /**
463 gezelter 1287 * add all special interaction pairs (including excluded
464     * interactions) in a molecule into the appropriate lists.
465 gezelter 507 */
466 gezelter 1287 void addInteractionPairs(Molecule* mol);
467 gezelter 246
468 gezelter 507 /**
469 gezelter 1287 * remove all special interaction pairs which belong to a molecule
470     * from the appropriate lists.
471 gezelter 507 */
472 gezelter 1287 void removeInteractionPairs(Molecule* mol);
473 gezelter 246
474 gezelter 1528 /** Returns the set of atom types present in this simulation */
475     set<AtomType*> getSimulatedAtomTypes();
476 tim 292
477 gezelter 1528 friend ostream& operator <<(ostream& o, SimInfo& info);
478 tim 326
479 tim 963 void getCutoff(RealType& rcut, RealType& rsw);
480 gezelter 246
481 gezelter 507 private:
482 gezelter 246
483 gezelter 1530 /** fill up the simtype struct and other simulation-related variables */
484     void setupSimVariables();
485 gezelter 246
486    
487 chrisfen 998 /** Determine if we need to accumulate the simulation box dipole */
488     void setupAccumulateBoxDipole();
489    
490 gezelter 507 /** Calculates the number of degress of freedom in the whole system */
491     void calcNdf();
492     void calcNdfRaw();
493     void calcNdfTrans();
494 gezelter 246
495 gezelter 507 /**
496 gezelter 1528 * Adds molecule stamp and the total number of the molecule with
497     * same molecule stamp in the whole system.
498 gezelter 507 */
499     void addMoleculeStamp(MoleculeStamp* molStamp, int nmol);
500 gezelter 246
501 gezelter 1528 // Other classes holdingn important information
502     ForceField* forceField_; /**< provides access to defined atom types, bond types, etc. */
503     Globals* simParams_; /**< provides access to simulation parameters set by user */
504 gezelter 246
505 gezelter 1528 /// Counts of local objects
506 gezelter 1277 int nAtoms_; /**< number of atoms in local processor */
507     int nBonds_; /**< number of bonds in local processor */
508     int nBends_; /**< number of bends in local processor */
509     int nTorsions_; /**< number of torsions in local processor */
510     int nInversions_; /**< number of inversions in local processor */
511     int nRigidBodies_; /**< number of rigid bodies in local processor */
512     int nIntegrableObjects_; /**< number of integrable objects in local processor */
513     int nCutoffGroups_; /**< number of cutoff groups in local processor */
514     int nConstraints_; /**< number of constraints in local processors */
515 gezelter 1528
516     /// Counts of global objects
517     int nGlobalMols_; /**< number of molecules in the system (GLOBAL) */
518     int nGlobalAtoms_; /**< number of atoms in the system (GLOBAL) */
519     int nGlobalCutoffGroups_; /**< number of cutoff groups in this system (GLOBAL) */
520     int nGlobalIntegrableObjects_; /**< number of integrable objects in this system */
521     int nGlobalRigidBodies_; /**< number of rigid bodies in this system (GLOBAL) */
522    
523     /// Degress of freedom
524     int ndf_; /**< number of degress of freedom (excludes constraints) (LOCAL) */
525     int fdf_local; /**< number of frozen degrees of freedom (LOCAL) */
526     int fdf_; /**< number of frozen degrees of freedom (GLOBAL) */
527     int ndfRaw_; /**< number of degress of freedom (includes constraints), (LOCAL) */
528     int ndfTrans_; /**< number of translation degress of freedom, (LOCAL) */
529     int nZconstraint_; /**< number of z-constraint molecules (GLOBAL) */
530 gezelter 246
531 gezelter 1528 /// logicals
532     bool usesPeriodicBoundaries_; /**< use periodic boundary conditions? */
533     bool usesDirectionalAtoms_; /**< are there atoms with position AND orientation? */
534     bool usesMetallicAtoms_; /**< are there transition metal atoms? */
535     bool usesElectrostaticAtoms_; /**< are there electrostatic atoms? */
536     bool usesAtomicVirial_; /**< are we computing atomic virials? */
537     bool requiresPrepair_; /**< does this simulation require a pre-pair loop? */
538     bool requiresSkipCorrection_; /**< does this simulation require a skip-correction? */
539     bool requiresSelfCorrection_; /**< does this simulation require a self-correction? */
540 gezelter 246
541 gezelter 1535 public:
542     bool usesElectrostaticAtoms() { return usesElectrostaticAtoms_; }
543     bool usesDirectionalAtoms() { return usesDirectionalAtoms_; }
544     bool usesMetallicAtoms() { return usesMetallicAtoms_; }
545 gezelter 1546 bool usesAtomicVirial() { return usesAtomicVirial_; }
546     bool requiresPrepair() { return requiresPrepair_; }
547     bool requiresSkipCorrection() { return requiresSkipCorrection_;}
548     bool requiresSelfCorrection() { return requiresSelfCorrection_;}
549 gezelter 1535
550     private:
551 gezelter 1528 /// Data structures holding primary simulation objects
552     map<int, Molecule*> molecules_; /**< map holding pointers to LOCAL molecules */
553 gezelter 1535
554 gezelter 1528 /// Stamps are templates for objects that are then used to create
555     /// groups of objects. For example, a molecule stamp contains
556     /// information on how to build that molecule (i.e. the topology,
557     /// the atoms, the bonds, etc.) Once the system is built, the
558     /// stamps are no longer useful.
559     vector<int> molStampIds_; /**< stamp id for molecules in the system */
560     vector<MoleculeStamp*> moleculeStamps_; /**< molecule stamps array */
561    
562     /**
563     * A vector that maps between the global index of an atom, and the
564     * global index of cutoff group the atom belong to. It is filled
565     * by SimCreator once and only once, since it never changed during
566     * the simulation. It should be nGlobalAtoms_ in size.
567     */
568     vector<int> globalGroupMembership_;
569 gezelter 1547 public:
570     vector<int> getGlobalGroupMembership() { return globalGroupMembership_; }
571     private:
572 gezelter 1528
573     /**
574     * A vector that maps between the global index of an atom and the
575     * global index of the molecule the atom belongs to. It is filled
576     * by SimCreator once and only once, since it is never changed
577     * during the simulation. It shoudl be nGlobalAtoms_ in size.
578     */
579 gezelter 1544 vector<int> globalMolMembership_;
580    
581     /**
582     * A vector that maps between the local index of an atom and the
583     * index of the AtomType.
584     */
585     vector<int> identArray_;
586 gezelter 1545 public:
587 gezelter 1544 vector<int> getIdentArray() { return identArray_; }
588 gezelter 1545 private:
589 gezelter 1528
590     /// lists to handle atoms needing special treatment in the non-bonded interactions
591     PairList excludedInteractions_; /**< atoms excluded from interacting with each other */
592     PairList oneTwoInteractions_; /**< atoms that are directly Bonded */
593     PairList oneThreeInteractions_; /**< atoms sharing a Bend */
594     PairList oneFourInteractions_; /**< atoms sharing a Torsion */
595    
596     PropertyMap properties_; /**< Generic Properties can be added */
597     SnapshotManager* sman_; /**< SnapshotManager (handles particle positions, etc.) */
598    
599 gezelter 507 /**
600 gezelter 1528 * The reason to have a local index manager is that when molecule
601     * is migrating to other processors, the atoms and the
602     * rigid-bodies will release their local indices to
603     * LocalIndexManager. Combining the information of molecule
604     * migrating to current processor, Migrator class can query the
605     * LocalIndexManager to make a efficient data moving plan.
606 gezelter 507 */
607     LocalIndexManager localIndexMan_;
608 gezelter 246
609 tim 1024 // unparsed MetaData block for storing in Dump and EOR files:
610 gezelter 1528 string rawMetaData_;
611 tim 1024
612 gezelter 1528 // file names
613     string finalConfigFileName_;
614     string dumpFileName_;
615     string statFileName_;
616     string restFileName_;
617 chrisfen 417
618 gezelter 246
619 gezelter 1467 bool fortranInitialized_; /** flag to indicate whether the fortran side is initialized */
620 gezelter 1126
621     bool calcBoxDipole_; /**< flag to indicate whether or not we calculate
622     the simulation box dipole moment */
623    
624     bool useAtomicVirial_; /**< flag to indicate whether or not we use
625     Atomic Virials to calculate the pressure */
626 gezelter 1528
627     public:
628     /**
629     * return an integral objects by its global index. In MPI
630     * version, if the StuntDouble with specified global index does
631     * not belong to local processor, a NULL will be return.
632 tim 1024 */
633 gezelter 1528 StuntDouble* getIOIndexToIntegrableObject(int index);
634     void setIOIndexToIntegrableObject(const vector<StuntDouble*>& v);
635 tim 1024
636 gezelter 1528 private:
637     vector<StuntDouble*> IOIndexToIntegrableObject;
638    
639 gezelter 507 public:
640 gezelter 246
641 gezelter 507 /**
642     * Finds the processor where a molecule resides
643     * @return the id of the processor which contains the molecule
644     * @param globalIndex global Index of the molecule
645     */
646     int getMolToProc(int globalIndex) {
647     //assert(globalIndex < molToProcMap_.size());
648     return molToProcMap_[globalIndex];
649     }
650 gezelter 1528
651 gezelter 507 /**
652     * Set MolToProcMap array
653     * @see #SimCreator::divideMolecules
654     */
655 gezelter 1528 void setMolToProcMap(const vector<int>& molToProcMap) {
656 gezelter 507 molToProcMap_ = molToProcMap;
657     }
658 gezelter 246
659 gezelter 507 private:
660 gezelter 246
661 gezelter 507 /**
662 gezelter 1241 * The size of molToProcMap_ is equal to total number of molecules
663     * in the system. It maps a molecule to the processor on which it
664     * resides. it is filled by SimCreator once and only once.
665 gezelter 507 */
666 gezelter 1528 vector<int> molToProcMap_;
667 tim 292
668 gezelter 507 };
669 gezelter 2
670 gezelter 1390 } //namespace OpenMD
671 gezelter 246 #endif //BRAINS_SIMMODEL_HPP
672 gezelter 2

Properties

Name Value
svn:keywords Author Id Revision Date