ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/brains/SimInfo.hpp
Revision: 2982
Committed: Wed Aug 30 18:42:29 2006 UTC (18 years ago) by tim
File size: 20186 byte(s)
Log Message:
Massive changes preparing for release of OOPSE-4: The main difference
is that the new MD file format (.md, .dump, .eor) now contains meta-data
information along with the configuration information.

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