ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/branches/new_design/OOPSE-3.0/src/brains/SimInfo.hpp
Revision: 1841
Committed: Fri Dec 3 17:59:45 2004 UTC (19 years, 9 months ago) by tim
File size: 18681 byte(s)
Log Message:
begin to fix bugs

File Contents

# User Rev Content
1 tim 1804 /*
2     * Copyright (C) 2000-2004 Object Oriented Parallel Simulation Engine (OOPSE) project
3     *
4     * Contact: oopse@oopse.org
5     *
6     * This program is free software; you can redistribute it and/or
7     * modify it under the terms of the GNU Lesser General Public License
8     * as published by the Free Software Foundation; either version 2.1
9     * of the License, or (at your option) any later version.
10     * All we ask is that proper credit is given for our work, which includes
11     * - but is not limited to - adding the above copyright notice to the beginning
12     * of your source code files, and to any copyright notice that you may distribute
13     * with programs based on this work.
14     *
15     * This program is distributed in the hope that it will be useful,
16     * but WITHOUT ANY WARRANTY; without even the implied warranty of
17     * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18     * GNU Lesser General Public License for more details.
19     *
20     * You should have received a copy of the GNU Lesser General Public License
21     * along with this program; if not, write to the Free Software
22     * Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
23     *
24     */
25    
26     /**
27     * @file SimInfo.hpp
28     * @author tlin
29     * @date 11/02/2004
30     * @version 1.0
31     */
32    
33     #ifndef BRAINS_SIMMODEL_HPP
34     #define BRAINS_SIMMODEL_HPP
35    
36     #include <iostream>
37     #include <set>
38     #include <utility>
39     #include <vector>
40    
41     #include "brains/Exclude.hpp"
42     #include "io/Globals.hpp"
43     #include "math/Vector3.hpp"
44     #include "types/MoleculeStamp.hpp"
45     #include "UseTheForce/ForceField.hpp"
46     #include "utils/PropertyMap.hpp"
47     #include "utils/LocalIndexManager.hpp"
48    
49     //another nonsense macro declaration
50     #define __C
51     #include "brains/fSimulation.h"
52    
53     namespace oopse{
54    
55     //forward decalration
56     class SnapshotManager;
57     class Molecule;
58    
59     /**
60     * @class SimInfo SimInfo.hpp "brains/SimInfo.hpp"
61     * @brief As one of the heavy weight class of OOPSE, SimInfo
62     * One of the major changes in SimInfo class is the data struct. It only maintains a list of molecules.
63     * And the Molecule class will maintain all of the concrete objects (atoms, bond, bend, torsions, rigid bodies,
64     * cutoff groups, constrains).
65     * Another major change is the index. No matter single version or parallel version, atoms and
66     * rigid bodies have both global index and local index. Local index is not important to molecule as well as
67     * cutoff group.
68     */
69     class SimInfo {
70     public:
71     typedef std::map<int, Molecule*>::iterator MoleculeIterator;
72    
73     /**
74     * Constructor of SimInfo
75     * @param molStampPairs MoleculeStamp Array. The first element of the pair is molecule stamp, the
76     * second element is the total number of molecules with the same molecule stamp in the system
77     * @param ff pointer of a concrete ForceField instance
78 tim 1841 * @param simParams
79 tim 1804 * @note
80     */
81 tim 1841 SimInfo(std::vector<std::pair<MoleculeStamp*, int> >& molStampPairs, ForceField* ff, Globals* simParams);
82 tim 1804 virtual ~SimInfo();
83    
84     /**
85     * Adds a molecule
86     * @return return true if adding successfully, return false if the molecule is already in SimInfo
87     * @param mol molecule to be added
88     */
89     bool addMolecule(Molecule* mol);
90    
91     /**
92     * Removes a molecule from SimInfo
93     * @return true if removing successfully, return false if molecule is not in this SimInfo
94     */
95     bool removeMolecule(Molecule* mol);
96    
97     /** Returns the total number of molecules in the system. */
98     int getNGlobalMolecules() {
99     return nGlobalMols_;
100     }
101    
102     /** Returns the total number of atoms in the system. */
103     int getNGlobalAtoms() {
104     return nGlobalAtoms_;
105     }
106    
107     /** Returns the total number of cutoff groups in the system. */
108     int getNGlobalCutoffGroups() {
109     return nGlobalCutoffGroups_;
110     }
111    
112     /**
113     * Returns the total number of integrable objects (total number of rigid bodies plus the total number
114     * of atoms which do not belong to the rigid bodies) in the system
115     */
116     int getNGlobalIntegrableObjects() {
117     return nGlobalIntegrableObjects_;
118     }
119 tim 1841
120     /**
121     * Returns the total number of integrable objects (total number of rigid bodies plus the total number
122     * of atoms which do not belong to the rigid bodies) in the system
123     */
124     int getNGlobalRigidBodies() {
125     return nGlobalRigidBodies_;
126     }
127 tim 1804
128     /**
129     * Returns the number of local molecules.
130     * @return the number of local molecules
131     */
132     int getNMolecules() {
133     return molecules_.size();
134     }
135    
136     /** Returns the number of local atoms */
137     unsigned int getNAtoms() {
138     return nAtoms_;
139     }
140    
141     /** Returns the number of local bonds */
142     unsigned int getNBonds(){
143     return nBonds_;
144     }
145    
146     /** Returns the number of local bends */
147     unsigned int getNBends() {
148     return nBends_;
149     }
150    
151     /** Returns the number of local torsions */
152     unsigned int getNTorsions() {
153     return nTorsions_;
154     }
155    
156     /** Returns the number of local rigid bodies */
157     unsigned int getNRigidBodies() {
158     return nRigidBodies_;
159     }
160    
161     /** Returns the number of local integrable objects */
162     unsigned int getNIntegrableObjects() {
163     return nIntegrableObjects_;
164     }
165    
166     /** Returns the number of local cutoff groups */
167     unsigned int getNCutoffGroups() {
168     return nCutoffGroups_;
169     }
170    
171     /** Returns the total number of constraints in this SimInfo */
172     unsigned int getNConstraints() {
173     return nConstraints_;
174     }
175    
176     /**
177     * Returns the first molecule in this SimInfo and intialize the iterator.
178     * @return the first molecule, return NULL if there is not molecule in this SimInfo
179     * @param i the iterator of molecule array (user shouldn't change it)
180     */
181     Molecule* beginMolecule(MoleculeIterator& i);
182    
183     /**
184     * Returns the next avaliable Molecule based on the iterator.
185     * @return the next avaliable molecule, return NULL if reaching the end of the array
186     * @param i the iterator of molecule array
187     */
188     Molecule* nextMolecule(MoleculeIterator& i);
189    
190     /** Returns the number of degrees of freedom */
191     int getNdf() {
192     return ndf_;
193     }
194    
195     /** Returns the number of raw degrees of freedom */
196     int getNdfRaw() {
197     return ndfRaw_;
198     }
199    
200     /** Returns the number of translational degrees of freedom */
201     int getNdfTrans() {
202     return ndfTrans_;
203     }
204    
205     //getNZconstraint and setNZconstraint ruin the coherent of SimInfo class, need refactorying
206    
207     /** Returns the total number of z-constraint molecules in the system */
208     int getNZconstraint() {
209     return nZconstraint_;
210     }
211    
212     /**
213     * Sets the number of z-constraint molecules in the system.
214     */
215     void setNZconstraint(int nZconstraint) {
216     nZconstraint_ = nZconstraint;
217     }
218    
219     /** Returns the snapshot manager. */
220     SnapshotManager* getSnapshotManager() {
221     return sman_;
222     }
223    
224     /** Sets the snapshot manager. */
225 tim 1807 void setSnapshotManager(SnapshotManager* sman);
226    
227 tim 1804 /** Returns the force field */
228     ForceField* getForceField() {
229     return forceField_;
230     }
231    
232 tim 1841 Globals* getSimParams() {
233     return simParams_;
234 tim 1804 }
235    
236     /** Returns the velocity of center of mass of the whole system.*/
237     Vector3d getComVel();
238    
239     /** Returns the center of the mass of the whole system.*/
240     Vector3d getCom();
241    
242     /** Returns the seed (used for random number generator) */
243     int getSeed() {
244     return seed_;
245     }
246    
247     /** Sets the seed*/
248     void setSeed(int seed) {
249     seed_ = seed;
250     }
251    
252     /** main driver function to interact with fortran during the initialization and molecule migration */
253     void update();
254    
255     /** Returns the local index manager */
256     LocalIndexManager* getLocalIndexManager() {
257     return &localIndexMan_;
258     }
259    
260     int getMoleculeStampId(int globalIndex) {
261     //assert(globalIndex < molStampIds_.size())
262     return molStampIds_[globalIndex];
263     }
264    
265     /** Returns the molecule stamp */
266     MoleculeStamp* getMoleculeStamp(int id) {
267     return moleculeStamps_[id];
268     }
269    
270     /**
271     * Finds a molecule with a specified global index
272     * @return a pointer point to found molecule
273     * @param index
274     */
275     Molecule* getMoleculeByGlobalIndex(int index) {
276     MoleculeIterator i;
277     i = molecules_.find(index);
278    
279     return i != molecules_.end() ? i->second : NULL;
280     }
281    
282     /** Calculate the maximum cutoff radius based on the atom types */
283     double calcMaxCutoffRadius();
284    
285     double getRcut() {
286     return rcut_;
287     }
288    
289     double getRsw() {
290     return rsw_;
291     }
292    
293     std::string getFinalConfigFileName() {
294     return finalConfigFileName_;
295     }
296    
297     void setFinalConfigFileName(const std::string& fileName) {
298     finalConfigFileName_ = fileName;
299     }
300    
301     std::string getDumpFileName() {
302     return dumpFileName_;
303     }
304    
305     void setDumpFileName(const std::string& fileName) {
306     dumpFileName_ = fileName;
307     }
308    
309     std::string getStatFileName() {
310     return statFileName_;
311     }
312    
313     void setStatFileName(const std::string& fileName) {
314     statFileName_ = fileName;
315     }
316    
317     /**
318     * Returns the pointer of internal globalGroupMembership_ array. This array will be filled by SimCreator class
319     * @see #SimCreator::setGlobalIndex
320     */
321     int* getGlobalGroupMembershipPointer() {
322     return &globalGroupMembership_[0];
323     }
324    
325     /**
326     * Returns the pointer of internal globalMolMembership_ array. This array will be filled by SimCreator class
327     * @see #SimCreator::setGlobalIndex
328     */
329     int* getGlobalMolMembershipPointer() {
330     return &globalMolMembership_[0];
331     }
332    
333    
334     bool isFortranInitialized() {
335     return fortranInitialized_;
336     }
337    
338     //below functions are just forward functions
339     //To compose or to inherit is always a hot debate. In general, is-a relation need subclassing, in the
340     //the other hand, has-a relation need composing.
341     /**
342     * Adds property into property map
343     * @param genData GenericData to be added into PropertyMap
344     */
345     void addProperty(GenericData* genData);
346    
347     /**
348     * Removes property from PropertyMap by name
349     * @param propName the name of property to be removed
350     */
351     void removeProperty(const std::string& propName);
352    
353     /**
354     * clear all of the properties
355     */
356     void clearProperties();
357    
358     /**
359     * Returns all names of properties
360     * @return all names of properties
361     */
362     std::vector<std::string> getPropertyNames();
363    
364     /**
365     * Returns all of the properties in PropertyMap
366     * @return all of the properties in PropertyMap
367     */
368     std::vector<GenericData*> getProperties();
369    
370     /**
371     * Returns property
372     * @param propName name of property
373     * @return a pointer point to property with propName. If no property named propName
374     * exists, return NULL
375     */
376     GenericData* getPropertyByName(const std::string& propName);
377    
378 tim 1832 friend std::ostream& operator <<(std::ostream& o, SimInfo& info);
379 tim 1804
380     private:
381    
382    
383     /** Returns the unique atom types of local processor in an array */
384     std::set<AtomType*> getUniqueAtomTypes();
385    
386     /** fill up the simtype struct*/
387     void setupSimType();
388    
389     /**
390     * Setup Fortran Simulation
391     * @see #setupFortranParallel
392     */
393     void setupFortranSim();
394    
395     /** Figure out the radius of cutoff, radius of switching function and pass them to fortran */
396     void setupCutoff();
397    
398     /** Calculates the number of degress of freedom in the whole system */
399     void calcNdf();
400     void calcNdfRaw();
401     void calcNdfTrans();
402    
403     void addExcludePairs(Molecule* mol);
404     void removeExcludePairs(Molecule* mol);
405    
406     /**
407     * Adds molecule stamp and the total number of the molecule with same molecule stamp in the whole
408     * system.
409     */
410     void addMoleculeStamp(MoleculeStamp* molStamp, int nmol);
411    
412 tim 1841 ForceField* forceField_;
413     Globals* simParams_;
414    
415 tim 1804 std::map<int, Molecule*> molecules_; /**< Molecule array */
416    
417     //degress of freedom
418     int ndf_; /**< number of degress of freedom (excludes constraints), ndf_ is local */
419     int ndfRaw_; /**< number of degress of freedom (includes constraints), ndfRaw_ is local */
420     int ndfTrans_; /**< number of translation degress of freedom, ndfTrans_ is local */
421     int nZconstraint_; /** number of z-constraint molecules, nZconstraint_ is global */
422    
423     //number of global objects
424     int nGlobalMols_; /**< number of molecules in the system */
425     int nGlobalAtoms_; /**< number of atoms in the system */
426     int nGlobalCutoffGroups_; /**< number of cutoff groups in this system */
427     int nGlobalIntegrableObjects_; /**< number of integrable objects in this system */
428 tim 1841 int nGlobalRigidBodies_; /**< number of rigid bodies in this system */
429 tim 1804 /**
430     * the size of globalGroupMembership_ is nGlobalAtoms. Its index is global index of an atom, and the
431     * corresponding content is the global index of cutoff group this atom belong to.
432     * It is filled by SimCreator once and only once, since it is never changed during the simulation.
433     */
434     std::vector<int> globalGroupMembership_;
435    
436     /**
437     * the size of globalGroupMembership_ is nGlobalAtoms. Its index is global index of an atom, and the
438     * corresponding content is the global index of molecule this atom belong to.
439     * It is filled by SimCreator once and only once, since it is never changed during the simulation.
440     */
441     std::vector<int> globalMolMembership_;
442    
443    
444     std::vector<int> molStampIds_; /**< stamp id array of all molecules in the system */
445     std::vector<MoleculeStamp*> moleculeStamps_; /**< molecule stamps array */
446    
447     //number of local objects
448     int nAtoms_; /**< number of atoms in local processor */
449     int nBonds_; /**< number of bonds in local processor */
450     int nBends_; /**< number of bends in local processor */
451     int nTorsions_; /**< number of torsions in local processor */
452     int nRigidBodies_; /**< number of rigid bodies in local processor */
453     int nIntegrableObjects_; /**< number of integrable objects in local processor */
454     int nCutoffGroups_; /**< number of cutoff groups in local processor */
455     int nConstraints_; /**< number of constraints in local processors */
456    
457     simtype fInfo_; /**< A dual struct shared by c++/fortran which indicates the atom types in simulation*/
458 tim 1841 Exclude exclude_;
459 tim 1804 PropertyMap properties_; /**< Generic Property */
460     SnapshotManager* sman_; /**< SnapshotManager */
461 tim 1841
462 tim 1804 int seed_; /**< seed for random number generator */
463    
464     /**
465     * The reason to have a local index manager is that when molecule is migrating to other processors,
466     * the atoms and the rigid-bodies will release their local indices to LocalIndexManager. Combining the
467     * information of molecule migrating to current processor, Migrator class can query the LocalIndexManager
468     * to make a efficient data moving plan.
469     */
470     LocalIndexManager localIndexMan_;
471    
472     //file names
473     std::string finalConfigFileName_;
474     std::string dumpFileName_;
475     std::string statFileName_;
476    
477     double rcut_; /**< cutoff radius*/
478     double rsw_; /**< radius of switching function*/
479    
480     bool fortranInitialized_; /**< flag indicate whether fortran side is initialized */
481    
482     #ifdef IS_MPI
483     //in Parallel version, we need MolToProc
484     public:
485    
486     /**
487     * Finds the processor where a molecule resides
488     * @return the id of the processor which contains the molecule
489     * @param globalIndex global Index of the molecule
490     */
491     int getMolToProc(int globalIndex) {
492     //assert(globalIndex < molToProcMap_.size());
493     return molToProcMap_[globalIndex];
494     }
495    
496     /**
497     * Returns the pointer of internal molToProcMap array. This array will be filled by SimCreator class
498     * @see #SimCreator::divideMolecules
499     */
500     int* getMolToProcMapPointer() {
501     return &molToProcMap_[0];
502     }
503    
504     private:
505    
506     void setupFortranParallel();
507    
508     /**
509     * The size of molToProcMap_ is equal to total number of molecules in the system.
510     * It maps a molecule to the processor on which it resides. it is filled by SimCreator once and only
511     * once.
512     */
513     std::vector<int> molToProcMap_;
514     #endif
515    
516     };
517    
518     } //namespace oopse
519     #endif //BRAINS_SIMMODEL_HPP