1 |
/* |
2 |
* Copyright (C) 2000-2004 Object Oriented Parallel Simulation Engine (OOPSE) project |
3 |
* |
4 |
* Contact: oopse@oopse.org |
5 |
* |
6 |
* This program is free software; you can redistribute it and/or |
7 |
* modify it under the terms of the GNU Lesser General Public License |
8 |
* as published by the Free Software Foundation; either version 2.1 |
9 |
* of the License, or (at your option) any later version. |
10 |
* All we ask is that proper credit is given for our work, which includes |
11 |
* - but is not limited to - adding the above copyright notice to the beginning |
12 |
* of your source code files, and to any copyright notice that you may distribute |
13 |
* with programs based on this work. |
14 |
* |
15 |
* This program is distributed in the hope that it will be useful, |
16 |
* but WITHOUT ANY WARRANTY; without even the implied warranty of |
17 |
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
18 |
* GNU Lesser General Public License for more details. |
19 |
* |
20 |
* You should have received a copy of the GNU Lesser General Public License |
21 |
* along with this program; if not, write to the Free Software |
22 |
* Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. |
23 |
* |
24 |
*/ |
25 |
|
26 |
/** |
27 |
* @file SimInfo.hpp |
28 |
* @author tlin |
29 |
* @date 11/02/2004 |
30 |
* @version 1.0 |
31 |
*/ |
32 |
|
33 |
#ifndef BRAINS_SIMMODEL_HPP |
34 |
#define BRAINS_SIMMODEL_HPP |
35 |
|
36 |
#include <iostream> |
37 |
#include <vector> |
38 |
#include <utility> |
39 |
|
40 |
#include "brains/fSimulation.h" |
41 |
#include "primitives/Molecule.hpp" |
42 |
#include "types/MoleculeStamp.hpp" |
43 |
#include "utils/PropertyMap.hpp" |
44 |
#include "io/Globals.hpp" |
45 |
|
46 |
namespace oopse{ |
47 |
|
48 |
/** |
49 |
* @class SimInfo SimInfo.hpp "brains/SimInfo.hpp" |
50 |
* @brief As one of the heavy weight class of OOPSE, SimInfo |
51 |
* One of the major changes in SimInfo class is the data struct. It only maintains a list of molecules. |
52 |
* And the Molecule class will maintain all of the concrete objects (atoms, bond, bend, torsions, rigid bodies, |
53 |
* cutoff groups, constrains). |
54 |
* Another major change is the index. No matter single version or parallel version, atoms and |
55 |
* rigid bodies have both global index and local index. Local index is not important to molecule as well as |
56 |
* cutoff group. |
57 |
*/ |
58 |
class SimInfo { |
59 |
public: |
60 |
typedef MoleculeIterator MoleculeIterator; |
61 |
|
62 |
/** |
63 |
* Constructor of SimInfo |
64 |
* @param molStampPairs MoleculeStamp Array. The first element of the pair is molecule stamp, the |
65 |
* second element is the total number of molecules with the same molecule stamp in the system |
66 |
* @param ff pointer of a concrete ForceField instance |
67 |
* @param globals |
68 |
* @note |
69 |
* <p> |
70 |
* The major change of SimInfo |
71 |
* </p> |
72 |
*/ |
73 |
SimInfo(const std::vector<std::pair<MoleculeStamp*, int> >& molStampPairs, ForceField* ff, Globals* globals); |
74 |
virtual ~SimInfo(); |
75 |
|
76 |
/** |
77 |
* Adds a molecule |
78 |
* @return return true if adding successfully, return false if the molecule is already in SimInfo |
79 |
* @param mol molecule to be added |
80 |
*/ |
81 |
bool addMolecule(Molecule* mol); |
82 |
|
83 |
/** |
84 |
* Removes a molecule from SimInfo |
85 |
* @return true if removing successfully, return false if molecule is not in this SimInfo |
86 |
*/ |
87 |
bool removeMolecule(Molecule* mol); |
88 |
|
89 |
/** Returns the total number of molecules in the system. */ |
90 |
int getNGlobalMolecules() { |
91 |
return nGlobalMols_; |
92 |
} |
93 |
|
94 |
/** Returns the total number of atoms in the system. */ |
95 |
int getNGlobalAtoms() { |
96 |
return nGlobalAtoms_; |
97 |
} |
98 |
|
99 |
/** Returns the total number of cutoff groups in the system. */ |
100 |
int getNGlobalCutoffGroups() { |
101 |
return nGlobalCutoffGroups_; |
102 |
} |
103 |
|
104 |
/** |
105 |
* Returns the number of local molecules. |
106 |
* @return the number of local molecules |
107 |
*/ |
108 |
int getNMolecules() { |
109 |
return molecules_.size(); |
110 |
} |
111 |
|
112 |
/** Returns the number of local atoms */ |
113 |
unsigned int getNAtoms() { |
114 |
return nAtoms_; |
115 |
} |
116 |
|
117 |
/** Returns the number of local bonds */ |
118 |
unsigned int getNBonds(){ |
119 |
return nBonds_; |
120 |
} |
121 |
|
122 |
/** Returns the number of local bends */ |
123 |
unsigned int getNBends() { |
124 |
return nBends_; |
125 |
} |
126 |
|
127 |
/** Returns the number of local torsions */ |
128 |
unsigned int getNTorsions() { |
129 |
return nTorsions_; |
130 |
} |
131 |
|
132 |
/** Returns the number of local rigid bodies */ |
133 |
unsigned int getNRigidBodies() { |
134 |
return nRigidBodies_; |
135 |
} |
136 |
|
137 |
/** Returns the number of local integrable objects */ |
138 |
unsigned int getNIntegrableObjects() { |
139 |
return nIntegrableObjects_; |
140 |
} |
141 |
|
142 |
/** Returns the number of local cutoff groups */ |
143 |
unsigned int getNCutoffGroups() { |
144 |
return nCutoffGroups_; |
145 |
} |
146 |
|
147 |
/** Returns the total number of constraints in this SimInfo */ |
148 |
unsigned int getNConstraints() { |
149 |
return nConstraints_; |
150 |
} |
151 |
|
152 |
/** |
153 |
* Returns the first molecule in this SimInfo and intialize the iterator. |
154 |
* @return the first molecule, return NULL if there is not molecule in this SimInfo |
155 |
* @param i the iterator of molecule array (user shouldn't change it) |
156 |
*/ |
157 |
Molecule* beginMolecule(MoleculeIterator& i); |
158 |
|
159 |
/** |
160 |
* Returns the next avaliable Molecule based on the iterator. |
161 |
* @return the next avaliable molecule, return NULL if reaching the end of the array |
162 |
* @param i the iterator of molecule array |
163 |
*/ |
164 |
Molecule* nextMolecule(MoleculeIterator& i); |
165 |
|
166 |
/** Returns the number of degrees of freedom */ |
167 |
int getNdf() { |
168 |
return ndf_; |
169 |
} |
170 |
|
171 |
/** Returns the number of raw degrees of freedom */ |
172 |
int getNdfRaw() { |
173 |
return ndfRaw_; |
174 |
} |
175 |
|
176 |
/** Returns the number of translational degrees of freedom */ |
177 |
int getNdfTrans() { |
178 |
return ndfTrans_; |
179 |
} |
180 |
|
181 |
//getNZconstraint and setNZconstraint ruin the coherent of SimInfo class, need refactorying |
182 |
|
183 |
/** Returns the total number of z-constraint molecules in the system */ |
184 |
int getNZconstraint() { |
185 |
return nZconstraint_; |
186 |
} |
187 |
|
188 |
/** |
189 |
* Sets the number of z-constraint molecules in the system. |
190 |
*/ |
191 |
int setNZconstraint(int nZconstraint) { |
192 |
nZconstraint_ = nZconstraint; |
193 |
} |
194 |
|
195 |
/** Returns the snapshot manager. */ |
196 |
SnapshotManager* getSnapshotManager() { |
197 |
return sman_; |
198 |
} |
199 |
|
200 |
/** Sets the snapshot manager. */ |
201 |
void setSnapshotManager(SnapshotManager* sman) { |
202 |
sman_ = sman; |
203 |
} |
204 |
|
205 |
/** Returns the force field */ |
206 |
ForceField* getForceField() { |
207 |
return forceField_; |
208 |
} |
209 |
|
210 |
Globals* getGlobals() { |
211 |
return globals_; |
212 |
} |
213 |
|
214 |
/** Returns the velocity of center of mass of the whole system.*/ |
215 |
Vector3d getComVel(); |
216 |
|
217 |
/** Returns the center of the mass of the whole system.*/ |
218 |
Vector3d getCom(); |
219 |
|
220 |
/** Returns the seed (used for random number generator) */ |
221 |
int getSeed() { |
222 |
return seed_; |
223 |
} |
224 |
|
225 |
/** Sets the seed*/ |
226 |
void setSeed(int seed) { |
227 |
seed_ = seed; |
228 |
} |
229 |
|
230 |
/** main driver function to interact with fortran during the initialization and molecule migration */ |
231 |
void update(); |
232 |
|
233 |
/** Returns the local index manager */ |
234 |
LocalIndexManager* getLocalIndexManager() { |
235 |
return &localIndexMan_; |
236 |
} |
237 |
|
238 |
int getMoleculeStampId(int globalIndex) { |
239 |
//assert(globalIndex < molStampIds_.size()) |
240 |
return molStampIds_[globalIndex]; |
241 |
} |
242 |
|
243 |
/** Returns the molecule stamp */ |
244 |
MoleculeStamp* getMoleculeStamp(int id) { |
245 |
return moleculeStamps_[id]; |
246 |
} |
247 |
|
248 |
/** |
249 |
* Finds a molecule with a specified global index |
250 |
* @return a pointer point to found molecule |
251 |
* @param index |
252 |
*/ |
253 |
Molecule* getMoleculeByGlobalIndex(int index) { |
254 |
std::map<int, Molecule*> i; |
255 |
i = molecules_.find(index); |
256 |
|
257 |
return i != molecules_.end() ? i->second : NULL; |
258 |
} |
259 |
|
260 |
/** Calculate the maximum cutoff radius based on the atom types */ |
261 |
double calcMaxCutoffRadius(); |
262 |
|
263 |
double getRcut() { |
264 |
return rcut_; |
265 |
} |
266 |
|
267 |
void setRcut(double rcut) { |
268 |
rcut_ = rcut; |
269 |
} |
270 |
|
271 |
double getRsw() { |
272 |
return rsw_; |
273 |
} |
274 |
void setRsw(double rsw) { |
275 |
rsw_ = rsw; |
276 |
} |
277 |
|
278 |
std::string getFinalConfigFileName() { |
279 |
return finalConfigFileName_; |
280 |
} |
281 |
|
282 |
void setFinalConfigFileName(const std::string& fileName) { |
283 |
finalConfigFileName_ = fileName; |
284 |
} |
285 |
|
286 |
|
287 |
std::string getDumpFileName() { |
288 |
return dumpFileName_; |
289 |
} |
290 |
|
291 |
void setDumpFileName(const std::string& fileName) { |
292 |
dumpFileName_ = fileName; |
293 |
} |
294 |
|
295 |
std::string getStatFileName() { |
296 |
return statFileName_; |
297 |
} |
298 |
|
299 |
void setStatFileName(const std::string& fileName) { |
300 |
statFileName_ = fileName; |
301 |
} |
302 |
|
303 |
/** |
304 |
* Returns the pointer of internal globalGroupMembership_ array. This array will be filled by SimCreator class |
305 |
* @see #SimCreator::setGlobalIndex |
306 |
*/ |
307 |
int* getGlobalGroupMembershipPointer() { |
308 |
return globalGroupMembership_[0]; |
309 |
} |
310 |
|
311 |
/** |
312 |
* Returns the pointer of internal globalMolMembership_ array. This array will be filled by SimCreator class |
313 |
* @see #SimCreator::setGlobalIndex |
314 |
*/ |
315 |
int* getGlobalMolMembershipPointer() { |
316 |
return globalMolMembership_[0]; |
317 |
} |
318 |
|
319 |
//below functions are just forward functions |
320 |
//To compose or to inherit is always a hot debate. In general, is-a relation need subclassing, in the |
321 |
//the other hand, has-a relation need composing. |
322 |
/** |
323 |
* Adds property into property map |
324 |
* @param genData GenericData to be added into PropertyMap |
325 |
*/ |
326 |
void addProperty(GenericData* genData); |
327 |
|
328 |
/** |
329 |
* Removes property from PropertyMap by name |
330 |
* @param propName the name of property to be removed |
331 |
*/ |
332 |
void removeProperty(const std::string& propName); |
333 |
|
334 |
/** |
335 |
* clear all of the properties |
336 |
*/ |
337 |
void clearProperties(); |
338 |
|
339 |
/** |
340 |
* Returns all names of properties |
341 |
* @return all names of properties |
342 |
*/ |
343 |
std::vector<std::string> getPropertyNames(); |
344 |
|
345 |
/** |
346 |
* Returns all of the properties in PropertyMap |
347 |
* @return all of the properties in PropertyMap |
348 |
*/ |
349 |
std::vector<GenericData*> getProperties(); |
350 |
|
351 |
/** |
352 |
* Returns property |
353 |
* @param propName name of property |
354 |
* @return a pointer point to property with propName. If no property named propName |
355 |
* exists, return NULL |
356 |
*/ |
357 |
GenericData* getPropertyByName(const std::string& propName); |
358 |
|
359 |
friend std::ostream& operator <<(ostream& o, SimInfo& info); |
360 |
|
361 |
private: |
362 |
|
363 |
|
364 |
/** Returns the unique atom types of local processor in an array */ |
365 |
std::set<AtomType*> SimInfo::getUniqueAtomTypes(); |
366 |
|
367 |
/** fill up the simtype struct*/ |
368 |
void setupSimType(); |
369 |
|
370 |
/** |
371 |
* Setup Fortran Simulation |
372 |
* @see #setupFortranParallel |
373 |
*/ |
374 |
void setupFortranSim(); |
375 |
|
376 |
/** Calculates the number of degress of freedom in the whole system */ |
377 |
void calcNdf(); |
378 |
void calcNdfRaw(); |
379 |
void calcNdfTrans(); |
380 |
|
381 |
void addExcludePairs(Molecule* mol); |
382 |
void removeExcludePairs(Molecule* mol); |
383 |
|
384 |
/** |
385 |
* Adds molecule stamp and the total number of the molecule with same molecule stamp in the whole |
386 |
* system. |
387 |
*/ |
388 |
void addMoleculeStamp(MoleculeStamp* molStamp, int nmol); |
389 |
|
390 |
std::map<int, Molecule*> molecules_; /**< Molecule array */ |
391 |
|
392 |
//degress of freedom |
393 |
int ndf_; /**< number of degress of freedom (excludes constraints), ndf_ is local */ |
394 |
int ndfRaw_; /**< number of degress of freedom (includes constraints), ndfRaw_ is local */ |
395 |
int ndfTrans_; /**< number of translation degress of freedom, ndfTrans_ is local */ |
396 |
int nZconstraint_; /** number of z-constraint molecules, nZconstraint_ is global */ |
397 |
|
398 |
//number of global objects |
399 |
int nGlobalMols_; /**< number of molecules in the system */ |
400 |
int nGlobalAtoms_; /**< number of atoms in the system */ |
401 |
int nGlobalCutoffGroups_; /**< number of cutoff groups in this system */ |
402 |
|
403 |
/** |
404 |
* the size of globalGroupMembership_ is nGlobalAtoms. Its index is global index of an atom, and the |
405 |
* corresponding content is the global index of cutoff group this atom belong to. |
406 |
* It is filled by SimCreator once and only once, since it is never changed during the simulation. |
407 |
*/ |
408 |
std::vector<int> globalGroupMembership_; |
409 |
|
410 |
/** |
411 |
* the size of globalGroupMembership_ is nGlobalAtoms. Its index is global index of an atom, and the |
412 |
* corresponding content is the global index of molecule this atom belong to. |
413 |
* It is filled by SimCreator once and only once, since it is never changed during the simulation. |
414 |
*/ |
415 |
std::vector<int> globalMolMembership_; |
416 |
|
417 |
|
418 |
std::vector<int> molStampIds_; /**< stamp id array of all molecules in the system */ |
419 |
std::vector<MoleculeStamp*> moleculeStamps_; /**< molecule stamps array */ |
420 |
|
421 |
//number of local objects |
422 |
int nAtoms_; /**< number of atoms in local processor */ |
423 |
int nBonds_; /**< number of bonds in local processor */ |
424 |
int nBends_; /**< number of bends in local processor */ |
425 |
int nTorsions_; /**< number of torsions in local processor */ |
426 |
int nRigidBodies_; /**< number of rigid bodies in local processor */ |
427 |
int nIntegrableObjects_; /**< number of integrable objects in local processor */ |
428 |
int nCutoffGroups_; /**< number of cutoff groups in local processor */ |
429 |
int nConstraints_; /**< number of constraints in local processors */ |
430 |
|
431 |
simtype fInfo_; /**< A dual struct shared by c++/fortran which indicates the atom types in simulation*/ |
432 |
Exclude exclude_; |
433 |
ForceField* forceField_; |
434 |
PropertyMap properties_; /**< Generic Property */ |
435 |
SnapshotManager* sman_; /**< SnapshotManager */ |
436 |
Globals* globals_; |
437 |
int seed_; /**< seed for random number generator */ |
438 |
|
439 |
/** |
440 |
* The reason to have a local index manager is that when molecule is migrating to other processors, |
441 |
* the atoms and the rigid-bodies will release their local indices to LocalIndexManager. Combining the |
442 |
* information of molecule migrating to current processor, Migrator class can query the LocalIndexManager |
443 |
* to make a efficient data moving plan. |
444 |
*/ |
445 |
LocalIndexManager localIndexMan_; |
446 |
|
447 |
//file names |
448 |
std::string finalConfigFileName_; |
449 |
std::string dumpFileName_; |
450 |
std::string statFileName_; |
451 |
|
452 |
double rcut_; /**< cutoff radius*/ |
453 |
double rsw_; /**< radius of switching function*/ |
454 |
|
455 |
#ifdef IS_MPI |
456 |
//in Parallel version, we need MolToProc |
457 |
public: |
458 |
|
459 |
/** |
460 |
* Finds the processor where a molecule resides |
461 |
* @return the id of the processor which contains the molecule |
462 |
* @param globalIndex global Index of the molecule |
463 |
*/ |
464 |
int getMolToProc(int globalIndex) { |
465 |
//assert(globalIndex < molToProcMap_.size()); |
466 |
return molToProcMap_[globalIndex]; |
467 |
} |
468 |
|
469 |
/** |
470 |
* Returns the pointer of internal molToProcMap array. This array will be filled by SimCreator class |
471 |
* @see #SimCreator::divideMolecules |
472 |
*/ |
473 |
int* getMolToProcMapPointer() { |
474 |
return &molToProcMap_[0]; |
475 |
} |
476 |
|
477 |
private: |
478 |
|
479 |
void setupFortranParallel(); |
480 |
|
481 |
/** |
482 |
* The size of molToProcMap_ is equal to total number of molecules in the system. |
483 |
* It maps a molecule to the processor on which it resides. it is filled by SimCreator once and only |
484 |
* once. |
485 |
*/ |
486 |
std::vector<int> molToProcMap_; |
487 |
#endif |
488 |
|
489 |
}; |
490 |
|
491 |
} //namespace oopse |
492 |
#endif //BRAINS_SIMMODEL_HPP |