ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.hpp
Revision: 555
Committed: Mon May 30 14:01:52 2005 UTC (19 years, 11 months ago) by chuckv
Original Path: trunk/src/brains/SimInfo.hpp
File size: 18444 byte(s)
Log Message:
Added method to remove system angular momentum to velocitizer and added method to calculate system angular momentum to siminfo.

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     * 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 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 tim 3 #include "brains/Exclude.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 tim 316
66 gezelter 246 //another nonsense macro declaration
67 gezelter 2 #define __C
68 tim 3 #include "brains/fSimulation.h"
69 gezelter 2
70 gezelter 246 namespace oopse{
71 gezelter 2
72 gezelter 507 //forward decalration
73     class SnapshotManager;
74     class Molecule;
75     class SelectionManager;
76     /**
77     * @class SimInfo SimInfo.hpp "brains/SimInfo.hpp"
78     * @brief As one of the heavy weight class of OOPSE, SimInfo
79     * One of the major changes in SimInfo class is the data struct. It only maintains a list of molecules.
80     * And the Molecule class will maintain all of the concrete objects (atoms, bond, bend, torsions, rigid bodies,
81     * cutoff groups, constrains).
82     * Another major change is the index. No matter single version or parallel version, atoms and
83     * rigid bodies have both global index and local index. Local index is not important to molecule as well as
84     * cutoff group.
85     */
86     class SimInfo {
87     public:
88     typedef std::map<int, Molecule*>::iterator MoleculeIterator;
89 gezelter 2
90 gezelter 507 /**
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     SimInfo(MakeStamps* stamps, std::vector<std::pair<MoleculeStamp*, int> >& molStampPairs, ForceField* ff, Globals* simParams);
99     virtual ~SimInfo();
100 gezelter 2
101 gezelter 507 /**
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 2
108 gezelter 507 /**
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 2
114 gezelter 507 /** Returns the total number of molecules in the system. */
115     int getNGlobalMolecules() {
116     return nGlobalMols_;
117     }
118 gezelter 2
119 gezelter 507 /** Returns the total number of atoms in the system. */
120     int getNGlobalAtoms() {
121     return nGlobalAtoms_;
122     }
123 gezelter 2
124 gezelter 507 /** Returns the total number of cutoff groups in the system. */
125     int getNGlobalCutoffGroups() {
126     return nGlobalCutoffGroups_;
127     }
128 gezelter 2
129 gezelter 507 /**
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 2
137 gezelter 507 /**
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 2
145 gezelter 507 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 2
154 gezelter 507 /** Returns the number of local atoms */
155     unsigned int getNAtoms() {
156     return nAtoms_;
157     }
158 gezelter 2
159 gezelter 507 /** Returns the number of local bonds */
160     unsigned int getNBonds(){
161     return nBonds_;
162     }
163 gezelter 2
164 gezelter 507 /** Returns the number of local bends */
165     unsigned int getNBends() {
166     return nBends_;
167     }
168 gezelter 2
169 gezelter 507 /** Returns the number of local torsions */
170     unsigned int getNTorsions() {
171     return nTorsions_;
172     }
173 gezelter 2
174 gezelter 507 /** Returns the number of local rigid bodies */
175     unsigned int getNRigidBodies() {
176     return nRigidBodies_;
177     }
178 gezelter 2
179 gezelter 507 /** Returns the number of local integrable objects */
180     unsigned int getNIntegrableObjects() {
181     return nIntegrableObjects_;
182     }
183 gezelter 2
184 gezelter 507 /** Returns the number of local cutoff groups */
185     unsigned int getNCutoffGroups() {
186     return nCutoffGroups_;
187     }
188 gezelter 2
189 gezelter 507 /** Returns the total number of constraints in this SimInfo */
190     unsigned int getNConstraints() {
191     return nConstraints_;
192     }
193 gezelter 246
194 gezelter 507 /**
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 2
201 gezelter 507 /**
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 2
208 gezelter 507 /** Returns the number of degrees of freedom */
209     int getNdf() {
210     return ndf_;
211     }
212 gezelter 2
213 gezelter 507 /** Returns the number of raw degrees of freedom */
214     int getNdfRaw() {
215     return ndfRaw_;
216     }
217 gezelter 2
218 gezelter 507 /** Returns the number of translational degrees of freedom */
219     int getNdfTrans() {
220     return ndfTrans_;
221     }
222 gezelter 2
223 gezelter 507 //getNZconstraint and setNZconstraint ruin the coherent of SimInfo class, need refactorying
224 gezelter 246
225 gezelter 507 /** Returns the total number of z-constraint molecules in the system */
226     int getNZconstraint() {
227     return nZconstraint_;
228     }
229 gezelter 2
230 gezelter 507 /**
231     * Sets the number of z-constraint molecules in the system.
232     */
233     void setNZconstraint(int nZconstraint) {
234     nZconstraint_ = nZconstraint;
235     }
236 gezelter 246
237 gezelter 507 /** Returns the snapshot manager. */
238     SnapshotManager* getSnapshotManager() {
239     return sman_;
240     }
241 gezelter 2
242 gezelter 507 /** Sets the snapshot manager. */
243     void setSnapshotManager(SnapshotManager* sman);
244 gezelter 246
245 gezelter 507 /** Returns the force field */
246     ForceField* getForceField() {
247     return forceField_;
248     }
249 gezelter 2
250 gezelter 507 Globals* getSimParams() {
251     return simParams_;
252     }
253 gezelter 2
254 gezelter 507 /** Returns the velocity of center of mass of the whole system.*/
255     Vector3d getComVel();
256 gezelter 2
257 gezelter 507 /** Returns the center of the mass of the whole system.*/
258     Vector3d getCom();
259 chuckv 555 /** Returns the center of the mass and Center of Mass velocity of the whole system.*/
260     void getComAll(Vector3d& com,Vector3d& comVel);
261 gezelter 2
262 chuckv 555 /** Returns intertia tensor for the entire system and system Angular Momentum.*/
263     void getInertiaTensor(Mat3x3d &intertiaTensor,Vector3d &angularMomentum);
264    
265     /** Returns system angular momentum */
266     Vector3d getAngularMomentum();
267    
268 gezelter 507 /** main driver function to interact with fortran during the initialization and molecule migration */
269     void update();
270 gezelter 2
271 gezelter 507 /** Returns the local index manager */
272     LocalIndexManager* getLocalIndexManager() {
273     return &localIndexMan_;
274     }
275 gezelter 2
276 gezelter 507 int getMoleculeStampId(int globalIndex) {
277     //assert(globalIndex < molStampIds_.size())
278     return molStampIds_[globalIndex];
279     }
280 gezelter 2
281 gezelter 507 /** Returns the molecule stamp */
282     MoleculeStamp* getMoleculeStamp(int id) {
283     return moleculeStamps_[id];
284     }
285 gezelter 2
286 gezelter 507 /** Return the total number of the molecule stamps */
287     int getNMoleculeStamp() {
288     return moleculeStamps_.size();
289     }
290     /**
291     * Finds a molecule with a specified global index
292     * @return a pointer point to found molecule
293     * @param index
294     */
295     Molecule* getMoleculeByGlobalIndex(int index) {
296     MoleculeIterator i;
297     i = molecules_.find(index);
298 gezelter 2
299 gezelter 507 return i != molecules_.end() ? i->second : NULL;
300     }
301 gezelter 2
302 gezelter 507 /** Calculate the maximum cutoff radius based on the atom types */
303     double calcMaxCutoffRadius();
304 gezelter 2
305 gezelter 507 double getRcut() {
306     return rcut_;
307     }
308 gezelter 2
309 gezelter 507 double getRsw() {
310     return rsw_;
311     }
312 gezelter 246
313 gezelter 507 std::string getFinalConfigFileName() {
314     return finalConfigFileName_;
315     }
316 gezelter 246
317 gezelter 507 void setFinalConfigFileName(const std::string& fileName) {
318     finalConfigFileName_ = fileName;
319     }
320 gezelter 2
321 gezelter 507 std::string getDumpFileName() {
322     return dumpFileName_;
323     }
324 gezelter 246
325 gezelter 507 void setDumpFileName(const std::string& fileName) {
326     dumpFileName_ = fileName;
327     }
328 gezelter 2
329 gezelter 507 std::string getStatFileName() {
330     return statFileName_;
331     }
332 gezelter 246
333 gezelter 507 void setStatFileName(const std::string& fileName) {
334     statFileName_ = fileName;
335     }
336 chrisfen 417
337 gezelter 507 std::string getRestFileName() {
338     return restFileName_;
339     }
340 chrisfen 417
341 gezelter 507 void setRestFileName(const std::string& fileName) {
342     restFileName_ = fileName;
343     }
344 gezelter 2
345 gezelter 507 /**
346     * Sets GlobalGroupMembership
347     * @see #SimCreator::setGlobalIndex
348     */
349     void setGlobalGroupMembership(const std::vector<int>& globalGroupMembership) {
350     assert(globalGroupMembership.size() == nGlobalAtoms_);
351     globalGroupMembership_ = globalGroupMembership;
352     }
353 gezelter 2
354 gezelter 507 /**
355     * Sets GlobalMolMembership
356     * @see #SimCreator::setGlobalIndex
357     */
358     void setGlobalMolMembership(const std::vector<int>& globalMolMembership) {
359     assert(globalMolMembership.size() == nGlobalAtoms_);
360     globalMolMembership_ = globalMolMembership;
361     }
362 gezelter 246
363    
364 gezelter 507 bool isFortranInitialized() {
365     return fortranInitialized_;
366     }
367 gezelter 246
368 gezelter 507 //below functions are just forward functions
369     //To compose or to inherit is always a hot debate. In general, is-a relation need subclassing, in the
370     //the other hand, has-a relation need composing.
371     /**
372     * Adds property into property map
373     * @param genData GenericData to be added into PropertyMap
374     */
375     void addProperty(GenericData* genData);
376 gezelter 246
377 gezelter 507 /**
378     * Removes property from PropertyMap by name
379     * @param propName the name of property to be removed
380     */
381     void removeProperty(const std::string& propName);
382 gezelter 246
383 gezelter 507 /**
384     * clear all of the properties
385     */
386     void clearProperties();
387 gezelter 246
388 gezelter 507 /**
389     * Returns all names of properties
390     * @return all names of properties
391     */
392     std::vector<std::string> getPropertyNames();
393 gezelter 246
394 gezelter 507 /**
395     * Returns all of the properties in PropertyMap
396     * @return all of the properties in PropertyMap
397     */
398     std::vector<GenericData*> getProperties();
399 gezelter 246
400 gezelter 507 /**
401     * Returns property
402     * @param propName name of property
403     * @return a pointer point to property with propName. If no property named propName
404     * exists, return NULL
405     */
406     GenericData* getPropertyByName(const std::string& propName);
407 gezelter 246
408 gezelter 507 /**
409     * add all exclude pairs of a molecule into exclude list.
410     */
411     void addExcludePairs(Molecule* mol);
412 gezelter 246
413 gezelter 507 /**
414     * remove all exclude pairs which belong to a molecule from exclude list
415     */
416 gezelter 246
417 gezelter 507 void removeExcludePairs(Molecule* mol);
418 tim 292
419    
420 gezelter 507 /** Returns the unique atom types of local processor in an array */
421     std::set<AtomType*> getUniqueAtomTypes();
422 tim 292
423 gezelter 507 friend std::ostream& operator <<(std::ostream& o, SimInfo& info);
424 tim 326
425 gezelter 507 void getCutoff(double& rcut, double& rsw);
426 gezelter 246
427 gezelter 507 private:
428 gezelter 246
429 gezelter 507 /** fill up the simtype struct*/
430     void setupSimType();
431 gezelter 246
432 gezelter 507 /**
433     * Setup Fortran Simulation
434     * @see #setupFortranParallel
435     */
436     void setupFortranSim();
437 gezelter 246
438 gezelter 507 /** Figure out the radius of cutoff, radius of switching function and pass them to fortran */
439     void setupCutoff();
440 gezelter 246
441 gezelter 507 /** Calculates the number of degress of freedom in the whole system */
442     void calcNdf();
443     void calcNdfRaw();
444     void calcNdfTrans();
445 gezelter 246
446 gezelter 507 /**
447     * Adds molecule stamp and the total number of the molecule with same molecule stamp in the whole
448     * system.
449     */
450     void addMoleculeStamp(MoleculeStamp* molStamp, int nmol);
451 gezelter 246
452 gezelter 507 MakeStamps* stamps_;
453     ForceField* forceField_;
454     Globals* simParams_;
455 gezelter 246
456 gezelter 507 std::map<int, Molecule*> molecules_; /**< Molecule array */
457 gezelter 246
458 gezelter 507 //degress of freedom
459     int ndf_; /**< number of degress of freedom (excludes constraints), ndf_ is local */
460     int ndfRaw_; /**< number of degress of freedom (includes constraints), ndfRaw_ is local */
461     int ndfTrans_; /**< number of translation degress of freedom, ndfTrans_ is local */
462     int nZconstraint_; /** number of z-constraint molecules, nZconstraint_ is global */
463 gezelter 246
464 gezelter 507 //number of global objects
465     int nGlobalMols_; /**< number of molecules in the system */
466     int nGlobalAtoms_; /**< number of atoms in the system */
467     int nGlobalCutoffGroups_; /**< number of cutoff groups in this system */
468     int nGlobalIntegrableObjects_; /**< number of integrable objects in this system */
469     int nGlobalRigidBodies_; /**< number of rigid bodies in this system */
470     /**
471     * the size of globalGroupMembership_ is nGlobalAtoms. Its index is global index of an atom, and the
472     * corresponding content is the global index of cutoff group this atom belong to.
473     * It is filled by SimCreator once and only once, since it never changed during the simulation.
474     */
475     std::vector<int> globalGroupMembership_;
476 gezelter 246
477 gezelter 507 /**
478     * the size of globalGroupMembership_ is nGlobalAtoms. Its index is global index of an atom, and the
479     * corresponding content is the global index of molecule this atom belong to.
480     * It is filled by SimCreator once and only once, since it is never changed during the simulation.
481     */
482     std::vector<int> globalMolMembership_;
483 gezelter 246
484    
485 gezelter 507 std::vector<int> molStampIds_; /**< stamp id array of all molecules in the system */
486     std::vector<MoleculeStamp*> moleculeStamps_; /**< molecule stamps array */
487 gezelter 246
488 gezelter 507 //number of local objects
489     int nAtoms_; /**< number of atoms in local processor */
490     int nBonds_; /**< number of bonds in local processor */
491     int nBends_; /**< number of bends in local processor */
492     int nTorsions_; /**< number of torsions in local processor */
493     int nRigidBodies_; /**< number of rigid bodies in local processor */
494     int nIntegrableObjects_; /**< number of integrable objects in local processor */
495     int nCutoffGroups_; /**< number of cutoff groups in local processor */
496     int nConstraints_; /**< number of constraints in local processors */
497 gezelter 246
498 gezelter 507 simtype fInfo_; /**< A dual struct shared by c++/fortran which indicates the atom types in simulation*/
499     Exclude exclude_;
500     PropertyMap properties_; /**< Generic Property */
501     SnapshotManager* sman_; /**< SnapshotManager */
502 gezelter 246
503 gezelter 507 /**
504     * The reason to have a local index manager is that when molecule is migrating to other processors,
505     * the atoms and the rigid-bodies will release their local indices to LocalIndexManager. Combining the
506     * information of molecule migrating to current processor, Migrator class can query the LocalIndexManager
507     * to make a efficient data moving plan.
508     */
509     LocalIndexManager localIndexMan_;
510 gezelter 246
511 gezelter 507 //file names
512     std::string finalConfigFileName_;
513     std::string dumpFileName_;
514     std::string statFileName_;
515     std::string restFileName_;
516 chrisfen 417
517 gezelter 507 double rcut_; /**< cutoff radius*/
518     double rsw_; /**< radius of switching function*/
519 gezelter 246
520 gezelter 507 bool fortranInitialized_; /**< flag indicate whether fortran side is initialized */
521 tim 292
522 gezelter 246 #ifdef IS_MPI
523     //in Parallel version, we need MolToProc
524 gezelter 507 public:
525 gezelter 246
526 gezelter 507 /**
527     * Finds the processor where a molecule resides
528     * @return the id of the processor which contains the molecule
529     * @param globalIndex global Index of the molecule
530     */
531     int getMolToProc(int globalIndex) {
532     //assert(globalIndex < molToProcMap_.size());
533     return molToProcMap_[globalIndex];
534     }
535 gezelter 246
536 gezelter 507 /**
537     * Set MolToProcMap array
538     * @see #SimCreator::divideMolecules
539     */
540     void setMolToProcMap(const std::vector<int>& molToProcMap) {
541     molToProcMap_ = molToProcMap;
542     }
543 gezelter 246
544 gezelter 507 private:
545 gezelter 246
546 gezelter 507 void setupFortranParallel();
547 gezelter 246
548 gezelter 507 /**
549     * The size of molToProcMap_ is equal to total number of molecules in the system.
550     * It maps a molecule to the processor on which it resides. it is filled by SimCreator once and only
551     * once.
552     */
553     std::vector<int> molToProcMap_;
554 tim 292
555 gezelter 246 #endif
556    
557 gezelter 507 };
558 gezelter 2
559 gezelter 246 } //namespace oopse
560     #endif //BRAINS_SIMMODEL_HPP
561 gezelter 2