ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimInfo.hpp
Revision: 1846
Committed: Thu Jan 31 15:55:47 2013 UTC (12 years, 3 months ago) by gezelter
File size: 23112 byte(s)
Log Message:
Bug fix in Tetrahedrality parameter Z exposed another bug in SimInfo.

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

Properties

Name Value
svn:keywords Author Id Revision Date