ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimInfo.cpp
Revision: 1969
Committed: Wed Feb 26 14:14:50 2014 UTC (11 years, 2 months ago) by gezelter
File size: 35811 byte(s)
Log Message:
Fixes to deal with deprecation of MPI C++ bindings.  We've reverted back to the
C calls.

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, 234107 (2008).
39 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40 * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 */
42
43 /**
44 * @file SimInfo.cpp
45 * @author tlin
46 * @date 11/02/2004
47 * @version 1.0
48 */
49
50 #ifdef IS_MPI
51 #include <mpi.h>
52 #endif
53 #include <algorithm>
54 #include <set>
55 #include <map>
56
57 #include "brains/SimInfo.hpp"
58 #include "math/Vector3.hpp"
59 #include "primitives/Molecule.hpp"
60 #include "primitives/StuntDouble.hpp"
61 #include "utils/MemoryUtils.hpp"
62 #include "utils/simError.h"
63 #include "selection/SelectionManager.hpp"
64 #include "io/ForceFieldOptions.hpp"
65 #include "brains/ForceField.hpp"
66 #include "nonbonded/SwitchingFunction.hpp"
67
68 using namespace std;
69 namespace OpenMD {
70
71 SimInfo::SimInfo(ForceField* ff, Globals* simParams) :
72 forceField_(ff), simParams_(simParams),
73 ndf_(0), fdf_local(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0),
74 nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0),
75 nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
76 nGlobalFluctuatingCharges_(0), nGlobalBonds_(0), nGlobalBends_(0),
77 nGlobalTorsions_(0), nGlobalInversions_(0), nAtoms_(0), nBonds_(0),
78 nBends_(0), nTorsions_(0), nInversions_(0), nRigidBodies_(0),
79 nIntegrableObjects_(0), nCutoffGroups_(0), nConstraints_(0),
80 nFluctuatingCharges_(0), sman_(NULL), topologyDone_(false),
81 calcBoxDipole_(false), useAtomicVirial_(true) {
82
83 MoleculeStamp* molStamp;
84 int nMolWithSameStamp;
85 int nCutoffAtoms = 0; // number of atoms belong to cutoff groups
86 int nGroups = 0; //total cutoff groups defined in meta-data file
87 CutoffGroupStamp* cgStamp;
88 RigidBodyStamp* rbStamp;
89 int nRigidAtoms = 0;
90
91 vector<Component*> components = simParams->getComponents();
92
93 for (vector<Component*>::iterator i = components.begin();
94 i !=components.end(); ++i) {
95 molStamp = (*i)->getMoleculeStamp();
96 if ( (*i)->haveRegion() ) {
97 molStamp->setRegion( (*i)->getRegion() );
98 } else {
99 // set the region to a disallowed value:
100 molStamp->setRegion( -1 );
101 }
102
103 nMolWithSameStamp = (*i)->getNMol();
104
105 addMoleculeStamp(molStamp, nMolWithSameStamp);
106
107 //calculate atoms in molecules
108 nGlobalAtoms_ += molStamp->getNAtoms() * nMolWithSameStamp;
109 nGlobalBonds_ += molStamp->getNBonds() * nMolWithSameStamp;
110 nGlobalBends_ += molStamp->getNBends() * nMolWithSameStamp;
111 nGlobalTorsions_ += molStamp->getNTorsions() * nMolWithSameStamp;
112 nGlobalInversions_ += molStamp->getNInversions() * nMolWithSameStamp;
113
114 //calculate atoms in cutoff groups
115 int nAtomsInGroups = 0;
116 int nCutoffGroupsInStamp = molStamp->getNCutoffGroups();
117
118 for (int j=0; j < nCutoffGroupsInStamp; j++) {
119 cgStamp = molStamp->getCutoffGroupStamp(j);
120 nAtomsInGroups += cgStamp->getNMembers();
121 }
122
123 nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
124
125 nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;
126
127 //calculate atoms in rigid bodies
128 int nAtomsInRigidBodies = 0;
129 int nRigidBodiesInStamp = molStamp->getNRigidBodies();
130
131 for (int j=0; j < nRigidBodiesInStamp; j++) {
132 rbStamp = molStamp->getRigidBodyStamp(j);
133 nAtomsInRigidBodies += rbStamp->getNMembers();
134 }
135
136 nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
137 nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;
138
139 }
140
141 //every free atom (atom does not belong to cutoff groups) is a cutoff
142 //group therefore the total number of cutoff groups in the system is
143 //equal to the total number of atoms minus number of atoms belong to
144 //cutoff group defined in meta-data file plus the number of cutoff
145 //groups defined in meta-data file
146
147 nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
148
149 //every free atom (atom does not belong to rigid bodies) is an
150 //integrable object therefore the total number of integrable objects
151 //in the system is equal to the total number of atoms minus number of
152 //atoms belong to rigid body defined in meta-data file plus the number
153 //of rigid bodies defined in meta-data file
154 nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms
155 + nGlobalRigidBodies_;
156
157 nGlobalMols_ = molStampIds_.size();
158 molToProcMap_.resize(nGlobalMols_);
159 }
160
161 SimInfo::~SimInfo() {
162 map<int, Molecule*>::iterator i;
163 for (i = molecules_.begin(); i != molecules_.end(); ++i) {
164 delete i->second;
165 }
166 molecules_.clear();
167
168 delete sman_;
169 delete simParams_;
170 delete forceField_;
171 }
172
173
174 bool SimInfo::addMolecule(Molecule* mol) {
175 MoleculeIterator i;
176
177 i = molecules_.find(mol->getGlobalIndex());
178 if (i == molecules_.end() ) {
179
180 molecules_.insert(make_pair(mol->getGlobalIndex(), mol));
181
182 nAtoms_ += mol->getNAtoms();
183 nBonds_ += mol->getNBonds();
184 nBends_ += mol->getNBends();
185 nTorsions_ += mol->getNTorsions();
186 nInversions_ += mol->getNInversions();
187 nRigidBodies_ += mol->getNRigidBodies();
188 nIntegrableObjects_ += mol->getNIntegrableObjects();
189 nCutoffGroups_ += mol->getNCutoffGroups();
190 nConstraints_ += mol->getNConstraintPairs();
191
192 addInteractionPairs(mol);
193
194 return true;
195 } else {
196 return false;
197 }
198 }
199
200 bool SimInfo::removeMolecule(Molecule* mol) {
201 MoleculeIterator i;
202 i = molecules_.find(mol->getGlobalIndex());
203
204 if (i != molecules_.end() ) {
205
206 assert(mol == i->second);
207
208 nAtoms_ -= mol->getNAtoms();
209 nBonds_ -= mol->getNBonds();
210 nBends_ -= mol->getNBends();
211 nTorsions_ -= mol->getNTorsions();
212 nInversions_ -= mol->getNInversions();
213 nRigidBodies_ -= mol->getNRigidBodies();
214 nIntegrableObjects_ -= mol->getNIntegrableObjects();
215 nCutoffGroups_ -= mol->getNCutoffGroups();
216 nConstraints_ -= mol->getNConstraintPairs();
217
218 removeInteractionPairs(mol);
219 molecules_.erase(mol->getGlobalIndex());
220
221 delete mol;
222
223 return true;
224 } else {
225 return false;
226 }
227 }
228
229
230 Molecule* SimInfo::beginMolecule(MoleculeIterator& i) {
231 i = molecules_.begin();
232 return i == molecules_.end() ? NULL : i->second;
233 }
234
235 Molecule* SimInfo::nextMolecule(MoleculeIterator& i) {
236 ++i;
237 return i == molecules_.end() ? NULL : i->second;
238 }
239
240
241 void SimInfo::calcNdf() {
242 int ndf_local, nfq_local;
243 MoleculeIterator i;
244 vector<StuntDouble*>::iterator j;
245 vector<Atom*>::iterator k;
246
247 Molecule* mol;
248 StuntDouble* sd;
249 Atom* atom;
250
251 ndf_local = 0;
252 nfq_local = 0;
253
254 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
255
256 for (sd = mol->beginIntegrableObject(j); sd != NULL;
257 sd = mol->nextIntegrableObject(j)) {
258
259 ndf_local += 3;
260
261 if (sd->isDirectional()) {
262 if (sd->isLinear()) {
263 ndf_local += 2;
264 } else {
265 ndf_local += 3;
266 }
267 }
268 }
269
270 for (atom = mol->beginFluctuatingCharge(k); atom != NULL;
271 atom = mol->nextFluctuatingCharge(k)) {
272 if (atom->isFluctuatingCharge()) {
273 nfq_local++;
274 }
275 }
276 }
277
278 ndfLocal_ = ndf_local;
279
280 // n_constraints is local, so subtract them on each processor
281 ndf_local -= nConstraints_;
282
283 #ifdef IS_MPI
284 MPI_Allreduce(&ndf_local, &ndf_, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
285 MPI_Allreduce(&nfq_local, &nGlobalFluctuatingCharges_, 1,
286 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
287 // MPI::COMM_WORLD.Allreduce(&ndf_local, &ndf_, 1, MPI::INT,MPI::SUM);
288 // MPI::COMM_WORLD.Allreduce(&nfq_local, &nGlobalFluctuatingCharges_, 1,
289 // MPI::INT, MPI::SUM);
290 #else
291 ndf_ = ndf_local;
292 nGlobalFluctuatingCharges_ = nfq_local;
293 #endif
294
295 // nZconstraints_ is global, as are the 3 COM translations for the
296 // entire system:
297 ndf_ = ndf_ - 3 - nZconstraint_;
298
299 }
300
301 int SimInfo::getFdf() {
302 #ifdef IS_MPI
303 MPI_Allreduce(&fdf_local, &fdf_, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
304 // MPI::COMM_WORLD.Allreduce(&fdf_local, &fdf_, 1, MPI::INT, MPI::SUM);
305 #else
306 fdf_ = fdf_local;
307 #endif
308 return fdf_;
309 }
310
311 unsigned int SimInfo::getNLocalCutoffGroups(){
312 int nLocalCutoffAtoms = 0;
313 Molecule* mol;
314 MoleculeIterator mi;
315 CutoffGroup* cg;
316 Molecule::CutoffGroupIterator ci;
317
318 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
319
320 for (cg = mol->beginCutoffGroup(ci); cg != NULL;
321 cg = mol->nextCutoffGroup(ci)) {
322 nLocalCutoffAtoms += cg->getNumAtom();
323
324 }
325 }
326
327 return nAtoms_ - nLocalCutoffAtoms + nCutoffGroups_;
328 }
329
330 void SimInfo::calcNdfRaw() {
331 int ndfRaw_local;
332
333 MoleculeIterator i;
334 vector<StuntDouble*>::iterator j;
335 Molecule* mol;
336 StuntDouble* sd;
337
338 // Raw degrees of freedom that we have to set
339 ndfRaw_local = 0;
340
341 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
342
343 for (sd = mol->beginIntegrableObject(j); sd != NULL;
344 sd = mol->nextIntegrableObject(j)) {
345
346 ndfRaw_local += 3;
347
348 if (sd->isDirectional()) {
349 if (sd->isLinear()) {
350 ndfRaw_local += 2;
351 } else {
352 ndfRaw_local += 3;
353 }
354 }
355
356 }
357 }
358
359 #ifdef IS_MPI
360 MPI_Allreduce(&ndfRaw_local, &ndfRaw_, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
361 // MPI::COMM_WORLD.Allreduce(&ndfRaw_local, &ndfRaw_, 1, MPI::INT, MPI::SUM);
362 #else
363 ndfRaw_ = ndfRaw_local;
364 #endif
365 }
366
367 void SimInfo::calcNdfTrans() {
368 int ndfTrans_local;
369
370 ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
371
372
373 #ifdef IS_MPI
374 MPI_Allreduce(&ndfTrans_local, &ndfTrans_, 1,
375 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
376 // MPI::COMM_WORLD.Allreduce(&ndfTrans_local, &ndfTrans_, 1,
377 // MPI::INT, MPI::SUM);
378 #else
379 ndfTrans_ = ndfTrans_local;
380 #endif
381
382 ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
383
384 }
385
386 void SimInfo::addInteractionPairs(Molecule* mol) {
387 ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
388 vector<Bond*>::iterator bondIter;
389 vector<Bend*>::iterator bendIter;
390 vector<Torsion*>::iterator torsionIter;
391 vector<Inversion*>::iterator inversionIter;
392 Bond* bond;
393 Bend* bend;
394 Torsion* torsion;
395 Inversion* inversion;
396 int a;
397 int b;
398 int c;
399 int d;
400
401 // atomGroups can be used to add special interaction maps between
402 // groups of atoms that are in two separate rigid bodies.
403 // However, most site-site interactions between two rigid bodies
404 // are probably not special, just the ones between the physically
405 // bonded atoms. Interactions *within* a single rigid body should
406 // always be excluded. These are done at the bottom of this
407 // function.
408
409 map<int, set<int> > atomGroups;
410 Molecule::RigidBodyIterator rbIter;
411 RigidBody* rb;
412 Molecule::IntegrableObjectIterator ii;
413 StuntDouble* sd;
414
415 for (sd = mol->beginIntegrableObject(ii); sd != NULL;
416 sd = mol->nextIntegrableObject(ii)) {
417
418 if (sd->isRigidBody()) {
419 rb = static_cast<RigidBody*>(sd);
420 vector<Atom*> atoms = rb->getAtoms();
421 set<int> rigidAtoms;
422 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
423 rigidAtoms.insert(atoms[i]->getGlobalIndex());
424 }
425 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
426 atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
427 }
428 } else {
429 set<int> oneAtomSet;
430 oneAtomSet.insert(sd->getGlobalIndex());
431 atomGroups.insert(map<int, set<int> >::value_type(sd->getGlobalIndex(), oneAtomSet));
432 }
433 }
434
435
436 for (bond= mol->beginBond(bondIter); bond != NULL;
437 bond = mol->nextBond(bondIter)) {
438
439 a = bond->getAtomA()->getGlobalIndex();
440 b = bond->getAtomB()->getGlobalIndex();
441
442 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
443 oneTwoInteractions_.addPair(a, b);
444 } else {
445 excludedInteractions_.addPair(a, b);
446 }
447 }
448
449 for (bend= mol->beginBend(bendIter); bend != NULL;
450 bend = mol->nextBend(bendIter)) {
451
452 a = bend->getAtomA()->getGlobalIndex();
453 b = bend->getAtomB()->getGlobalIndex();
454 c = bend->getAtomC()->getGlobalIndex();
455
456 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
457 oneTwoInteractions_.addPair(a, b);
458 oneTwoInteractions_.addPair(b, c);
459 } else {
460 excludedInteractions_.addPair(a, b);
461 excludedInteractions_.addPair(b, c);
462 }
463
464 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
465 oneThreeInteractions_.addPair(a, c);
466 } else {
467 excludedInteractions_.addPair(a, c);
468 }
469 }
470
471 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
472 torsion = mol->nextTorsion(torsionIter)) {
473
474 a = torsion->getAtomA()->getGlobalIndex();
475 b = torsion->getAtomB()->getGlobalIndex();
476 c = torsion->getAtomC()->getGlobalIndex();
477 d = torsion->getAtomD()->getGlobalIndex();
478
479 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
480 oneTwoInteractions_.addPair(a, b);
481 oneTwoInteractions_.addPair(b, c);
482 oneTwoInteractions_.addPair(c, d);
483 } else {
484 excludedInteractions_.addPair(a, b);
485 excludedInteractions_.addPair(b, c);
486 excludedInteractions_.addPair(c, d);
487 }
488
489 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
490 oneThreeInteractions_.addPair(a, c);
491 oneThreeInteractions_.addPair(b, d);
492 } else {
493 excludedInteractions_.addPair(a, c);
494 excludedInteractions_.addPair(b, d);
495 }
496
497 if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
498 oneFourInteractions_.addPair(a, d);
499 } else {
500 excludedInteractions_.addPair(a, d);
501 }
502 }
503
504 for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
505 inversion = mol->nextInversion(inversionIter)) {
506
507 a = inversion->getAtomA()->getGlobalIndex();
508 b = inversion->getAtomB()->getGlobalIndex();
509 c = inversion->getAtomC()->getGlobalIndex();
510 d = inversion->getAtomD()->getGlobalIndex();
511
512 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
513 oneTwoInteractions_.addPair(a, b);
514 oneTwoInteractions_.addPair(a, c);
515 oneTwoInteractions_.addPair(a, d);
516 } else {
517 excludedInteractions_.addPair(a, b);
518 excludedInteractions_.addPair(a, c);
519 excludedInteractions_.addPair(a, d);
520 }
521
522 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
523 oneThreeInteractions_.addPair(b, c);
524 oneThreeInteractions_.addPair(b, d);
525 oneThreeInteractions_.addPair(c, d);
526 } else {
527 excludedInteractions_.addPair(b, c);
528 excludedInteractions_.addPair(b, d);
529 excludedInteractions_.addPair(c, d);
530 }
531 }
532
533 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
534 rb = mol->nextRigidBody(rbIter)) {
535 vector<Atom*> atoms = rb->getAtoms();
536 for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
537 for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
538 a = atoms[i]->getGlobalIndex();
539 b = atoms[j]->getGlobalIndex();
540 excludedInteractions_.addPair(a, b);
541 }
542 }
543 }
544
545 }
546
547 void SimInfo::removeInteractionPairs(Molecule* mol) {
548 ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
549 vector<Bond*>::iterator bondIter;
550 vector<Bend*>::iterator bendIter;
551 vector<Torsion*>::iterator torsionIter;
552 vector<Inversion*>::iterator inversionIter;
553 Bond* bond;
554 Bend* bend;
555 Torsion* torsion;
556 Inversion* inversion;
557 int a;
558 int b;
559 int c;
560 int d;
561
562 map<int, set<int> > atomGroups;
563 Molecule::RigidBodyIterator rbIter;
564 RigidBody* rb;
565 Molecule::IntegrableObjectIterator ii;
566 StuntDouble* sd;
567
568 for (sd = mol->beginIntegrableObject(ii); sd != NULL;
569 sd = mol->nextIntegrableObject(ii)) {
570
571 if (sd->isRigidBody()) {
572 rb = static_cast<RigidBody*>(sd);
573 vector<Atom*> atoms = rb->getAtoms();
574 set<int> rigidAtoms;
575 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
576 rigidAtoms.insert(atoms[i]->getGlobalIndex());
577 }
578 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
579 atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
580 }
581 } else {
582 set<int> oneAtomSet;
583 oneAtomSet.insert(sd->getGlobalIndex());
584 atomGroups.insert(map<int, set<int> >::value_type(sd->getGlobalIndex(), oneAtomSet));
585 }
586 }
587
588 for (bond= mol->beginBond(bondIter); bond != NULL;
589 bond = mol->nextBond(bondIter)) {
590
591 a = bond->getAtomA()->getGlobalIndex();
592 b = bond->getAtomB()->getGlobalIndex();
593
594 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
595 oneTwoInteractions_.removePair(a, b);
596 } else {
597 excludedInteractions_.removePair(a, b);
598 }
599 }
600
601 for (bend= mol->beginBend(bendIter); bend != NULL;
602 bend = mol->nextBend(bendIter)) {
603
604 a = bend->getAtomA()->getGlobalIndex();
605 b = bend->getAtomB()->getGlobalIndex();
606 c = bend->getAtomC()->getGlobalIndex();
607
608 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
609 oneTwoInteractions_.removePair(a, b);
610 oneTwoInteractions_.removePair(b, c);
611 } else {
612 excludedInteractions_.removePair(a, b);
613 excludedInteractions_.removePair(b, c);
614 }
615
616 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
617 oneThreeInteractions_.removePair(a, c);
618 } else {
619 excludedInteractions_.removePair(a, c);
620 }
621 }
622
623 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
624 torsion = mol->nextTorsion(torsionIter)) {
625
626 a = torsion->getAtomA()->getGlobalIndex();
627 b = torsion->getAtomB()->getGlobalIndex();
628 c = torsion->getAtomC()->getGlobalIndex();
629 d = torsion->getAtomD()->getGlobalIndex();
630
631 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
632 oneTwoInteractions_.removePair(a, b);
633 oneTwoInteractions_.removePair(b, c);
634 oneTwoInteractions_.removePair(c, d);
635 } else {
636 excludedInteractions_.removePair(a, b);
637 excludedInteractions_.removePair(b, c);
638 excludedInteractions_.removePair(c, d);
639 }
640
641 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
642 oneThreeInteractions_.removePair(a, c);
643 oneThreeInteractions_.removePair(b, d);
644 } else {
645 excludedInteractions_.removePair(a, c);
646 excludedInteractions_.removePair(b, d);
647 }
648
649 if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
650 oneFourInteractions_.removePair(a, d);
651 } else {
652 excludedInteractions_.removePair(a, d);
653 }
654 }
655
656 for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
657 inversion = mol->nextInversion(inversionIter)) {
658
659 a = inversion->getAtomA()->getGlobalIndex();
660 b = inversion->getAtomB()->getGlobalIndex();
661 c = inversion->getAtomC()->getGlobalIndex();
662 d = inversion->getAtomD()->getGlobalIndex();
663
664 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
665 oneTwoInteractions_.removePair(a, b);
666 oneTwoInteractions_.removePair(a, c);
667 oneTwoInteractions_.removePair(a, d);
668 } else {
669 excludedInteractions_.removePair(a, b);
670 excludedInteractions_.removePair(a, c);
671 excludedInteractions_.removePair(a, d);
672 }
673
674 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
675 oneThreeInteractions_.removePair(b, c);
676 oneThreeInteractions_.removePair(b, d);
677 oneThreeInteractions_.removePair(c, d);
678 } else {
679 excludedInteractions_.removePair(b, c);
680 excludedInteractions_.removePair(b, d);
681 excludedInteractions_.removePair(c, d);
682 }
683 }
684
685 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
686 rb = mol->nextRigidBody(rbIter)) {
687 vector<Atom*> atoms = rb->getAtoms();
688 for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
689 for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
690 a = atoms[i]->getGlobalIndex();
691 b = atoms[j]->getGlobalIndex();
692 excludedInteractions_.removePair(a, b);
693 }
694 }
695 }
696
697 }
698
699
700 void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) {
701 int curStampId;
702
703 //index from 0
704 curStampId = moleculeStamps_.size();
705
706 moleculeStamps_.push_back(molStamp);
707 molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
708 }
709
710
711 /**
712 * update
713 *
714 * Performs the global checks and variable settings after the
715 * objects have been created.
716 *
717 */
718 void SimInfo::update() {
719 setupSimVariables();
720 calcNdf();
721 calcNdfRaw();
722 calcNdfTrans();
723 }
724
725 /**
726 * getSimulatedAtomTypes
727 *
728 * Returns an STL set of AtomType* that are actually present in this
729 * simulation. Must query all processors to assemble this information.
730 *
731 */
732 set<AtomType*> SimInfo::getSimulatedAtomTypes() {
733 SimInfo::MoleculeIterator mi;
734 Molecule* mol;
735 Molecule::AtomIterator ai;
736 Atom* atom;
737 set<AtomType*> atomTypes;
738
739 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
740 for(atom = mol->beginAtom(ai); atom != NULL;
741 atom = mol->nextAtom(ai)) {
742 atomTypes.insert(atom->getAtomType());
743 }
744 }
745
746 #ifdef IS_MPI
747
748 // loop over the found atom types on this processor, and add their
749 // numerical idents to a vector:
750
751 vector<int> foundTypes;
752 set<AtomType*>::iterator i;
753 for (i = atomTypes.begin(); i != atomTypes.end(); ++i)
754 foundTypes.push_back( (*i)->getIdent() );
755
756 // count_local holds the number of found types on this processor
757 int count_local = foundTypes.size();
758
759 int nproc;
760 MPI_Comm_size( MPI_COMM_WORLD, &nproc);
761 // int nproc = MPI::COMM_WORLD.Get_size();
762
763 // we need arrays to hold the counts and displacement vectors for
764 // all processors
765 vector<int> counts(nproc, 0);
766 vector<int> disps(nproc, 0);
767
768 // fill the counts array
769 MPI_Allgather(&count_local, 1, MPI_INT, &counts[0],
770 1, MPI_INT, MPI_COMM_WORLD);
771 // MPI::COMM_WORLD.Allgather(&count_local, 1, MPI::INT, &counts[0],
772 // 1, MPI::INT);
773
774 // use the processor counts to compute the displacement array
775 disps[0] = 0;
776 int totalCount = counts[0];
777 for (int iproc = 1; iproc < nproc; iproc++) {
778 disps[iproc] = disps[iproc-1] + counts[iproc-1];
779 totalCount += counts[iproc];
780 }
781
782 // we need a (possibly redundant) set of all found types:
783 vector<int> ftGlobal(totalCount);
784
785 // now spray out the foundTypes to all the other processors:
786 MPI_Allgatherv(&foundTypes[0], count_local, MPI_INT,
787 &ftGlobal[0], &counts[0], &disps[0],
788 MPI_INT, MPI_COMM_WORLD);
789 // MPI::COMM_WORLD.Allgatherv(&foundTypes[0], count_local, MPI::INT,
790 // &ftGlobal[0], &counts[0], &disps[0],
791 // MPI::INT);
792
793 vector<int>::iterator j;
794
795 // foundIdents is a stl set, so inserting an already found ident
796 // will have no effect.
797 set<int> foundIdents;
798
799 for (j = ftGlobal.begin(); j != ftGlobal.end(); ++j)
800 foundIdents.insert((*j));
801
802 // now iterate over the foundIdents and get the actual atom types
803 // that correspond to these:
804 set<int>::iterator it;
805 for (it = foundIdents.begin(); it != foundIdents.end(); ++it)
806 atomTypes.insert( forceField_->getAtomType((*it)) );
807
808 #endif
809
810 return atomTypes;
811 }
812
813
814 int getGlobalCountOfType(AtomType* atype) {
815 /*
816 set<AtomType*> atypes = getSimulatedAtomTypes();
817 map<AtomType*, int> counts_;
818
819 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
820 for(atom = mol->beginAtom(ai); atom != NULL;
821 atom = mol->nextAtom(ai)) {
822 atom->getAtomType();
823 }
824 }
825 */
826 return 0;
827 }
828
829 void SimInfo::setupSimVariables() {
830 useAtomicVirial_ = simParams_->getUseAtomicVirial();
831 // we only call setAccumulateBoxDipole if the accumulateBoxDipole
832 // parameter is true
833 calcBoxDipole_ = false;
834 if ( simParams_->haveAccumulateBoxDipole() )
835 if ( simParams_->getAccumulateBoxDipole() ) {
836 calcBoxDipole_ = true;
837 }
838
839 set<AtomType*>::iterator i;
840 set<AtomType*> atomTypes;
841 atomTypes = getSimulatedAtomTypes();
842 bool usesElectrostatic = false;
843 bool usesMetallic = false;
844 bool usesDirectional = false;
845 bool usesFluctuatingCharges = false;
846 //loop over all of the atom types
847 for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
848 usesElectrostatic |= (*i)->isElectrostatic();
849 usesMetallic |= (*i)->isMetal();
850 usesDirectional |= (*i)->isDirectional();
851 usesFluctuatingCharges |= (*i)->isFluctuatingCharge();
852 }
853
854 #ifdef IS_MPI
855 int temp;
856
857 temp = usesDirectional;
858 MPI_Allreduce(MPI_IN_PLACE, &temp, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
859 usesDirectionalAtoms_ = (temp == 0) ? false : true;
860
861 // MPI::COMM_WORLD.Allreduce(&temp, &usesDirectionalAtoms_, 1, MPI::BOOL,
862 // MPI::LOR);
863
864 temp = usesMetallic;
865 MPI_Allreduce(MPI_IN_PLACE, &temp, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
866 usesMetallicAtoms_ = (temp == 0) ? false : true;
867
868 // MPI::COMM_WORLD.Allreduce(&temp, &usesMetallicAtoms_, 1, MPI::BOOL,
869 // MPI::LOR);
870
871 temp = usesElectrostatic;
872 MPI_Allreduce(MPI_IN_PLACE, &temp, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
873 usesElectrostaticAtoms_ = (temp == 0) ? false : true;
874
875 // MPI::COMM_WORLD.Allreduce(&temp, &usesElectrostaticAtoms_, 1, MPI::BOOL,
876 // MPI::LOR);
877
878 temp = usesFluctuatingCharges;
879 MPI_Allreduce(MPI_IN_PLACE, &temp, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
880 usesFluctuatingCharges_ = (temp == 0) ? false : true;
881
882 // MPI::COMM_WORLD.Allreduce(&temp, &usesFluctuatingCharges_, 1, MPI::BOOL,
883 // MPI::LOR);
884
885 #else
886
887 usesDirectionalAtoms_ = usesDirectional;
888 usesMetallicAtoms_ = usesMetallic;
889 usesElectrostaticAtoms_ = usesElectrostatic;
890 usesFluctuatingCharges_ = usesFluctuatingCharges;
891
892 #endif
893
894 requiresPrepair_ = usesMetallicAtoms_ ? true : false;
895 requiresSkipCorrection_ = usesElectrostaticAtoms_ ? true : false;
896 requiresSelfCorrection_ = usesElectrostaticAtoms_ ? true : false;
897 }
898
899
900 vector<int> SimInfo::getGlobalAtomIndices() {
901 SimInfo::MoleculeIterator mi;
902 Molecule* mol;
903 Molecule::AtomIterator ai;
904 Atom* atom;
905
906 vector<int> GlobalAtomIndices(getNAtoms(), 0);
907
908 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
909
910 for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
911 GlobalAtomIndices[atom->getLocalIndex()] = atom->getGlobalIndex();
912 }
913 }
914 return GlobalAtomIndices;
915 }
916
917
918 vector<int> SimInfo::getGlobalGroupIndices() {
919 SimInfo::MoleculeIterator mi;
920 Molecule* mol;
921 Molecule::CutoffGroupIterator ci;
922 CutoffGroup* cg;
923
924 vector<int> GlobalGroupIndices;
925
926 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
927
928 //local index of cutoff group is trivial, it only depends on the
929 //order of travesing
930 for (cg = mol->beginCutoffGroup(ci); cg != NULL;
931 cg = mol->nextCutoffGroup(ci)) {
932 GlobalGroupIndices.push_back(cg->getGlobalIndex());
933 }
934 }
935 return GlobalGroupIndices;
936 }
937
938
939 void SimInfo::prepareTopology() {
940
941 //calculate mass ratio of cutoff group
942 SimInfo::MoleculeIterator mi;
943 Molecule* mol;
944 Molecule::CutoffGroupIterator ci;
945 CutoffGroup* cg;
946 Molecule::AtomIterator ai;
947 Atom* atom;
948 RealType totalMass;
949
950 /**
951 * The mass factor is the relative mass of an atom to the total
952 * mass of the cutoff group it belongs to. By default, all atoms
953 * are their own cutoff groups, and therefore have mass factors of
954 * 1. We need some special handling for massless atoms, which
955 * will be treated as carrying the entire mass of the cutoff
956 * group.
957 */
958 massFactors_.clear();
959 massFactors_.resize(getNAtoms(), 1.0);
960
961 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
962 for (cg = mol->beginCutoffGroup(ci); cg != NULL;
963 cg = mol->nextCutoffGroup(ci)) {
964
965 totalMass = cg->getMass();
966 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
967 // Check for massless groups - set mfact to 1 if true
968 if (totalMass != 0)
969 massFactors_[atom->getLocalIndex()] = atom->getMass()/totalMass;
970 else
971 massFactors_[atom->getLocalIndex()] = 1.0;
972 }
973 }
974 }
975
976 // Build the identArray_ and regions_
977
978 identArray_.clear();
979 identArray_.reserve(getNAtoms());
980 regions_.clear();
981 regions_.reserve(getNAtoms());
982
983 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
984 int reg = mol->getRegion();
985 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
986 identArray_.push_back(atom->getIdent());
987 regions_.push_back(reg);
988 }
989 }
990
991 topologyDone_ = true;
992 }
993
994 void SimInfo::addProperty(GenericData* genData) {
995 properties_.addProperty(genData);
996 }
997
998 void SimInfo::removeProperty(const string& propName) {
999 properties_.removeProperty(propName);
1000 }
1001
1002 void SimInfo::clearProperties() {
1003 properties_.clearProperties();
1004 }
1005
1006 vector<string> SimInfo::getPropertyNames() {
1007 return properties_.getPropertyNames();
1008 }
1009
1010 vector<GenericData*> SimInfo::getProperties() {
1011 return properties_.getProperties();
1012 }
1013
1014 GenericData* SimInfo::getPropertyByName(const string& propName) {
1015 return properties_.getPropertyByName(propName);
1016 }
1017
1018 void SimInfo::setSnapshotManager(SnapshotManager* sman) {
1019 if (sman_ == sman) {
1020 return;
1021 }
1022 delete sman_;
1023 sman_ = sman;
1024
1025 SimInfo::MoleculeIterator mi;
1026 Molecule::AtomIterator ai;
1027 Molecule::RigidBodyIterator rbIter;
1028 Molecule::CutoffGroupIterator cgIter;
1029 Molecule::BondIterator bondIter;
1030 Molecule::BendIterator bendIter;
1031 Molecule::TorsionIterator torsionIter;
1032 Molecule::InversionIterator inversionIter;
1033
1034 Molecule* mol;
1035 Atom* atom;
1036 RigidBody* rb;
1037 CutoffGroup* cg;
1038 Bond* bond;
1039 Bend* bend;
1040 Torsion* torsion;
1041 Inversion* inversion;
1042
1043 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
1044
1045 for (atom = mol->beginAtom(ai); atom != NULL;
1046 atom = mol->nextAtom(ai)) {
1047 atom->setSnapshotManager(sman_);
1048 }
1049 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
1050 rb = mol->nextRigidBody(rbIter)) {
1051 rb->setSnapshotManager(sman_);
1052 }
1053 for (cg = mol->beginCutoffGroup(cgIter); cg != NULL;
1054 cg = mol->nextCutoffGroup(cgIter)) {
1055 cg->setSnapshotManager(sman_);
1056 }
1057 for (bond = mol->beginBond(bondIter); bond != NULL;
1058 bond = mol->nextBond(bondIter)) {
1059 bond->setSnapshotManager(sman_);
1060 }
1061 for (bend = mol->beginBend(bendIter); bend != NULL;
1062 bend = mol->nextBend(bendIter)) {
1063 bend->setSnapshotManager(sman_);
1064 }
1065 for (torsion = mol->beginTorsion(torsionIter); torsion != NULL;
1066 torsion = mol->nextTorsion(torsionIter)) {
1067 torsion->setSnapshotManager(sman_);
1068 }
1069 for (inversion = mol->beginInversion(inversionIter); inversion != NULL;
1070 inversion = mol->nextInversion(inversionIter)) {
1071 inversion->setSnapshotManager(sman_);
1072 }
1073 }
1074 }
1075
1076
1077 ostream& operator <<(ostream& o, SimInfo& info) {
1078
1079 return o;
1080 }
1081
1082
1083 StuntDouble* SimInfo::getIOIndexToIntegrableObject(int index) {
1084 if (index >= int(IOIndexToIntegrableObject.size())) {
1085 sprintf(painCave.errMsg,
1086 "SimInfo::getIOIndexToIntegrableObject Error: Integrable Object\n"
1087 "\tindex exceeds number of known objects!\n");
1088 painCave.isFatal = 1;
1089 simError();
1090 return NULL;
1091 } else
1092 return IOIndexToIntegrableObject.at(index);
1093 }
1094
1095 void SimInfo::setIOIndexToIntegrableObject(const vector<StuntDouble*>& v) {
1096 IOIndexToIntegrableObject= v;
1097 }
1098
1099 int SimInfo::getNGlobalConstraints() {
1100 int nGlobalConstraints;
1101 #ifdef IS_MPI
1102 MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1,
1103 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
1104 // MPI::COMM_WORLD.Allreduce(&nConstraints_, &nGlobalConstraints, 1,
1105 // MPI::INT, MPI::SUM);
1106 #else
1107 nGlobalConstraints = nConstraints_;
1108 #endif
1109 return nGlobalConstraints;
1110 }
1111
1112 }//end namespace OpenMD
1113

Properties

Name Value
svn:keywords Author Id Revision Date