ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimInfo.cpp
Revision: 2022
Committed: Fri Sep 26 22:22:28 2014 UTC (10 years, 7 months ago) by gezelter
File size: 34673 byte(s)
Log Message:
Added support for accumulateBoxQuadrupole flag

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

Properties

Name Value
svn:keywords Author Id Revision Date