ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/brains/SimInfo.cpp
Revision: 3054
Committed: Wed Oct 18 21:58:48 2006 UTC (18 years, 6 months ago) by gezelter
File size: 46109 byte(s)
Log Message:
fixing a wrapVector problem in staticProps, also making Shifted force
and electrostatic damping the default behavior

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