ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/brains/SimInfo.hpp
Revision: 3126
Committed: Fri Apr 6 21:53:43 2007 UTC (17 years, 5 months ago) by gezelter
File size: 20779 byte(s)
Log Message:
Massive update to do virials (both atomic and cutoff-group) correctly.
The rigid body constraint contributions had been missing and this was
masked by the use of cutoff groups...

File Contents

# User Rev Content
1 gezelter 2204 /*
2 gezelter 1930 * 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     * 1. Acknowledgement of the program authors must be made in any
10     * publication of scientific results based in part on use of the
11     * program. An acceptable form of acknowledgement is citation of
12     * the article in which the program was described (Matthew
13     * A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14     * J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15     * Parallel Simulation Engine for Molecular Dynamics,"
16     * J. Comput. Chem. 26, pp. 252-271 (2005))
17     *
18     * 2. Redistributions of source code must retain the above copyright
19     * notice, this list of conditions and the following disclaimer.
20     *
21     * 3. Redistributions in binary form must reproduce the above copyright
22     * notice, this list of conditions and the following disclaimer in the
23     * documentation and/or other materials provided with the
24     * distribution.
25     *
26     * This software is provided "AS IS," without a warranty of any
27     * kind. All express or implied conditions, representations and
28     * warranties, including any implied warranty of merchantability,
29     * fitness for a particular purpose or non-infringement, are hereby
30     * excluded. The University of Notre Dame and its licensors shall not
31     * be liable for any damages suffered by licensee as a result of
32     * using, modifying or distributing the software or its
33     * derivatives. In no event will the University of Notre Dame or its
34     * licensors be liable for any lost revenue, profit or data, or for
35     * direct, indirect, special, consequential, incidental or punitive
36     * damages, however caused and regardless of the theory of liability,
37     * arising out of the use of or inability to use software, even if the
38     * University of Notre Dame has been advised of the possibility of
39     * such damages.
40     */
41    
42     /**
43     * @file SimInfo.hpp
44     * @author tlin
45     * @date 11/02/2004
46     * @version 1.0
47     */
48 gezelter 1490
49 gezelter 1930 #ifndef BRAINS_SIMMODEL_HPP
50     #define BRAINS_SIMMODEL_HPP
51    
52     #include <iostream>
53     #include <set>
54     #include <utility>
55 gezelter 1490 #include <vector>
56    
57 tim 1492 #include "brains/Exclude.hpp"
58 gezelter 1930 #include "io/Globals.hpp"
59     #include "math/Vector3.hpp"
60 chuckv 2252 #include "math/SquareMatrix3.hpp"
61 gezelter 1930 #include "types/MoleculeStamp.hpp"
62     #include "UseTheForce/ForceField.hpp"
63     #include "utils/PropertyMap.hpp"
64     #include "utils/LocalIndexManager.hpp"
65 tim 2000
66 gezelter 1930 //another nonsense macro declaration
67 gezelter 1490 #define __C
68 tim 1492 #include "brains/fSimulation.h"
69 gezelter 1490
70 gezelter 1930 namespace oopse{
71 gezelter 1490
72 gezelter 2204 //forward decalration
73     class SnapshotManager;
74     class Molecule;
75     class SelectionManager;
76 tim 2982 class StuntDouble;
77 gezelter 2204 /**
78     * @class SimInfo SimInfo.hpp "brains/SimInfo.hpp"
79 gezelter 2733 * @brief One of the heavy weight classes of OOPSE, SimInfo maintains a list of molecules.
80     * The Molecule class maintains all of the concrete objects
81     * (atoms, bond, bend, torsions, rigid bodies, cutoff groups, constrains).
82     * In both the single and parallel versions, atoms and
83     * rigid bodies have both global and local indices. The local index is
84     * not relevant to molecules or cutoff groups.
85 gezelter 2204 */
86     class SimInfo {
87     public:
88     typedef std::map<int, Molecule*>::iterator MoleculeIterator;
89 gezelter 1490
90 gezelter 2204 /**
91     * Constructor of SimInfo
92     * @param molStampPairs MoleculeStamp Array. The first element of the pair is molecule stamp, the
93     * second element is the total number of molecules with the same molecule stamp in the system
94     * @param ff pointer of a concrete ForceField instance
95     * @param simParams
96     * @note
97     */
98 tim 2469 SimInfo(ForceField* ff, Globals* simParams);
99 gezelter 2204 virtual ~SimInfo();
100 gezelter 1490
101 gezelter 2204 /**
102     * Adds a molecule
103     * @return return true if adding successfully, return false if the molecule is already in SimInfo
104     * @param mol molecule to be added
105     */
106     bool addMolecule(Molecule* mol);
107 gezelter 1490
108 gezelter 2204 /**
109     * Removes a molecule from SimInfo
110     * @return true if removing successfully, return false if molecule is not in this SimInfo
111     */
112     bool removeMolecule(Molecule* mol);
113 gezelter 1490
114 gezelter 2204 /** Returns the total number of molecules in the system. */
115     int getNGlobalMolecules() {
116     return nGlobalMols_;
117     }
118 gezelter 1490
119 gezelter 2204 /** Returns the total number of atoms in the system. */
120     int getNGlobalAtoms() {
121     return nGlobalAtoms_;
122     }
123 gezelter 1490
124 gezelter 2204 /** Returns the total number of cutoff groups in the system. */
125     int getNGlobalCutoffGroups() {
126     return nGlobalCutoffGroups_;
127     }
128 gezelter 1490
129 gezelter 2204 /**
130     * Returns the total number of integrable objects (total number of rigid bodies plus the total number
131     * of atoms which do not belong to the rigid bodies) in the system
132     */
133     int getNGlobalIntegrableObjects() {
134     return nGlobalIntegrableObjects_;
135     }
136 gezelter 1490
137 gezelter 2204 /**
138     * Returns the total number of integrable objects (total number of rigid bodies plus the total number
139     * of atoms which do not belong to the rigid bodies) in the system
140     */
141     int getNGlobalRigidBodies() {
142     return nGlobalRigidBodies_;
143     }
144 gezelter 1490
145 gezelter 2204 int getNGlobalConstraints();
146     /**
147     * Returns the number of local molecules.
148     * @return the number of local molecules
149     */
150     int getNMolecules() {
151     return molecules_.size();
152     }
153 gezelter 1490
154 gezelter 2204 /** Returns the number of local atoms */
155     unsigned int getNAtoms() {
156     return nAtoms_;
157     }
158 gezelter 1490
159 gezelter 2204 /** Returns the number of local bonds */
160     unsigned int getNBonds(){
161     return nBonds_;
162     }
163 gezelter 1490
164 gezelter 2204 /** Returns the number of local bends */
165     unsigned int getNBends() {
166     return nBends_;
167     }
168 gezelter 1490
169 gezelter 2204 /** Returns the number of local torsions */
170     unsigned int getNTorsions() {
171     return nTorsions_;
172     }
173 gezelter 1490
174 gezelter 2204 /** Returns the number of local rigid bodies */
175     unsigned int getNRigidBodies() {
176     return nRigidBodies_;
177     }
178 gezelter 1490
179 gezelter 2204 /** Returns the number of local integrable objects */
180     unsigned int getNIntegrableObjects() {
181     return nIntegrableObjects_;
182     }
183 gezelter 1490
184 gezelter 2204 /** Returns the number of local cutoff groups */
185     unsigned int getNCutoffGroups() {
186     return nCutoffGroups_;
187     }
188 gezelter 1490
189 gezelter 2204 /** Returns the total number of constraints in this SimInfo */
190     unsigned int getNConstraints() {
191     return nConstraints_;
192     }
193 gezelter 1930
194 gezelter 2204 /**
195     * Returns the first molecule in this SimInfo and intialize the iterator.
196     * @return the first molecule, return NULL if there is not molecule in this SimInfo
197     * @param i the iterator of molecule array (user shouldn't change it)
198     */
199     Molecule* beginMolecule(MoleculeIterator& i);
200 gezelter 1490
201 gezelter 2204 /**
202     * Returns the next avaliable Molecule based on the iterator.
203     * @return the next avaliable molecule, return NULL if reaching the end of the array
204     * @param i the iterator of molecule array
205     */
206     Molecule* nextMolecule(MoleculeIterator& i);
207 gezelter 1490
208 gezelter 2204 /** Returns the number of degrees of freedom */
209     int getNdf() {
210 gezelter 2733 return ndf_ - getFdf();
211 gezelter 2204 }
212 gezelter 1490
213 gezelter 2204 /** Returns the number of raw degrees of freedom */
214     int getNdfRaw() {
215     return ndfRaw_;
216     }
217 gezelter 1490
218 gezelter 2204 /** Returns the number of translational degrees of freedom */
219     int getNdfTrans() {
220     return ndfTrans_;
221     }
222 gezelter 1490
223 gezelter 2733 /** sets the current number of frozen degrees of freedom */
224     void setFdf(int fdf) {
225     fdf_local = fdf;
226     }
227    
228     int getFdf();
229    
230 gezelter 2204 //getNZconstraint and setNZconstraint ruin the coherent of SimInfo class, need refactorying
231 gezelter 1930
232 gezelter 2204 /** Returns the total number of z-constraint molecules in the system */
233     int getNZconstraint() {
234     return nZconstraint_;
235     }
236 gezelter 1490
237 gezelter 2204 /**
238     * Sets the number of z-constraint molecules in the system.
239     */
240     void setNZconstraint(int nZconstraint) {
241     nZconstraint_ = nZconstraint;
242     }
243 gezelter 1930
244 gezelter 2204 /** Returns the snapshot manager. */
245     SnapshotManager* getSnapshotManager() {
246     return sman_;
247     }
248 gezelter 1490
249 gezelter 2204 /** Sets the snapshot manager. */
250     void setSnapshotManager(SnapshotManager* sman);
251 gezelter 1930
252 gezelter 2204 /** Returns the force field */
253     ForceField* getForceField() {
254     return forceField_;
255     }
256 gezelter 1490
257 gezelter 2204 Globals* getSimParams() {
258     return simParams_;
259     }
260 gezelter 1490
261 gezelter 2204 /** Returns the velocity of center of mass of the whole system.*/
262     Vector3d getComVel();
263 gezelter 1490
264 gezelter 2204 /** Returns the center of the mass of the whole system.*/
265     Vector3d getCom();
266 chuckv 2252 /** Returns the center of the mass and Center of Mass velocity of the whole system.*/
267     void getComAll(Vector3d& com,Vector3d& comVel);
268 gezelter 1490
269 chuckv 2252 /** Returns intertia tensor for the entire system and system Angular Momentum.*/
270     void getInertiaTensor(Mat3x3d &intertiaTensor,Vector3d &angularMomentum);
271    
272     /** Returns system angular momentum */
273     Vector3d getAngularMomentum();
274    
275 chuckv 3100 /** Returns volume of system as estimated by an ellipsoid defined by the radii of gyration*/
276     void getGyrationalVolume(RealType &vol);
277     /** Overloaded version of gyrational volume that also returns det(I) so dV/dr can be calculated*/
278     void getGyrationalVolume(RealType &vol, RealType &detI);
279 gezelter 2204 /** main driver function to interact with fortran during the initialization and molecule migration */
280     void update();
281 gezelter 1490
282 gezelter 2204 /** Returns the local index manager */
283     LocalIndexManager* getLocalIndexManager() {
284     return &localIndexMan_;
285     }
286 gezelter 1490
287 gezelter 2204 int getMoleculeStampId(int globalIndex) {
288     //assert(globalIndex < molStampIds_.size())
289     return molStampIds_[globalIndex];
290     }
291 gezelter 1490
292 gezelter 2204 /** Returns the molecule stamp */
293     MoleculeStamp* getMoleculeStamp(int id) {
294     return moleculeStamps_[id];
295     }
296 gezelter 1490
297 gezelter 2204 /** Return the total number of the molecule stamps */
298     int getNMoleculeStamp() {
299     return moleculeStamps_.size();
300     }
301     /**
302     * Finds a molecule with a specified global index
303     * @return a pointer point to found molecule
304     * @param index
305     */
306     Molecule* getMoleculeByGlobalIndex(int index) {
307     MoleculeIterator i;
308     i = molecules_.find(index);
309 gezelter 1490
310 gezelter 2204 return i != molecules_.end() ? i->second : NULL;
311     }
312 gezelter 1490
313 tim 2759 RealType getRcut() {
314 gezelter 2204 return rcut_;
315     }
316 gezelter 1490
317 tim 2759 RealType getRsw() {
318 gezelter 2204 return rsw_;
319     }
320 gezelter 2463
321 tim 2759 RealType getList() {
322 gezelter 2463 return rlist_;
323     }
324 gezelter 1930
325 gezelter 2204 std::string getFinalConfigFileName() {
326     return finalConfigFileName_;
327     }
328 tim 2982
329 gezelter 2204 void setFinalConfigFileName(const std::string& fileName) {
330     finalConfigFileName_ = fileName;
331     }
332 gezelter 1490
333 tim 2982 std::string getRawMetaData() {
334     return rawMetaData_;
335     }
336     void setRawMetaData(const std::string& rawMetaData) {
337     rawMetaData_ = rawMetaData;
338     }
339    
340 gezelter 2204 std::string getDumpFileName() {
341     return dumpFileName_;
342     }
343 gezelter 1930
344 gezelter 2204 void setDumpFileName(const std::string& fileName) {
345     dumpFileName_ = fileName;
346     }
347 gezelter 1490
348 gezelter 2204 std::string getStatFileName() {
349     return statFileName_;
350     }
351 gezelter 1930
352 gezelter 2204 void setStatFileName(const std::string& fileName) {
353     statFileName_ = fileName;
354     }
355 chrisfen 2101
356 gezelter 2204 std::string getRestFileName() {
357     return restFileName_;
358     }
359 chrisfen 2101
360 gezelter 2204 void setRestFileName(const std::string& fileName) {
361     restFileName_ = fileName;
362     }
363 gezelter 1490
364 gezelter 2204 /**
365     * Sets GlobalGroupMembership
366     * @see #SimCreator::setGlobalIndex
367     */
368     void setGlobalGroupMembership(const std::vector<int>& globalGroupMembership) {
369     assert(globalGroupMembership.size() == nGlobalAtoms_);
370     globalGroupMembership_ = globalGroupMembership;
371     }
372 gezelter 1490
373 gezelter 2204 /**
374     * Sets GlobalMolMembership
375     * @see #SimCreator::setGlobalIndex
376     */
377     void setGlobalMolMembership(const std::vector<int>& globalMolMembership) {
378     assert(globalMolMembership.size() == nGlobalAtoms_);
379     globalMolMembership_ = globalMolMembership;
380     }
381 gezelter 1930
382    
383 gezelter 2204 bool isFortranInitialized() {
384     return fortranInitialized_;
385     }
386 gezelter 1930
387 chrisfen 2917 bool getCalcBoxDipole() {
388     return calcBoxDipole_;
389     }
390    
391 gezelter 3126 bool getUseAtomicVirial() {
392     return useAtomicVirial_;
393     }
394    
395 gezelter 2204 //below functions are just forward functions
396     //To compose or to inherit is always a hot debate. In general, is-a relation need subclassing, in the
397     //the other hand, has-a relation need composing.
398     /**
399     * Adds property into property map
400     * @param genData GenericData to be added into PropertyMap
401     */
402     void addProperty(GenericData* genData);
403 gezelter 1930
404 gezelter 2204 /**
405     * Removes property from PropertyMap by name
406     * @param propName the name of property to be removed
407     */
408     void removeProperty(const std::string& propName);
409 gezelter 1930
410 gezelter 2204 /**
411     * clear all of the properties
412     */
413     void clearProperties();
414 gezelter 1930
415 gezelter 2204 /**
416     * Returns all names of properties
417     * @return all names of properties
418     */
419     std::vector<std::string> getPropertyNames();
420 gezelter 1930
421 gezelter 2204 /**
422     * Returns all of the properties in PropertyMap
423     * @return all of the properties in PropertyMap
424     */
425     std::vector<GenericData*> getProperties();
426 gezelter 1930
427 gezelter 2204 /**
428     * Returns property
429     * @param propName name of property
430     * @return a pointer point to property with propName. If no property named propName
431     * exists, return NULL
432     */
433     GenericData* getPropertyByName(const std::string& propName);
434 gezelter 1930
435 gezelter 2204 /**
436     * add all exclude pairs of a molecule into exclude list.
437     */
438     void addExcludePairs(Molecule* mol);
439 gezelter 1930
440 gezelter 2204 /**
441     * remove all exclude pairs which belong to a molecule from exclude list
442     */
443 gezelter 1930
444 gezelter 2204 void removeExcludePairs(Molecule* mol);
445 tim 1976
446    
447 gezelter 2204 /** Returns the unique atom types of local processor in an array */
448     std::set<AtomType*> getUniqueAtomTypes();
449 tim 1976
450 gezelter 2204 friend std::ostream& operator <<(std::ostream& o, SimInfo& info);
451 tim 2010
452 tim 2759 void getCutoff(RealType& rcut, RealType& rsw);
453 gezelter 1930
454 gezelter 2204 private:
455 gezelter 1930
456 gezelter 2204 /** fill up the simtype struct*/
457     void setupSimType();
458 gezelter 1930
459 gezelter 2204 /**
460     * Setup Fortran Simulation
461     * @see #setupFortranParallel
462     */
463     void setupFortranSim();
464 gezelter 1930
465 gezelter 2204 /** Figure out the radius of cutoff, radius of switching function and pass them to fortran */
466     void setupCutoff();
467 gezelter 1930
468 chrisfen 2297 /** Figure out which coulombic correction method to use and pass to fortran */
469 chrisfen 2303 void setupElectrostaticSummationMethod( int isError );
470 chrisfen 2297
471 chrisfen 2425 /** Figure out which polynomial type to use for the switching function */
472     void setupSwitchingFunction();
473    
474 chrisfen 2917 /** Determine if we need to accumulate the simulation box dipole */
475     void setupAccumulateBoxDipole();
476    
477 gezelter 2204 /** Calculates the number of degress of freedom in the whole system */
478     void calcNdf();
479     void calcNdfRaw();
480     void calcNdfTrans();
481 gezelter 1930
482 tim 2469 ForceField* forceField_;
483     Globals* simParams_;
484    
485     std::map<int, Molecule*> molecules_; /**< Molecule array */
486    
487 gezelter 2204 /**
488     * Adds molecule stamp and the total number of the molecule with same molecule stamp in the whole
489     * system.
490     */
491     void addMoleculeStamp(MoleculeStamp* molStamp, int nmol);
492 gezelter 1930
493 gezelter 2204 //degress of freedom
494     int ndf_; /**< number of degress of freedom (excludes constraints), ndf_ is local */
495 gezelter 2733 int fdf_local; /**< number of frozen degrees of freedom */
496     int fdf_; /**< number of frozen degrees of freedom */
497 gezelter 2204 int ndfRaw_; /**< number of degress of freedom (includes constraints), ndfRaw_ is local */
498     int ndfTrans_; /**< number of translation degress of freedom, ndfTrans_ is local */
499     int nZconstraint_; /** number of z-constraint molecules, nZconstraint_ is global */
500 gezelter 1930
501 gezelter 2204 //number of global objects
502     int nGlobalMols_; /**< number of molecules in the system */
503     int nGlobalAtoms_; /**< number of atoms in the system */
504     int nGlobalCutoffGroups_; /**< number of cutoff groups in this system */
505     int nGlobalIntegrableObjects_; /**< number of integrable objects in this system */
506     int nGlobalRigidBodies_; /**< number of rigid bodies in this system */
507     /**
508     * the size of globalGroupMembership_ is nGlobalAtoms. Its index is global index of an atom, and the
509     * corresponding content is the global index of cutoff group this atom belong to.
510     * It is filled by SimCreator once and only once, since it never changed during the simulation.
511     */
512     std::vector<int> globalGroupMembership_;
513 gezelter 1930
514 gezelter 2204 /**
515     * the size of globalGroupMembership_ is nGlobalAtoms. Its index is global index of an atom, and the
516     * corresponding content is the global index of molecule this atom belong to.
517     * It is filled by SimCreator once and only once, since it is never changed during the simulation.
518     */
519     std::vector<int> globalMolMembership_;
520 gezelter 1930
521    
522 gezelter 2204 std::vector<int> molStampIds_; /**< stamp id array of all molecules in the system */
523     std::vector<MoleculeStamp*> moleculeStamps_; /**< molecule stamps array */
524 gezelter 1930
525 gezelter 2204 //number of local objects
526     int nAtoms_; /**< number of atoms in local processor */
527     int nBonds_; /**< number of bonds in local processor */
528     int nBends_; /**< number of bends in local processor */
529     int nTorsions_; /**< number of torsions in local processor */
530     int nRigidBodies_; /**< number of rigid bodies in local processor */
531     int nIntegrableObjects_; /**< number of integrable objects in local processor */
532     int nCutoffGroups_; /**< number of cutoff groups in local processor */
533     int nConstraints_; /**< number of constraints in local processors */
534 gezelter 1930
535 gezelter 2204 simtype fInfo_; /**< A dual struct shared by c++/fortran which indicates the atom types in simulation*/
536     Exclude exclude_;
537     PropertyMap properties_; /**< Generic Property */
538     SnapshotManager* sman_; /**< SnapshotManager */
539 gezelter 1930
540 gezelter 2204 /**
541     * The reason to have a local index manager is that when molecule is migrating to other processors,
542     * the atoms and the rigid-bodies will release their local indices to LocalIndexManager. Combining the
543     * information of molecule migrating to current processor, Migrator class can query the LocalIndexManager
544     * to make a efficient data moving plan.
545     */
546     LocalIndexManager localIndexMan_;
547 gezelter 1930
548 tim 2982 // unparsed MetaData block for storing in Dump and EOR files:
549     std::string rawMetaData_;
550    
551 gezelter 2204 //file names
552     std::string finalConfigFileName_;
553     std::string dumpFileName_;
554     std::string statFileName_;
555     std::string restFileName_;
556 chrisfen 2101
557 tim 2759 RealType rcut_; /**< cutoff radius*/
558     RealType rsw_; /**< radius of switching function*/
559     RealType rlist_; /**< neighbor list radius */
560 gezelter 1930
561 gezelter 3126 bool fortranInitialized_; /**< flag indicate whether fortran side
562     is initialized */
563    
564     bool calcBoxDipole_; /**< flag to indicate whether or not we calculate
565     the simulation box dipole moment */
566    
567     bool useAtomicVirial_; /**< flag to indicate whether or not we use
568     Atomic Virials to calculate the pressure */
569 tim 1976
570 tim 2982 public:
571     /**
572     * return an integral objects by its global index. In MPI version, if the StuntDouble with specified
573     * global index does not belong to local processor, a NULL will be return.
574     */
575     StuntDouble* getIOIndexToIntegrableObject(int index);
576     void setIOIndexToIntegrableObject(const std::vector<StuntDouble*>& v);
577     private:
578     std::vector<StuntDouble*> IOIndexToIntegrableObject;
579     //public:
580     //void setStuntDoubleFromGlobalIndex(std::vector<StuntDouble*> v);
581     /**
582     * return a StuntDouble by its global index. In MPI version, if the StuntDouble with specified
583     * global index does not belong to local processor, a NULL will be return.
584     */
585     //StuntDouble* getStuntDoubleFromGlobalIndex(int index);
586     //private:
587     //std::vector<StuntDouble*> sdByGlobalIndex_;
588    
589 gezelter 1930 #ifdef IS_MPI
590     //in Parallel version, we need MolToProc
591 gezelter 2204 public:
592 gezelter 1930
593 gezelter 2204 /**
594     * Finds the processor where a molecule resides
595     * @return the id of the processor which contains the molecule
596     * @param globalIndex global Index of the molecule
597     */
598     int getMolToProc(int globalIndex) {
599     //assert(globalIndex < molToProcMap_.size());
600     return molToProcMap_[globalIndex];
601     }
602 gezelter 1930
603 gezelter 2204 /**
604     * Set MolToProcMap array
605     * @see #SimCreator::divideMolecules
606     */
607     void setMolToProcMap(const std::vector<int>& molToProcMap) {
608     molToProcMap_ = molToProcMap;
609     }
610 tim 2982
611    
612 gezelter 1930
613 gezelter 2204 private:
614 gezelter 1930
615 gezelter 2204 void setupFortranParallel();
616 gezelter 1930
617 gezelter 2204 /**
618     * The size of molToProcMap_ is equal to total number of molecules in the system.
619     * It maps a molecule to the processor on which it resides. it is filled by SimCreator once and only
620     * once.
621     */
622     std::vector<int> molToProcMap_;
623 tim 1976
624 gezelter 1930 #endif
625    
626 gezelter 2204 };
627 gezelter 1490
628 gezelter 1930 } //namespace oopse
629     #endif //BRAINS_SIMMODEL_HPP
630 gezelter 1490