ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1095
Committed: Tue Dec 5 00:17:24 2006 UTC (18 years, 4 months ago) by chuckv
Original Path: trunk/src/brains/SimInfo.cpp
File size: 46415 byte(s)
Log Message:
Added interface to change number of neighbors in calculating neighbor list.

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.cpp
44 * @author tlin
45 * @date 11/02/2004
46 * @version 1.0
47 */
48
49 #include <algorithm>
50 #include <set>
51 #include <map>
52
53 #include "brains/SimInfo.hpp"
54 #include "math/Vector3.hpp"
55 #include "primitives/Molecule.hpp"
56 #include "primitives/StuntDouble.hpp"
57 #include "UseTheForce/fCutoffPolicy.h"
58 #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
59 #include "UseTheForce/DarkSide/fElectrostaticScreeningMethod.h"
60 #include "UseTheForce/DarkSide/fSwitchingFunctionType.h"
61 #include "UseTheForce/doForces_interface.h"
62 #include "UseTheForce/DarkSide/neighborLists_interface.h"
63 #include "UseTheForce/DarkSide/electrostatic_interface.h"
64 #include "UseTheForce/DarkSide/switcheroo_interface.h"
65 #include "utils/MemoryUtils.hpp"
66 #include "utils/simError.h"
67 #include "selection/SelectionManager.hpp"
68 #include "io/ForceFieldOptions.hpp"
69 #include "UseTheForce/ForceField.hpp"
70
71
72 #ifdef IS_MPI
73 #include "UseTheForce/mpiComponentPlan.h"
74 #include "UseTheForce/DarkSide/simParallel_interface.h"
75 #endif
76
77 namespace oopse {
78 std::set<int> getRigidSet(int index, std::map<int, std::set<int> >& container) {
79 std::map<int, std::set<int> >::iterator i = container.find(index);
80 std::set<int> result;
81 if (i != container.end()) {
82 result = i->second;
83 }
84
85 return result;
86 }
87
88 SimInfo::SimInfo(ForceField* ff, Globals* simParams) :
89 forceField_(ff), simParams_(simParams),
90 ndf_(0), fdf_local(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0),
91 nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0),
92 nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
93 nAtoms_(0), nBonds_(0), nBends_(0), nTorsions_(0), nRigidBodies_(0),
94 nIntegrableObjects_(0), nCutoffGroups_(0), nConstraints_(0),
95 sman_(NULL), fortranInitialized_(false), calcBoxDipole_(false) {
96
97 MoleculeStamp* molStamp;
98 int nMolWithSameStamp;
99 int nCutoffAtoms = 0; // number of atoms belong to cutoff groups
100 int nGroups = 0; //total cutoff groups defined in meta-data file
101 CutoffGroupStamp* cgStamp;
102 RigidBodyStamp* rbStamp;
103 int nRigidAtoms = 0;
104 std::vector<Component*> components = simParams->getComponents();
105
106 for (std::vector<Component*>::iterator i = components.begin(); i !=components.end(); ++i) {
107 molStamp = (*i)->getMoleculeStamp();
108 nMolWithSameStamp = (*i)->getNMol();
109
110 addMoleculeStamp(molStamp, nMolWithSameStamp);
111
112 //calculate atoms in molecules
113 nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp;
114
115 //calculate atoms in cutoff groups
116 int nAtomsInGroups = 0;
117 int nCutoffGroupsInStamp = molStamp->getNCutoffGroups();
118
119 for (int j=0; j < nCutoffGroupsInStamp; j++) {
120 cgStamp = molStamp->getCutoffGroupStamp(j);
121 nAtomsInGroups += cgStamp->getNMembers();
122 }
123
124 nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
125
126 nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;
127
128 //calculate atoms in rigid bodies
129 int nAtomsInRigidBodies = 0;
130 int nRigidBodiesInStamp = molStamp->getNRigidBodies();
131
132 for (int j=0; j < nRigidBodiesInStamp; j++) {
133 rbStamp = molStamp->getRigidBodyStamp(j);
134 nAtomsInRigidBodies += rbStamp->getNMembers();
135 }
136
137 nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
138 nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;
139
140 }
141
142 //every free atom (atom does not belong to cutoff groups) is a cutoff
143 //group therefore the total number of cutoff groups in the system is
144 //equal to the total number of atoms minus number of atoms belong to
145 //cutoff group defined in meta-data file plus the number of cutoff
146 //groups defined in meta-data file
147 nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
148
149 //every free atom (atom does not belong to rigid bodies) is an
150 //integrable object therefore the total number of integrable objects
151 //in the system is equal to the total number of atoms minus number of
152 //atoms belong to rigid body defined in meta-data file plus the number
153 //of rigid bodies defined in meta-data file
154 nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms
155 + nGlobalRigidBodies_;
156
157 nGlobalMols_ = molStampIds_.size();
158
159 #ifdef IS_MPI
160 molToProcMap_.resize(nGlobalMols_);
161 #endif
162
163 }
164
165 SimInfo::~SimInfo() {
166 std::map<int, Molecule*>::iterator i;
167 for (i = molecules_.begin(); i != molecules_.end(); ++i) {
168 delete i->second;
169 }
170 molecules_.clear();
171
172 delete sman_;
173 delete simParams_;
174 delete forceField_;
175 }
176
177 int SimInfo::getNGlobalConstraints() {
178 int nGlobalConstraints;
179 #ifdef IS_MPI
180 MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM,
181 MPI_COMM_WORLD);
182 #else
183 nGlobalConstraints = nConstraints_;
184 #endif
185 return nGlobalConstraints;
186 }
187
188 bool SimInfo::addMolecule(Molecule* mol) {
189 MoleculeIterator i;
190
191 i = molecules_.find(mol->getGlobalIndex());
192 if (i == molecules_.end() ) {
193
194 molecules_.insert(std::make_pair(mol->getGlobalIndex(), mol));
195
196 nAtoms_ += mol->getNAtoms();
197 nBonds_ += mol->getNBonds();
198 nBends_ += mol->getNBends();
199 nTorsions_ += mol->getNTorsions();
200 nRigidBodies_ += mol->getNRigidBodies();
201 nIntegrableObjects_ += mol->getNIntegrableObjects();
202 nCutoffGroups_ += mol->getNCutoffGroups();
203 nConstraints_ += mol->getNConstraintPairs();
204
205 addExcludePairs(mol);
206
207 return true;
208 } else {
209 return false;
210 }
211 }
212
213 bool SimInfo::removeMolecule(Molecule* mol) {
214 MoleculeIterator i;
215 i = molecules_.find(mol->getGlobalIndex());
216
217 if (i != molecules_.end() ) {
218
219 assert(mol == i->second);
220
221 nAtoms_ -= mol->getNAtoms();
222 nBonds_ -= mol->getNBonds();
223 nBends_ -= mol->getNBends();
224 nTorsions_ -= mol->getNTorsions();
225 nRigidBodies_ -= mol->getNRigidBodies();
226 nIntegrableObjects_ -= mol->getNIntegrableObjects();
227 nCutoffGroups_ -= mol->getNCutoffGroups();
228 nConstraints_ -= mol->getNConstraintPairs();
229
230 removeExcludePairs(mol);
231 molecules_.erase(mol->getGlobalIndex());
232
233 delete mol;
234
235 return true;
236 } else {
237 return false;
238 }
239
240
241 }
242
243
244 Molecule* SimInfo::beginMolecule(MoleculeIterator& i) {
245 i = molecules_.begin();
246 return i == molecules_.end() ? NULL : i->second;
247 }
248
249 Molecule* SimInfo::nextMolecule(MoleculeIterator& i) {
250 ++i;
251 return i == molecules_.end() ? NULL : i->second;
252 }
253
254
255 void SimInfo::calcNdf() {
256 int ndf_local;
257 MoleculeIterator i;
258 std::vector<StuntDouble*>::iterator j;
259 Molecule* mol;
260 StuntDouble* integrableObject;
261
262 ndf_local = 0;
263
264 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
265 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
266 integrableObject = mol->nextIntegrableObject(j)) {
267
268 ndf_local += 3;
269
270 if (integrableObject->isDirectional()) {
271 if (integrableObject->isLinear()) {
272 ndf_local += 2;
273 } else {
274 ndf_local += 3;
275 }
276 }
277
278 }
279 }
280
281 // n_constraints is local, so subtract them on each processor
282 ndf_local -= nConstraints_;
283
284 #ifdef IS_MPI
285 MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
286 #else
287 ndf_ = ndf_local;
288 #endif
289
290 // nZconstraints_ is global, as are the 3 COM translations for the
291 // entire system:
292 ndf_ = ndf_ - 3 - nZconstraint_;
293
294 }
295
296 int SimInfo::getFdf() {
297 #ifdef IS_MPI
298 MPI_Allreduce(&fdf_local,&fdf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
299 #else
300 fdf_ = fdf_local;
301 #endif
302 return fdf_;
303 }
304
305 void SimInfo::calcNdfRaw() {
306 int ndfRaw_local;
307
308 MoleculeIterator i;
309 std::vector<StuntDouble*>::iterator j;
310 Molecule* mol;
311 StuntDouble* integrableObject;
312
313 // Raw degrees of freedom that we have to set
314 ndfRaw_local = 0;
315
316 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
317 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
318 integrableObject = mol->nextIntegrableObject(j)) {
319
320 ndfRaw_local += 3;
321
322 if (integrableObject->isDirectional()) {
323 if (integrableObject->isLinear()) {
324 ndfRaw_local += 2;
325 } else {
326 ndfRaw_local += 3;
327 }
328 }
329
330 }
331 }
332
333 #ifdef IS_MPI
334 MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
335 #else
336 ndfRaw_ = ndfRaw_local;
337 #endif
338 }
339
340 void SimInfo::calcNdfTrans() {
341 int ndfTrans_local;
342
343 ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
344
345
346 #ifdef IS_MPI
347 MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
348 #else
349 ndfTrans_ = ndfTrans_local;
350 #endif
351
352 ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
353
354 }
355
356 void SimInfo::addExcludePairs(Molecule* mol) {
357 std::vector<Bond*>::iterator bondIter;
358 std::vector<Bend*>::iterator bendIter;
359 std::vector<Torsion*>::iterator torsionIter;
360 Bond* bond;
361 Bend* bend;
362 Torsion* torsion;
363 int a;
364 int b;
365 int c;
366 int d;
367
368 std::map<int, std::set<int> > atomGroups;
369
370 Molecule::RigidBodyIterator rbIter;
371 RigidBody* rb;
372 Molecule::IntegrableObjectIterator ii;
373 StuntDouble* integrableObject;
374
375 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
376 integrableObject = mol->nextIntegrableObject(ii)) {
377
378 if (integrableObject->isRigidBody()) {
379 rb = static_cast<RigidBody*>(integrableObject);
380 std::vector<Atom*> atoms = rb->getAtoms();
381 std::set<int> rigidAtoms;
382 for (int i = 0; i < atoms.size(); ++i) {
383 rigidAtoms.insert(atoms[i]->getGlobalIndex());
384 }
385 for (int i = 0; i < atoms.size(); ++i) {
386 atomGroups.insert(std::map<int, std::set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
387 }
388 } else {
389 std::set<int> oneAtomSet;
390 oneAtomSet.insert(integrableObject->getGlobalIndex());
391 atomGroups.insert(std::map<int, std::set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));
392 }
393 }
394
395
396
397 for (bond= mol->beginBond(bondIter); bond != NULL; bond = mol->nextBond(bondIter)) {
398 a = bond->getAtomA()->getGlobalIndex();
399 b = bond->getAtomB()->getGlobalIndex();
400 exclude_.addPair(a, b);
401 }
402
403 for (bend= mol->beginBend(bendIter); bend != NULL; bend = mol->nextBend(bendIter)) {
404 a = bend->getAtomA()->getGlobalIndex();
405 b = bend->getAtomB()->getGlobalIndex();
406 c = bend->getAtomC()->getGlobalIndex();
407 std::set<int> rigidSetA = getRigidSet(a, atomGroups);
408 std::set<int> rigidSetB = getRigidSet(b, atomGroups);
409 std::set<int> rigidSetC = getRigidSet(c, atomGroups);
410
411 exclude_.addPairs(rigidSetA, rigidSetB);
412 exclude_.addPairs(rigidSetA, rigidSetC);
413 exclude_.addPairs(rigidSetB, rigidSetC);
414
415 //exclude_.addPair(a, b);
416 //exclude_.addPair(a, c);
417 //exclude_.addPair(b, c);
418 }
419
420 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; torsion = mol->nextTorsion(torsionIter)) {
421 a = torsion->getAtomA()->getGlobalIndex();
422 b = torsion->getAtomB()->getGlobalIndex();
423 c = torsion->getAtomC()->getGlobalIndex();
424 d = torsion->getAtomD()->getGlobalIndex();
425 std::set<int> rigidSetA = getRigidSet(a, atomGroups);
426 std::set<int> rigidSetB = getRigidSet(b, atomGroups);
427 std::set<int> rigidSetC = getRigidSet(c, atomGroups);
428 std::set<int> rigidSetD = getRigidSet(d, atomGroups);
429
430 exclude_.addPairs(rigidSetA, rigidSetB);
431 exclude_.addPairs(rigidSetA, rigidSetC);
432 exclude_.addPairs(rigidSetA, rigidSetD);
433 exclude_.addPairs(rigidSetB, rigidSetC);
434 exclude_.addPairs(rigidSetB, rigidSetD);
435 exclude_.addPairs(rigidSetC, rigidSetD);
436
437 /*
438 exclude_.addPairs(rigidSetA.begin(), rigidSetA.end(), rigidSetB.begin(), rigidSetB.end());
439 exclude_.addPairs(rigidSetA.begin(), rigidSetA.end(), rigidSetC.begin(), rigidSetC.end());
440 exclude_.addPairs(rigidSetA.begin(), rigidSetA.end(), rigidSetD.begin(), rigidSetD.end());
441 exclude_.addPairs(rigidSetB.begin(), rigidSetB.end(), rigidSetC.begin(), rigidSetC.end());
442 exclude_.addPairs(rigidSetB.begin(), rigidSetB.end(), rigidSetD.begin(), rigidSetD.end());
443 exclude_.addPairs(rigidSetC.begin(), rigidSetC.end(), rigidSetD.begin(), rigidSetD.end());
444
445
446 exclude_.addPair(a, b);
447 exclude_.addPair(a, c);
448 exclude_.addPair(a, d);
449 exclude_.addPair(b, c);
450 exclude_.addPair(b, d);
451 exclude_.addPair(c, d);
452 */
453 }
454
455 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
456 std::vector<Atom*> atoms = rb->getAtoms();
457 for (int i = 0; i < atoms.size() -1 ; ++i) {
458 for (int j = i + 1; j < atoms.size(); ++j) {
459 a = atoms[i]->getGlobalIndex();
460 b = atoms[j]->getGlobalIndex();
461 exclude_.addPair(a, b);
462 }
463 }
464 }
465
466 }
467
468 void SimInfo::removeExcludePairs(Molecule* mol) {
469 std::vector<Bond*>::iterator bondIter;
470 std::vector<Bend*>::iterator bendIter;
471 std::vector<Torsion*>::iterator torsionIter;
472 Bond* bond;
473 Bend* bend;
474 Torsion* torsion;
475 int a;
476 int b;
477 int c;
478 int d;
479
480 std::map<int, std::set<int> > atomGroups;
481
482 Molecule::RigidBodyIterator rbIter;
483 RigidBody* rb;
484 Molecule::IntegrableObjectIterator ii;
485 StuntDouble* integrableObject;
486
487 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
488 integrableObject = mol->nextIntegrableObject(ii)) {
489
490 if (integrableObject->isRigidBody()) {
491 rb = static_cast<RigidBody*>(integrableObject);
492 std::vector<Atom*> atoms = rb->getAtoms();
493 std::set<int> rigidAtoms;
494 for (int i = 0; i < atoms.size(); ++i) {
495 rigidAtoms.insert(atoms[i]->getGlobalIndex());
496 }
497 for (int i = 0; i < atoms.size(); ++i) {
498 atomGroups.insert(std::map<int, std::set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
499 }
500 } else {
501 std::set<int> oneAtomSet;
502 oneAtomSet.insert(integrableObject->getGlobalIndex());
503 atomGroups.insert(std::map<int, std::set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));
504 }
505 }
506
507
508 for (bond= mol->beginBond(bondIter); bond != NULL; bond = mol->nextBond(bondIter)) {
509 a = bond->getAtomA()->getGlobalIndex();
510 b = bond->getAtomB()->getGlobalIndex();
511 exclude_.removePair(a, b);
512 }
513
514 for (bend= mol->beginBend(bendIter); bend != NULL; bend = mol->nextBend(bendIter)) {
515 a = bend->getAtomA()->getGlobalIndex();
516 b = bend->getAtomB()->getGlobalIndex();
517 c = bend->getAtomC()->getGlobalIndex();
518
519 std::set<int> rigidSetA = getRigidSet(a, atomGroups);
520 std::set<int> rigidSetB = getRigidSet(b, atomGroups);
521 std::set<int> rigidSetC = getRigidSet(c, atomGroups);
522
523 exclude_.removePairs(rigidSetA, rigidSetB);
524 exclude_.removePairs(rigidSetA, rigidSetC);
525 exclude_.removePairs(rigidSetB, rigidSetC);
526
527 //exclude_.removePair(a, b);
528 //exclude_.removePair(a, c);
529 //exclude_.removePair(b, c);
530 }
531
532 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; torsion = mol->nextTorsion(torsionIter)) {
533 a = torsion->getAtomA()->getGlobalIndex();
534 b = torsion->getAtomB()->getGlobalIndex();
535 c = torsion->getAtomC()->getGlobalIndex();
536 d = torsion->getAtomD()->getGlobalIndex();
537
538 std::set<int> rigidSetA = getRigidSet(a, atomGroups);
539 std::set<int> rigidSetB = getRigidSet(b, atomGroups);
540 std::set<int> rigidSetC = getRigidSet(c, atomGroups);
541 std::set<int> rigidSetD = getRigidSet(d, atomGroups);
542
543 exclude_.removePairs(rigidSetA, rigidSetB);
544 exclude_.removePairs(rigidSetA, rigidSetC);
545 exclude_.removePairs(rigidSetA, rigidSetD);
546 exclude_.removePairs(rigidSetB, rigidSetC);
547 exclude_.removePairs(rigidSetB, rigidSetD);
548 exclude_.removePairs(rigidSetC, rigidSetD);
549
550 /*
551 exclude_.removePairs(rigidSetA.begin(), rigidSetA.end(), rigidSetB.begin(), rigidSetB.end());
552 exclude_.removePairs(rigidSetA.begin(), rigidSetA.end(), rigidSetC.begin(), rigidSetC.end());
553 exclude_.removePairs(rigidSetA.begin(), rigidSetA.end(), rigidSetD.begin(), rigidSetD.end());
554 exclude_.removePairs(rigidSetB.begin(), rigidSetB.end(), rigidSetC.begin(), rigidSetC.end());
555 exclude_.removePairs(rigidSetB.begin(), rigidSetB.end(), rigidSetD.begin(), rigidSetD.end());
556 exclude_.removePairs(rigidSetC.begin(), rigidSetC.end(), rigidSetD.begin(), rigidSetD.end());
557
558
559 exclude_.removePair(a, b);
560 exclude_.removePair(a, c);
561 exclude_.removePair(a, d);
562 exclude_.removePair(b, c);
563 exclude_.removePair(b, d);
564 exclude_.removePair(c, d);
565 */
566 }
567
568 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
569 std::vector<Atom*> atoms = rb->getAtoms();
570 for (int i = 0; i < atoms.size() -1 ; ++i) {
571 for (int j = i + 1; j < atoms.size(); ++j) {
572 a = atoms[i]->getGlobalIndex();
573 b = atoms[j]->getGlobalIndex();
574 exclude_.removePair(a, b);
575 }
576 }
577 }
578
579 }
580
581
582 void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) {
583 int curStampId;
584
585 //index from 0
586 curStampId = moleculeStamps_.size();
587
588 moleculeStamps_.push_back(molStamp);
589 molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
590 }
591
592 void SimInfo::update() {
593
594 setupSimType();
595
596 #ifdef IS_MPI
597 setupFortranParallel();
598 #endif
599
600 setupFortranSim();
601
602 //setup fortran force field
603 /** @deprecate */
604 int isError = 0;
605
606 setupCutoff();
607
608 setupElectrostaticSummationMethod( isError );
609 setupSwitchingFunction();
610 setupAccumulateBoxDipole();
611
612 if(isError){
613 sprintf( painCave.errMsg,
614 "ForceField error: There was an error initializing the forceField in fortran.\n" );
615 painCave.isFatal = 1;
616 simError();
617 }
618
619 calcNdf();
620 calcNdfRaw();
621 calcNdfTrans();
622
623 fortranInitialized_ = true;
624 }
625
626 std::set<AtomType*> SimInfo::getUniqueAtomTypes() {
627 SimInfo::MoleculeIterator mi;
628 Molecule* mol;
629 Molecule::AtomIterator ai;
630 Atom* atom;
631 std::set<AtomType*> atomTypes;
632
633 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
634
635 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
636 atomTypes.insert(atom->getAtomType());
637 }
638
639 }
640
641 return atomTypes;
642 }
643
644 void SimInfo::setupSimType() {
645 std::set<AtomType*>::iterator i;
646 std::set<AtomType*> atomTypes;
647 atomTypes = getUniqueAtomTypes();
648
649 int useLennardJones = 0;
650 int useElectrostatic = 0;
651 int useEAM = 0;
652 int useSC = 0;
653 int useCharge = 0;
654 int useDirectional = 0;
655 int useDipole = 0;
656 int useGayBerne = 0;
657 int useSticky = 0;
658 int useStickyPower = 0;
659 int useShape = 0;
660 int useFLARB = 0; //it is not in AtomType yet
661 int useDirectionalAtom = 0;
662 int useElectrostatics = 0;
663 //usePBC and useRF are from simParams
664 int usePBC = simParams_->getUsePeriodicBoundaryConditions();
665 int useRF;
666 int useSF;
667 int useSP;
668 int useBoxDipole;
669 std::string myMethod;
670
671 // set the useRF logical
672 useRF = 0;
673 useSF = 0;
674 useSP = 0;
675
676
677 if (simParams_->haveElectrostaticSummationMethod()) {
678 std::string myMethod = simParams_->getElectrostaticSummationMethod();
679 toUpper(myMethod);
680 if (myMethod == "REACTION_FIELD"){
681 useRF = 1;
682 } else if (myMethod == "SHIFTED_FORCE"){
683 useSF = 1;
684 } else if (myMethod == "SHIFTED_POTENTIAL"){
685 useSP = 1;
686 }
687 }
688
689 if (simParams_->haveAccumulateBoxDipole())
690 if (simParams_->getAccumulateBoxDipole())
691 useBoxDipole = 1;
692
693 //loop over all of the atom types
694 for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
695 useLennardJones |= (*i)->isLennardJones();
696 useElectrostatic |= (*i)->isElectrostatic();
697 useEAM |= (*i)->isEAM();
698 useSC |= (*i)->isSC();
699 useCharge |= (*i)->isCharge();
700 useDirectional |= (*i)->isDirectional();
701 useDipole |= (*i)->isDipole();
702 useGayBerne |= (*i)->isGayBerne();
703 useSticky |= (*i)->isSticky();
704 useStickyPower |= (*i)->isStickyPower();
705 useShape |= (*i)->isShape();
706 }
707
708 if (useSticky || useStickyPower || useDipole || useGayBerne || useShape) {
709 useDirectionalAtom = 1;
710 }
711
712 if (useCharge || useDipole) {
713 useElectrostatics = 1;
714 }
715
716 #ifdef IS_MPI
717 int temp;
718
719 temp = usePBC;
720 MPI_Allreduce(&temp, &usePBC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
721
722 temp = useDirectionalAtom;
723 MPI_Allreduce(&temp, &useDirectionalAtom, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
724
725 temp = useLennardJones;
726 MPI_Allreduce(&temp, &useLennardJones, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
727
728 temp = useElectrostatics;
729 MPI_Allreduce(&temp, &useElectrostatics, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
730
731 temp = useCharge;
732 MPI_Allreduce(&temp, &useCharge, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
733
734 temp = useDipole;
735 MPI_Allreduce(&temp, &useDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
736
737 temp = useSticky;
738 MPI_Allreduce(&temp, &useSticky, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
739
740 temp = useStickyPower;
741 MPI_Allreduce(&temp, &useStickyPower, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
742
743 temp = useGayBerne;
744 MPI_Allreduce(&temp, &useGayBerne, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
745
746 temp = useEAM;
747 MPI_Allreduce(&temp, &useEAM, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
748
749 temp = useSC;
750 MPI_Allreduce(&temp, &useSC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
751
752 temp = useShape;
753 MPI_Allreduce(&temp, &useShape, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
754
755 temp = useFLARB;
756 MPI_Allreduce(&temp, &useFLARB, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
757
758 temp = useRF;
759 MPI_Allreduce(&temp, &useRF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
760
761 temp = useSF;
762 MPI_Allreduce(&temp, &useSF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
763
764 temp = useSP;
765 MPI_Allreduce(&temp, &useSP, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
766
767 temp = useBoxDipole;
768 MPI_Allreduce(&temp, &useBoxDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
769
770 #endif
771
772 fInfo_.SIM_uses_PBC = usePBC;
773 fInfo_.SIM_uses_DirectionalAtoms = useDirectionalAtom;
774 fInfo_.SIM_uses_LennardJones = useLennardJones;
775 fInfo_.SIM_uses_Electrostatics = useElectrostatics;
776 fInfo_.SIM_uses_Charges = useCharge;
777 fInfo_.SIM_uses_Dipoles = useDipole;
778 fInfo_.SIM_uses_Sticky = useSticky;
779 fInfo_.SIM_uses_StickyPower = useStickyPower;
780 fInfo_.SIM_uses_GayBerne = useGayBerne;
781 fInfo_.SIM_uses_EAM = useEAM;
782 fInfo_.SIM_uses_SC = useSC;
783 fInfo_.SIM_uses_Shapes = useShape;
784 fInfo_.SIM_uses_FLARB = useFLARB;
785 fInfo_.SIM_uses_RF = useRF;
786 fInfo_.SIM_uses_SF = useSF;
787 fInfo_.SIM_uses_SP = useSP;
788 fInfo_.SIM_uses_BoxDipole = useBoxDipole;
789 }
790
791 void SimInfo::setupFortranSim() {
792 int isError;
793 int nExclude;
794 std::vector<int> fortranGlobalGroupMembership;
795
796 nExclude = exclude_.getSize();
797 isError = 0;
798
799 //globalGroupMembership_ is filled by SimCreator
800 for (int i = 0; i < nGlobalAtoms_; i++) {
801 fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
802 }
803
804 //calculate mass ratio of cutoff group
805 std::vector<RealType> mfact;
806 SimInfo::MoleculeIterator mi;
807 Molecule* mol;
808 Molecule::CutoffGroupIterator ci;
809 CutoffGroup* cg;
810 Molecule::AtomIterator ai;
811 Atom* atom;
812 RealType totalMass;
813
814 //to avoid memory reallocation, reserve enough space for mfact
815 mfact.reserve(getNCutoffGroups());
816
817 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
818 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
819
820 totalMass = cg->getMass();
821 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
822 // Check for massless groups - set mfact to 1 if true
823 if (totalMass != 0)
824 mfact.push_back(atom->getMass()/totalMass);
825 else
826 mfact.push_back( 1.0 );
827 }
828
829 }
830 }
831
832 //fill ident array of local atoms (it is actually ident of AtomType, it is so confusing !!!)
833 std::vector<int> identArray;
834
835 //to avoid memory reallocation, reserve enough space identArray
836 identArray.reserve(getNAtoms());
837
838 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
839 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
840 identArray.push_back(atom->getIdent());
841 }
842 }
843
844 //fill molMembershipArray
845 //molMembershipArray is filled by SimCreator
846 std::vector<int> molMembershipArray(nGlobalAtoms_);
847 for (int i = 0; i < nGlobalAtoms_; i++) {
848 molMembershipArray[i] = globalMolMembership_[i] + 1;
849 }
850
851 //setup fortran simulation
852 int nGlobalExcludes = 0;
853 int* globalExcludes = NULL;
854 int* excludeList = exclude_.getExcludeList();
855 setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray[0], &nExclude, excludeList ,
856 &nGlobalExcludes, globalExcludes, &molMembershipArray[0],
857 &mfact[0], &nCutoffGroups_, &fortranGlobalGroupMembership[0], &isError);
858
859 if( isError ){
860
861 sprintf( painCave.errMsg,
862 "There was an error setting the simulation information in fortran.\n" );
863 painCave.isFatal = 1;
864 painCave.severity = OOPSE_ERROR;
865 simError();
866 }
867
868 #ifdef IS_MPI
869 sprintf( checkPointMsg,
870 "succesfully sent the simulation information to fortran.\n");
871 MPIcheckPoint();
872 #endif // is_mpi
873
874 // Setup number of neighbors in neighbor list if present
875 if (simParams_->haveNeighborListNeighbors()) {
876 setNeighbors(simParams_->getNeighborListNeighbors());
877 }
878
879
880 }
881
882
883 #ifdef IS_MPI
884 void SimInfo::setupFortranParallel() {
885
886 //SimInfo is responsible for creating localToGlobalAtomIndex and localToGlobalGroupIndex
887 std::vector<int> localToGlobalAtomIndex(getNAtoms(), 0);
888 std::vector<int> localToGlobalCutoffGroupIndex;
889 SimInfo::MoleculeIterator mi;
890 Molecule::AtomIterator ai;
891 Molecule::CutoffGroupIterator ci;
892 Molecule* mol;
893 Atom* atom;
894 CutoffGroup* cg;
895 mpiSimData parallelData;
896 int isError;
897
898 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
899
900 //local index(index in DataStorge) of atom is important
901 for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
902 localToGlobalAtomIndex[atom->getLocalIndex()] = atom->getGlobalIndex() + 1;
903 }
904
905 //local index of cutoff group is trivial, it only depends on the order of travesing
906 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
907 localToGlobalCutoffGroupIndex.push_back(cg->getGlobalIndex() + 1);
908 }
909
910 }
911
912 //fill up mpiSimData struct
913 parallelData.nMolGlobal = getNGlobalMolecules();
914 parallelData.nMolLocal = getNMolecules();
915 parallelData.nAtomsGlobal = getNGlobalAtoms();
916 parallelData.nAtomsLocal = getNAtoms();
917 parallelData.nGroupsGlobal = getNGlobalCutoffGroups();
918 parallelData.nGroupsLocal = getNCutoffGroups();
919 parallelData.myNode = worldRank;
920 MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors));
921
922 //pass mpiSimData struct and index arrays to fortran
923 setFsimParallel(&parallelData, &(parallelData.nAtomsLocal),
924 &localToGlobalAtomIndex[0], &(parallelData.nGroupsLocal),
925 &localToGlobalCutoffGroupIndex[0], &isError);
926
927 if (isError) {
928 sprintf(painCave.errMsg,
929 "mpiRefresh errror: fortran didn't like something we gave it.\n");
930 painCave.isFatal = 1;
931 simError();
932 }
933
934 sprintf(checkPointMsg, " mpiRefresh successful.\n");
935 MPIcheckPoint();
936
937
938 }
939
940 #endif
941
942 void SimInfo::setupCutoff() {
943
944 ForceFieldOptions& forceFieldOptions_ = forceField_->getForceFieldOptions();
945
946 // Check the cutoff policy
947 int cp = TRADITIONAL_CUTOFF_POLICY; // Set to traditional by default
948
949 std::string myPolicy;
950 if (forceFieldOptions_.haveCutoffPolicy()){
951 myPolicy = forceFieldOptions_.getCutoffPolicy();
952 }else if (simParams_->haveCutoffPolicy()) {
953 myPolicy = simParams_->getCutoffPolicy();
954 }
955
956 if (!myPolicy.empty()){
957 toUpper(myPolicy);
958 if (myPolicy == "MIX") {
959 cp = MIX_CUTOFF_POLICY;
960 } else {
961 if (myPolicy == "MAX") {
962 cp = MAX_CUTOFF_POLICY;
963 } else {
964 if (myPolicy == "TRADITIONAL") {
965 cp = TRADITIONAL_CUTOFF_POLICY;
966 } else {
967 // throw error
968 sprintf( painCave.errMsg,
969 "SimInfo error: Unknown cutoffPolicy. (Input file specified %s .)\n\tcutoffPolicy must be one of: \"Mix\", \"Max\", or \"Traditional\".", myPolicy.c_str() );
970 painCave.isFatal = 1;
971 simError();
972 }
973 }
974 }
975 }
976 notifyFortranCutoffPolicy(&cp);
977
978 // Check the Skin Thickness for neighborlists
979 RealType skin;
980 if (simParams_->haveSkinThickness()) {
981 skin = simParams_->getSkinThickness();
982 notifyFortranSkinThickness(&skin);
983 }
984
985 // Check if the cutoff was set explicitly:
986 if (simParams_->haveCutoffRadius()) {
987 rcut_ = simParams_->getCutoffRadius();
988 if (simParams_->haveSwitchingRadius()) {
989 rsw_ = simParams_->getSwitchingRadius();
990 } else {
991 if (fInfo_.SIM_uses_Charges |
992 fInfo_.SIM_uses_Dipoles |
993 fInfo_.SIM_uses_RF) {
994
995 rsw_ = 0.85 * rcut_;
996 sprintf(painCave.errMsg,
997 "SimCreator Warning: No value was set for the switchingRadius.\n"
998 "\tOOPSE will use a default value of 85 percent of the cutoffRadius.\n"
999 "\tswitchingRadius = %f. for this simulation\n", rsw_);
1000 painCave.isFatal = 0;
1001 simError();
1002 } else {
1003 rsw_ = rcut_;
1004 sprintf(painCave.errMsg,
1005 "SimCreator Warning: No value was set for the switchingRadius.\n"
1006 "\tOOPSE will use the same value as the cutoffRadius.\n"
1007 "\tswitchingRadius = %f. for this simulation\n", rsw_);
1008 painCave.isFatal = 0;
1009 simError();
1010 }
1011 }
1012
1013 notifyFortranCutoffs(&rcut_, &rsw_);
1014
1015 } else {
1016
1017 // For electrostatic atoms, we'll assume a large safe value:
1018 if (fInfo_.SIM_uses_Charges | fInfo_.SIM_uses_Dipoles | fInfo_.SIM_uses_RF) {
1019 sprintf(painCave.errMsg,
1020 "SimCreator Warning: No value was set for the cutoffRadius.\n"
1021 "\tOOPSE will use a default value of 15.0 angstroms"
1022 "\tfor the cutoffRadius.\n");
1023 painCave.isFatal = 0;
1024 simError();
1025 rcut_ = 15.0;
1026
1027 if (simParams_->haveElectrostaticSummationMethod()) {
1028 std::string myMethod = simParams_->getElectrostaticSummationMethod();
1029 toUpper(myMethod);
1030 if (myMethod == "SHIFTED_POTENTIAL" || myMethod == "SHIFTED_FORCE") {
1031 if (simParams_->haveSwitchingRadius()){
1032 sprintf(painCave.errMsg,
1033 "SimInfo Warning: A value was set for the switchingRadius\n"
1034 "\teven though the electrostaticSummationMethod was\n"
1035 "\tset to %s\n", myMethod.c_str());
1036 painCave.isFatal = 1;
1037 simError();
1038 }
1039 }
1040 }
1041
1042 if (simParams_->haveSwitchingRadius()){
1043 rsw_ = simParams_->getSwitchingRadius();
1044 } else {
1045 sprintf(painCave.errMsg,
1046 "SimCreator Warning: No value was set for switchingRadius.\n"
1047 "\tOOPSE will use a default value of\n"
1048 "\t0.85 * cutoffRadius for the switchingRadius\n");
1049 painCave.isFatal = 0;
1050 simError();
1051 rsw_ = 0.85 * rcut_;
1052 }
1053 notifyFortranCutoffs(&rcut_, &rsw_);
1054 } else {
1055 // We didn't set rcut explicitly, and we don't have electrostatic atoms, so
1056 // We'll punt and let fortran figure out the cutoffs later.
1057
1058 notifyFortranYouAreOnYourOwn();
1059
1060 }
1061 }
1062 }
1063
1064 void SimInfo::setupElectrostaticSummationMethod( int isError ) {
1065
1066 int errorOut;
1067 int esm = NONE;
1068 int sm = UNDAMPED;
1069 RealType alphaVal;
1070 RealType dielectric;
1071
1072 errorOut = isError;
1073
1074 if (simParams_->haveElectrostaticSummationMethod()) {
1075 std::string myMethod = simParams_->getElectrostaticSummationMethod();
1076 toUpper(myMethod);
1077 if (myMethod == "NONE") {
1078 esm = NONE;
1079 } else {
1080 if (myMethod == "SWITCHING_FUNCTION") {
1081 esm = SWITCHING_FUNCTION;
1082 } else {
1083 if (myMethod == "SHIFTED_POTENTIAL") {
1084 esm = SHIFTED_POTENTIAL;
1085 } else {
1086 if (myMethod == "SHIFTED_FORCE") {
1087 esm = SHIFTED_FORCE;
1088 } else {
1089 if (myMethod == "REACTION_FIELD") {
1090 esm = REACTION_FIELD;
1091 dielectric = simParams_->getDielectric();
1092 if (!simParams_->haveDielectric()) {
1093 // throw warning
1094 sprintf( painCave.errMsg,
1095 "SimInfo warning: dielectric was not specified in the input file\n\tfor the reaction field correction method.\n"
1096 "\tA default value of %f will be used for the dielectric.\n", dielectric);
1097 painCave.isFatal = 0;
1098 simError();
1099 }
1100 } else {
1101 // throw error
1102 sprintf( painCave.errMsg,
1103 "SimInfo error: Unknown electrostaticSummationMethod.\n"
1104 "\t(Input file specified %s .)\n"
1105 "\telectrostaticSummationMethod must be one of: \"none\",\n"
1106 "\t\"shifted_potential\", \"shifted_force\", or \n"
1107 "\t\"reaction_field\".\n", myMethod.c_str() );
1108 painCave.isFatal = 1;
1109 simError();
1110 }
1111 }
1112 }
1113 }
1114 }
1115 }
1116
1117 if (simParams_->haveElectrostaticScreeningMethod()) {
1118 std::string myScreen = simParams_->getElectrostaticScreeningMethod();
1119 toUpper(myScreen);
1120 if (myScreen == "UNDAMPED") {
1121 sm = UNDAMPED;
1122 } else {
1123 if (myScreen == "DAMPED") {
1124 sm = DAMPED;
1125 if (!simParams_->haveDampingAlpha()) {
1126 // first set a cutoff dependent alpha value
1127 // we assume alpha depends linearly with rcut from 0 to 20.5 ang
1128 alphaVal = 0.5125 - rcut_* 0.025;
1129 // for values rcut > 20.5, alpha is zero
1130 if (alphaVal < 0) alphaVal = 0;
1131
1132 // throw warning
1133 sprintf( painCave.errMsg,
1134 "SimInfo warning: dampingAlpha was not specified in the input file.\n"
1135 "\tA default value of %f (1/ang) will be used for the cutoff of\n\t%f (ang).\n", alphaVal, rcut_);
1136 painCave.isFatal = 0;
1137 simError();
1138 } else {
1139 alphaVal = simParams_->getDampingAlpha();
1140 }
1141
1142 } else {
1143 // throw error
1144 sprintf( painCave.errMsg,
1145 "SimInfo error: Unknown electrostaticScreeningMethod.\n"
1146 "\t(Input file specified %s .)\n"
1147 "\telectrostaticScreeningMethod must be one of: \"undamped\"\n"
1148 "or \"damped\".\n", myScreen.c_str() );
1149 painCave.isFatal = 1;
1150 simError();
1151 }
1152 }
1153 }
1154
1155 // let's pass some summation method variables to fortran
1156 setElectrostaticSummationMethod( &esm );
1157 setFortranElectrostaticMethod( &esm );
1158 setScreeningMethod( &sm );
1159 setDampingAlpha( &alphaVal );
1160 setReactionFieldDielectric( &dielectric );
1161 initFortranFF( &errorOut );
1162 }
1163
1164 void SimInfo::setupSwitchingFunction() {
1165 int ft = CUBIC;
1166
1167 if (simParams_->haveSwitchingFunctionType()) {
1168 std::string funcType = simParams_->getSwitchingFunctionType();
1169 toUpper(funcType);
1170 if (funcType == "CUBIC") {
1171 ft = CUBIC;
1172 } else {
1173 if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
1174 ft = FIFTH_ORDER_POLY;
1175 } else {
1176 // throw error
1177 sprintf( painCave.errMsg,
1178 "SimInfo error: Unknown switchingFunctionType. (Input file specified %s .)\n\tswitchingFunctionType must be one of: \"cubic\" or \"fifth_order_polynomial\".", funcType.c_str() );
1179 painCave.isFatal = 1;
1180 simError();
1181 }
1182 }
1183 }
1184
1185 // send switching function notification to switcheroo
1186 setFunctionType(&ft);
1187
1188 }
1189
1190 void SimInfo::setupAccumulateBoxDipole() {
1191
1192 // we only call setAccumulateBoxDipole if the accumulateBoxDipole parameter is true
1193 if ( simParams_->haveAccumulateBoxDipole() )
1194 if ( simParams_->getAccumulateBoxDipole() ) {
1195 setAccumulateBoxDipole();
1196 calcBoxDipole_ = true;
1197 }
1198
1199 }
1200
1201 void SimInfo::addProperty(GenericData* genData) {
1202 properties_.addProperty(genData);
1203 }
1204
1205 void SimInfo::removeProperty(const std::string& propName) {
1206 properties_.removeProperty(propName);
1207 }
1208
1209 void SimInfo::clearProperties() {
1210 properties_.clearProperties();
1211 }
1212
1213 std::vector<std::string> SimInfo::getPropertyNames() {
1214 return properties_.getPropertyNames();
1215 }
1216
1217 std::vector<GenericData*> SimInfo::getProperties() {
1218 return properties_.getProperties();
1219 }
1220
1221 GenericData* SimInfo::getPropertyByName(const std::string& propName) {
1222 return properties_.getPropertyByName(propName);
1223 }
1224
1225 void SimInfo::setSnapshotManager(SnapshotManager* sman) {
1226 if (sman_ == sman) {
1227 return;
1228 }
1229 delete sman_;
1230 sman_ = sman;
1231
1232 Molecule* mol;
1233 RigidBody* rb;
1234 Atom* atom;
1235 SimInfo::MoleculeIterator mi;
1236 Molecule::RigidBodyIterator rbIter;
1237 Molecule::AtomIterator atomIter;;
1238
1239 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
1240
1241 for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) {
1242 atom->setSnapshotManager(sman_);
1243 }
1244
1245 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
1246 rb->setSnapshotManager(sman_);
1247 }
1248 }
1249
1250 }
1251
1252 Vector3d SimInfo::getComVel(){
1253 SimInfo::MoleculeIterator i;
1254 Molecule* mol;
1255
1256 Vector3d comVel(0.0);
1257 RealType totalMass = 0.0;
1258
1259
1260 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1261 RealType mass = mol->getMass();
1262 totalMass += mass;
1263 comVel += mass * mol->getComVel();
1264 }
1265
1266 #ifdef IS_MPI
1267 RealType tmpMass = totalMass;
1268 Vector3d tmpComVel(comVel);
1269 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1270 MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1271 #endif
1272
1273 comVel /= totalMass;
1274
1275 return comVel;
1276 }
1277
1278 Vector3d SimInfo::getCom(){
1279 SimInfo::MoleculeIterator i;
1280 Molecule* mol;
1281
1282 Vector3d com(0.0);
1283 RealType totalMass = 0.0;
1284
1285 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1286 RealType mass = mol->getMass();
1287 totalMass += mass;
1288 com += mass * mol->getCom();
1289 }
1290
1291 #ifdef IS_MPI
1292 RealType tmpMass = totalMass;
1293 Vector3d tmpCom(com);
1294 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1295 MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1296 #endif
1297
1298 com /= totalMass;
1299
1300 return com;
1301
1302 }
1303
1304 std::ostream& operator <<(std::ostream& o, SimInfo& info) {
1305
1306 return o;
1307 }
1308
1309
1310 /*
1311 Returns center of mass and center of mass velocity in one function call.
1312 */
1313
1314 void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){
1315 SimInfo::MoleculeIterator i;
1316 Molecule* mol;
1317
1318
1319 RealType totalMass = 0.0;
1320
1321
1322 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1323 RealType mass = mol->getMass();
1324 totalMass += mass;
1325 com += mass * mol->getCom();
1326 comVel += mass * mol->getComVel();
1327 }
1328
1329 #ifdef IS_MPI
1330 RealType tmpMass = totalMass;
1331 Vector3d tmpCom(com);
1332 Vector3d tmpComVel(comVel);
1333 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1334 MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1335 MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1336 #endif
1337
1338 com /= totalMass;
1339 comVel /= totalMass;
1340 }
1341
1342 /*
1343 Return intertia tensor for entire system and angular momentum Vector.
1344
1345
1346 [ Ixx -Ixy -Ixz ]
1347 J =| -Iyx Iyy -Iyz |
1348 [ -Izx -Iyz Izz ]
1349 */
1350
1351 void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){
1352
1353
1354 RealType xx = 0.0;
1355 RealType yy = 0.0;
1356 RealType zz = 0.0;
1357 RealType xy = 0.0;
1358 RealType xz = 0.0;
1359 RealType yz = 0.0;
1360 Vector3d com(0.0);
1361 Vector3d comVel(0.0);
1362
1363 getComAll(com, comVel);
1364
1365 SimInfo::MoleculeIterator i;
1366 Molecule* mol;
1367
1368 Vector3d thisq(0.0);
1369 Vector3d thisv(0.0);
1370
1371 RealType thisMass = 0.0;
1372
1373
1374
1375
1376 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1377
1378 thisq = mol->getCom()-com;
1379 thisv = mol->getComVel()-comVel;
1380 thisMass = mol->getMass();
1381 // Compute moment of intertia coefficients.
1382 xx += thisq[0]*thisq[0]*thisMass;
1383 yy += thisq[1]*thisq[1]*thisMass;
1384 zz += thisq[2]*thisq[2]*thisMass;
1385
1386 // compute products of intertia
1387 xy += thisq[0]*thisq[1]*thisMass;
1388 xz += thisq[0]*thisq[2]*thisMass;
1389 yz += thisq[1]*thisq[2]*thisMass;
1390
1391 angularMomentum += cross( thisq, thisv ) * thisMass;
1392
1393 }
1394
1395
1396 inertiaTensor(0,0) = yy + zz;
1397 inertiaTensor(0,1) = -xy;
1398 inertiaTensor(0,2) = -xz;
1399 inertiaTensor(1,0) = -xy;
1400 inertiaTensor(1,1) = xx + zz;
1401 inertiaTensor(1,2) = -yz;
1402 inertiaTensor(2,0) = -xz;
1403 inertiaTensor(2,1) = -yz;
1404 inertiaTensor(2,2) = xx + yy;
1405
1406 #ifdef IS_MPI
1407 Mat3x3d tmpI(inertiaTensor);
1408 Vector3d tmpAngMom;
1409 MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1410 MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1411 #endif
1412
1413 return;
1414 }
1415
1416 //Returns the angular momentum of the system
1417 Vector3d SimInfo::getAngularMomentum(){
1418
1419 Vector3d com(0.0);
1420 Vector3d comVel(0.0);
1421 Vector3d angularMomentum(0.0);
1422
1423 getComAll(com,comVel);
1424
1425 SimInfo::MoleculeIterator i;
1426 Molecule* mol;
1427
1428 Vector3d thisr(0.0);
1429 Vector3d thisp(0.0);
1430
1431 RealType thisMass;
1432
1433 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1434 thisMass = mol->getMass();
1435 thisr = mol->getCom()-com;
1436 thisp = (mol->getComVel()-comVel)*thisMass;
1437
1438 angularMomentum += cross( thisr, thisp );
1439
1440 }
1441
1442 #ifdef IS_MPI
1443 Vector3d tmpAngMom;
1444 MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1445 #endif
1446
1447 return angularMomentum;
1448 }
1449
1450 StuntDouble* SimInfo::getIOIndexToIntegrableObject(int index) {
1451 return IOIndexToIntegrableObject.at(index);
1452 }
1453
1454 void SimInfo::setIOIndexToIntegrableObject(const std::vector<StuntDouble*>& v) {
1455 IOIndexToIntegrableObject= v;
1456 }
1457
1458 /*
1459 void SimInfo::setStuntDoubleFromGlobalIndex(std::vector<StuntDouble*> v) {
1460 assert( v.size() == nAtoms_ + nRigidBodies_);
1461 sdByGlobalIndex_ = v;
1462 }
1463
1464 StuntDouble* SimInfo::getStuntDoubleFromGlobalIndex(int index) {
1465 //assert(index < nAtoms_ + nRigidBodies_);
1466 return sdByGlobalIndex_.at(index);
1467 }
1468 */
1469 }//end namespace oopse
1470