ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/brains/SimInfo.cpp
Revision: 2552
Committed: Thu Jan 12 16:47:25 2006 UTC (18 years, 5 months ago) by chrisfen
File size: 43521 byte(s)
Log Message:
unifying function name in electrostatics

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