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, 2 months ago) by tim
File size: 18776 byte(s)
Log Message:
Adding single precision capabilities to c++ side

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 "math/SquareMatrix3.hpp"
61 #include "types/MoleculeStamp.hpp"
62 #include "UseTheForce/ForceField.hpp"
63 #include "utils/PropertyMap.hpp"
64 #include "utils/LocalIndexManager.hpp"
65
66 //another nonsense macro declaration
67 #define __C
68 #include "brains/fSimulation.h"
69
70 namespace oopse{
71
72 //forward decalration
73 class SnapshotManager;
74 class Molecule;
75 class SelectionManager;
76 /**
77 * @class SimInfo SimInfo.hpp "brains/SimInfo.hpp"
78 * @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 */
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(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_ - getFdf();
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 /** 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 //getNZconstraint and setNZconstraint ruin the coherent of SimInfo class, need refactorying
230
231 /** Returns the total number of z-constraint molecules in the system */
232 int getNZconstraint() {
233 return nZconstraint_;
234 }
235
236 /**
237 * Sets the number of z-constraint molecules in the system.
238 */
239 void setNZconstraint(int nZconstraint) {
240 nZconstraint_ = nZconstraint;
241 }
242
243 /** Returns the snapshot manager. */
244 SnapshotManager* getSnapshotManager() {
245 return sman_;
246 }
247
248 /** Sets the snapshot manager. */
249 void setSnapshotManager(SnapshotManager* sman);
250
251 /** Returns the force field */
252 ForceField* getForceField() {
253 return forceField_;
254 }
255
256 Globals* getSimParams() {
257 return simParams_;
258 }
259
260 /** Returns the velocity of center of mass of the whole system.*/
261 Vector3d getComVel();
262
263 /** Returns the center of the mass of the whole system.*/
264 Vector3d getCom();
265 /** Returns the center of the mass and Center of Mass velocity of the whole system.*/
266 void getComAll(Vector3d& com,Vector3d& comVel);
267
268 /** 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 /** main driver function to interact with fortran during the initialization and molecule migration */
275 void update();
276
277 /** Returns the local index manager */
278 LocalIndexManager* getLocalIndexManager() {
279 return &localIndexMan_;
280 }
281
282 int getMoleculeStampId(int globalIndex) {
283 //assert(globalIndex < molStampIds_.size())
284 return molStampIds_[globalIndex];
285 }
286
287 /** Returns the molecule stamp */
288 MoleculeStamp* getMoleculeStamp(int id) {
289 return moleculeStamps_[id];
290 }
291
292 /** 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
305 return i != molecules_.end() ? i->second : NULL;
306 }
307
308 RealType getRcut() {
309 return rcut_;
310 }
311
312 RealType getRsw() {
313 return rsw_;
314 }
315
316 RealType getList() {
317 return rlist_;
318 }
319
320 std::string getFinalConfigFileName() {
321 return finalConfigFileName_;
322 }
323
324 void setFinalConfigFileName(const std::string& fileName) {
325 finalConfigFileName_ = fileName;
326 }
327
328 std::string getDumpFileName() {
329 return dumpFileName_;
330 }
331
332 void setDumpFileName(const std::string& fileName) {
333 dumpFileName_ = fileName;
334 }
335
336 std::string getStatFileName() {
337 return statFileName_;
338 }
339
340 void setStatFileName(const std::string& fileName) {
341 statFileName_ = fileName;
342 }
343
344 std::string getRestFileName() {
345 return restFileName_;
346 }
347
348 void setRestFileName(const std::string& fileName) {
349 restFileName_ = fileName;
350 }
351
352 /**
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
361 /**
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
370
371 bool isFortranInitialized() {
372 return fortranInitialized_;
373 }
374
375 //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
384 /**
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
390 /**
391 * clear all of the properties
392 */
393 void clearProperties();
394
395 /**
396 * Returns all names of properties
397 * @return all names of properties
398 */
399 std::vector<std::string> getPropertyNames();
400
401 /**
402 * Returns all of the properties in PropertyMap
403 * @return all of the properties in PropertyMap
404 */
405 std::vector<GenericData*> getProperties();
406
407 /**
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
415 /**
416 * add all exclude pairs of a molecule into exclude list.
417 */
418 void addExcludePairs(Molecule* mol);
419
420 /**
421 * remove all exclude pairs which belong to a molecule from exclude list
422 */
423
424 void removeExcludePairs(Molecule* mol);
425
426
427 /** Returns the unique atom types of local processor in an array */
428 std::set<AtomType*> getUniqueAtomTypes();
429
430 friend std::ostream& operator <<(std::ostream& o, SimInfo& info);
431
432 void getCutoff(RealType& rcut, RealType& rsw);
433
434 private:
435
436 /** fill up the simtype struct*/
437 void setupSimType();
438
439 /**
440 * Setup Fortran Simulation
441 * @see #setupFortranParallel
442 */
443 void setupFortranSim();
444
445 /** Figure out the radius of cutoff, radius of switching function and pass them to fortran */
446 void setupCutoff();
447
448 /** Figure out which coulombic correction method to use and pass to fortran */
449 void setupElectrostaticSummationMethod( int isError );
450
451 /** Figure out which polynomial type to use for the switching function */
452 void setupSwitchingFunction();
453
454 /** Calculates the number of degress of freedom in the whole system */
455 void calcNdf();
456 void calcNdfRaw();
457 void calcNdfTrans();
458
459 ForceField* forceField_;
460 Globals* simParams_;
461
462 std::map<int, Molecule*> molecules_; /**< Molecule array */
463
464 /**
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
470 //degress of freedom
471 int ndf_; /**< number of degress of freedom (excludes constraints), ndf_ is local */
472 int fdf_local; /**< number of frozen degrees of freedom */
473 int fdf_; /**< number of frozen degrees of freedom */
474 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
478 //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
491 /**
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
498
499 std::vector<int> molStampIds_; /**< stamp id array of all molecules in the system */
500 std::vector<MoleculeStamp*> moleculeStamps_; /**< molecule stamps array */
501
502 //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
512 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
517 /**
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
525 //file names
526 std::string finalConfigFileName_;
527 std::string dumpFileName_;
528 std::string statFileName_;
529 std::string restFileName_;
530
531 RealType rcut_; /**< cutoff radius*/
532 RealType rsw_; /**< radius of switching function*/
533 RealType rlist_; /**< neighbor list radius */
534
535 bool fortranInitialized_; /**< flag indicate whether fortran side is initialized */
536
537 #ifdef IS_MPI
538 //in Parallel version, we need MolToProc
539 public:
540
541 /**
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
551 /**
552 * Set MolToProcMap array
553 * @see #SimCreator::divideMolecules
554 */
555 void setMolToProcMap(const std::vector<int>& molToProcMap) {
556 molToProcMap_ = molToProcMap;
557 }
558
559 private:
560
561 void setupFortranParallel();
562
563 /**
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
570 #endif
571
572 };
573
574 } //namespace oopse
575 #endif //BRAINS_SIMMODEL_HPP
576