OpenMD 3.0
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
SimInfo.cpp
Go to the documentation of this file.
1/*
2 * Copyright (c) 2004-present, The University of Notre Dame. All rights
3 * reserved.
4 *
5 * Redistribution and use in source and binary forms, with or without
6 * modification, are permitted provided that the following conditions are met:
7 *
8 * 1. Redistributions of source code must retain the above copyright notice,
9 * this list of conditions and the following disclaimer.
10 *
11 * 2. Redistributions in binary form must reproduce the above copyright notice,
12 * this list of conditions and the following disclaimer in the documentation
13 * and/or other materials provided with the distribution.
14 *
15 * 3. Neither the name of the copyright holder nor the names of its
16 * contributors may be used to endorse or promote products derived from
17 * this software without specific prior written permission.
18 *
19 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
20 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
21 * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
22 * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
23 * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
24 * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
25 * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
26 * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
27 * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
28 * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
29 * POSSIBILITY OF SUCH DAMAGE.
30 *
31 * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
32 * research, please cite the appropriate papers when you publish your
33 * work. Good starting points are:
34 *
35 * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
36 * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
37 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
38 * [4] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
39 * [5] Kuang & Gezelter, Mol. Phys., 110, 691-701 (2012).
40 * [6] Lamichhane, Gezelter & Newman, J. Chem. Phys. 141, 134109 (2014).
41 * [7] Lamichhane, Newman & Gezelter, J. Chem. Phys. 141, 134110 (2014).
42 * [8] Bhattarai, Newman & Gezelter, Phys. Rev. B 99, 094106 (2019).
43 */
44
45/**
46 * @file SimInfo.cpp
47 * @author tlin
48 * @date 11/02/2004
49 * @version 1.0
50 */
51
52#include "brains/SimInfo.hpp"
53
54#include <algorithm>
55#include <cstdint>
56#include <map>
57#include <memory>
58#include <random>
59#include <set>
60
61#ifdef IS_MPI
62#include <mpi.h>
63#endif
64
65#include "brains/ForceField.hpp"
66#include "io/ForceFieldOptions.hpp"
67#include "math/Vector3.hpp"
68#include "nonbonded/SwitchingFunction.hpp"
71#include "selection/SelectionManager.hpp"
72#include "utils/MemoryUtils.hpp"
73#include "utils/RandNumGen.hpp"
74#include "utils/simError.h"
75
76using namespace std;
77namespace OpenMD {
78
80 forceField_(ff), simParams_(simParams), randNumGen_ {nullptr}, nAtoms_(0),
81 nBonds_(0), nBends_(0), nTorsions_(0), nInversions_(0), nRigidBodies_(0),
82 nIntegrableObjects_(0), nCutoffGroups_(0), nConstraints_(0),
83 nFluctuatingCharges_(0), nGlobalMols_(0), nGlobalAtoms_(0),
84 nGlobalCutoffGroups_(0), nGlobalIntegrableObjects_(0),
85 nGlobalRigidBodies_(0), nGlobalFluctuatingCharges_(0), nGlobalBonds_(0),
86 nGlobalBends_(0), nGlobalTorsions_(0), nGlobalInversions_(0),
87 nGlobalConstraints_(0), hasNGlobalConstraints_(false), ndf_(0),
88 fdf_local(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0), sman_(NULL),
89 topologyDone_(false), calcBoxDipole_(false), calcBoxQuadrupole_(false),
90 useAtomicVirial_(true) {
91 MoleculeStamp* molStamp;
92 int nMolWithSameStamp;
93 int nCutoffAtoms = 0; // number of atoms belong to cutoff groups
94 int nGroups = 0; // total cutoff groups defined in meta-data file
95 CutoffGroupStamp* cgStamp;
96 RigidBodyStamp* rbStamp;
97 int nRigidAtoms = 0;
98
99 vector<Component*> components = simParams->getComponents();
100
101 for (vector<Component*>::iterator i = components.begin();
102 i != components.end(); ++i) {
103 molStamp = (*i)->getMoleculeStamp();
104 if ((*i)->haveRegion()) {
105 molStamp->setRegion((*i)->getRegion());
106 } else {
107 // set the region to a disallowed value:
108 molStamp->setRegion(-1);
109 }
110
111 nMolWithSameStamp = (*i)->getNMol();
112
113 addMoleculeStamp(molStamp, nMolWithSameStamp);
114
115 // calculate atoms in molecules
116 nGlobalAtoms_ += molStamp->getNAtoms() * nMolWithSameStamp;
117 nGlobalBonds_ += molStamp->getNBonds() * nMolWithSameStamp;
118 nGlobalBends_ += molStamp->getNBends() * nMolWithSameStamp;
119 nGlobalTorsions_ += molStamp->getNTorsions() * nMolWithSameStamp;
120 nGlobalInversions_ += molStamp->getNInversions() * nMolWithSameStamp;
121
122 // calculate atoms in cutoff groups
123 int nAtomsInGroups = 0;
124 int nCutoffGroupsInStamp = molStamp->getNCutoffGroups();
125
126 for (int j = 0; j < nCutoffGroupsInStamp; j++) {
127 cgStamp = molStamp->getCutoffGroupStamp(j);
128 nAtomsInGroups += cgStamp->getNMembers();
129 }
130
131 nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
132
133 nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;
134
135 // calculate atoms in rigid bodies
136 int nAtomsInRigidBodies = 0;
137 int nRigidBodiesInStamp = molStamp->getNRigidBodies();
138
139 for (int j = 0; j < nRigidBodiesInStamp; j++) {
140 rbStamp = molStamp->getRigidBodyStamp(j);
141 nAtomsInRigidBodies += rbStamp->getNMembers();
142 }
143
144 nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
145 nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;
146 }
147
148 // every free atom (atom does not belong to cutoff groups) is a cutoff
149 // group therefore the total number of cutoff groups in the system is
150 // equal to the total number of atoms minus number of atoms belong to
151 // cutoff group defined in meta-data file plus the number of cutoff
152 // groups defined in meta-data file
153
154 nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
155
156 // every free atom (atom does not belong to rigid bodies) is an
157 // integrable object therefore the total number of integrable objects
158 // in the system is equal to the total number of atoms minus number of
159 // atoms belong to rigid body defined in meta-data file plus the number
160 // of rigid bodies defined in meta-data file
161 nGlobalIntegrableObjects_ =
162 nGlobalAtoms_ - nRigidAtoms + nGlobalRigidBodies_;
163
164 nGlobalMols_ = molStampIds_.size();
165 molToProcMap_.resize(nGlobalMols_);
166
167 // Initialize the random number generator on each processing
168 // element using either the user-supplied seed or a default seed
169 std::uint_fast32_t seed;
170
171 if (simParams_->haveSeed())
172 seed = static_cast<std::uint_fast32_t>(simParams_->getSeed());
173 else
174 seed = std::mt19937::default_seed;
175
176 randNumGen_ = std::make_shared<Utils::RandNumGen>(seed);
177 }
178
179 SimInfo::~SimInfo() {
180 Utils::deletePointers(molecules_);
181
182 delete sman_;
183 delete simParams_;
184 delete forceField_;
185 }
186
188 MoleculeIterator i;
189
190 i = molecules_.find(mol->getGlobalIndex());
191 if (i == molecules_.end()) {
192 molecules_.insert(make_pair(mol->getGlobalIndex(), mol));
193
194 nAtoms_ += mol->getNAtoms();
195 nBonds_ += mol->getNBonds();
196 nBends_ += mol->getNBends();
197 nTorsions_ += mol->getNTorsions();
198 nInversions_ += mol->getNInversions();
199 nRigidBodies_ += mol->getNRigidBodies();
200 nIntegrableObjects_ += mol->getNIntegrableObjects();
201 nCutoffGroups_ += mol->getNCutoffGroups();
202 nConstraints_ += mol->getNConstraintPairs();
203
205
206 return true;
207 } else {
208 return false;
209 }
210 }
211
213 MoleculeIterator i;
214 i = molecules_.find(mol->getGlobalIndex());
215
216 if (i != molecules_.end()) {
217 assert(mol == i->second);
218
219 nAtoms_ -= mol->getNAtoms();
220 nBonds_ -= mol->getNBonds();
221 nBends_ -= mol->getNBends();
222 nTorsions_ -= mol->getNTorsions();
223 nInversions_ -= mol->getNInversions();
224 nRigidBodies_ -= mol->getNRigidBodies();
225 nIntegrableObjects_ -= mol->getNIntegrableObjects();
226 nCutoffGroups_ -= mol->getNCutoffGroups();
227 nConstraints_ -= mol->getNConstraintPairs();
228
230 molecules_.erase(mol->getGlobalIndex());
231
232 delete mol;
233
234 return true;
235 } else {
236 return false;
237 }
238 }
239
240 Molecule* SimInfo::beginMolecule(MoleculeIterator& i) {
241 i = molecules_.begin();
242 return i == molecules_.end() ? NULL : i->second;
243 }
244
245 Molecule* SimInfo::nextMolecule(MoleculeIterator& i) {
246 ++i;
247 return i == molecules_.end() ? NULL : i->second;
248 }
249
250 void SimInfo::calcNdf() {
251 int ndf_local, nfq_local;
252 MoleculeIterator i;
253 vector<StuntDouble*>::iterator j;
254 vector<Atom*>::iterator k;
255
256 Molecule* mol;
257 StuntDouble* sd;
258 Atom* atom;
259
260 ndf_local = 0;
261 nfq_local = 0;
262
263 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
264 for (sd = mol->beginIntegrableObject(j); sd != NULL;
265 sd = mol->nextIntegrableObject(j)) {
266 ndf_local += 3;
267
268 if (sd->isDirectional()) {
269 if (sd->isLinear()) {
270 ndf_local += 2;
271 } else {
272 ndf_local += 3;
273 }
274 }
275 }
276
277 for (atom = mol->beginFluctuatingCharge(k); atom != NULL;
278 atom = mol->nextFluctuatingCharge(k)) {
279 if (atom->isFluctuatingCharge()) { nfq_local++; }
280 }
281 }
282
283 ndfLocal_ = ndf_local;
284
285 // n_constraints is local, so subtract them on each processor
286 ndf_local -= nConstraints_;
287
288#ifdef IS_MPI
289 MPI_Allreduce(&ndf_local, &ndf_, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
290 MPI_Allreduce(&nfq_local, &nGlobalFluctuatingCharges_, 1, MPI_INT, MPI_SUM,
291 MPI_COMM_WORLD);
292#else
293 ndf_ = ndf_local;
294 nGlobalFluctuatingCharges_ = nfq_local;
295#endif
296
297 // nZconstraints_ is global, as are the 3 COM translations for the
298 // entire system:
299 ndf_ = ndf_ - 3 - nZconstraint_;
300 }
301
302 int SimInfo::getFdf() {
303#ifdef IS_MPI
304 MPI_Allreduce(&fdf_local, &fdf_, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
305#else
306 fdf_ = fdf_local;
307#endif
308 return fdf_;
309 }
310
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 for (cg = mol->beginCutoffGroup(ci); cg != NULL;
320 cg = mol->nextCutoffGroup(ci)) {
321 nLocalCutoffAtoms += cg->getNumAtom();
322 }
323 }
324
325 return nAtoms_ - nLocalCutoffAtoms + nCutoffGroups_;
326 }
327
328 void SimInfo::calcNdfRaw() {
329 int ndfRaw_local;
330
331 MoleculeIterator i;
332 vector<StuntDouble*>::iterator j;
333 Molecule* mol;
334 StuntDouble* sd;
335
336 // Raw degrees of freedom that we have to set
337 ndfRaw_local = 0;
338
339 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
340 for (sd = mol->beginIntegrableObject(j); sd != NULL;
341 sd = mol->nextIntegrableObject(j)) {
342 ndfRaw_local += 3;
343
344 if (sd->isDirectional()) {
345 if (sd->isLinear()) {
346 ndfRaw_local += 2;
347 } else {
348 ndfRaw_local += 3;
349 }
350 }
351 }
352 }
353
354#ifdef IS_MPI
355 MPI_Allreduce(&ndfRaw_local, &ndfRaw_, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
356#else
357 ndfRaw_ = ndfRaw_local;
358#endif
359 }
360
361 void SimInfo::calcNdfTrans() {
362 int ndfTrans_local;
363
364 ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
365
366#ifdef IS_MPI
367 MPI_Allreduce(&ndfTrans_local, &ndfTrans_, 1, MPI_INT, MPI_SUM,
368 MPI_COMM_WORLD);
369#else
370 ndfTrans_ = ndfTrans_local;
371#endif
372
373 ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
374 }
375
377 ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
378 vector<Bond*>::iterator bondIter;
379 vector<Bend*>::iterator bendIter;
380 vector<Torsion*>::iterator torsionIter;
381 vector<Inversion*>::iterator inversionIter;
382 Bond* bond;
383 Bend* bend;
384 Torsion* torsion;
385 Inversion* inversion;
386 int a;
387 int b;
388 int c;
389 int d;
390
391 // atomGroups can be used to add special interaction maps between
392 // groups of atoms that are in two separate rigid bodies.
393 // However, most site-site interactions between two rigid bodies
394 // are probably not special, just the ones between the physically
395 // bonded atoms. Interactions *within* a single rigid body should
396 // always be excluded. These are done at the bottom of this
397 // function.
398
399 map<int, set<int>> atomGroups;
400 Molecule::RigidBodyIterator rbIter;
401 RigidBody* rb;
402 Molecule::IntegrableObjectIterator ii;
403 StuntDouble* sd;
404
405 for (sd = mol->beginIntegrableObject(ii); sd != NULL;
406 sd = mol->nextIntegrableObject(ii)) {
407 if (sd->isRigidBody()) {
408 rb = static_cast<RigidBody*>(sd);
409 vector<Atom*> atoms = rb->getAtoms();
410 set<int> rigidAtoms;
411 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
412 rigidAtoms.insert(atoms[i]->getGlobalIndex());
413 }
414 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
415 atomGroups.insert(map<int, set<int>>::value_type(
416 atoms[i]->getGlobalIndex(), rigidAtoms));
417 }
418 } else {
419 set<int> oneAtomSet;
420 oneAtomSet.insert(sd->getGlobalIndex());
421 atomGroups.insert(
422 map<int, set<int>>::value_type(sd->getGlobalIndex(), oneAtomSet));
423 }
424 }
425
426 for (bond = mol->beginBond(bondIter); bond != NULL;
427 bond = mol->nextBond(bondIter)) {
428 a = bond->getAtomA()->getGlobalIndex();
429 b = bond->getAtomB()->getGlobalIndex();
430
431 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
432 oneTwoInteractions_.addPair(a, b);
433 } else {
434 excludedInteractions_.addPair(a, b);
435 }
436 }
437
438 for (bend = mol->beginBend(bendIter); bend != NULL;
439 bend = mol->nextBend(bendIter)) {
440 a = bend->getAtomA()->getGlobalIndex();
441 b = bend->getAtomB()->getGlobalIndex();
442 c = bend->getAtomC()->getGlobalIndex();
443
444 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
445 oneTwoInteractions_.addPair(a, b);
446 oneTwoInteractions_.addPair(b, c);
447 } else {
448 excludedInteractions_.addPair(a, b);
449 excludedInteractions_.addPair(b, c);
450 }
451
452 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
453 oneThreeInteractions_.addPair(a, c);
454 } else {
455 excludedInteractions_.addPair(a, c);
456 }
457 }
458
459 for (torsion = mol->beginTorsion(torsionIter); torsion != NULL;
460 torsion = mol->nextTorsion(torsionIter)) {
461 a = torsion->getAtomA()->getGlobalIndex();
462 b = torsion->getAtomB()->getGlobalIndex();
463 c = torsion->getAtomC()->getGlobalIndex();
464 d = torsion->getAtomD()->getGlobalIndex();
465
466 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
467 oneTwoInteractions_.addPair(a, b);
468 oneTwoInteractions_.addPair(b, c);
469 oneTwoInteractions_.addPair(c, d);
470 } else {
471 excludedInteractions_.addPair(a, b);
472 excludedInteractions_.addPair(b, c);
473 excludedInteractions_.addPair(c, d);
474 }
475
476 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
477 oneThreeInteractions_.addPair(a, c);
478 oneThreeInteractions_.addPair(b, d);
479 } else {
480 excludedInteractions_.addPair(a, c);
481 excludedInteractions_.addPair(b, d);
482 }
483
484 if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
485 oneFourInteractions_.addPair(a, d);
486 } else {
487 excludedInteractions_.addPair(a, d);
488 }
489 }
490
491 for (inversion = mol->beginInversion(inversionIter); inversion != NULL;
492 inversion = mol->nextInversion(inversionIter)) {
493 a = inversion->getAtomA()->getGlobalIndex();
494 b = inversion->getAtomB()->getGlobalIndex();
495 c = inversion->getAtomC()->getGlobalIndex();
496 d = inversion->getAtomD()->getGlobalIndex();
497
498 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
499 oneTwoInteractions_.addPair(a, b);
500 oneTwoInteractions_.addPair(a, c);
501 oneTwoInteractions_.addPair(a, d);
502 } else {
503 excludedInteractions_.addPair(a, b);
504 excludedInteractions_.addPair(a, c);
505 excludedInteractions_.addPair(a, d);
506 }
507
508 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
509 oneThreeInteractions_.addPair(b, c);
510 oneThreeInteractions_.addPair(b, d);
511 oneThreeInteractions_.addPair(c, d);
512 } else {
513 excludedInteractions_.addPair(b, c);
514 excludedInteractions_.addPair(b, d);
515 excludedInteractions_.addPair(c, d);
516 }
517 }
518
519 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
520 rb = mol->nextRigidBody(rbIter)) {
521 vector<Atom*> atoms = rb->getAtoms();
522 for (int i = 0; i < static_cast<int>(atoms.size()) - 1; ++i) {
523 for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
524 a = atoms[i]->getGlobalIndex();
525 b = atoms[j]->getGlobalIndex();
526 excludedInteractions_.addPair(a, b);
527 }
528 }
529 }
530 }
531
533 ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
534 vector<Bond*>::iterator bondIter;
535 vector<Bend*>::iterator bendIter;
536 vector<Torsion*>::iterator torsionIter;
537 vector<Inversion*>::iterator inversionIter;
538 Bond* bond;
539 Bend* bend;
540 Torsion* torsion;
541 Inversion* inversion;
542 int a;
543 int b;
544 int c;
545 int d;
546
547 map<int, set<int>> atomGroups;
548 Molecule::RigidBodyIterator rbIter;
549 RigidBody* rb;
550 Molecule::IntegrableObjectIterator ii;
551 StuntDouble* sd;
552
553 for (sd = mol->beginIntegrableObject(ii); sd != NULL;
554 sd = mol->nextIntegrableObject(ii)) {
555 if (sd->isRigidBody()) {
556 rb = static_cast<RigidBody*>(sd);
557 vector<Atom*> atoms = rb->getAtoms();
558 set<int> rigidAtoms;
559 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
560 rigidAtoms.insert(atoms[i]->getGlobalIndex());
561 }
562 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
563 atomGroups.insert(map<int, set<int>>::value_type(
564 atoms[i]->getGlobalIndex(), rigidAtoms));
565 }
566 } else {
567 set<int> oneAtomSet;
568 oneAtomSet.insert(sd->getGlobalIndex());
569 atomGroups.insert(
570 map<int, set<int>>::value_type(sd->getGlobalIndex(), oneAtomSet));
571 }
572 }
573
574 for (bond = mol->beginBond(bondIter); bond != NULL;
575 bond = mol->nextBond(bondIter)) {
576 a = bond->getAtomA()->getGlobalIndex();
577 b = bond->getAtomB()->getGlobalIndex();
578
579 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
580 oneTwoInteractions_.removePair(a, b);
581 } else {
582 excludedInteractions_.removePair(a, b);
583 }
584 }
585
586 for (bend = mol->beginBend(bendIter); bend != NULL;
587 bend = mol->nextBend(bendIter)) {
588 a = bend->getAtomA()->getGlobalIndex();
589 b = bend->getAtomB()->getGlobalIndex();
590 c = bend->getAtomC()->getGlobalIndex();
591
592 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
593 oneTwoInteractions_.removePair(a, b);
594 oneTwoInteractions_.removePair(b, c);
595 } else {
596 excludedInteractions_.removePair(a, b);
597 excludedInteractions_.removePair(b, c);
598 }
599
600 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
601 oneThreeInteractions_.removePair(a, c);
602 } else {
603 excludedInteractions_.removePair(a, c);
604 }
605 }
606
607 for (torsion = mol->beginTorsion(torsionIter); torsion != NULL;
608 torsion = mol->nextTorsion(torsionIter)) {
609 a = torsion->getAtomA()->getGlobalIndex();
610 b = torsion->getAtomB()->getGlobalIndex();
611 c = torsion->getAtomC()->getGlobalIndex();
612 d = torsion->getAtomD()->getGlobalIndex();
613
614 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
615 oneTwoInteractions_.removePair(a, b);
616 oneTwoInteractions_.removePair(b, c);
617 oneTwoInteractions_.removePair(c, d);
618 } else {
619 excludedInteractions_.removePair(a, b);
620 excludedInteractions_.removePair(b, c);
621 excludedInteractions_.removePair(c, d);
622 }
623
624 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
625 oneThreeInteractions_.removePair(a, c);
626 oneThreeInteractions_.removePair(b, d);
627 } else {
628 excludedInteractions_.removePair(a, c);
629 excludedInteractions_.removePair(b, d);
630 }
631
632 if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
633 oneFourInteractions_.removePair(a, d);
634 } else {
635 excludedInteractions_.removePair(a, d);
636 }
637 }
638
639 for (inversion = mol->beginInversion(inversionIter); inversion != NULL;
640 inversion = mol->nextInversion(inversionIter)) {
641 a = inversion->getAtomA()->getGlobalIndex();
642 b = inversion->getAtomB()->getGlobalIndex();
643 c = inversion->getAtomC()->getGlobalIndex();
644 d = inversion->getAtomD()->getGlobalIndex();
645
646 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
647 oneTwoInteractions_.removePair(a, b);
648 oneTwoInteractions_.removePair(a, c);
649 oneTwoInteractions_.removePair(a, d);
650 } else {
651 excludedInteractions_.removePair(a, b);
652 excludedInteractions_.removePair(a, c);
653 excludedInteractions_.removePair(a, d);
654 }
655
656 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
657 oneThreeInteractions_.removePair(b, c);
658 oneThreeInteractions_.removePair(b, d);
659 oneThreeInteractions_.removePair(c, d);
660 } else {
661 excludedInteractions_.removePair(b, c);
662 excludedInteractions_.removePair(b, d);
663 excludedInteractions_.removePair(c, d);
664 }
665 }
666
667 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
668 rb = mol->nextRigidBody(rbIter)) {
669 vector<Atom*> atoms = rb->getAtoms();
670 for (int i = 0; i < static_cast<int>(atoms.size()) - 1; ++i) {
671 for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
672 a = atoms[i]->getGlobalIndex();
673 b = atoms[j]->getGlobalIndex();
674 excludedInteractions_.removePair(a, b);
675 }
676 }
677 }
678 }
679
680 void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) {
681 int curStampId;
682
683 // index from 0
684 curStampId = moleculeStamps_.size();
685
686 moleculeStamps_.push_back(molStamp);
687 molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
688 }
689
690 /**
691 * update
692 *
693 * Performs the global checks and variable settings after the
694 * objects have been created.
695 *
696 */
698 setupSimVariables();
699 calcNConstraints();
700 calcNdf();
701 calcNdfRaw();
702 calcNdfTrans();
703 }
704
705 /**
706 * getSimulatedAtomTypes
707 *
708 * Returns an STL set of AtomType* that are actually present in this
709 * simulation. Must query all processors to assemble this information.
710 *
711 */
713 SimInfo::MoleculeIterator mi;
714 Molecule* mol;
715 Molecule::AtomIterator ai;
716 Atom* atom;
717 AtomTypeSet atomTypes;
718
719 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
720 for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
721 atomTypes.insert(atom->getAtomType());
722 }
723 }
724
725#ifdef IS_MPI
726
727 // loop over the found atom types on this processor, and add their
728 // numerical idents to a vector:
729
730 vector<int> foundTypes;
731 AtomTypeSet::iterator i;
732 for (i = atomTypes.begin(); i != atomTypes.end(); ++i)
733 foundTypes.push_back((*i)->getIdent());
734
735 // count_local holds the number of found types on this processor
736 int count_local = foundTypes.size();
737
738 int nproc;
739 MPI_Comm_size(MPI_COMM_WORLD, &nproc);
740
741 // we need arrays to hold the counts and displacement vectors for
742 // all processors
743 vector<int> counts(nproc, 0);
744 vector<int> disps(nproc, 0);
745
746 // fill the counts array
747 MPI_Allgather(&count_local, 1, MPI_INT, &counts[0], 1, MPI_INT,
748 MPI_COMM_WORLD);
749
750 // use the processor counts to compute the displacement array
751 disps[0] = 0;
752 int totalCount = counts[0];
753 for (int iproc = 1; iproc < nproc; iproc++) {
754 disps[iproc] = disps[iproc - 1] + counts[iproc - 1];
755 totalCount += counts[iproc];
756 }
757
758 // we need a (possibly redundant) set of all found types:
759 vector<int> ftGlobal(totalCount);
760
761 // now spray out the foundTypes to all the other processors:
762 MPI_Allgatherv(&foundTypes[0], count_local, MPI_INT, &ftGlobal[0],
763 &counts[0], &disps[0], MPI_INT, MPI_COMM_WORLD);
764
765 vector<int>::iterator j;
766
767 // foundIdents is a stl set, so inserting an already found ident
768 // will have no effect.
769 set<int> foundIdents;
770
771 for (j = ftGlobal.begin(); j != ftGlobal.end(); ++j)
772 foundIdents.insert((*j));
773
774 // now iterate over the foundIdents and get the actual atom types
775 // that correspond to these:
776 set<int>::iterator it;
777 for (it = foundIdents.begin(); it != foundIdents.end(); ++it)
778 atomTypes.insert(forceField_->getAtomType((*it)));
779
780#endif
781
782 return atomTypes;
783 }
784
785 int getGlobalCountOfType(AtomType*) {
786 /*
787 AtomTypeSet atypes = getSimulatedAtomTypes();
788 map<AtomType*, int> counts_;
789
790 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
791 for(atom = mol->beginAtom(ai); atom != NULL;
792 atom = mol->nextAtom(ai)) {
793 atom->getAtomType();
794 }
795 }
796 */
797 return 0;
798 }
799
800 void SimInfo::setupSimVariables() {
801 useAtomicVirial_ = simParams_->getUseAtomicVirial();
802 // we only call setAccumulateBoxDipole if the accumulateBoxDipole
803 // parameter is true
804 calcBoxDipole_ = false;
805 if (simParams_->haveAccumulateBoxDipole())
806 if (simParams_->getAccumulateBoxDipole()) { calcBoxDipole_ = true; }
807 // we only call setAccumulateBoxQuadrupole if the accumulateBoxQuadrupole
808 // parameter is true
809 calcBoxQuadrupole_ = false;
810 if (simParams_->haveAccumulateBoxQuadrupole())
811 if (simParams_->getAccumulateBoxQuadrupole()) {
812 calcBoxQuadrupole_ = true;
813 }
814
815 AtomTypeSet::iterator i;
816 AtomTypeSet atomTypes;
817 atomTypes = getSimulatedAtomTypes();
818 bool usesElectrostatic = false;
819 bool usesMetallic = false;
820 bool usesDirectional = false;
821 bool usesFluctuatingCharges = false;
822 // loop over all of the atom types
823 for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
824 usesElectrostatic |= (*i)->isElectrostatic();
825 usesMetallic |= (*i)->isMetal();
826 usesDirectional |= (*i)->isDirectional();
827 usesFluctuatingCharges |= (*i)->isFluctuatingCharge();
828 }
829
830#ifdef IS_MPI
831 int temp;
832
833 temp = usesDirectional;
834 MPI_Allreduce(MPI_IN_PLACE, &temp, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
835 usesDirectionalAtoms_ = (temp == 0) ? false : true;
836
837 temp = usesMetallic;
838 MPI_Allreduce(MPI_IN_PLACE, &temp, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
839 usesMetallicAtoms_ = (temp == 0) ? false : true;
840
841 temp = usesElectrostatic;
842 MPI_Allreduce(MPI_IN_PLACE, &temp, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
843 usesElectrostaticAtoms_ = (temp == 0) ? false : true;
844
845 temp = usesFluctuatingCharges;
846 MPI_Allreduce(MPI_IN_PLACE, &temp, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
847 usesFluctuatingCharges_ = (temp == 0) ? false : true;
848#else
849
850 usesDirectionalAtoms_ = usesDirectional;
851 usesMetallicAtoms_ = usesMetallic;
852 usesElectrostaticAtoms_ = usesElectrostatic;
853 usesFluctuatingCharges_ = usesFluctuatingCharges;
854
855#endif
856
857 requiresPrepair_ = usesMetallicAtoms_ ? true : false;
858 requiresSkipCorrection_ = usesElectrostaticAtoms_ ? true : false;
859 requiresSelfCorrection_ = usesElectrostaticAtoms_ ? true : false;
860 }
861
863 SimInfo::MoleculeIterator mi;
864 Molecule* mol;
865 Molecule::AtomIterator ai;
866 Atom* atom;
867
868 vector<int> GlobalAtomIndices(getNAtoms(), 0);
869
870 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
871 for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
872 GlobalAtomIndices[atom->getLocalIndex()] = atom->getGlobalIndex();
873 }
874 }
875 return GlobalAtomIndices;
876 }
877
879 SimInfo::MoleculeIterator mi;
880 Molecule* mol;
881 Molecule::CutoffGroupIterator ci;
882 CutoffGroup* cg;
883
884 vector<int> GlobalGroupIndices;
885
886 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
887 // local index of cutoff group is trivial, it only depends on the
888 // order of travesing
889 for (cg = mol->beginCutoffGroup(ci); cg != NULL;
890 cg = mol->nextCutoffGroup(ci)) {
891 GlobalGroupIndices.push_back(cg->getGlobalIndex());
892 }
893 }
894 return GlobalGroupIndices;
895 }
896
898 // calculate mass ratio of cutoff group
899 SimInfo::MoleculeIterator mi;
900 Molecule* mol;
901 Molecule::CutoffGroupIterator ci;
902 CutoffGroup* cg;
903 Molecule::AtomIterator ai;
904 Atom* atom;
905 RealType totalMass;
906
907 /**
908 * The mass factor is the relative mass of an atom to the total
909 * mass of the cutoff group it belongs to. By default, all atoms
910 * are their own cutoff groups, and therefore have mass factors of
911 * 1. We need some special handling for massless atoms, which
912 * will be treated as carrying the entire mass of the cutoff
913 * group.
914 */
915 massFactors_.clear();
916 massFactors_.resize(getNAtoms(), 1.0);
917
918 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
919 for (cg = mol->beginCutoffGroup(ci); cg != NULL;
920 cg = mol->nextCutoffGroup(ci)) {
921 totalMass = cg->getMass();
922 for (atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
923 // Check for massless groups - set mfact to 1 if true
924 if (totalMass != 0)
925 massFactors_[atom->getLocalIndex()] = atom->getMass() / totalMass;
926 else
927 massFactors_[atom->getLocalIndex()] = 1.0;
928 }
929 }
930 }
931
932 // Build the identArray_ and regions_
933
934 identArray_.clear();
935 identArray_.reserve(getNAtoms());
936 regions_.clear();
937 regions_.reserve(getNAtoms());
938
939 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
940 int reg = mol->getRegion();
941 for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
942 identArray_.push_back(atom->getIdent());
943 regions_.push_back(reg);
944 }
945 }
946
947 topologyDone_ = true;
948 }
949
950 void SimInfo::addProperty(std::shared_ptr<GenericData> genData) {
951 properties_.addProperty(genData);
952 }
953
954 void SimInfo::removeProperty(const string& propName) {
955 properties_.removeProperty(propName);
956 }
957
958 std::vector<string> SimInfo::getPropertyNames() {
959 return properties_.getPropertyNames();
960 }
961
962 std::vector<std::shared_ptr<GenericData>> SimInfo::getProperties() {
963 return properties_.getProperties();
964 }
965
966 std::shared_ptr<GenericData> SimInfo::getPropertyByName(
967 const string& propName) {
968 return properties_.getPropertyByName(propName);
969 }
970
972 if (sman_ == sman) { return; }
973 delete sman_;
974 sman_ = sman;
975
976 SimInfo::MoleculeIterator mi;
977 Molecule::AtomIterator ai;
978 Molecule::RigidBodyIterator rbIter;
979 Molecule::CutoffGroupIterator cgIter;
980 Molecule::BondIterator bondIter;
981 Molecule::BendIterator bendIter;
982 Molecule::TorsionIterator torsionIter;
983 Molecule::InversionIterator inversionIter;
984
985 Molecule* mol;
986 Atom* atom;
987 RigidBody* rb;
988 CutoffGroup* cg;
989 Bond* bond;
990 Bend* bend;
991 Torsion* torsion;
992 Inversion* inversion;
993
994 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
995 for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
996 atom->setSnapshotManager(sman_);
997 }
998 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
999 rb = mol->nextRigidBody(rbIter)) {
1000 rb->setSnapshotManager(sman_);
1001 }
1002 for (cg = mol->beginCutoffGroup(cgIter); cg != NULL;
1003 cg = mol->nextCutoffGroup(cgIter)) {
1004 cg->setSnapshotManager(sman_);
1005 }
1006 for (bond = mol->beginBond(bondIter); bond != NULL;
1007 bond = mol->nextBond(bondIter)) {
1008 bond->setSnapshotManager(sman_);
1009 }
1010 for (bend = mol->beginBend(bendIter); bend != NULL;
1011 bend = mol->nextBend(bendIter)) {
1012 bend->setSnapshotManager(sman_);
1013 }
1014 for (torsion = mol->beginTorsion(torsionIter); torsion != NULL;
1015 torsion = mol->nextTorsion(torsionIter)) {
1016 torsion->setSnapshotManager(sman_);
1017 }
1018 for (inversion = mol->beginInversion(inversionIter); inversion != NULL;
1019 inversion = mol->nextInversion(inversionIter)) {
1020 inversion->setSnapshotManager(sman_);
1021 }
1022 }
1023 }
1024
1025 ostream& operator<<(ostream& o, SimInfo&) { return o; }
1026
1028 if (index >= int(IOIndexToIntegrableObject.size())) {
1029 snprintf(
1030 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
1031 "SimInfo::getIOIndexToIntegrableObject Error: Integrable Object\n"
1032 "\tindex exceeds number of known objects!\n");
1033 painCave.isFatal = 1;
1034 simError();
1035 return NULL;
1036 } else
1037 return IOIndexToIntegrableObject.at(index);
1038 }
1039
1040 void SimInfo::setIOIndexToIntegrableObject(const vector<StuntDouble*>& v) {
1041 IOIndexToIntegrableObject = v;
1042 }
1043
1044 void SimInfo::calcNConstraints() {
1045#ifdef IS_MPI
1046 MPI_Allreduce(&nConstraints_, &nGlobalConstraints_, 1, MPI_INT, MPI_SUM,
1047 MPI_COMM_WORLD);
1048#else
1049 nGlobalConstraints_ = nConstraints_;
1050#endif
1051 }
1052} // namespace OpenMD
AtomType is what OpenMD looks to for unchanging data about an atom.
Definition AtomType.hpp:66
AtomType * getAtomType(const std::string &at)
getAtomType by string
size_t getNIntegrableObjects()
Returns the total number of integrable objects in this molecule.
Definition Molecule.hpp:179
size_t getNInversions()
Returns the total number of improper torsions in this molecule.
Definition Molecule.hpp:173
int getGlobalIndex()
Returns the global index of this molecule.
Definition Molecule.hpp:107
size_t getNBends()
Returns the total number of bends in this molecule.
Definition Molecule.hpp:167
size_t getNConstraintPairs()
Returns the total number of constraints in this molecule.
Definition Molecule.hpp:185
size_t getNAtoms()
Returns the total number of atoms in this molecule.
Definition Molecule.hpp:161
size_t getNRigidBodies()
Returns the total number of rigid bodies in this molecule.
Definition Molecule.hpp:176
size_t getNBonds()
Returns the total number of bonds in this molecule.
Definition Molecule.hpp:164
size_t getNCutoffGroups()
Returns the total number of cutoff groups in this molecule.
Definition Molecule.hpp:182
size_t getNTorsions()
Returns the total number of torsions in this molecule.
Definition Molecule.hpp:170
void addPair(int i, int j)
Adds a pair into this PairList class.
Definition PairList.cpp:68
void removePair(int i, int j)
Remove a pair from PairList class.
Definition PairList.cpp:104
std::shared_ptr< GenericData > getPropertyByName(const std::string &propName)
Returns property.
std::vector< std::shared_ptr< GenericData > > getProperties()
Returns all of the properties in PropertyMap.
void addProperty(std::shared_ptr< GenericData > genData)
Adds property into property map.
bool removeProperty(const std::string &propName)
Removes property from PropertyMap by name.
std::vector< std::string > getPropertyNames()
Returns all names of properties.
std::vector< Atom * > getAtoms()
Returns the atoms of this rigid body.
void setSnapshotManager(SnapshotManager *sman)
Sets the Snapshot Manager of this ShortRangeInteraction.
One of the heavy-weight classes of OpenMD, SimInfo maintains objects and variables relating to the cu...
Definition SimInfo.hpp:93
bool removeMolecule(Molecule *mol)
Removes a molecule from SimInfo.
Definition SimInfo.cpp:212
std::vector< std::shared_ptr< GenericData > > getProperties()
Returns all of the properties in PropertyMap.
Definition SimInfo.cpp:962
std::vector< int > getGlobalGroupIndices()
returns a vector which maps the local cutoff group index on this processor to the global cutoff group...
Definition SimInfo.cpp:878
void removeInteractionPairs(Molecule *mol)
remove all special interaction pairs which belong to a molecule from the appropriate lists.
Definition SimInfo.cpp:532
Molecule * beginMolecule(MoleculeIterator &i)
Returns the first molecule in this SimInfo and intialize the iterator.
Definition SimInfo.cpp:240
unsigned int getNAtoms()
Returns the number of local atoms.
Definition SimInfo.hpp:172
void prepareTopology()
Do final bookkeeping before Force managers need their data.
Definition SimInfo.cpp:897
std::shared_ptr< GenericData > getPropertyByName(const std::string &propName)
Returns property.
Definition SimInfo.cpp:966
void setSnapshotManager(SnapshotManager *sman)
Sets the snapshot manager.
Definition SimInfo.cpp:971
bool addMolecule(Molecule *mol)
Adds a molecule.
Definition SimInfo.cpp:187
std::vector< std::string > getPropertyNames()
Returns all names of properties.
Definition SimInfo.cpp:958
void addInteractionPairs(Molecule *mol)
add all special interaction pairs (including excluded interactions) in a molecule into the appropriat...
Definition SimInfo.cpp:376
void update()
update
Definition SimInfo.cpp:697
unsigned int getNLocalCutoffGroups()
Returns the number of effective cutoff groups on local processor.
Definition SimInfo.cpp:311
Molecule * nextMolecule(MoleculeIterator &i)
Returns the next avaliable Molecule based on the iterator.
Definition SimInfo.cpp:245
StuntDouble * getIOIndexToIntegrableObject(int index)
return an integral objects by its global index.
Definition SimInfo.cpp:1027
void addProperty(std::shared_ptr< GenericData > genData)
Adds property into property map.
Definition SimInfo.cpp:950
std::vector< int > getGlobalAtomIndices()
returns a vector which maps the local atom index on this processor to the global atom index.
Definition SimInfo.cpp:862
void removeProperty(const std::string &propName)
Removes property from PropertyMap by name.
Definition SimInfo.cpp:954
AtomTypeSet getSimulatedAtomTypes()
Returns the set of atom types present in this simulation.
Definition SimInfo.cpp:712
SimInfo(ForceField *ff, Globals *simParams)
Constructor of SimInfo.
Definition SimInfo.cpp:79
SnapshotManager class is an abstract class which maintains a series of snapshots.
"Don't move, or you're dead! Stand up! Captain, we've got them!"
bool isLinear()
Tests the if this stuntDouble is a linear rigidbody.
bool isRigidBody()
Tests if this stuntDouble is a rigid body.
int getGlobalIndex()
Returns the global index of this stuntDouble.
bool isDirectional()
Tests if this stuntDouble is a directional one.
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.