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: 1907
Committed: Thu Jan 6 22:31:07 2005 UTC (19 years, 8 months ago) by tim
File size: 18553 byte(s)
Log Message:
constraint is almost working

File Contents

# User Rev Content
1 tim 1842 /*
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     * @param simParams
79     * @note
80     */
81     SimInfo(std::vector<std::pair<MoleculeStamp*, int> >& molStampPairs, ForceField* ff, Globals* simParams);
82     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    
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 1907
128     int getNGlobalConstraints();
129 tim 1842 /**
130     * Returns the number of local molecules.
131     * @return the number of local molecules
132     */
133     int getNMolecules() {
134     return molecules_.size();
135     }
136    
137     /** Returns the number of local atoms */
138     unsigned int getNAtoms() {
139     return nAtoms_;
140     }
141    
142     /** Returns the number of local bonds */
143     unsigned int getNBonds(){
144     return nBonds_;
145     }
146    
147     /** Returns the number of local bends */
148     unsigned int getNBends() {
149     return nBends_;
150     }
151    
152     /** Returns the number of local torsions */
153     unsigned int getNTorsions() {
154     return nTorsions_;
155     }
156    
157     /** Returns the number of local rigid bodies */
158     unsigned int getNRigidBodies() {
159     return nRigidBodies_;
160     }
161    
162     /** Returns the number of local integrable objects */
163     unsigned int getNIntegrableObjects() {
164     return nIntegrableObjects_;
165     }
166    
167     /** Returns the number of local cutoff groups */
168     unsigned int getNCutoffGroups() {
169     return nCutoffGroups_;
170     }
171    
172     /** Returns the total number of constraints in this SimInfo */
173     unsigned int getNConstraints() {
174     return nConstraints_;
175     }
176    
177     /**
178     * Returns the first molecule in this SimInfo and intialize the iterator.
179     * @return the first molecule, return NULL if there is not molecule in this SimInfo
180     * @param i the iterator of molecule array (user shouldn't change it)
181     */
182     Molecule* beginMolecule(MoleculeIterator& i);
183    
184     /**
185     * Returns the next avaliable Molecule based on the iterator.
186     * @return the next avaliable molecule, return NULL if reaching the end of the array
187     * @param i the iterator of molecule array
188     */
189     Molecule* nextMolecule(MoleculeIterator& i);
190    
191     /** Returns the number of degrees of freedom */
192     int getNdf() {
193     return ndf_;
194     }
195    
196     /** Returns the number of raw degrees of freedom */
197     int getNdfRaw() {
198     return ndfRaw_;
199     }
200    
201     /** Returns the number of translational degrees of freedom */
202     int getNdfTrans() {
203     return ndfTrans_;
204     }
205    
206     //getNZconstraint and setNZconstraint ruin the coherent of SimInfo class, need refactorying
207    
208     /** Returns the total number of z-constraint molecules in the system */
209     int getNZconstraint() {
210     return nZconstraint_;
211     }
212    
213     /**
214     * Sets the number of z-constraint molecules in the system.
215     */
216     void setNZconstraint(int nZconstraint) {
217     nZconstraint_ = nZconstraint;
218     }
219    
220     /** Returns the snapshot manager. */
221     SnapshotManager* getSnapshotManager() {
222     return sman_;
223     }
224    
225     /** Sets the snapshot manager. */
226     void setSnapshotManager(SnapshotManager* sman);
227    
228     /** Returns the force field */
229     ForceField* getForceField() {
230     return forceField_;
231     }
232    
233     Globals* getSimParams() {
234     return simParams_;
235     }
236    
237     /** Returns the velocity of center of mass of the whole system.*/
238     Vector3d getComVel();
239    
240     /** Returns the center of the mass of the whole system.*/
241     Vector3d getCom();
242    
243     /** Returns the seed (used for random number generator) */
244     int getSeed() {
245     return seed_;
246     }
247    
248     /** Sets the seed*/
249     void setSeed(int seed) {
250     seed_ = seed;
251     }
252    
253     /** main driver function to interact with fortran during the initialization and molecule migration */
254     void update();
255    
256     /** Returns the local index manager */
257     LocalIndexManager* getLocalIndexManager() {
258     return &localIndexMan_;
259     }
260    
261     int getMoleculeStampId(int globalIndex) {
262     //assert(globalIndex < molStampIds_.size())
263     return molStampIds_[globalIndex];
264     }
265    
266     /** Returns the molecule stamp */
267     MoleculeStamp* getMoleculeStamp(int id) {
268     return moleculeStamps_[id];
269     }
270 tim 1903
271     /** Return the total number of the molecule stamps */
272     int getNMoleculeStamp() {
273     return moleculeStamps_.size();
274     }
275 tim 1842 /**
276     * Finds a molecule with a specified global index
277     * @return a pointer point to found molecule
278     * @param index
279     */
280     Molecule* getMoleculeByGlobalIndex(int index) {
281     MoleculeIterator i;
282     i = molecules_.find(index);
283    
284     return i != molecules_.end() ? i->second : NULL;
285     }
286    
287     /** Calculate the maximum cutoff radius based on the atom types */
288     double calcMaxCutoffRadius();
289    
290     double getRcut() {
291     return rcut_;
292     }
293    
294     double getRsw() {
295     return rsw_;
296     }
297    
298     std::string getFinalConfigFileName() {
299     return finalConfigFileName_;
300     }
301    
302     void setFinalConfigFileName(const std::string& fileName) {
303     finalConfigFileName_ = fileName;
304     }
305    
306     std::string getDumpFileName() {
307     return dumpFileName_;
308     }
309    
310     void setDumpFileName(const std::string& fileName) {
311     dumpFileName_ = fileName;
312     }
313    
314     std::string getStatFileName() {
315     return statFileName_;
316     }
317    
318     void setStatFileName(const std::string& fileName) {
319     statFileName_ = fileName;
320     }
321    
322     /**
323     * Sets GlobalGroupMembership
324     * @see #SimCreator::setGlobalIndex
325     */
326     void setGlobalGroupMembership(const std::vector<int>& globalGroupMembership) {
327     assert(globalGroupMembership.size() == nGlobalAtoms_);
328 tim 1843 globalGroupMembership_ = globalGroupMembership;
329 tim 1842 }
330    
331     /**
332     * Sets GlobalMolMembership
333     * @see #SimCreator::setGlobalIndex
334     */
335     void setGlobalMolMembership(const std::vector<int>& globalMolMembership) {
336     assert(globalMolMembership.size() == nGlobalAtoms_);
337     globalMolMembership_ = globalMolMembership;
338     }
339    
340    
341     bool isFortranInitialized() {
342     return fortranInitialized_;
343     }
344    
345     //below functions are just forward functions
346     //To compose or to inherit is always a hot debate. In general, is-a relation need subclassing, in the
347     //the other hand, has-a relation need composing.
348     /**
349     * Adds property into property map
350     * @param genData GenericData to be added into PropertyMap
351     */
352     void addProperty(GenericData* genData);
353    
354     /**
355     * Removes property from PropertyMap by name
356     * @param propName the name of property to be removed
357     */
358     void removeProperty(const std::string& propName);
359    
360     /**
361     * clear all of the properties
362     */
363     void clearProperties();
364    
365     /**
366     * Returns all names of properties
367     * @return all names of properties
368     */
369     std::vector<std::string> getPropertyNames();
370    
371     /**
372     * Returns all of the properties in PropertyMap
373     * @return all of the properties in PropertyMap
374     */
375     std::vector<GenericData*> getProperties();
376    
377     /**
378     * Returns property
379     * @param propName name of property
380     * @return a pointer point to property with propName. If no property named propName
381     * exists, return NULL
382     */
383     GenericData* getPropertyByName(const std::string& propName);
384 tim 1856
385     /**
386     * add all exclude pairs of a molecule into exclude list.
387     */
388     void addExcludePairs(Molecule* mol);
389    
390     /**
391     * remove all exclude pairs which belong to a molecule from exclude list
392     */
393    
394     void removeExcludePairs(Molecule* mol);
395 tim 1842
396     friend std::ostream& operator <<(std::ostream& o, SimInfo& info);
397    
398     private:
399    
400    
401     /** Returns the unique atom types of local processor in an array */
402     std::set<AtomType*> getUniqueAtomTypes();
403    
404     /** fill up the simtype struct*/
405     void setupSimType();
406    
407     /**
408     * Setup Fortran Simulation
409     * @see #setupFortranParallel
410     */
411     void setupFortranSim();
412    
413     /** Figure out the radius of cutoff, radius of switching function and pass them to fortran */
414     void setupCutoff();
415    
416     /** Calculates the number of degress of freedom in the whole system */
417     void calcNdf();
418     void calcNdfRaw();
419     void calcNdfTrans();
420    
421     /**
422     * Adds molecule stamp and the total number of the molecule with same molecule stamp in the whole
423     * system.
424     */
425     void addMoleculeStamp(MoleculeStamp* molStamp, int nmol);
426    
427     ForceField* forceField_;
428     Globals* simParams_;
429    
430     std::map<int, Molecule*> molecules_; /**< Molecule array */
431    
432     //degress of freedom
433     int ndf_; /**< number of degress of freedom (excludes constraints), ndf_ is local */
434     int ndfRaw_; /**< number of degress of freedom (includes constraints), ndfRaw_ is local */
435     int ndfTrans_; /**< number of translation degress of freedom, ndfTrans_ is local */
436     int nZconstraint_; /** number of z-constraint molecules, nZconstraint_ is global */
437    
438     //number of global objects
439     int nGlobalMols_; /**< number of molecules in the system */
440     int nGlobalAtoms_; /**< number of atoms in the system */
441     int nGlobalCutoffGroups_; /**< number of cutoff groups in this system */
442     int nGlobalIntegrableObjects_; /**< number of integrable objects in this system */
443     int nGlobalRigidBodies_; /**< number of rigid bodies in this system */
444     /**
445     * the size of globalGroupMembership_ is nGlobalAtoms. Its index is global index of an atom, and the
446     * corresponding content is the global index of cutoff group this atom belong to.
447 tim 1867 * It is filled by SimCreator once and only once, since it never changed during the simulation.
448 tim 1842 */
449     std::vector<int> globalGroupMembership_;
450    
451     /**
452     * the size of globalGroupMembership_ is nGlobalAtoms. Its index is global index of an atom, and the
453     * corresponding content is the global index of molecule this atom belong to.
454     * It is filled by SimCreator once and only once, since it is never changed during the simulation.
455     */
456     std::vector<int> globalMolMembership_;
457    
458    
459     std::vector<int> molStampIds_; /**< stamp id array of all molecules in the system */
460     std::vector<MoleculeStamp*> moleculeStamps_; /**< molecule stamps array */
461    
462     //number of local objects
463     int nAtoms_; /**< number of atoms in local processor */
464     int nBonds_; /**< number of bonds in local processor */
465     int nBends_; /**< number of bends in local processor */
466     int nTorsions_; /**< number of torsions in local processor */
467     int nRigidBodies_; /**< number of rigid bodies in local processor */
468     int nIntegrableObjects_; /**< number of integrable objects in local processor */
469     int nCutoffGroups_; /**< number of cutoff groups in local processor */
470     int nConstraints_; /**< number of constraints in local processors */
471    
472     simtype fInfo_; /**< A dual struct shared by c++/fortran which indicates the atom types in simulation*/
473     Exclude exclude_;
474     PropertyMap properties_; /**< Generic Property */
475     SnapshotManager* sman_; /**< SnapshotManager */
476    
477     int seed_; /**< seed for random number generator */
478    
479     /**
480     * The reason to have a local index manager is that when molecule is migrating to other processors,
481     * the atoms and the rigid-bodies will release their local indices to LocalIndexManager. Combining the
482     * information of molecule migrating to current processor, Migrator class can query the LocalIndexManager
483     * to make a efficient data moving plan.
484     */
485     LocalIndexManager localIndexMan_;
486    
487     //file names
488     std::string finalConfigFileName_;
489     std::string dumpFileName_;
490     std::string statFileName_;
491    
492     double rcut_; /**< cutoff radius*/
493     double rsw_; /**< radius of switching function*/
494    
495     bool fortranInitialized_; /**< flag indicate whether fortran side is initialized */
496    
497     #ifdef IS_MPI
498     //in Parallel version, we need MolToProc
499     public:
500    
501     /**
502     * Finds the processor where a molecule resides
503     * @return the id of the processor which contains the molecule
504     * @param globalIndex global Index of the molecule
505     */
506     int getMolToProc(int globalIndex) {
507     //assert(globalIndex < molToProcMap_.size());
508     return molToProcMap_[globalIndex];
509     }
510    
511     /**
512     * Set MolToProcMap array
513     * @see #SimCreator::divideMolecules
514     */
515     void setMolToProcMap(const std::vector<int>& molToProcMap) {
516     molToProcMap_ = molToProcMap;
517     }
518    
519     private:
520    
521     void setupFortranParallel();
522    
523     /**
524     * The size of molToProcMap_ is equal to total number of molecules in the system.
525     * It maps a molecule to the processor on which it resides. it is filled by SimCreator once and only
526     * once.
527     */
528     std::vector<int> molToProcMap_;
529     #endif
530    
531     };
532    
533     } //namespace oopse
534     #endif //BRAINS_SIMMODEL_HPP
535