ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/brains/SimInfo.hpp
Revision: 2759
Committed: Wed May 17 21:51:42 2006 UTC (18 years, 4 months ago) by tim
File size: 18776 byte(s)
Log Message:
Adding single precision capabilities to c++ side

File Contents

# User Rev Content
1 gezelter 2204 /*
2 gezelter 1930 * 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     * 1. Acknowledgement of the program authors must be made in any
10     * publication of scientific results based in part on use of the
11     * program. An acceptable form of acknowledgement is citation of
12     * the article in which the program was described (Matthew
13     * A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14     * J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15     * Parallel Simulation Engine for Molecular Dynamics,"
16     * J. Comput. Chem. 26, pp. 252-271 (2005))
17     *
18     * 2. Redistributions of source code must retain the above copyright
19     * notice, this list of conditions and the following disclaimer.
20     *
21     * 3. Redistributions in binary form must reproduce the above copyright
22     * notice, this list of conditions and the following disclaimer in the
23     * documentation and/or other materials provided with the
24     * distribution.
25     *
26     * This software is provided "AS IS," without a warranty of any
27     * kind. All express or implied conditions, representations and
28     * warranties, including any implied warranty of merchantability,
29     * fitness for a particular purpose or non-infringement, are hereby
30     * excluded. The University of Notre Dame and its licensors shall not
31     * be liable for any damages suffered by licensee as a result of
32     * using, modifying or distributing the software or its
33     * derivatives. In no event will the University of Notre Dame or its
34     * licensors be liable for any lost revenue, profit or data, or for
35     * direct, indirect, special, consequential, incidental or punitive
36     * damages, however caused and regardless of the theory of liability,
37     * arising out of the use of or inability to use software, even if the
38     * University of Notre Dame has been advised of the possibility of
39     * such damages.
40     */
41    
42     /**
43     * @file SimInfo.hpp
44     * @author tlin
45     * @date 11/02/2004
46     * @version 1.0
47     */
48 gezelter 1490
49 gezelter 1930 #ifndef BRAINS_SIMMODEL_HPP
50     #define BRAINS_SIMMODEL_HPP
51    
52     #include <iostream>
53     #include <set>
54     #include <utility>
55 gezelter 1490 #include <vector>
56    
57 tim 1492 #include "brains/Exclude.hpp"
58 gezelter 1930 #include "io/Globals.hpp"
59     #include "math/Vector3.hpp"
60 chuckv 2252 #include "math/SquareMatrix3.hpp"
61 gezelter 1930 #include "types/MoleculeStamp.hpp"
62     #include "UseTheForce/ForceField.hpp"
63     #include "utils/PropertyMap.hpp"
64     #include "utils/LocalIndexManager.hpp"
65 tim 2000
66 gezelter 1930 //another nonsense macro declaration
67 gezelter 1490 #define __C
68 tim 1492 #include "brains/fSimulation.h"
69 gezelter 1490
70 gezelter 1930 namespace oopse{
71 gezelter 1490
72 gezelter 2204 //forward decalration
73     class SnapshotManager;
74     class Molecule;
75     class SelectionManager;
76     /**
77     * @class SimInfo SimInfo.hpp "brains/SimInfo.hpp"
78 gezelter 2733 * @brief One of the heavy weight classes of OOPSE, SimInfo maintains a list of molecules.
79     * The Molecule class maintains all of the concrete objects
80     * (atoms, bond, bend, torsions, rigid bodies, cutoff groups, constrains).
81     * In both the single and parallel versions, atoms and
82     * rigid bodies have both global and local indices. The local index is
83     * not relevant to molecules or cutoff groups.
84 gezelter 2204 */
85     class SimInfo {
86     public:
87     typedef std::map<int, Molecule*>::iterator MoleculeIterator;
88 gezelter 1490
89 gezelter 2204 /**
90     * Constructor of SimInfo
91     * @param molStampPairs MoleculeStamp Array. The first element of the pair is molecule stamp, the
92     * second element is the total number of molecules with the same molecule stamp in the system
93     * @param ff pointer of a concrete ForceField instance
94     * @param simParams
95     * @note
96     */
97 tim 2469 SimInfo(ForceField* ff, Globals* simParams);
98 gezelter 2204 virtual ~SimInfo();
99 gezelter 1490
100 gezelter 2204 /**
101     * Adds a molecule
102     * @return return true if adding successfully, return false if the molecule is already in SimInfo
103     * @param mol molecule to be added
104     */
105     bool addMolecule(Molecule* mol);
106 gezelter 1490
107 gezelter 2204 /**
108     * Removes a molecule from SimInfo
109     * @return true if removing successfully, return false if molecule is not in this SimInfo
110     */
111     bool removeMolecule(Molecule* mol);
112 gezelter 1490
113 gezelter 2204 /** Returns the total number of molecules in the system. */
114     int getNGlobalMolecules() {
115     return nGlobalMols_;
116     }
117 gezelter 1490
118 gezelter 2204 /** Returns the total number of atoms in the system. */
119     int getNGlobalAtoms() {
120     return nGlobalAtoms_;
121     }
122 gezelter 1490
123 gezelter 2204 /** Returns the total number of cutoff groups in the system. */
124     int getNGlobalCutoffGroups() {
125     return nGlobalCutoffGroups_;
126     }
127 gezelter 1490
128 gezelter 2204 /**
129     * Returns the total number of integrable objects (total number of rigid bodies plus the total number
130     * of atoms which do not belong to the rigid bodies) in the system
131     */
132     int getNGlobalIntegrableObjects() {
133     return nGlobalIntegrableObjects_;
134     }
135 gezelter 1490
136 gezelter 2204 /**
137     * Returns the total number of integrable objects (total number of rigid bodies plus the total number
138     * of atoms which do not belong to the rigid bodies) in the system
139     */
140     int getNGlobalRigidBodies() {
141     return nGlobalRigidBodies_;
142     }
143 gezelter 1490
144 gezelter 2204 int getNGlobalConstraints();
145     /**
146     * Returns the number of local molecules.
147     * @return the number of local molecules
148     */
149     int getNMolecules() {
150     return molecules_.size();
151     }
152 gezelter 1490
153 gezelter 2204 /** Returns the number of local atoms */
154     unsigned int getNAtoms() {
155     return nAtoms_;
156     }
157 gezelter 1490
158 gezelter 2204 /** Returns the number of local bonds */
159     unsigned int getNBonds(){
160     return nBonds_;
161     }
162 gezelter 1490
163 gezelter 2204 /** Returns the number of local bends */
164     unsigned int getNBends() {
165     return nBends_;
166     }
167 gezelter 1490
168 gezelter 2204 /** Returns the number of local torsions */
169     unsigned int getNTorsions() {
170     return nTorsions_;
171     }
172 gezelter 1490
173 gezelter 2204 /** Returns the number of local rigid bodies */
174     unsigned int getNRigidBodies() {
175     return nRigidBodies_;
176     }
177 gezelter 1490
178 gezelter 2204 /** Returns the number of local integrable objects */
179     unsigned int getNIntegrableObjects() {
180     return nIntegrableObjects_;
181     }
182 gezelter 1490
183 gezelter 2204 /** Returns the number of local cutoff groups */
184     unsigned int getNCutoffGroups() {
185     return nCutoffGroups_;
186     }
187 gezelter 1490
188 gezelter 2204 /** Returns the total number of constraints in this SimInfo */
189     unsigned int getNConstraints() {
190     return nConstraints_;
191     }
192 gezelter 1930
193 gezelter 2204 /**
194     * Returns the first molecule in this SimInfo and intialize the iterator.
195     * @return the first molecule, return NULL if there is not molecule in this SimInfo
196     * @param i the iterator of molecule array (user shouldn't change it)
197     */
198     Molecule* beginMolecule(MoleculeIterator& i);
199 gezelter 1490
200 gezelter 2204 /**
201     * Returns the next avaliable Molecule based on the iterator.
202     * @return the next avaliable molecule, return NULL if reaching the end of the array
203     * @param i the iterator of molecule array
204     */
205     Molecule* nextMolecule(MoleculeIterator& i);
206 gezelter 1490
207 gezelter 2204 /** Returns the number of degrees of freedom */
208     int getNdf() {
209 gezelter 2733 return ndf_ - getFdf();
210 gezelter 2204 }
211 gezelter 1490
212 gezelter 2204 /** Returns the number of raw degrees of freedom */
213     int getNdfRaw() {
214     return ndfRaw_;
215     }
216 gezelter 1490
217 gezelter 2204 /** Returns the number of translational degrees of freedom */
218     int getNdfTrans() {
219     return ndfTrans_;
220     }
221 gezelter 1490
222 gezelter 2733 /** sets the current number of frozen degrees of freedom */
223     void setFdf(int fdf) {
224     fdf_local = fdf;
225     }
226    
227     int getFdf();
228    
229 gezelter 2204 //getNZconstraint and setNZconstraint ruin the coherent of SimInfo class, need refactorying
230 gezelter 1930
231 gezelter 2204 /** Returns the total number of z-constraint molecules in the system */
232     int getNZconstraint() {
233     return nZconstraint_;
234     }
235 gezelter 1490
236 gezelter 2204 /**
237     * Sets the number of z-constraint molecules in the system.
238     */
239     void setNZconstraint(int nZconstraint) {
240     nZconstraint_ = nZconstraint;
241     }
242 gezelter 1930
243 gezelter 2204 /** Returns the snapshot manager. */
244     SnapshotManager* getSnapshotManager() {
245     return sman_;
246     }
247 gezelter 1490
248 gezelter 2204 /** Sets the snapshot manager. */
249     void setSnapshotManager(SnapshotManager* sman);
250 gezelter 1930
251 gezelter 2204 /** Returns the force field */
252     ForceField* getForceField() {
253     return forceField_;
254     }
255 gezelter 1490
256 gezelter 2204 Globals* getSimParams() {
257     return simParams_;
258     }
259 gezelter 1490
260 gezelter 2204 /** Returns the velocity of center of mass of the whole system.*/
261     Vector3d getComVel();
262 gezelter 1490
263 gezelter 2204 /** Returns the center of the mass of the whole system.*/
264     Vector3d getCom();
265 chuckv 2252 /** Returns the center of the mass and Center of Mass velocity of the whole system.*/
266     void getComAll(Vector3d& com,Vector3d& comVel);
267 gezelter 1490
268 chuckv 2252 /** Returns intertia tensor for the entire system and system Angular Momentum.*/
269     void getInertiaTensor(Mat3x3d &intertiaTensor,Vector3d &angularMomentum);
270    
271     /** Returns system angular momentum */
272     Vector3d getAngularMomentum();
273    
274 gezelter 2204 /** main driver function to interact with fortran during the initialization and molecule migration */
275     void update();
276 gezelter 1490
277 gezelter 2204 /** Returns the local index manager */
278     LocalIndexManager* getLocalIndexManager() {
279     return &localIndexMan_;
280     }
281 gezelter 1490
282 gezelter 2204 int getMoleculeStampId(int globalIndex) {
283     //assert(globalIndex < molStampIds_.size())
284     return molStampIds_[globalIndex];
285     }
286 gezelter 1490
287 gezelter 2204 /** Returns the molecule stamp */
288     MoleculeStamp* getMoleculeStamp(int id) {
289     return moleculeStamps_[id];
290     }
291 gezelter 1490
292 gezelter 2204 /** Return the total number of the molecule stamps */
293     int getNMoleculeStamp() {
294     return moleculeStamps_.size();
295     }
296     /**
297     * Finds a molecule with a specified global index
298     * @return a pointer point to found molecule
299     * @param index
300     */
301     Molecule* getMoleculeByGlobalIndex(int index) {
302     MoleculeIterator i;
303     i = molecules_.find(index);
304 gezelter 1490
305 gezelter 2204 return i != molecules_.end() ? i->second : NULL;
306     }
307 gezelter 1490
308 tim 2759 RealType getRcut() {
309 gezelter 2204 return rcut_;
310     }
311 gezelter 1490
312 tim 2759 RealType getRsw() {
313 gezelter 2204 return rsw_;
314     }
315 gezelter 2463
316 tim 2759 RealType getList() {
317 gezelter 2463 return rlist_;
318     }
319 gezelter 1930
320 gezelter 2204 std::string getFinalConfigFileName() {
321     return finalConfigFileName_;
322     }
323 gezelter 1930
324 gezelter 2204 void setFinalConfigFileName(const std::string& fileName) {
325     finalConfigFileName_ = fileName;
326     }
327 gezelter 1490
328 gezelter 2204 std::string getDumpFileName() {
329     return dumpFileName_;
330     }
331 gezelter 1930
332 gezelter 2204 void setDumpFileName(const std::string& fileName) {
333     dumpFileName_ = fileName;
334     }
335 gezelter 1490
336 gezelter 2204 std::string getStatFileName() {
337     return statFileName_;
338     }
339 gezelter 1930
340 gezelter 2204 void setStatFileName(const std::string& fileName) {
341     statFileName_ = fileName;
342     }
343 chrisfen 2101
344 gezelter 2204 std::string getRestFileName() {
345     return restFileName_;
346     }
347 chrisfen 2101
348 gezelter 2204 void setRestFileName(const std::string& fileName) {
349     restFileName_ = fileName;
350     }
351 gezelter 1490
352 gezelter 2204 /**
353     * Sets GlobalGroupMembership
354     * @see #SimCreator::setGlobalIndex
355     */
356     void setGlobalGroupMembership(const std::vector<int>& globalGroupMembership) {
357     assert(globalGroupMembership.size() == nGlobalAtoms_);
358     globalGroupMembership_ = globalGroupMembership;
359     }
360 gezelter 1490
361 gezelter 2204 /**
362     * Sets GlobalMolMembership
363     * @see #SimCreator::setGlobalIndex
364     */
365     void setGlobalMolMembership(const std::vector<int>& globalMolMembership) {
366     assert(globalMolMembership.size() == nGlobalAtoms_);
367     globalMolMembership_ = globalMolMembership;
368     }
369 gezelter 1930
370    
371 gezelter 2204 bool isFortranInitialized() {
372     return fortranInitialized_;
373     }
374 gezelter 1930
375 gezelter 2204 //below functions are just forward functions
376     //To compose or to inherit is always a hot debate. In general, is-a relation need subclassing, in the
377     //the other hand, has-a relation need composing.
378     /**
379     * Adds property into property map
380     * @param genData GenericData to be added into PropertyMap
381     */
382     void addProperty(GenericData* genData);
383 gezelter 1930
384 gezelter 2204 /**
385     * Removes property from PropertyMap by name
386     * @param propName the name of property to be removed
387     */
388     void removeProperty(const std::string& propName);
389 gezelter 1930
390 gezelter 2204 /**
391     * clear all of the properties
392     */
393     void clearProperties();
394 gezelter 1930
395 gezelter 2204 /**
396     * Returns all names of properties
397     * @return all names of properties
398     */
399     std::vector<std::string> getPropertyNames();
400 gezelter 1930
401 gezelter 2204 /**
402     * Returns all of the properties in PropertyMap
403     * @return all of the properties in PropertyMap
404     */
405     std::vector<GenericData*> getProperties();
406 gezelter 1930
407 gezelter 2204 /**
408     * Returns property
409     * @param propName name of property
410     * @return a pointer point to property with propName. If no property named propName
411     * exists, return NULL
412     */
413     GenericData* getPropertyByName(const std::string& propName);
414 gezelter 1930
415 gezelter 2204 /**
416     * add all exclude pairs of a molecule into exclude list.
417     */
418     void addExcludePairs(Molecule* mol);
419 gezelter 1930
420 gezelter 2204 /**
421     * remove all exclude pairs which belong to a molecule from exclude list
422     */
423 gezelter 1930
424 gezelter 2204 void removeExcludePairs(Molecule* mol);
425 tim 1976
426    
427 gezelter 2204 /** Returns the unique atom types of local processor in an array */
428     std::set<AtomType*> getUniqueAtomTypes();
429 tim 1976
430 gezelter 2204 friend std::ostream& operator <<(std::ostream& o, SimInfo& info);
431 tim 2010
432 tim 2759 void getCutoff(RealType& rcut, RealType& rsw);
433 gezelter 1930
434 gezelter 2204 private:
435 gezelter 1930
436 gezelter 2204 /** fill up the simtype struct*/
437     void setupSimType();
438 gezelter 1930
439 gezelter 2204 /**
440     * Setup Fortran Simulation
441     * @see #setupFortranParallel
442     */
443     void setupFortranSim();
444 gezelter 1930
445 gezelter 2204 /** Figure out the radius of cutoff, radius of switching function and pass them to fortran */
446     void setupCutoff();
447 gezelter 1930
448 chrisfen 2297 /** Figure out which coulombic correction method to use and pass to fortran */
449 chrisfen 2303 void setupElectrostaticSummationMethod( int isError );
450 chrisfen 2297
451 chrisfen 2425 /** Figure out which polynomial type to use for the switching function */
452     void setupSwitchingFunction();
453    
454 gezelter 2204 /** Calculates the number of degress of freedom in the whole system */
455     void calcNdf();
456     void calcNdfRaw();
457     void calcNdfTrans();
458 gezelter 1930
459 tim 2469 ForceField* forceField_;
460     Globals* simParams_;
461    
462     std::map<int, Molecule*> molecules_; /**< Molecule array */
463    
464 gezelter 2204 /**
465     * Adds molecule stamp and the total number of the molecule with same molecule stamp in the whole
466     * system.
467     */
468     void addMoleculeStamp(MoleculeStamp* molStamp, int nmol);
469 gezelter 1930
470 gezelter 2204 //degress of freedom
471     int ndf_; /**< number of degress of freedom (excludes constraints), ndf_ is local */
472 gezelter 2733 int fdf_local; /**< number of frozen degrees of freedom */
473     int fdf_; /**< number of frozen degrees of freedom */
474 gezelter 2204 int ndfRaw_; /**< number of degress of freedom (includes constraints), ndfRaw_ is local */
475     int ndfTrans_; /**< number of translation degress of freedom, ndfTrans_ is local */
476     int nZconstraint_; /** number of z-constraint molecules, nZconstraint_ is global */
477 gezelter 1930
478 gezelter 2204 //number of global objects
479     int nGlobalMols_; /**< number of molecules in the system */
480     int nGlobalAtoms_; /**< number of atoms in the system */
481     int nGlobalCutoffGroups_; /**< number of cutoff groups in this system */
482     int nGlobalIntegrableObjects_; /**< number of integrable objects in this system */
483     int nGlobalRigidBodies_; /**< number of rigid bodies in this system */
484     /**
485     * the size of globalGroupMembership_ is nGlobalAtoms. Its index is global index of an atom, and the
486     * corresponding content is the global index of cutoff group this atom belong to.
487     * It is filled by SimCreator once and only once, since it never changed during the simulation.
488     */
489     std::vector<int> globalGroupMembership_;
490 gezelter 1930
491 gezelter 2204 /**
492     * the size of globalGroupMembership_ is nGlobalAtoms. Its index is global index of an atom, and the
493     * corresponding content is the global index of molecule this atom belong to.
494     * It is filled by SimCreator once and only once, since it is never changed during the simulation.
495     */
496     std::vector<int> globalMolMembership_;
497 gezelter 1930
498    
499 gezelter 2204 std::vector<int> molStampIds_; /**< stamp id array of all molecules in the system */
500     std::vector<MoleculeStamp*> moleculeStamps_; /**< molecule stamps array */
501 gezelter 1930
502 gezelter 2204 //number of local objects
503     int nAtoms_; /**< number of atoms in local processor */
504     int nBonds_; /**< number of bonds in local processor */
505     int nBends_; /**< number of bends in local processor */
506     int nTorsions_; /**< number of torsions in local processor */
507     int nRigidBodies_; /**< number of rigid bodies in local processor */
508     int nIntegrableObjects_; /**< number of integrable objects in local processor */
509     int nCutoffGroups_; /**< number of cutoff groups in local processor */
510     int nConstraints_; /**< number of constraints in local processors */
511 gezelter 1930
512 gezelter 2204 simtype fInfo_; /**< A dual struct shared by c++/fortran which indicates the atom types in simulation*/
513     Exclude exclude_;
514     PropertyMap properties_; /**< Generic Property */
515     SnapshotManager* sman_; /**< SnapshotManager */
516 gezelter 1930
517 gezelter 2204 /**
518     * The reason to have a local index manager is that when molecule is migrating to other processors,
519     * the atoms and the rigid-bodies will release their local indices to LocalIndexManager. Combining the
520     * information of molecule migrating to current processor, Migrator class can query the LocalIndexManager
521     * to make a efficient data moving plan.
522     */
523     LocalIndexManager localIndexMan_;
524 gezelter 1930
525 gezelter 2204 //file names
526     std::string finalConfigFileName_;
527     std::string dumpFileName_;
528     std::string statFileName_;
529     std::string restFileName_;
530 chrisfen 2101
531 tim 2759 RealType rcut_; /**< cutoff radius*/
532     RealType rsw_; /**< radius of switching function*/
533     RealType rlist_; /**< neighbor list radius */
534 gezelter 1930
535 gezelter 2204 bool fortranInitialized_; /**< flag indicate whether fortran side is initialized */
536 tim 1976
537 gezelter 1930 #ifdef IS_MPI
538     //in Parallel version, we need MolToProc
539 gezelter 2204 public:
540 gezelter 1930
541 gezelter 2204 /**
542     * Finds the processor where a molecule resides
543     * @return the id of the processor which contains the molecule
544     * @param globalIndex global Index of the molecule
545     */
546     int getMolToProc(int globalIndex) {
547     //assert(globalIndex < molToProcMap_.size());
548     return molToProcMap_[globalIndex];
549     }
550 gezelter 1930
551 gezelter 2204 /**
552     * Set MolToProcMap array
553     * @see #SimCreator::divideMolecules
554     */
555     void setMolToProcMap(const std::vector<int>& molToProcMap) {
556     molToProcMap_ = molToProcMap;
557     }
558 gezelter 1930
559 gezelter 2204 private:
560 gezelter 1930
561 gezelter 2204 void setupFortranParallel();
562 gezelter 1930
563 gezelter 2204 /**
564     * The size of molToProcMap_ is equal to total number of molecules in the system.
565     * It maps a molecule to the processor on which it resides. it is filled by SimCreator once and only
566     * once.
567     */
568     std::vector<int> molToProcMap_;
569 tim 1976
570 gezelter 1930 #endif
571    
572 gezelter 2204 };
573 gezelter 1490
574 gezelter 1930 } //namespace oopse
575     #endif //BRAINS_SIMMODEL_HPP
576 gezelter 1490