ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1581
Committed: Mon Jun 13 22:13:12 2011 UTC (13 years, 10 months ago) by gezelter
File size: 38107 byte(s)
Log Message:
bug fixes

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. Redistributions of source code must retain the above copyright
10 * notice, this list of conditions and the following disclaimer.
11 *
12 * 2. Redistributions in binary form must reproduce the above copyright
13 * notice, this list of conditions and the following disclaimer in the
14 * documentation and/or other materials provided with the
15 * distribution.
16 *
17 * This software is provided "AS IS," without a warranty of any
18 * kind. All express or implied conditions, representations and
19 * warranties, including any implied warranty of merchantability,
20 * fitness for a particular purpose or non-infringement, are hereby
21 * excluded. The University of Notre Dame and its licensors shall not
22 * be liable for any damages suffered by licensee as a result of
23 * using, modifying or distributing the software or its
24 * derivatives. In no event will the University of Notre Dame or its
25 * licensors be liable for any lost revenue, profit or data, or for
26 * direct, indirect, special, consequential, incidental or punitive
27 * damages, however caused and regardless of the theory of liability,
28 * arising out of the use of or inability to use software, even if the
29 * University of Notre Dame has been advised of the possibility of
30 * such damages.
31 *
32 * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33 * research, please cite the appropriate papers when you publish your
34 * work. Good starting points are:
35 *
36 * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37 * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39 * [4] Vardeman & Gezelter, in progress (2009).
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 "utils/MemoryUtils.hpp"
58 #include "utils/simError.h"
59 #include "selection/SelectionManager.hpp"
60 #include "io/ForceFieldOptions.hpp"
61 #include "UseTheForce/ForceField.hpp"
62 #include "nonbonded/SwitchingFunction.hpp"
63
64 using namespace std;
65 namespace OpenMD {
66
67 SimInfo::SimInfo(ForceField* ff, Globals* simParams) :
68 forceField_(ff), simParams_(simParams),
69 ndf_(0), fdf_local(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0),
70 nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0),
71 nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
72 nAtoms_(0), nBonds_(0), nBends_(0), nTorsions_(0), nInversions_(0),
73 nRigidBodies_(0), nIntegrableObjects_(0), nCutoffGroups_(0),
74 nConstraints_(0), sman_(NULL), topologyDone_(false),
75 calcBoxDipole_(false), useAtomicVirial_(true) {
76
77 MoleculeStamp* molStamp;
78 int nMolWithSameStamp;
79 int nCutoffAtoms = 0; // number of atoms belong to cutoff groups
80 int nGroups = 0; //total cutoff groups defined in meta-data file
81 CutoffGroupStamp* cgStamp;
82 RigidBodyStamp* rbStamp;
83 int nRigidAtoms = 0;
84
85 vector<Component*> components = simParams->getComponents();
86
87 for (vector<Component*>::iterator i = components.begin(); i !=components.end(); ++i) {
88 molStamp = (*i)->getMoleculeStamp();
89 nMolWithSameStamp = (*i)->getNMol();
90
91 addMoleculeStamp(molStamp, nMolWithSameStamp);
92
93 //calculate atoms in molecules
94 nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp;
95
96 //calculate atoms in cutoff groups
97 int nAtomsInGroups = 0;
98 int nCutoffGroupsInStamp = molStamp->getNCutoffGroups();
99
100 for (int j=0; j < nCutoffGroupsInStamp; j++) {
101 cgStamp = molStamp->getCutoffGroupStamp(j);
102 nAtomsInGroups += cgStamp->getNMembers();
103 }
104
105 nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
106
107 nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;
108
109 //calculate atoms in rigid bodies
110 int nAtomsInRigidBodies = 0;
111 int nRigidBodiesInStamp = molStamp->getNRigidBodies();
112
113 for (int j=0; j < nRigidBodiesInStamp; j++) {
114 rbStamp = molStamp->getRigidBodyStamp(j);
115 nAtomsInRigidBodies += rbStamp->getNMembers();
116 }
117
118 nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
119 nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;
120
121 }
122
123 //every free atom (atom does not belong to cutoff groups) is a cutoff
124 //group therefore the total number of cutoff groups in the system is
125 //equal to the total number of atoms minus number of atoms belong to
126 //cutoff group defined in meta-data file plus the number of cutoff
127 //groups defined in meta-data file
128
129 nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
130
131 //every free atom (atom does not belong to rigid bodies) is an
132 //integrable object therefore the total number of integrable objects
133 //in the system is equal to the total number of atoms minus number of
134 //atoms belong to rigid body defined in meta-data file plus the number
135 //of rigid bodies defined in meta-data file
136 nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms
137 + nGlobalRigidBodies_;
138
139 nGlobalMols_ = molStampIds_.size();
140 molToProcMap_.resize(nGlobalMols_);
141 }
142
143 SimInfo::~SimInfo() {
144 map<int, Molecule*>::iterator i;
145 for (i = molecules_.begin(); i != molecules_.end(); ++i) {
146 delete i->second;
147 }
148 molecules_.clear();
149
150 delete sman_;
151 delete simParams_;
152 delete forceField_;
153 }
154
155
156 bool SimInfo::addMolecule(Molecule* mol) {
157 MoleculeIterator i;
158
159 i = molecules_.find(mol->getGlobalIndex());
160 if (i == molecules_.end() ) {
161
162 molecules_.insert(make_pair(mol->getGlobalIndex(), mol));
163
164 nAtoms_ += mol->getNAtoms();
165 nBonds_ += mol->getNBonds();
166 nBends_ += mol->getNBends();
167 nTorsions_ += mol->getNTorsions();
168 nInversions_ += mol->getNInversions();
169 nRigidBodies_ += mol->getNRigidBodies();
170 nIntegrableObjects_ += mol->getNIntegrableObjects();
171 nCutoffGroups_ += mol->getNCutoffGroups();
172 nConstraints_ += mol->getNConstraintPairs();
173
174 addInteractionPairs(mol);
175
176 return true;
177 } else {
178 return false;
179 }
180 }
181
182 bool SimInfo::removeMolecule(Molecule* mol) {
183 MoleculeIterator i;
184 i = molecules_.find(mol->getGlobalIndex());
185
186 if (i != molecules_.end() ) {
187
188 assert(mol == i->second);
189
190 nAtoms_ -= mol->getNAtoms();
191 nBonds_ -= mol->getNBonds();
192 nBends_ -= mol->getNBends();
193 nTorsions_ -= mol->getNTorsions();
194 nInversions_ -= mol->getNInversions();
195 nRigidBodies_ -= mol->getNRigidBodies();
196 nIntegrableObjects_ -= mol->getNIntegrableObjects();
197 nCutoffGroups_ -= mol->getNCutoffGroups();
198 nConstraints_ -= mol->getNConstraintPairs();
199
200 removeInteractionPairs(mol);
201 molecules_.erase(mol->getGlobalIndex());
202
203 delete mol;
204
205 return true;
206 } else {
207 return false;
208 }
209 }
210
211
212 Molecule* SimInfo::beginMolecule(MoleculeIterator& i) {
213 i = molecules_.begin();
214 return i == molecules_.end() ? NULL : i->second;
215 }
216
217 Molecule* SimInfo::nextMolecule(MoleculeIterator& i) {
218 ++i;
219 return i == molecules_.end() ? NULL : i->second;
220 }
221
222
223 void SimInfo::calcNdf() {
224 int ndf_local;
225 MoleculeIterator i;
226 vector<StuntDouble*>::iterator j;
227 Molecule* mol;
228 StuntDouble* integrableObject;
229
230 ndf_local = 0;
231
232 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
233 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
234 integrableObject = mol->nextIntegrableObject(j)) {
235
236 ndf_local += 3;
237
238 if (integrableObject->isDirectional()) {
239 if (integrableObject->isLinear()) {
240 ndf_local += 2;
241 } else {
242 ndf_local += 3;
243 }
244 }
245
246 }
247 }
248
249 // n_constraints is local, so subtract them on each processor
250 ndf_local -= nConstraints_;
251
252 #ifdef IS_MPI
253 MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
254 #else
255 ndf_ = ndf_local;
256 #endif
257
258 // nZconstraints_ is global, as are the 3 COM translations for the
259 // entire system:
260 ndf_ = ndf_ - 3 - nZconstraint_;
261
262 }
263
264 int SimInfo::getFdf() {
265 #ifdef IS_MPI
266 MPI_Allreduce(&fdf_local,&fdf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
267 #else
268 fdf_ = fdf_local;
269 #endif
270 return fdf_;
271 }
272
273 unsigned int SimInfo::getNLocalCutoffGroups(){
274 int nLocalCutoffAtoms = 0;
275 Molecule* mol;
276 MoleculeIterator mi;
277 CutoffGroup* cg;
278 Molecule::CutoffGroupIterator ci;
279
280 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
281
282 for (cg = mol->beginCutoffGroup(ci); cg != NULL;
283 cg = mol->nextCutoffGroup(ci)) {
284 nLocalCutoffAtoms += cg->getNumAtom();
285
286 }
287 }
288
289 return nAtoms_ - nLocalCutoffAtoms + nCutoffGroups_;
290 }
291
292 void SimInfo::calcNdfRaw() {
293 int ndfRaw_local;
294
295 MoleculeIterator i;
296 vector<StuntDouble*>::iterator j;
297 Molecule* mol;
298 StuntDouble* integrableObject;
299
300 // Raw degrees of freedom that we have to set
301 ndfRaw_local = 0;
302
303 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
304 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
305 integrableObject = mol->nextIntegrableObject(j)) {
306
307 ndfRaw_local += 3;
308
309 if (integrableObject->isDirectional()) {
310 if (integrableObject->isLinear()) {
311 ndfRaw_local += 2;
312 } else {
313 ndfRaw_local += 3;
314 }
315 }
316
317 }
318 }
319
320 #ifdef IS_MPI
321 MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
322 #else
323 ndfRaw_ = ndfRaw_local;
324 #endif
325 }
326
327 void SimInfo::calcNdfTrans() {
328 int ndfTrans_local;
329
330 ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
331
332
333 #ifdef IS_MPI
334 MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
335 #else
336 ndfTrans_ = ndfTrans_local;
337 #endif
338
339 ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
340
341 }
342
343 void SimInfo::addInteractionPairs(Molecule* mol) {
344 ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
345 vector<Bond*>::iterator bondIter;
346 vector<Bend*>::iterator bendIter;
347 vector<Torsion*>::iterator torsionIter;
348 vector<Inversion*>::iterator inversionIter;
349 Bond* bond;
350 Bend* bend;
351 Torsion* torsion;
352 Inversion* inversion;
353 int a;
354 int b;
355 int c;
356 int d;
357
358 // atomGroups can be used to add special interaction maps between
359 // groups of atoms that are in two separate rigid bodies.
360 // However, most site-site interactions between two rigid bodies
361 // are probably not special, just the ones between the physically
362 // bonded atoms. Interactions *within* a single rigid body should
363 // always be excluded. These are done at the bottom of this
364 // function.
365
366 map<int, set<int> > atomGroups;
367 Molecule::RigidBodyIterator rbIter;
368 RigidBody* rb;
369 Molecule::IntegrableObjectIterator ii;
370 StuntDouble* integrableObject;
371
372 for (integrableObject = mol->beginIntegrableObject(ii);
373 integrableObject != NULL;
374 integrableObject = mol->nextIntegrableObject(ii)) {
375
376 if (integrableObject->isRigidBody()) {
377 rb = static_cast<RigidBody*>(integrableObject);
378 vector<Atom*> atoms = rb->getAtoms();
379 set<int> rigidAtoms;
380 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
381 rigidAtoms.insert(atoms[i]->getGlobalIndex());
382 }
383 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
384 atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
385 }
386 } else {
387 set<int> oneAtomSet;
388 oneAtomSet.insert(integrableObject->getGlobalIndex());
389 atomGroups.insert(map<int, set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));
390 }
391 }
392
393 for (bond= mol->beginBond(bondIter); bond != NULL;
394 bond = mol->nextBond(bondIter)) {
395
396 a = bond->getAtomA()->getGlobalIndex();
397 b = bond->getAtomB()->getGlobalIndex();
398
399 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
400 oneTwoInteractions_.addPair(a, b);
401 } else {
402 excludedInteractions_.addPair(a, b);
403 }
404 }
405
406 for (bend= mol->beginBend(bendIter); bend != NULL;
407 bend = mol->nextBend(bendIter)) {
408
409 a = bend->getAtomA()->getGlobalIndex();
410 b = bend->getAtomB()->getGlobalIndex();
411 c = bend->getAtomC()->getGlobalIndex();
412
413 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
414 oneTwoInteractions_.addPair(a, b);
415 oneTwoInteractions_.addPair(b, c);
416 } else {
417 excludedInteractions_.addPair(a, b);
418 excludedInteractions_.addPair(b, c);
419 }
420
421 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
422 oneThreeInteractions_.addPair(a, c);
423 } else {
424 excludedInteractions_.addPair(a, c);
425 }
426 }
427
428 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
429 torsion = mol->nextTorsion(torsionIter)) {
430
431 a = torsion->getAtomA()->getGlobalIndex();
432 b = torsion->getAtomB()->getGlobalIndex();
433 c = torsion->getAtomC()->getGlobalIndex();
434 d = torsion->getAtomD()->getGlobalIndex();
435
436 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
437 oneTwoInteractions_.addPair(a, b);
438 oneTwoInteractions_.addPair(b, c);
439 oneTwoInteractions_.addPair(c, d);
440 } else {
441 excludedInteractions_.addPair(a, b);
442 excludedInteractions_.addPair(b, c);
443 excludedInteractions_.addPair(c, d);
444 }
445
446 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
447 oneThreeInteractions_.addPair(a, c);
448 oneThreeInteractions_.addPair(b, d);
449 } else {
450 excludedInteractions_.addPair(a, c);
451 excludedInteractions_.addPair(b, d);
452 }
453
454 if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
455 oneFourInteractions_.addPair(a, d);
456 } else {
457 excludedInteractions_.addPair(a, d);
458 }
459 }
460
461 for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
462 inversion = mol->nextInversion(inversionIter)) {
463
464 a = inversion->getAtomA()->getGlobalIndex();
465 b = inversion->getAtomB()->getGlobalIndex();
466 c = inversion->getAtomC()->getGlobalIndex();
467 d = inversion->getAtomD()->getGlobalIndex();
468
469 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
470 oneTwoInteractions_.addPair(a, b);
471 oneTwoInteractions_.addPair(a, c);
472 oneTwoInteractions_.addPair(a, d);
473 } else {
474 excludedInteractions_.addPair(a, b);
475 excludedInteractions_.addPair(a, c);
476 excludedInteractions_.addPair(a, d);
477 }
478
479 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
480 oneThreeInteractions_.addPair(b, c);
481 oneThreeInteractions_.addPair(b, d);
482 oneThreeInteractions_.addPair(c, d);
483 } else {
484 excludedInteractions_.addPair(b, c);
485 excludedInteractions_.addPair(b, d);
486 excludedInteractions_.addPair(c, d);
487 }
488 }
489
490 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
491 rb = mol->nextRigidBody(rbIter)) {
492 vector<Atom*> atoms = rb->getAtoms();
493 for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
494 for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
495 a = atoms[i]->getGlobalIndex();
496 b = atoms[j]->getGlobalIndex();
497 excludedInteractions_.addPair(a, b);
498 }
499 }
500 }
501
502 }
503
504 void SimInfo::removeInteractionPairs(Molecule* mol) {
505 ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
506 vector<Bond*>::iterator bondIter;
507 vector<Bend*>::iterator bendIter;
508 vector<Torsion*>::iterator torsionIter;
509 vector<Inversion*>::iterator inversionIter;
510 Bond* bond;
511 Bend* bend;
512 Torsion* torsion;
513 Inversion* inversion;
514 int a;
515 int b;
516 int c;
517 int d;
518
519 map<int, set<int> > atomGroups;
520 Molecule::RigidBodyIterator rbIter;
521 RigidBody* rb;
522 Molecule::IntegrableObjectIterator ii;
523 StuntDouble* integrableObject;
524
525 for (integrableObject = mol->beginIntegrableObject(ii);
526 integrableObject != NULL;
527 integrableObject = mol->nextIntegrableObject(ii)) {
528
529 if (integrableObject->isRigidBody()) {
530 rb = static_cast<RigidBody*>(integrableObject);
531 vector<Atom*> atoms = rb->getAtoms();
532 set<int> rigidAtoms;
533 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
534 rigidAtoms.insert(atoms[i]->getGlobalIndex());
535 }
536 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
537 atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
538 }
539 } else {
540 set<int> oneAtomSet;
541 oneAtomSet.insert(integrableObject->getGlobalIndex());
542 atomGroups.insert(map<int, set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));
543 }
544 }
545
546 for (bond= mol->beginBond(bondIter); bond != NULL;
547 bond = mol->nextBond(bondIter)) {
548
549 a = bond->getAtomA()->getGlobalIndex();
550 b = bond->getAtomB()->getGlobalIndex();
551
552 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
553 oneTwoInteractions_.removePair(a, b);
554 } else {
555 excludedInteractions_.removePair(a, b);
556 }
557 }
558
559 for (bend= mol->beginBend(bendIter); bend != NULL;
560 bend = mol->nextBend(bendIter)) {
561
562 a = bend->getAtomA()->getGlobalIndex();
563 b = bend->getAtomB()->getGlobalIndex();
564 c = bend->getAtomC()->getGlobalIndex();
565
566 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
567 oneTwoInteractions_.removePair(a, b);
568 oneTwoInteractions_.removePair(b, c);
569 } else {
570 excludedInteractions_.removePair(a, b);
571 excludedInteractions_.removePair(b, c);
572 }
573
574 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
575 oneThreeInteractions_.removePair(a, c);
576 } else {
577 excludedInteractions_.removePair(a, c);
578 }
579 }
580
581 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
582 torsion = mol->nextTorsion(torsionIter)) {
583
584 a = torsion->getAtomA()->getGlobalIndex();
585 b = torsion->getAtomB()->getGlobalIndex();
586 c = torsion->getAtomC()->getGlobalIndex();
587 d = torsion->getAtomD()->getGlobalIndex();
588
589 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
590 oneTwoInteractions_.removePair(a, b);
591 oneTwoInteractions_.removePair(b, c);
592 oneTwoInteractions_.removePair(c, d);
593 } else {
594 excludedInteractions_.removePair(a, b);
595 excludedInteractions_.removePair(b, c);
596 excludedInteractions_.removePair(c, d);
597 }
598
599 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
600 oneThreeInteractions_.removePair(a, c);
601 oneThreeInteractions_.removePair(b, d);
602 } else {
603 excludedInteractions_.removePair(a, c);
604 excludedInteractions_.removePair(b, d);
605 }
606
607 if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
608 oneFourInteractions_.removePair(a, d);
609 } else {
610 excludedInteractions_.removePair(a, d);
611 }
612 }
613
614 for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
615 inversion = mol->nextInversion(inversionIter)) {
616
617 a = inversion->getAtomA()->getGlobalIndex();
618 b = inversion->getAtomB()->getGlobalIndex();
619 c = inversion->getAtomC()->getGlobalIndex();
620 d = inversion->getAtomD()->getGlobalIndex();
621
622 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
623 oneTwoInteractions_.removePair(a, b);
624 oneTwoInteractions_.removePair(a, c);
625 oneTwoInteractions_.removePair(a, d);
626 } else {
627 excludedInteractions_.removePair(a, b);
628 excludedInteractions_.removePair(a, c);
629 excludedInteractions_.removePair(a, d);
630 }
631
632 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
633 oneThreeInteractions_.removePair(b, c);
634 oneThreeInteractions_.removePair(b, d);
635 oneThreeInteractions_.removePair(c, d);
636 } else {
637 excludedInteractions_.removePair(b, c);
638 excludedInteractions_.removePair(b, d);
639 excludedInteractions_.removePair(c, d);
640 }
641 }
642
643 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
644 rb = mol->nextRigidBody(rbIter)) {
645 vector<Atom*> atoms = rb->getAtoms();
646 for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
647 for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
648 a = atoms[i]->getGlobalIndex();
649 b = atoms[j]->getGlobalIndex();
650 excludedInteractions_.removePair(a, b);
651 }
652 }
653 }
654
655 }
656
657
658 void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) {
659 int curStampId;
660
661 //index from 0
662 curStampId = moleculeStamps_.size();
663
664 moleculeStamps_.push_back(molStamp);
665 molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
666 }
667
668
669 /**
670 * update
671 *
672 * Performs the global checks and variable settings after the
673 * objects have been created.
674 *
675 */
676 void SimInfo::update() {
677 setupSimVariables();
678 calcNdf();
679 calcNdfRaw();
680 calcNdfTrans();
681 }
682
683 /**
684 * getSimulatedAtomTypes
685 *
686 * Returns an STL set of AtomType* that are actually present in this
687 * simulation. Must query all processors to assemble this information.
688 *
689 */
690 set<AtomType*> SimInfo::getSimulatedAtomTypes() {
691 SimInfo::MoleculeIterator mi;
692 Molecule* mol;
693 Molecule::AtomIterator ai;
694 Atom* atom;
695 set<AtomType*> atomTypes;
696
697 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
698 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
699 atomTypes.insert(atom->getAtomType());
700 }
701 }
702
703 #ifdef IS_MPI
704
705 // loop over the found atom types on this processor, and add their
706 // numerical idents to a vector:
707
708 vector<int> foundTypes;
709 set<AtomType*>::iterator i;
710 for (i = atomTypes.begin(); i != atomTypes.end(); ++i)
711 foundTypes.push_back( (*i)->getIdent() );
712
713 // count_local holds the number of found types on this processor
714 int count_local = foundTypes.size();
715
716 // count holds the total number of found types on all processors
717 // (some will be redundant with the ones found locally):
718 int count;
719 MPI::COMM_WORLD.Allreduce(&count_local, &count, 1, MPI::INT, MPI::SUM);
720
721 // create a vector to hold the globally found types, and resize it:
722 vector<int> ftGlobal;
723 ftGlobal.resize(count);
724 vector<int> counts;
725
726 int nproc = MPI::COMM_WORLD.Get_size();
727 counts.resize(nproc);
728 vector<int> disps;
729 disps.resize(nproc);
730
731 // now spray out the foundTypes to all the other processors:
732
733 MPI::COMM_WORLD.Allgatherv(&foundTypes[0], count_local, MPI::INT,
734 &ftGlobal[0], &counts[0], &disps[0], MPI::INT);
735
736 // foundIdents is a stl set, so inserting an already found ident
737 // will have no effect.
738 set<int> foundIdents;
739 vector<int>::iterator j;
740 for (j = ftGlobal.begin(); j != ftGlobal.end(); ++j)
741 foundIdents.insert((*j));
742
743 // now iterate over the foundIdents and get the actual atom types
744 // that correspond to these:
745 set<int>::iterator it;
746 for (it = foundIdents.begin(); it != foundIdents.end(); ++it)
747 atomTypes.insert( forceField_->getAtomType((*it)) );
748
749 #endif
750
751 return atomTypes;
752 }
753
754 void SimInfo::setupSimVariables() {
755 useAtomicVirial_ = simParams_->getUseAtomicVirial();
756 // we only call setAccumulateBoxDipole if the accumulateBoxDipole parameter is true
757 calcBoxDipole_ = false;
758 if ( simParams_->haveAccumulateBoxDipole() )
759 if ( simParams_->getAccumulateBoxDipole() ) {
760 calcBoxDipole_ = true;
761 }
762
763 set<AtomType*>::iterator i;
764 set<AtomType*> atomTypes;
765 atomTypes = getSimulatedAtomTypes();
766 int usesElectrostatic = 0;
767 int usesMetallic = 0;
768 int usesDirectional = 0;
769 //loop over all of the atom types
770 for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
771 usesElectrostatic |= (*i)->isElectrostatic();
772 usesMetallic |= (*i)->isMetal();
773 usesDirectional |= (*i)->isDirectional();
774 }
775
776 #ifdef IS_MPI
777 int temp;
778 temp = usesDirectional;
779 MPI_Allreduce(&temp, &usesDirectionalAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
780
781 temp = usesMetallic;
782 MPI_Allreduce(&temp, &usesMetallicAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
783
784 temp = usesElectrostatic;
785 MPI_Allreduce(&temp, &usesElectrostaticAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
786 #endif
787 }
788
789
790 vector<int> SimInfo::getGlobalAtomIndices() {
791 SimInfo::MoleculeIterator mi;
792 Molecule* mol;
793 Molecule::AtomIterator ai;
794 Atom* atom;
795
796 vector<int> GlobalAtomIndices(getNAtoms(), 0);
797
798 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
799
800 for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
801 GlobalAtomIndices[atom->getLocalIndex()] = atom->getGlobalIndex();
802 }
803 }
804 return GlobalAtomIndices;
805 }
806
807
808 vector<int> SimInfo::getGlobalGroupIndices() {
809 SimInfo::MoleculeIterator mi;
810 Molecule* mol;
811 Molecule::CutoffGroupIterator ci;
812 CutoffGroup* cg;
813
814 vector<int> GlobalGroupIndices;
815
816 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
817
818 //local index of cutoff group is trivial, it only depends on the
819 //order of travesing
820 for (cg = mol->beginCutoffGroup(ci); cg != NULL;
821 cg = mol->nextCutoffGroup(ci)) {
822 GlobalGroupIndices.push_back(cg->getGlobalIndex());
823 }
824 }
825 return GlobalGroupIndices;
826 }
827
828
829 void SimInfo::prepareTopology() {
830 int nExclude, nOneTwo, nOneThree, nOneFour;
831
832 //calculate mass ratio of cutoff group
833 SimInfo::MoleculeIterator mi;
834 Molecule* mol;
835 Molecule::CutoffGroupIterator ci;
836 CutoffGroup* cg;
837 Molecule::AtomIterator ai;
838 Atom* atom;
839 RealType totalMass;
840
841 /**
842 * The mass factor is the relative mass of an atom to the total
843 * mass of the cutoff group it belongs to. By default, all atoms
844 * are their own cutoff groups, and therefore have mass factors of
845 * 1. We need some special handling for massless atoms, which
846 * will be treated as carrying the entire mass of the cutoff
847 * group.
848 */
849 massFactors_.clear();
850 massFactors_.resize(getNAtoms(), 1.0);
851
852 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
853 for (cg = mol->beginCutoffGroup(ci); cg != NULL;
854 cg = mol->nextCutoffGroup(ci)) {
855
856 totalMass = cg->getMass();
857 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
858 // Check for massless groups - set mfact to 1 if true
859 if (totalMass != 0)
860 massFactors_[atom->getLocalIndex()] = atom->getMass()/totalMass;
861 else
862 massFactors_[atom->getLocalIndex()] = 1.0;
863 }
864 }
865 }
866
867 // Build the identArray_
868
869 identArray_.clear();
870 identArray_.reserve(getNAtoms());
871 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
872 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
873 identArray_.push_back(atom->getIdent());
874 }
875 }
876
877 //scan topology
878
879 nExclude = excludedInteractions_.getSize();
880 nOneTwo = oneTwoInteractions_.getSize();
881 nOneThree = oneThreeInteractions_.getSize();
882 nOneFour = oneFourInteractions_.getSize();
883
884 int* excludeList = excludedInteractions_.getPairList();
885 int* oneTwoList = oneTwoInteractions_.getPairList();
886 int* oneThreeList = oneThreeInteractions_.getPairList();
887 int* oneFourList = oneFourInteractions_.getPairList();
888
889 topologyDone_ = true;
890 }
891
892 void SimInfo::addProperty(GenericData* genData) {
893 properties_.addProperty(genData);
894 }
895
896 void SimInfo::removeProperty(const string& propName) {
897 properties_.removeProperty(propName);
898 }
899
900 void SimInfo::clearProperties() {
901 properties_.clearProperties();
902 }
903
904 vector<string> SimInfo::getPropertyNames() {
905 return properties_.getPropertyNames();
906 }
907
908 vector<GenericData*> SimInfo::getProperties() {
909 return properties_.getProperties();
910 }
911
912 GenericData* SimInfo::getPropertyByName(const string& propName) {
913 return properties_.getPropertyByName(propName);
914 }
915
916 void SimInfo::setSnapshotManager(SnapshotManager* sman) {
917 if (sman_ == sman) {
918 return;
919 }
920 delete sman_;
921 sman_ = sman;
922
923 Molecule* mol;
924 RigidBody* rb;
925 Atom* atom;
926 CutoffGroup* cg;
927 SimInfo::MoleculeIterator mi;
928 Molecule::RigidBodyIterator rbIter;
929 Molecule::AtomIterator atomIter;
930 Molecule::CutoffGroupIterator cgIter;
931
932 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
933
934 for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) {
935 atom->setSnapshotManager(sman_);
936 }
937
938 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
939 rb->setSnapshotManager(sman_);
940 }
941
942 for (cg = mol->beginCutoffGroup(cgIter); cg != NULL; cg = mol->nextCutoffGroup(cgIter)) {
943 cg->setSnapshotManager(sman_);
944 }
945 }
946
947 }
948
949 Vector3d SimInfo::getComVel(){
950 SimInfo::MoleculeIterator i;
951 Molecule* mol;
952
953 Vector3d comVel(0.0);
954 RealType totalMass = 0.0;
955
956
957 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
958 RealType mass = mol->getMass();
959 totalMass += mass;
960 comVel += mass * mol->getComVel();
961 }
962
963 #ifdef IS_MPI
964 RealType tmpMass = totalMass;
965 Vector3d tmpComVel(comVel);
966 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
967 MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
968 #endif
969
970 comVel /= totalMass;
971
972 return comVel;
973 }
974
975 Vector3d SimInfo::getCom(){
976 SimInfo::MoleculeIterator i;
977 Molecule* mol;
978
979 Vector3d com(0.0);
980 RealType totalMass = 0.0;
981
982 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
983 RealType mass = mol->getMass();
984 totalMass += mass;
985 com += mass * mol->getCom();
986 }
987
988 #ifdef IS_MPI
989 RealType tmpMass = totalMass;
990 Vector3d tmpCom(com);
991 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
992 MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
993 #endif
994
995 com /= totalMass;
996
997 return com;
998
999 }
1000
1001 ostream& operator <<(ostream& o, SimInfo& info) {
1002
1003 return o;
1004 }
1005
1006
1007 /*
1008 Returns center of mass and center of mass velocity in one function call.
1009 */
1010
1011 void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){
1012 SimInfo::MoleculeIterator i;
1013 Molecule* mol;
1014
1015
1016 RealType totalMass = 0.0;
1017
1018
1019 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1020 RealType mass = mol->getMass();
1021 totalMass += mass;
1022 com += mass * mol->getCom();
1023 comVel += mass * mol->getComVel();
1024 }
1025
1026 #ifdef IS_MPI
1027 RealType tmpMass = totalMass;
1028 Vector3d tmpCom(com);
1029 Vector3d tmpComVel(comVel);
1030 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1031 MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1032 MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1033 #endif
1034
1035 com /= totalMass;
1036 comVel /= totalMass;
1037 }
1038
1039 /*
1040 Return intertia tensor for entire system and angular momentum Vector.
1041
1042
1043 [ Ixx -Ixy -Ixz ]
1044 J =| -Iyx Iyy -Iyz |
1045 [ -Izx -Iyz Izz ]
1046 */
1047
1048 void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){
1049
1050
1051 RealType xx = 0.0;
1052 RealType yy = 0.0;
1053 RealType zz = 0.0;
1054 RealType xy = 0.0;
1055 RealType xz = 0.0;
1056 RealType yz = 0.0;
1057 Vector3d com(0.0);
1058 Vector3d comVel(0.0);
1059
1060 getComAll(com, comVel);
1061
1062 SimInfo::MoleculeIterator i;
1063 Molecule* mol;
1064
1065 Vector3d thisq(0.0);
1066 Vector3d thisv(0.0);
1067
1068 RealType thisMass = 0.0;
1069
1070
1071
1072
1073 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1074
1075 thisq = mol->getCom()-com;
1076 thisv = mol->getComVel()-comVel;
1077 thisMass = mol->getMass();
1078 // Compute moment of intertia coefficients.
1079 xx += thisq[0]*thisq[0]*thisMass;
1080 yy += thisq[1]*thisq[1]*thisMass;
1081 zz += thisq[2]*thisq[2]*thisMass;
1082
1083 // compute products of intertia
1084 xy += thisq[0]*thisq[1]*thisMass;
1085 xz += thisq[0]*thisq[2]*thisMass;
1086 yz += thisq[1]*thisq[2]*thisMass;
1087
1088 angularMomentum += cross( thisq, thisv ) * thisMass;
1089
1090 }
1091
1092
1093 inertiaTensor(0,0) = yy + zz;
1094 inertiaTensor(0,1) = -xy;
1095 inertiaTensor(0,2) = -xz;
1096 inertiaTensor(1,0) = -xy;
1097 inertiaTensor(1,1) = xx + zz;
1098 inertiaTensor(1,2) = -yz;
1099 inertiaTensor(2,0) = -xz;
1100 inertiaTensor(2,1) = -yz;
1101 inertiaTensor(2,2) = xx + yy;
1102
1103 #ifdef IS_MPI
1104 Mat3x3d tmpI(inertiaTensor);
1105 Vector3d tmpAngMom;
1106 MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1107 MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1108 #endif
1109
1110 return;
1111 }
1112
1113 //Returns the angular momentum of the system
1114 Vector3d SimInfo::getAngularMomentum(){
1115
1116 Vector3d com(0.0);
1117 Vector3d comVel(0.0);
1118 Vector3d angularMomentum(0.0);
1119
1120 getComAll(com,comVel);
1121
1122 SimInfo::MoleculeIterator i;
1123 Molecule* mol;
1124
1125 Vector3d thisr(0.0);
1126 Vector3d thisp(0.0);
1127
1128 RealType thisMass;
1129
1130 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1131 thisMass = mol->getMass();
1132 thisr = mol->getCom()-com;
1133 thisp = (mol->getComVel()-comVel)*thisMass;
1134
1135 angularMomentum += cross( thisr, thisp );
1136
1137 }
1138
1139 #ifdef IS_MPI
1140 Vector3d tmpAngMom;
1141 MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1142 #endif
1143
1144 return angularMomentum;
1145 }
1146
1147 StuntDouble* SimInfo::getIOIndexToIntegrableObject(int index) {
1148 return IOIndexToIntegrableObject.at(index);
1149 }
1150
1151 void SimInfo::setIOIndexToIntegrableObject(const vector<StuntDouble*>& v) {
1152 IOIndexToIntegrableObject= v;
1153 }
1154
1155 /* Returns the Volume of the simulation based on a ellipsoid with semi-axes
1156 based on the radius of gyration V=4/3*Pi*R_1*R_2*R_3
1157 where R_i are related to the principle inertia moments R_i = sqrt(C*I_i/N), this reduces to
1158 V = 4/3*Pi*(C/N)^3/2*sqrt(det(I)). See S.E. Baltazar et. al. Comp. Mat. Sci. 37 (2006) 526-536.
1159 */
1160 void SimInfo::getGyrationalVolume(RealType &volume){
1161 Mat3x3d intTensor;
1162 RealType det;
1163 Vector3d dummyAngMom;
1164 RealType sysconstants;
1165 RealType geomCnst;
1166
1167 geomCnst = 3.0/2.0;
1168 /* Get the inertial tensor and angular momentum for free*/
1169 getInertiaTensor(intTensor,dummyAngMom);
1170
1171 det = intTensor.determinant();
1172 sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1173 volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(det);
1174 return;
1175 }
1176
1177 void SimInfo::getGyrationalVolume(RealType &volume, RealType &detI){
1178 Mat3x3d intTensor;
1179 Vector3d dummyAngMom;
1180 RealType sysconstants;
1181 RealType geomCnst;
1182
1183 geomCnst = 3.0/2.0;
1184 /* Get the inertial tensor and angular momentum for free*/
1185 getInertiaTensor(intTensor,dummyAngMom);
1186
1187 detI = intTensor.determinant();
1188 sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1189 volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(detI);
1190 return;
1191 }
1192 /*
1193 void SimInfo::setStuntDoubleFromGlobalIndex(vector<StuntDouble*> v) {
1194 assert( v.size() == nAtoms_ + nRigidBodies_);
1195 sdByGlobalIndex_ = v;
1196 }
1197
1198 StuntDouble* SimInfo::getStuntDoubleFromGlobalIndex(int index) {
1199 //assert(index < nAtoms_ + nRigidBodies_);
1200 return sdByGlobalIndex_.at(index);
1201 }
1202 */
1203 int SimInfo::getNGlobalConstraints() {
1204 int nGlobalConstraints;
1205 #ifdef IS_MPI
1206 MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM,
1207 MPI_COMM_WORLD);
1208 #else
1209 nGlobalConstraints = nConstraints_;
1210 #endif
1211 return nGlobalConstraints;
1212 }
1213
1214 }//end namespace OpenMD
1215

Properties

Name Value
svn:keywords Author Id Revision Date