ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.hpp
Revision: 1390
Committed: Wed Nov 25 20:02:06 2009 UTC (15 years, 5 months ago) by gezelter
Original Path: trunk/src/brains/SimInfo.hpp
File size: 21426 byte(s)
Log Message:
Almost all of the changes necessary to create OpenMD out of our old
project (OOPSE-4)

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 tim 316
66 gezelter 246 //another nonsense macro declaration
67 gezelter 1390 #define __OPENMD_C
68 tim 3 #include "brains/fSimulation.h"
69 gezelter 2
70 gezelter 1390 namespace OpenMD{
71 gezelter 2
72 gezelter 507 //forward decalration
73     class SnapshotManager;
74     class Molecule;
75     class SelectionManager;
76 tim 1024 class StuntDouble;
77 gezelter 507 /**
78     * @class SimInfo SimInfo.hpp "brains/SimInfo.hpp"
79 gezelter 1390 * @brief One of the heavy weight classes of OpenMD, SimInfo maintains a list of molecules.
80 gezelter 1277 * The Molecule class maintains all of the concrete objects
81     * (atoms, bond, bend, torsions, inversions, rigid bodies, cutoff groups,
82     * constraints). 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     */
86 gezelter 507 class SimInfo {
87     public:
88     typedef std::map<int, Molecule*>::iterator MoleculeIterator;
89 gezelter 2
90 gezelter 507 /**
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 770 SimInfo(ForceField* ff, Globals* simParams);
99 gezelter 507 virtual ~SimInfo();
100 gezelter 2
101 gezelter 507 /**
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 2
108 gezelter 507 /**
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 2
114 gezelter 507 /** Returns the total number of molecules in the system. */
115     int getNGlobalMolecules() {
116     return nGlobalMols_;
117     }
118 gezelter 2
119 gezelter 507 /** Returns the total number of atoms in the system. */
120     int getNGlobalAtoms() {
121     return nGlobalAtoms_;
122     }
123 gezelter 2
124 gezelter 507 /** Returns the total number of cutoff groups in the system. */
125     int getNGlobalCutoffGroups() {
126     return nGlobalCutoffGroups_;
127     }
128 gezelter 2
129 gezelter 507 /**
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 2
137 gezelter 507 /**
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 2
145 gezelter 507 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 2
154 gezelter 507 /** Returns the number of local atoms */
155     unsigned int getNAtoms() {
156     return nAtoms_;
157     }
158 gezelter 2
159 gezelter 507 /** Returns the number of local bonds */
160     unsigned int getNBonds(){
161     return nBonds_;
162     }
163 gezelter 2
164 gezelter 507 /** Returns the number of local bends */
165     unsigned int getNBends() {
166     return nBends_;
167     }
168 gezelter 2
169 gezelter 507 /** Returns the number of local torsions */
170     unsigned int getNTorsions() {
171     return nTorsions_;
172     }
173 gezelter 2
174 gezelter 1277 /** Returns the number of local torsions */
175     unsigned int getNInversions() {
176     return nInversions_;
177     }
178 gezelter 507 /** Returns the number of local rigid bodies */
179     unsigned int getNRigidBodies() {
180     return nRigidBodies_;
181     }
182 gezelter 2
183 gezelter 507 /** Returns the number of local integrable objects */
184     unsigned int getNIntegrableObjects() {
185     return nIntegrableObjects_;
186     }
187 gezelter 2
188 gezelter 507 /** Returns the number of local cutoff groups */
189     unsigned int getNCutoffGroups() {
190     return nCutoffGroups_;
191     }
192 gezelter 2
193 gezelter 507 /** Returns the total number of constraints in this SimInfo */
194     unsigned int getNConstraints() {
195     return nConstraints_;
196     }
197 gezelter 246
198 gezelter 507 /**
199     * Returns the first molecule in this SimInfo and intialize the iterator.
200     * @return the first molecule, return NULL if there is not molecule in this SimInfo
201     * @param i the iterator of molecule array (user shouldn't change it)
202     */
203     Molecule* beginMolecule(MoleculeIterator& i);
204 gezelter 2
205 gezelter 507 /**
206     * Returns the next avaliable Molecule based on the iterator.
207     * @return the next avaliable molecule, return NULL if reaching the end of the array
208     * @param i the iterator of molecule array
209     */
210     Molecule* nextMolecule(MoleculeIterator& i);
211 gezelter 2
212 gezelter 507 /** Returns the number of degrees of freedom */
213     int getNdf() {
214 gezelter 945 return ndf_ - getFdf();
215 gezelter 507 }
216 gezelter 2
217 gezelter 507 /** Returns the number of raw degrees of freedom */
218     int getNdfRaw() {
219     return ndfRaw_;
220     }
221 gezelter 2
222 gezelter 507 /** Returns the number of translational degrees of freedom */
223     int getNdfTrans() {
224     return ndfTrans_;
225     }
226 gezelter 2
227 gezelter 945 /** sets the current number of frozen degrees of freedom */
228     void setFdf(int fdf) {
229     fdf_local = fdf;
230     }
231    
232     int getFdf();
233    
234 gezelter 507 //getNZconstraint and setNZconstraint ruin the coherent of SimInfo class, need refactorying
235 gezelter 246
236 gezelter 507 /** Returns the total number of z-constraint molecules in the system */
237     int getNZconstraint() {
238     return nZconstraint_;
239     }
240 gezelter 2
241 gezelter 507 /**
242     * Sets the number of z-constraint molecules in the system.
243     */
244     void setNZconstraint(int nZconstraint) {
245     nZconstraint_ = nZconstraint;
246     }
247 gezelter 246
248 gezelter 507 /** Returns the snapshot manager. */
249     SnapshotManager* getSnapshotManager() {
250     return sman_;
251     }
252 gezelter 2
253 gezelter 507 /** Sets the snapshot manager. */
254     void setSnapshotManager(SnapshotManager* sman);
255 gezelter 246
256 gezelter 507 /** Returns the force field */
257     ForceField* getForceField() {
258     return forceField_;
259     }
260 gezelter 2
261 gezelter 507 Globals* getSimParams() {
262     return simParams_;
263     }
264 gezelter 2
265 gezelter 507 /** Returns the velocity of center of mass of the whole system.*/
266     Vector3d getComVel();
267 gezelter 2
268 gezelter 507 /** Returns the center of the mass of the whole system.*/
269     Vector3d getCom();
270 chuckv 555 /** Returns the center of the mass and Center of Mass velocity of the whole system.*/
271     void getComAll(Vector3d& com,Vector3d& comVel);
272 gezelter 2
273 chuckv 555 /** Returns intertia tensor for the entire system and system Angular Momentum.*/
274     void getInertiaTensor(Mat3x3d &intertiaTensor,Vector3d &angularMomentum);
275    
276     /** Returns system angular momentum */
277     Vector3d getAngularMomentum();
278    
279 chuckv 1103 /** Returns volume of system as estimated by an ellipsoid defined by the radii of gyration*/
280     void getGyrationalVolume(RealType &vol);
281     /** Overloaded version of gyrational volume that also returns det(I) so dV/dr can be calculated*/
282     void getGyrationalVolume(RealType &vol, RealType &detI);
283 gezelter 507 /** main driver function to interact with fortran during the initialization and molecule migration */
284     void update();
285 gezelter 2
286 gezelter 507 /** Returns the local index manager */
287     LocalIndexManager* getLocalIndexManager() {
288     return &localIndexMan_;
289     }
290 gezelter 2
291 gezelter 507 int getMoleculeStampId(int globalIndex) {
292     //assert(globalIndex < molStampIds_.size())
293     return molStampIds_[globalIndex];
294     }
295 gezelter 2
296 gezelter 507 /** Returns the molecule stamp */
297     MoleculeStamp* getMoleculeStamp(int id) {
298     return moleculeStamps_[id];
299     }
300 gezelter 2
301 gezelter 507 /** Return the total number of the molecule stamps */
302     int getNMoleculeStamp() {
303     return moleculeStamps_.size();
304     }
305     /**
306     * Finds a molecule with a specified global index
307     * @return a pointer point to found molecule
308     * @param index
309     */
310     Molecule* getMoleculeByGlobalIndex(int index) {
311     MoleculeIterator i;
312     i = molecules_.find(index);
313 gezelter 2
314 gezelter 507 return i != molecules_.end() ? i->second : NULL;
315     }
316 gezelter 2
317 chuckv 1292 int getGlobalMolMembership(int id){
318     return globalMolMembership_[id];
319     }
320    
321 tim 963 RealType getRcut() {
322 gezelter 507 return rcut_;
323     }
324 gezelter 2
325 tim 963 RealType getRsw() {
326 gezelter 507 return rsw_;
327     }
328 gezelter 764
329 tim 963 RealType getList() {
330 gezelter 764 return rlist_;
331     }
332 gezelter 246
333 gezelter 507 std::string getFinalConfigFileName() {
334     return finalConfigFileName_;
335     }
336 tim 1024
337 gezelter 507 void setFinalConfigFileName(const std::string& fileName) {
338     finalConfigFileName_ = fileName;
339     }
340 gezelter 2
341 tim 1024 std::string getRawMetaData() {
342     return rawMetaData_;
343     }
344     void setRawMetaData(const std::string& rawMetaData) {
345     rawMetaData_ = rawMetaData;
346     }
347    
348 gezelter 507 std::string getDumpFileName() {
349     return dumpFileName_;
350     }
351 gezelter 246
352 gezelter 507 void setDumpFileName(const std::string& fileName) {
353     dumpFileName_ = fileName;
354     }
355 gezelter 2
356 gezelter 507 std::string getStatFileName() {
357     return statFileName_;
358     }
359 gezelter 246
360 gezelter 507 void setStatFileName(const std::string& fileName) {
361     statFileName_ = fileName;
362     }
363 chrisfen 417
364 gezelter 507 std::string getRestFileName() {
365     return restFileName_;
366     }
367 chrisfen 417
368 gezelter 507 void setRestFileName(const std::string& fileName) {
369     restFileName_ = fileName;
370     }
371 gezelter 2
372 gezelter 507 /**
373     * Sets GlobalGroupMembership
374     * @see #SimCreator::setGlobalIndex
375     */
376     void setGlobalGroupMembership(const std::vector<int>& globalGroupMembership) {
377 gezelter 1287 assert(globalGroupMembership.size() == static_cast<size_t>(nGlobalAtoms_));
378 gezelter 507 globalGroupMembership_ = globalGroupMembership;
379     }
380 gezelter 2
381 gezelter 507 /**
382     * Sets GlobalMolMembership
383     * @see #SimCreator::setGlobalIndex
384     */
385     void setGlobalMolMembership(const std::vector<int>& globalMolMembership) {
386 gezelter 1287 assert(globalMolMembership.size() == static_cast<size_t>(nGlobalAtoms_));
387 gezelter 507 globalMolMembership_ = globalMolMembership;
388     }
389 gezelter 246
390    
391 gezelter 507 bool isFortranInitialized() {
392     return fortranInitialized_;
393     }
394 gezelter 246
395 chrisfen 998 bool getCalcBoxDipole() {
396     return calcBoxDipole_;
397     }
398    
399 gezelter 1126 bool getUseAtomicVirial() {
400     return useAtomicVirial_;
401     }
402    
403 gezelter 507 //below functions are just forward functions
404     //To compose or to inherit is always a hot debate. In general, is-a relation need subclassing, in the
405     //the other hand, has-a relation need composing.
406     /**
407     * Adds property into property map
408     * @param genData GenericData to be added into PropertyMap
409     */
410     void addProperty(GenericData* genData);
411 gezelter 246
412 gezelter 507 /**
413     * Removes property from PropertyMap by name
414     * @param propName the name of property to be removed
415     */
416     void removeProperty(const std::string& propName);
417 gezelter 246
418 gezelter 507 /**
419     * clear all of the properties
420     */
421     void clearProperties();
422 gezelter 246
423 gezelter 507 /**
424     * Returns all names of properties
425     * @return all names of properties
426     */
427     std::vector<std::string> getPropertyNames();
428 gezelter 246
429 gezelter 507 /**
430     * Returns all of the properties in PropertyMap
431     * @return all of the properties in PropertyMap
432     */
433     std::vector<GenericData*> getProperties();
434 gezelter 246
435 gezelter 507 /**
436     * Returns property
437     * @param propName name of property
438     * @return a pointer point to property with propName. If no property named propName
439     * exists, return NULL
440     */
441     GenericData* getPropertyByName(const std::string& propName);
442 gezelter 246
443 gezelter 507 /**
444 gezelter 1287 * add all special interaction pairs (including excluded
445     * interactions) in a molecule into the appropriate lists.
446 gezelter 507 */
447 gezelter 1287 void addInteractionPairs(Molecule* mol);
448 gezelter 246
449 gezelter 507 /**
450 gezelter 1287 * remove all special interaction pairs which belong to a molecule
451     * from the appropriate lists.
452 gezelter 507 */
453 gezelter 1287 void removeInteractionPairs(Molecule* mol);
454 gezelter 246
455 tim 292
456 gezelter 507 /** Returns the unique atom types of local processor in an array */
457     std::set<AtomType*> getUniqueAtomTypes();
458 tim 292
459 gezelter 507 friend std::ostream& operator <<(std::ostream& o, SimInfo& info);
460 tim 326
461 tim 963 void getCutoff(RealType& rcut, RealType& rsw);
462 gezelter 246
463 gezelter 507 private:
464 gezelter 246
465 gezelter 507 /** fill up the simtype struct*/
466     void setupSimType();
467 gezelter 246
468 gezelter 507 /**
469     * Setup Fortran Simulation
470     * @see #setupFortranParallel
471     */
472     void setupFortranSim();
473 gezelter 246
474 gezelter 507 /** Figure out the radius of cutoff, radius of switching function and pass them to fortran */
475     void setupCutoff();
476 gezelter 246
477 chrisfen 598 /** Figure out which coulombic correction method to use and pass to fortran */
478 chrisfen 604 void setupElectrostaticSummationMethod( int isError );
479 chrisfen 598
480 chrisfen 726 /** Figure out which polynomial type to use for the switching function */
481     void setupSwitchingFunction();
482    
483 chrisfen 998 /** Determine if we need to accumulate the simulation box dipole */
484     void setupAccumulateBoxDipole();
485    
486 gezelter 507 /** Calculates the number of degress of freedom in the whole system */
487     void calcNdf();
488     void calcNdfRaw();
489     void calcNdfTrans();
490 gezelter 246
491 tim 770 ForceField* forceField_;
492     Globals* simParams_;
493    
494     std::map<int, Molecule*> molecules_; /**< Molecule array */
495    
496 gezelter 507 /**
497     * Adds molecule stamp and the total number of the molecule with same molecule stamp in the whole
498     * system.
499     */
500     void addMoleculeStamp(MoleculeStamp* molStamp, int nmol);
501 gezelter 246
502 gezelter 507 //degress of freedom
503     int ndf_; /**< number of degress of freedom (excludes constraints), ndf_ is local */
504 gezelter 945 int fdf_local; /**< number of frozen degrees of freedom */
505     int fdf_; /**< number of frozen degrees of freedom */
506 gezelter 507 int ndfRaw_; /**< number of degress of freedom (includes constraints), ndfRaw_ is local */
507     int ndfTrans_; /**< number of translation degress of freedom, ndfTrans_ is local */
508     int nZconstraint_; /** number of z-constraint molecules, nZconstraint_ is global */
509 gezelter 246
510 gezelter 507 //number of global objects
511     int nGlobalMols_; /**< number of molecules in the system */
512     int nGlobalAtoms_; /**< number of atoms in the system */
513     int nGlobalCutoffGroups_; /**< number of cutoff groups in this system */
514     int nGlobalIntegrableObjects_; /**< number of integrable objects in this system */
515     int nGlobalRigidBodies_; /**< number of rigid bodies in this system */
516     /**
517     * the size of globalGroupMembership_ is nGlobalAtoms. Its index is global index of an atom, and the
518     * corresponding content is the global index of cutoff group this atom belong to.
519     * It is filled by SimCreator once and only once, since it never changed during the simulation.
520     */
521     std::vector<int> globalGroupMembership_;
522 gezelter 246
523 gezelter 507 /**
524 gezelter 1313 * the size of globalMolMembership_ is nGlobalAtoms. Its index is global index of an atom, and the
525 gezelter 507 * corresponding content is the global index of molecule this atom belong to.
526     * It is filled by SimCreator once and only once, since it is never changed during the simulation.
527     */
528     std::vector<int> globalMolMembership_;
529 gezelter 246
530    
531 gezelter 507 std::vector<int> molStampIds_; /**< stamp id array of all molecules in the system */
532     std::vector<MoleculeStamp*> moleculeStamps_; /**< molecule stamps array */
533 gezelter 246
534 gezelter 507 //number of local objects
535 gezelter 1277 int nAtoms_; /**< number of atoms in local processor */
536     int nBonds_; /**< number of bonds in local processor */
537     int nBends_; /**< number of bends in local processor */
538     int nTorsions_; /**< number of torsions in local processor */
539     int nInversions_; /**< number of inversions in local processor */
540     int nRigidBodies_; /**< number of rigid bodies in local processor */
541     int nIntegrableObjects_; /**< number of integrable objects in local processor */
542     int nCutoffGroups_; /**< number of cutoff groups in local processor */
543     int nConstraints_; /**< number of constraints in local processors */
544 gezelter 246
545 gezelter 507 simtype fInfo_; /**< A dual struct shared by c++/fortran which indicates the atom types in simulation*/
546 gezelter 1287 PairList excludedInteractions_;
547     PairList oneTwoInteractions_;
548     PairList oneThreeInteractions_;
549     PairList oneFourInteractions_;
550 gezelter 507 PropertyMap properties_; /**< Generic Property */
551     SnapshotManager* sman_; /**< SnapshotManager */
552 gezelter 246
553 gezelter 507 /**
554     * The reason to have a local index manager is that when molecule is migrating to other processors,
555     * the atoms and the rigid-bodies will release their local indices to LocalIndexManager. Combining the
556     * information of molecule migrating to current processor, Migrator class can query the LocalIndexManager
557     * to make a efficient data moving plan.
558     */
559     LocalIndexManager localIndexMan_;
560 gezelter 246
561 tim 1024 // unparsed MetaData block for storing in Dump and EOR files:
562     std::string rawMetaData_;
563    
564 gezelter 507 //file names
565     std::string finalConfigFileName_;
566     std::string dumpFileName_;
567     std::string statFileName_;
568     std::string restFileName_;
569 chrisfen 417
570 tim 963 RealType rcut_; /**< cutoff radius*/
571     RealType rsw_; /**< radius of switching function*/
572     RealType rlist_; /**< neighbor list radius */
573 gezelter 246
574 gezelter 1386 int ljsp_; /**< use shifted potential for LJ*/
575     int ljsf_; /**< use shifted force for LJ*/
576 chrisfen 1129
577 gezelter 1126 bool fortranInitialized_; /**< flag indicate whether fortran side
578     is initialized */
579    
580     bool calcBoxDipole_; /**< flag to indicate whether or not we calculate
581     the simulation box dipole moment */
582    
583     bool useAtomicVirial_; /**< flag to indicate whether or not we use
584     Atomic Virials to calculate the pressure */
585 tim 292
586 tim 1024 public:
587     /**
588     * return an integral objects by its global index. In MPI version, if the StuntDouble with specified
589     * global index does not belong to local processor, a NULL will be return.
590     */
591     StuntDouble* getIOIndexToIntegrableObject(int index);
592     void setIOIndexToIntegrableObject(const std::vector<StuntDouble*>& v);
593     private:
594     std::vector<StuntDouble*> IOIndexToIntegrableObject;
595     //public:
596     //void setStuntDoubleFromGlobalIndex(std::vector<StuntDouble*> v);
597     /**
598     * return a StuntDouble by its global index. In MPI version, if the StuntDouble with specified
599     * global index does not belong to local processor, a NULL will be return.
600     */
601     //StuntDouble* getStuntDoubleFromGlobalIndex(int index);
602     //private:
603     //std::vector<StuntDouble*> sdByGlobalIndex_;
604    
605 gezelter 246 //in Parallel version, we need MolToProc
606 gezelter 507 public:
607 gezelter 246
608 gezelter 507 /**
609     * Finds the processor where a molecule resides
610     * @return the id of the processor which contains the molecule
611     * @param globalIndex global Index of the molecule
612     */
613     int getMolToProc(int globalIndex) {
614     //assert(globalIndex < molToProcMap_.size());
615     return molToProcMap_[globalIndex];
616     }
617 gezelter 246
618 gezelter 507 /**
619     * Set MolToProcMap array
620     * @see #SimCreator::divideMolecules
621     */
622     void setMolToProcMap(const std::vector<int>& molToProcMap) {
623     molToProcMap_ = molToProcMap;
624     }
625 gezelter 246
626 gezelter 507 private:
627 gezelter 246
628 gezelter 507 void setupFortranParallel();
629 gezelter 246
630 gezelter 507 /**
631 gezelter 1241 * The size of molToProcMap_ is equal to total number of molecules
632     * in the system. It maps a molecule to the processor on which it
633     * resides. it is filled by SimCreator once and only once.
634 gezelter 507 */
635     std::vector<int> molToProcMap_;
636 tim 292
637 gezelter 246
638 gezelter 507 };
639 gezelter 2
640 gezelter 1390 } //namespace OpenMD
641 gezelter 246 #endif //BRAINS_SIMMODEL_HPP
642 gezelter 2