ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.hpp
Revision: 1503
Committed: Sat Oct 2 19:54:41 2010 UTC (14 years, 7 months ago) by gezelter
File size: 21460 byte(s)
Log Message:
Changes to remove more of the low level stuff from the fortran side.

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

Properties

Name Value
svn:keywords Author Id Revision Date