ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/brains/SimInfo.hpp
Revision: 2204
Committed: Fri Apr 15 22:04:00 2005 UTC (19 years, 2 months ago) by gezelter
File size: 18021 byte(s)
Log Message:
xemacs has been drafted to perform our indentation services

File Contents

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