ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.hpp
Revision: 1530
Committed: Tue Dec 28 21:47:55 2010 UTC (14 years, 4 months ago) by gezelter
File size: 22964 byte(s)
Log Message:
Moved switching functions and force options over to the C++ side, and
removed them from Fortran.

File Contents

# User Rev Content
1 gezelter 507 /*
2 gezelter 246 * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3     *
4     * The University of Notre Dame grants you ("Licensee") a
5     * non-exclusive, royalty free, license to use, modify and
6     * redistribute this software in source and binary code form, provided
7     * that the following conditions are met:
8     *
9 gezelter 1390 * 1. Redistributions of source code must retain the above copyright
10 gezelter 246 * notice, this list of conditions and the following disclaimer.
11     *
12 gezelter 1390 * 2. Redistributions in binary form must reproduce the above copyright
13 gezelter 246 * notice, this list of conditions and the following disclaimer in the
14     * documentation and/or other materials provided with the
15     * distribution.
16     *
17     * This software is provided "AS IS," without a warranty of any
18     * kind. All express or implied conditions, representations and
19     * warranties, including any implied warranty of merchantability,
20     * fitness for a particular purpose or non-infringement, are hereby
21     * excluded. The University of Notre Dame and its licensors shall not
22     * be liable for any damages suffered by licensee as a result of
23     * using, modifying or distributing the software or its
24     * derivatives. In no event will the University of Notre Dame or its
25     * licensors be liable for any lost revenue, profit or data, or for
26     * direct, indirect, special, consequential, incidental or punitive
27     * damages, however caused and regardless of the theory of liability,
28     * arising out of the use of or inability to use software, even if the
29     * University of Notre Dame has been advised of the possibility of
30     * such damages.
31 gezelter 1390 *
32     * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33     * research, please cite the appropriate papers when you publish your
34     * work. Good starting points are:
35     *
36     * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37     * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38     * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39     * [4] Vardeman & Gezelter, in progress (2009).
40 gezelter 246 */
41    
42     /**
43     * @file SimInfo.hpp
44     * @author tlin
45     * @date 11/02/2004
46     * @version 1.0
47     */
48 gezelter 2
49 gezelter 246 #ifndef BRAINS_SIMMODEL_HPP
50     #define BRAINS_SIMMODEL_HPP
51    
52     #include <iostream>
53     #include <set>
54     #include <utility>
55 gezelter 2 #include <vector>
56    
57 gezelter 1287 #include "brains/PairList.hpp"
58 gezelter 246 #include "io/Globals.hpp"
59     #include "math/Vector3.hpp"
60 chuckv 555 #include "math/SquareMatrix3.hpp"
61 gezelter 246 #include "types/MoleculeStamp.hpp"
62     #include "UseTheForce/ForceField.hpp"
63     #include "utils/PropertyMap.hpp"
64     #include "utils/LocalIndexManager.hpp"
65 gezelter 1530 #include "nonbonded/SwitchingFunction.hpp"
66 tim 316
67 gezelter 246 //another nonsense macro declaration
68 gezelter 1390 #define __OPENMD_C
69 tim 3 #include "brains/fSimulation.h"
70 gezelter 2
71 gezelter 1528 using namespace std;
72 gezelter 1390 namespace OpenMD{
73 gezelter 1528 enum CutoffMethod {
74     HARD,
75     SWITCHING_FUNCTION,
76     SHIFTED_POTENTIAL,
77     SHIFTED_FORCE
78     };
79    
80 gezelter 507 //forward decalration
81     class SnapshotManager;
82     class Molecule;
83     class SelectionManager;
84 tim 1024 class StuntDouble;
85 gezelter 1528
86 gezelter 507 /**
87 gezelter 1528 * @class SimInfo SimInfo.hpp "brains/SimInfo.hpp"
88     *
89     * @brief One of the heavy-weight classes of OpenMD, SimInfo
90     * maintains objects and variables relating to the current
91     * simulation. This includes the master list of Molecules. The
92     * Molecule class maintains all of the concrete objects (Atoms,
93     * Bond, Bend, Torsions, Inversions, RigidBodies, CutoffGroups,
94     * Constraints). In both the single and parallel versions, Atoms and
95     * RigidBodies have both global and local indices.
96     */
97 gezelter 507 class SimInfo {
98     public:
99 gezelter 1528 typedef map<int, Molecule*>::iterator MoleculeIterator;
100    
101 gezelter 507 /**
102     * Constructor of SimInfo
103 gezelter 1528 *
104     * @param molStampPairs MoleculeStamp Array. The first element of
105     * the pair is molecule stamp, the second element is the total
106     * number of molecules with the same molecule stamp in the system
107     *
108 gezelter 507 * @param ff pointer of a concrete ForceField instance
109 gezelter 1528 *
110 gezelter 507 * @param simParams
111     */
112 tim 770 SimInfo(ForceField* ff, Globals* simParams);
113 gezelter 507 virtual ~SimInfo();
114 gezelter 2
115 gezelter 507 /**
116     * Adds a molecule
117 gezelter 1528 *
118     * @return return true if adding successfully, return false if the
119     * molecule is already in SimInfo
120     *
121 gezelter 507 * @param mol molecule to be added
122     */
123     bool addMolecule(Molecule* mol);
124 gezelter 2
125 gezelter 507 /**
126     * Removes a molecule from SimInfo
127 gezelter 1528 *
128     * @return true if removing successfully, return false if molecule
129     * is not in this SimInfo
130 gezelter 507 */
131     bool removeMolecule(Molecule* mol);
132 gezelter 2
133 gezelter 507 /** Returns the total number of molecules in the system. */
134     int getNGlobalMolecules() {
135     return nGlobalMols_;
136     }
137 gezelter 2
138 gezelter 507 /** Returns the total number of atoms in the system. */
139     int getNGlobalAtoms() {
140     return nGlobalAtoms_;
141     }
142 gezelter 2
143 gezelter 507 /** Returns the total number of cutoff groups in the system. */
144     int getNGlobalCutoffGroups() {
145     return nGlobalCutoffGroups_;
146     }
147 gezelter 2
148 gezelter 507 /**
149 gezelter 1528 * Returns the total number of integrable objects (total number of
150     * rigid bodies plus the total number of atoms which do not belong
151     * to the rigid bodies) in the system
152 gezelter 507 */
153     int getNGlobalIntegrableObjects() {
154     return nGlobalIntegrableObjects_;
155     }
156 gezelter 2
157 gezelter 507 /**
158 gezelter 1528 * Returns the total number of integrable objects (total number of
159     * rigid bodies plus the total number of atoms which do not belong
160     * to the rigid bodies) in the system
161 gezelter 507 */
162     int getNGlobalRigidBodies() {
163     return nGlobalRigidBodies_;
164     }
165 gezelter 2
166 gezelter 507 int getNGlobalConstraints();
167     /**
168     * Returns the number of local molecules.
169     * @return the number of local molecules
170     */
171     int getNMolecules() {
172     return molecules_.size();
173     }
174 gezelter 2
175 gezelter 507 /** Returns the number of local atoms */
176     unsigned int getNAtoms() {
177     return nAtoms_;
178     }
179 gezelter 2
180 gezelter 507 /** Returns the number of local bonds */
181     unsigned int getNBonds(){
182     return nBonds_;
183     }
184 gezelter 2
185 gezelter 507 /** Returns the number of local bends */
186     unsigned int getNBends() {
187     return nBends_;
188     }
189 gezelter 2
190 gezelter 507 /** Returns the number of local torsions */
191     unsigned int getNTorsions() {
192     return nTorsions_;
193     }
194 gezelter 2
195 gezelter 1277 /** Returns the number of local torsions */
196     unsigned int getNInversions() {
197     return nInversions_;
198     }
199 gezelter 507 /** Returns the number of local rigid bodies */
200     unsigned int getNRigidBodies() {
201     return nRigidBodies_;
202     }
203 gezelter 2
204 gezelter 507 /** Returns the number of local integrable objects */
205     unsigned int getNIntegrableObjects() {
206     return nIntegrableObjects_;
207     }
208 gezelter 2
209 gezelter 507 /** Returns the number of local cutoff groups */
210     unsigned int getNCutoffGroups() {
211     return nCutoffGroups_;
212     }
213 gezelter 2
214 gezelter 507 /** Returns the total number of constraints in this SimInfo */
215     unsigned int getNConstraints() {
216     return nConstraints_;
217     }
218 gezelter 246
219 gezelter 507 /**
220     * Returns the first molecule in this SimInfo and intialize the iterator.
221     * @return the first molecule, return NULL if there is not molecule in this SimInfo
222     * @param i the iterator of molecule array (user shouldn't change it)
223     */
224     Molecule* beginMolecule(MoleculeIterator& i);
225 gezelter 2
226 gezelter 507 /**
227     * Returns the next avaliable Molecule based on the iterator.
228     * @return the next avaliable molecule, return NULL if reaching the end of the array
229     * @param i the iterator of molecule array
230     */
231     Molecule* nextMolecule(MoleculeIterator& i);
232 gezelter 2
233 gezelter 507 /** Returns the number of degrees of freedom */
234     int getNdf() {
235 gezelter 945 return ndf_ - getFdf();
236 gezelter 507 }
237 gezelter 2
238 gezelter 507 /** Returns the number of raw degrees of freedom */
239     int getNdfRaw() {
240     return ndfRaw_;
241     }
242 gezelter 2
243 gezelter 507 /** Returns the number of translational degrees of freedom */
244     int getNdfTrans() {
245     return ndfTrans_;
246     }
247 gezelter 2
248 gezelter 945 /** sets the current number of frozen degrees of freedom */
249     void setFdf(int fdf) {
250     fdf_local = fdf;
251     }
252    
253     int getFdf();
254    
255 gezelter 1528 //getNZconstraint and setNZconstraint ruin the coherence of
256     //SimInfo class, need refactoring
257 gezelter 246
258 gezelter 507 /** Returns the total number of z-constraint molecules in the system */
259     int getNZconstraint() {
260     return nZconstraint_;
261     }
262 gezelter 2
263 gezelter 507 /**
264     * Sets the number of z-constraint molecules in the system.
265     */
266     void setNZconstraint(int nZconstraint) {
267     nZconstraint_ = nZconstraint;
268     }
269 gezelter 246
270 gezelter 507 /** Returns the snapshot manager. */
271     SnapshotManager* getSnapshotManager() {
272     return sman_;
273     }
274 gezelter 2
275 gezelter 507 /** Sets the snapshot manager. */
276     void setSnapshotManager(SnapshotManager* sman);
277 gezelter 246
278 gezelter 507 /** Returns the force field */
279     ForceField* getForceField() {
280     return forceField_;
281     }
282 gezelter 2
283 gezelter 507 Globals* getSimParams() {
284     return simParams_;
285     }
286 gezelter 2
287 gezelter 507 /** Returns the velocity of center of mass of the whole system.*/
288     Vector3d getComVel();
289 gezelter 2
290 gezelter 507 /** Returns the center of the mass of the whole system.*/
291     Vector3d getCom();
292 gezelter 1528 /** Returns the center of the mass and Center of Mass velocity of
293     the whole system.*/
294 chuckv 555 void getComAll(Vector3d& com,Vector3d& comVel);
295 gezelter 2
296 gezelter 1528 /** Returns intertia tensor for the entire system and system
297     Angular Momentum.*/
298 chuckv 555 void getInertiaTensor(Mat3x3d &intertiaTensor,Vector3d &angularMomentum);
299    
300     /** Returns system angular momentum */
301     Vector3d getAngularMomentum();
302    
303 gezelter 1528 /** Returns volume of system as estimated by an ellipsoid defined
304     by the radii of gyration*/
305 chuckv 1103 void getGyrationalVolume(RealType &vol);
306 gezelter 1528 /** Overloaded version of gyrational volume that also returns
307     det(I) so dV/dr can be calculated*/
308 chuckv 1103 void getGyrationalVolume(RealType &vol, RealType &detI);
309 gezelter 1528 /** main driver function to interact with fortran during the
310     initialization and molecule migration */
311 gezelter 507 void update();
312 gezelter 2
313 gezelter 507 /** Returns the local index manager */
314     LocalIndexManager* getLocalIndexManager() {
315     return &localIndexMan_;
316     }
317 gezelter 2
318 gezelter 507 int getMoleculeStampId(int globalIndex) {
319     //assert(globalIndex < molStampIds_.size())
320     return molStampIds_[globalIndex];
321     }
322 gezelter 2
323 gezelter 507 /** Returns the molecule stamp */
324     MoleculeStamp* getMoleculeStamp(int id) {
325     return moleculeStamps_[id];
326     }
327 gezelter 2
328 gezelter 507 /** Return the total number of the molecule stamps */
329     int getNMoleculeStamp() {
330     return moleculeStamps_.size();
331     }
332     /**
333     * Finds a molecule with a specified global index
334     * @return a pointer point to found molecule
335     * @param index
336     */
337     Molecule* getMoleculeByGlobalIndex(int index) {
338     MoleculeIterator i;
339     i = molecules_.find(index);
340 gezelter 2
341 gezelter 507 return i != molecules_.end() ? i->second : NULL;
342     }
343 gezelter 2
344 chuckv 1292 int getGlobalMolMembership(int id){
345     return globalMolMembership_[id];
346     }
347    
348 gezelter 1528 RealType getCutoffRadius() {
349     return cutoffRadius_;
350 gezelter 507 }
351 gezelter 2
352 gezelter 1528 RealType getSwitchingRadius() {
353     return switchingRadius_;
354 gezelter 507 }
355 gezelter 764
356 gezelter 1528 RealType getListRadius() {
357     return listRadius_;
358 gezelter 764 }
359 gezelter 246
360 gezelter 1528 string getFinalConfigFileName() {
361 gezelter 507 return finalConfigFileName_;
362     }
363 tim 1024
364 gezelter 1528 void setFinalConfigFileName(const string& fileName) {
365 gezelter 507 finalConfigFileName_ = fileName;
366     }
367 gezelter 2
368 gezelter 1528 string getRawMetaData() {
369 tim 1024 return rawMetaData_;
370     }
371 gezelter 1528 void setRawMetaData(const string& rawMetaData) {
372 tim 1024 rawMetaData_ = rawMetaData;
373     }
374    
375 gezelter 1528 string getDumpFileName() {
376 gezelter 507 return dumpFileName_;
377     }
378 gezelter 246
379 gezelter 1528 void setDumpFileName(const string& fileName) {
380 gezelter 507 dumpFileName_ = fileName;
381     }
382 gezelter 2
383 gezelter 1528 string getStatFileName() {
384 gezelter 507 return statFileName_;
385     }
386 gezelter 246
387 gezelter 1528 void setStatFileName(const string& fileName) {
388 gezelter 507 statFileName_ = fileName;
389     }
390 chrisfen 417
391 gezelter 1528 string getRestFileName() {
392 gezelter 507 return restFileName_;
393     }
394 chrisfen 417
395 gezelter 1528 void setRestFileName(const string& fileName) {
396 gezelter 507 restFileName_ = fileName;
397     }
398 gezelter 2
399 gezelter 507 /**
400     * Sets GlobalGroupMembership
401     * @see #SimCreator::setGlobalIndex
402     */
403 gezelter 1528 void setGlobalGroupMembership(const vector<int>& globalGroupMembership) {
404 gezelter 1287 assert(globalGroupMembership.size() == static_cast<size_t>(nGlobalAtoms_));
405 gezelter 507 globalGroupMembership_ = globalGroupMembership;
406     }
407 gezelter 2
408 gezelter 507 /**
409     * Sets GlobalMolMembership
410     * @see #SimCreator::setGlobalIndex
411     */
412 gezelter 1528 void setGlobalMolMembership(const vector<int>& globalMolMembership) {
413 gezelter 1287 assert(globalMolMembership.size() == static_cast<size_t>(nGlobalAtoms_));
414 gezelter 507 globalMolMembership_ = globalMolMembership;
415     }
416 gezelter 246
417    
418 gezelter 507 bool isFortranInitialized() {
419     return fortranInitialized_;
420     }
421 gezelter 246
422 chrisfen 998 bool getCalcBoxDipole() {
423     return calcBoxDipole_;
424     }
425    
426 gezelter 1126 bool getUseAtomicVirial() {
427     return useAtomicVirial_;
428     }
429    
430 gezelter 507 /**
431     * Adds property into property map
432     * @param genData GenericData to be added into PropertyMap
433     */
434     void addProperty(GenericData* genData);
435 gezelter 246
436 gezelter 507 /**
437     * Removes property from PropertyMap by name
438     * @param propName the name of property to be removed
439     */
440 gezelter 1528 void removeProperty(const string& propName);
441 gezelter 246
442 gezelter 507 /**
443     * clear all of the properties
444     */
445     void clearProperties();
446 gezelter 246
447 gezelter 507 /**
448     * Returns all names of properties
449     * @return all names of properties
450     */
451 gezelter 1528 vector<string> getPropertyNames();
452 gezelter 246
453 gezelter 507 /**
454     * Returns all of the properties in PropertyMap
455     * @return all of the properties in PropertyMap
456     */
457 gezelter 1528 vector<GenericData*> getProperties();
458 gezelter 246
459 gezelter 507 /**
460     * Returns property
461     * @param propName name of property
462     * @return a pointer point to property with propName. If no property named propName
463     * exists, return NULL
464     */
465 gezelter 1528 GenericData* getPropertyByName(const string& propName);
466 gezelter 246
467 gezelter 507 /**
468 gezelter 1287 * add all special interaction pairs (including excluded
469     * interactions) in a molecule into the appropriate lists.
470 gezelter 507 */
471 gezelter 1287 void addInteractionPairs(Molecule* mol);
472 gezelter 246
473 gezelter 507 /**
474 gezelter 1287 * remove all special interaction pairs which belong to a molecule
475     * from the appropriate lists.
476 gezelter 507 */
477 gezelter 1287 void removeInteractionPairs(Molecule* mol);
478 gezelter 246
479 tim 292
480 gezelter 507 /** Returns the unique atom types of local processor in an array */
481 gezelter 1528 set<AtomType*> getUniqueAtomTypes();
482    
483     /** Returns the set of atom types present in this simulation */
484     set<AtomType*> getSimulatedAtomTypes();
485 tim 292
486 gezelter 1528 friend ostream& operator <<(ostream& o, SimInfo& info);
487 tim 326
488 tim 963 void getCutoff(RealType& rcut, RealType& rsw);
489 gezelter 246
490 gezelter 507 private:
491 gezelter 246
492 gezelter 1530 /** fill up the simtype struct and other simulation-related variables */
493     void setupSimVariables();
494 gezelter 246
495 gezelter 507 /**
496     * Setup Fortran Simulation
497     * @see #setupFortranParallel
498     */
499     void setupFortranSim();
500 gezelter 246
501 gezelter 1530 /** Figure out the cutoff radius and cutoff method */
502     void setupCutoffs();
503     /** Figure out the switching radius and polynomial type for the switching function */
504     void setupSwitching();
505     /** Figure out the simulation variables associated with electrostatics */
506     void setupElectrostatics();
507 gezelter 1528 /** Figure out the neighbor list skin thickness */
508 gezelter 1530 void setupNeighborlists();
509 chrisfen 726
510 chrisfen 998 /** Determine if we need to accumulate the simulation box dipole */
511     void setupAccumulateBoxDipole();
512    
513 gezelter 507 /** Calculates the number of degress of freedom in the whole system */
514     void calcNdf();
515     void calcNdfRaw();
516     void calcNdfTrans();
517 gezelter 246
518 gezelter 507 /**
519 gezelter 1528 * Adds molecule stamp and the total number of the molecule with
520     * same molecule stamp in the whole system.
521 gezelter 507 */
522     void addMoleculeStamp(MoleculeStamp* molStamp, int nmol);
523 gezelter 246
524 gezelter 1528 // Other classes holdingn important information
525     ForceField* forceField_; /**< provides access to defined atom types, bond types, etc. */
526     Globals* simParams_; /**< provides access to simulation parameters set by user */
527 gezelter 246
528 gezelter 1528 /// Counts of local objects
529 gezelter 1277 int nAtoms_; /**< number of atoms in local processor */
530     int nBonds_; /**< number of bonds in local processor */
531     int nBends_; /**< number of bends in local processor */
532     int nTorsions_; /**< number of torsions in local processor */
533     int nInversions_; /**< number of inversions in local processor */
534     int nRigidBodies_; /**< number of rigid bodies in local processor */
535     int nIntegrableObjects_; /**< number of integrable objects in local processor */
536     int nCutoffGroups_; /**< number of cutoff groups in local processor */
537     int nConstraints_; /**< number of constraints in local processors */
538 gezelter 1528
539     /// Counts of global objects
540     int nGlobalMols_; /**< number of molecules in the system (GLOBAL) */
541     int nGlobalAtoms_; /**< number of atoms in the system (GLOBAL) */
542     int nGlobalCutoffGroups_; /**< number of cutoff groups in this system (GLOBAL) */
543     int nGlobalIntegrableObjects_; /**< number of integrable objects in this system */
544     int nGlobalRigidBodies_; /**< number of rigid bodies in this system (GLOBAL) */
545    
546     /// Degress of freedom
547     int ndf_; /**< number of degress of freedom (excludes constraints) (LOCAL) */
548     int fdf_local; /**< number of frozen degrees of freedom (LOCAL) */
549     int fdf_; /**< number of frozen degrees of freedom (GLOBAL) */
550     int ndfRaw_; /**< number of degress of freedom (includes constraints), (LOCAL) */
551     int ndfTrans_; /**< number of translation degress of freedom, (LOCAL) */
552     int nZconstraint_; /**< number of z-constraint molecules (GLOBAL) */
553 gezelter 246
554 gezelter 1528 /// logicals
555     bool usesPeriodicBoundaries_; /**< use periodic boundary conditions? */
556     bool usesDirectionalAtoms_; /**< are there atoms with position AND orientation? */
557     bool usesMetallicAtoms_; /**< are there transition metal atoms? */
558     bool usesElectrostaticAtoms_; /**< are there electrostatic atoms? */
559     bool usesAtomicVirial_; /**< are we computing atomic virials? */
560     bool requiresPrepair_; /**< does this simulation require a pre-pair loop? */
561     bool requiresSkipCorrection_; /**< does this simulation require a skip-correction? */
562     bool requiresSelfCorrection_; /**< does this simulation require a self-correction? */
563 gezelter 246
564 gezelter 1528 /// Data structures holding primary simulation objects
565     map<int, Molecule*> molecules_; /**< map holding pointers to LOCAL molecules */
566     simtype fInfo_; /**< A dual struct shared by C++
567     and Fortran to pass
568     information about what types
569     of calculation are
570     required */
571    
572     /// Stamps are templates for objects that are then used to create
573     /// groups of objects. For example, a molecule stamp contains
574     /// information on how to build that molecule (i.e. the topology,
575     /// the atoms, the bonds, etc.) Once the system is built, the
576     /// stamps are no longer useful.
577     vector<int> molStampIds_; /**< stamp id for molecules in the system */
578     vector<MoleculeStamp*> moleculeStamps_; /**< molecule stamps array */
579    
580     /**
581     * A vector that maps between the global index of an atom, and the
582     * global index of cutoff group the atom belong to. It is filled
583     * by SimCreator once and only once, since it never changed during
584     * the simulation. It should be nGlobalAtoms_ in size.
585     */
586     vector<int> globalGroupMembership_;
587    
588     /**
589     * A vector that maps between the global index of an atom and the
590     * global index of the molecule the atom belongs to. It is filled
591     * by SimCreator once and only once, since it is never changed
592     * during the simulation. It shoudl be nGlobalAtoms_ in size.
593     */
594     vector<int> globalMolMembership_;
595    
596     /// lists to handle atoms needing special treatment in the non-bonded interactions
597     PairList excludedInteractions_; /**< atoms excluded from interacting with each other */
598     PairList oneTwoInteractions_; /**< atoms that are directly Bonded */
599     PairList oneThreeInteractions_; /**< atoms sharing a Bend */
600     PairList oneFourInteractions_; /**< atoms sharing a Torsion */
601    
602     PropertyMap properties_; /**< Generic Properties can be added */
603     SnapshotManager* sman_; /**< SnapshotManager (handles particle positions, etc.) */
604    
605 gezelter 507 /**
606 gezelter 1528 * The reason to have a local index manager is that when molecule
607     * is migrating to other processors, the atoms and the
608     * rigid-bodies will release their local indices to
609     * LocalIndexManager. Combining the information of molecule
610     * migrating to current processor, Migrator class can query the
611     * LocalIndexManager to make a efficient data moving plan.
612 gezelter 507 */
613     LocalIndexManager localIndexMan_;
614 gezelter 246
615 tim 1024 // unparsed MetaData block for storing in Dump and EOR files:
616 gezelter 1528 string rawMetaData_;
617 tim 1024
618 gezelter 1528 // file names
619     string finalConfigFileName_;
620     string dumpFileName_;
621     string statFileName_;
622     string restFileName_;
623 chrisfen 417
624 gezelter 1528 RealType cutoffRadius_; /**< cutoff radius for non-bonded interactions */
625     RealType switchingRadius_; /**< inner radius of switching function */
626     RealType listRadius_; /**< Verlet neighbor list radius */
627     RealType skinThickness_; /**< Verlet neighbor list skin thickness */
628     CutoffMethod cutoffMethod_; /**< Cutoff Method for most non-bonded interactions */
629 gezelter 246
630 gezelter 1467 bool fortranInitialized_; /** flag to indicate whether the fortran side is initialized */
631 gezelter 1126
632     bool calcBoxDipole_; /**< flag to indicate whether or not we calculate
633     the simulation box dipole moment */
634    
635     bool useAtomicVirial_; /**< flag to indicate whether or not we use
636     Atomic Virials to calculate the pressure */
637 gezelter 1528
638     public:
639     /**
640     * return an integral objects by its global index. In MPI
641     * version, if the StuntDouble with specified global index does
642     * not belong to local processor, a NULL will be return.
643 tim 1024 */
644 gezelter 1528 StuntDouble* getIOIndexToIntegrableObject(int index);
645     void setIOIndexToIntegrableObject(const vector<StuntDouble*>& v);
646 tim 1024
647 gezelter 1528 private:
648     vector<StuntDouble*> IOIndexToIntegrableObject;
649    
650 gezelter 507 public:
651 gezelter 246
652 gezelter 507 /**
653     * Finds the processor where a molecule resides
654     * @return the id of the processor which contains the molecule
655     * @param globalIndex global Index of the molecule
656     */
657     int getMolToProc(int globalIndex) {
658     //assert(globalIndex < molToProcMap_.size());
659     return molToProcMap_[globalIndex];
660     }
661 gezelter 1528
662 gezelter 507 /**
663     * Set MolToProcMap array
664     * @see #SimCreator::divideMolecules
665     */
666 gezelter 1528 void setMolToProcMap(const vector<int>& molToProcMap) {
667 gezelter 507 molToProcMap_ = molToProcMap;
668     }
669 gezelter 246
670 gezelter 507 private:
671 gezelter 246
672 gezelter 507 void setupFortranParallel();
673 gezelter 246
674 gezelter 507 /**
675 gezelter 1241 * The size of molToProcMap_ is equal to total number of molecules
676     * in the system. It maps a molecule to the processor on which it
677     * resides. it is filled by SimCreator once and only once.
678 gezelter 507 */
679 gezelter 1528 vector<int> molToProcMap_;
680 tim 292
681 gezelter 507 };
682 gezelter 2
683 gezelter 1390 } //namespace OpenMD
684 gezelter 246 #endif //BRAINS_SIMMODEL_HPP
685 gezelter 2

Properties

Name Value
svn:keywords Author Id Revision Date