OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
InteractionManager.cpp
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#include "nonbonded/InteractionManager.hpp"
46
47#include "types/InversePowerSeriesInteractionType.hpp"
48#include "types/LennardJonesInteractionType.hpp"
49#include "types/MAWInteractionType.hpp"
50#include "types/MieInteractionType.hpp"
51#include "types/MorseInteractionType.hpp"
52#include "types/RepulsivePowerInteractionType.hpp"
53
54namespace OpenMD {
55
56 InteractionManager::InteractionManager() {
57 initialized_ = false;
58
59 lj_ = std::make_shared<LJ>();
60 gb_ = std::make_shared<GB>();
61 sticky_ = std::make_shared<Sticky>();
62 morse_ = std::make_shared<Morse>();
63 repulsivePower_ = std::make_shared<RepulsivePower>();
64 mie_ = std::make_shared<Mie>();
65 eam_ = std::make_shared<EAM>();
66 sc_ = std::make_shared<SC>();
67 electrostatic_ = std::make_shared<Electrostatic>();
68 maw_ = std::make_shared<MAW>();
69 inversePowerSeries_ = std::make_shared<InversePowerSeries>();
70 }
71
72 void InteractionManager::initialize() {
73 if (initialized_) return;
74
75 ForceField* forceField_ = info_->getForceField();
76
77 lj_->setForceField(forceField_);
78 gb_->setForceField(forceField_);
79 sticky_->setForceField(forceField_);
80 eam_->setForceField(forceField_);
81 eam_->setElectrostatic(electrostatic_.get());
82 sc_->setForceField(forceField_);
83 morse_->setForceField(forceField_);
84 electrostatic_->setSimInfo(info_);
85 electrostatic_->setForceField(forceField_);
86 maw_->setForceField(forceField_);
87 repulsivePower_->setForceField(forceField_);
88 mie_->setForceField(forceField_);
89 inversePowerSeries_->setForceField(forceField_);
90
91 ForceField::AtomTypeContainer* atomTypes = forceField_->getAtomTypes();
92 int nTypes = atomTypes->size();
93 sHash_.resize(nTypes);
94 iHash_.resize(nTypes);
95 interactions_.resize(nTypes);
96 ForceField::AtomTypeContainer::MapTypeIterator i1, i2;
97 AtomType* atype1;
98 AtomType* atype2;
99 int atid1, atid2;
100
101 // We only need to worry about the types that are actually in the
102 // simulation:
103
104 AtomTypeSet atypes = info_->getSimulatedAtomTypes();
105
106 lj_->setSimulatedAtomTypes(atypes);
107 gb_->setSimulatedAtomTypes(atypes);
108 sticky_->setSimulatedAtomTypes(atypes);
109 eam_->setSimulatedAtomTypes(atypes);
110 sc_->setSimulatedAtomTypes(atypes);
111 morse_->setSimulatedAtomTypes(atypes);
112 electrostatic_->setSimInfo(info_);
113 electrostatic_->setSimulatedAtomTypes(atypes);
114 maw_->setSimulatedAtomTypes(atypes);
115 repulsivePower_->setSimulatedAtomTypes(atypes);
116 mie_->setSimulatedAtomTypes(atypes);
117 inversePowerSeries_->setSimulatedAtomTypes(atypes);
118
119 AtomTypeSet::iterator at;
120 set<NonBondedInteractionPtr>::iterator it;
121
122 for (at = atypes.begin(); at != atypes.end(); ++at) {
123 atype1 = *at;
124 atid1 = atype1->getIdent();
125 iHash_[atid1].resize(nTypes);
126 interactions_[atid1].resize(nTypes);
127
128 // add it to the map:
129 pair<map<int, AtomType*>::iterator, bool> ret;
130 ret = typeMap_.insert(pair<int, AtomType*>(atid1, atype1));
131 if (ret.second == false) {
132 snprintf(
133 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
134 "InteractionManager already had a previous entry with ident %d\n",
135 atype1->getIdent());
136 painCave.severity = OPENMD_INFO;
137 painCave.isFatal = 0;
138 simError();
139 }
140
141 if (atype1->isLennardJones()) { sHash_[atid1] |= LJ_INTERACTION; }
142 if (atype1->isElectrostatic()) {
143 sHash_[atid1] |= ELECTROSTATIC_INTERACTION;
144 }
145 if (atype1->isSticky()) { sHash_[atid1] |= STICKY_INTERACTION; }
146 if (atype1->isStickyPower()) { sHash_[atid1] |= STICKY_INTERACTION; }
147 if (atype1->isEAM()) { sHash_[atid1] |= EAM_INTERACTION; }
148 if (atype1->isSC()) { sHash_[atid1] |= SC_INTERACTION; }
149 if (atype1->isGayBerne()) { sHash_[atid1] |= GB_INTERACTION; }
150 }
151 // Now, iterate over all known types and add to the interaction map:
152
153 map<int, AtomType*>::iterator it1, it2;
154 for (it1 = typeMap_.begin(); it1 != typeMap_.end(); ++it1) {
155 atype1 = (*it1).second;
156 atid1 = atype1->getIdent();
157
158 for (it2 = typeMap_.begin(); it2 != typeMap_.end(); ++it2) {
159 atype2 = (*it2).second;
160 atid2 = atype2->getIdent();
161
162 iHash_[atid1][atid2] = 0;
163
164 if (atype1->isLennardJones() && atype2->isLennardJones()) {
165 interactions_[atid1][atid2].insert(lj_);
166 iHash_[atid1][atid2] |= LJ_INTERACTION;
167 }
168 if (atype1->isElectrostatic() && atype2->isElectrostatic()) {
169 interactions_[atid1][atid2].insert(electrostatic_);
170 iHash_[atid1][atid2] |= ELECTROSTATIC_INTERACTION;
171 }
172
173 // A special case for calculating local fields:
174 if (info_->getSimParams()->getOutputElectricField()) {
175 if (atype1->isElectrostatic() || atype2->isElectrostatic()) {
176 interactions_[atid1][atid2].insert(electrostatic_);
177 iHash_[atid1][atid2] |= ELECTROSTATIC_INTERACTION;
178 }
179 }
180
181 if (atype1->isSticky() && atype2->isSticky()) {
182 interactions_[atid1][atid2].insert(sticky_);
183 iHash_[atid1][atid2] |= STICKY_INTERACTION;
184 }
185 if (atype1->isStickyPower() && atype2->isStickyPower()) {
186 interactions_[atid1][atid2].insert(sticky_);
187 iHash_[atid1][atid2] |= STICKY_INTERACTION;
188 }
189 if (atype1->isEAM() && atype2->isEAM()) {
190 interactions_[atid1][atid2].insert(eam_);
191 iHash_[atid1][atid2] |= EAM_INTERACTION;
192 }
193 if (atype1->isSC() && atype2->isSC()) {
194 interactions_[atid1][atid2].insert(sc_);
195 iHash_[atid1][atid2] |= SC_INTERACTION;
196 }
197 if (atype1->isGayBerne() && atype2->isGayBerne()) {
198 interactions_[atid1][atid2].insert(gb_);
199 iHash_[atid1][atid2] |= GB_INTERACTION;
200 }
201 if ((atype1->isGayBerne() && atype2->isLennardJones()) ||
202 (atype1->isLennardJones() && atype2->isGayBerne())) {
203 interactions_[atid1][atid2].insert(gb_);
204 iHash_[atid1][atid2] |= GB_INTERACTION;
205 }
206
207 // look for an explicitly-set non-bonded interaction type using the
208 // two atom types.
209 NonBondedInteractionType* nbiType =
210 forceField_->getNonBondedInteractionType(atype1->getName(),
211 atype2->getName());
212
213 if (nbiType != NULL) {
214 bool vdwExplicit = false;
215 bool metExplicit = false;
216 // bool hbExplicit = false;
217
218 if (nbiType->isLennardJones()) {
219 // We found an explicit Lennard-Jones interaction.
220 // override all other vdw entries for this pair of atom types:
221 for (it = interactions_[atid1][atid2].begin();
222 it != interactions_[atid1][atid2].end();) {
223 InteractionFamily ifam = (*it)->getFamily();
224 if (ifam == VANDERWAALS_FAMILY) {
225 iHash_[atid1][atid2] ^= (*it)->getHash();
226 interactions_[atid1][atid2].erase(it++);
227 } else {
228 ++it;
229 }
230 }
231 interactions_[atid1][atid2].insert(lj_);
232 iHash_[atid1][atid2] |= LJ_INTERACTION;
234 dynamic_cast<LennardJonesInteractionType*>(nbiType);
235 lj_->addExplicitInteraction(atype1, atype2, ljit->getSigma(),
236 ljit->getEpsilon());
237 vdwExplicit = true;
238 }
239
240 if (nbiType->isMorse()) {
241 if (vdwExplicit) {
242 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
243 "InteractionManager::initialize found more than one "
244 "explicit \n"
245 "\tvan der Waals interaction for atom types %s - %s\n",
246 atype1->getName().c_str(), atype2->getName().c_str());
247 painCave.severity = OPENMD_ERROR;
248 painCave.isFatal = 1;
249 simError();
250 }
251 // We found an explicit Morse interaction.
252 // override all other vdw entries for this pair of atom types:
253 for (it = interactions_[atid1][atid2].begin();
254 it != interactions_[atid1][atid2].end();) {
255 InteractionFamily ifam = (*it)->getFamily();
256 if (ifam == VANDERWAALS_FAMILY) {
257 iHash_[atid1][atid2] ^= (*it)->getHash();
258 interactions_[atid1][atid2].erase(it++);
259 } else {
260 ++it;
261 }
262 }
263 interactions_[atid1][atid2].insert(morse_);
264 iHash_[atid1][atid2] |= MORSE_INTERACTION;
266 dynamic_cast<MorseInteractionType*>(nbiType);
267 morse_->addExplicitInteraction(atype1, atype2, mit->getD(),
268 mit->getR(), mit->getBeta(),
269 mit->getInteractionType());
270 vdwExplicit = true;
271 }
272
273 if (nbiType->isRepulsivePower()) {
274 if (vdwExplicit) {
275 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
276 "InteractionManager::initialize found more than one "
277 "explicit \n"
278 "\tvan der Waals interaction for atom types %s - %s\n",
279 atype1->getName().c_str(), atype2->getName().c_str());
280 painCave.severity = OPENMD_ERROR;
281 painCave.isFatal = 1;
282 simError();
283 }
284 // We found an explicit RepulsivePower interaction.
285 // override all other vdw entries for this pair of atom types:
286 for (it = interactions_[atid1][atid2].begin();
287 it != interactions_[atid1][atid2].end();) {
288 InteractionFamily ifam = (*it)->getFamily();
289 if (ifam == VANDERWAALS_FAMILY) {
290 iHash_[atid1][atid2] ^= (*it)->getHash();
291 interactions_[atid1][atid2].erase(it++);
292 } else {
293 ++it;
294 }
295 }
296 interactions_[atid1][atid2].insert(repulsivePower_);
297 iHash_[atid1][atid2] |= REPULSIVEPOWER_INTERACTION;
299 dynamic_cast<RepulsivePowerInteractionType*>(nbiType);
300
301 repulsivePower_->addExplicitInteraction(
302 atype1, atype2, rpit->getSigma(), rpit->getEpsilon(),
303 rpit->getNrep());
304
305 vdwExplicit = true;
306 }
307
308 if (nbiType->isMie()) {
309 if (vdwExplicit) {
310 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
311 "InteractionManager::initialize found more than one "
312 "explicit \n"
313 "\tvan der Waals interaction for atom types %s - %s\n",
314 atype1->getName().c_str(), atype2->getName().c_str());
315 painCave.severity = OPENMD_ERROR;
316 painCave.isFatal = 1;
317 simError();
318 }
319 // We found an explicit Mie interaction.
320 // override all other vdw entries for this pair of atom types:
321 for (it = interactions_[atid1][atid2].begin();
322 it != interactions_[atid1][atid2].end();) {
323 InteractionFamily ifam = (*it)->getFamily();
324 if (ifam == VANDERWAALS_FAMILY) {
325 iHash_[atid1][atid2] ^= (*it)->getHash();
326 interactions_[atid1][atid2].erase(it++);
327 } else {
328 ++it;
329 }
330 }
331 interactions_[atid1][atid2].insert(mie_);
332 iHash_[atid1][atid2] |= MIE_INTERACTION;
333 MieInteractionType* mit =
334 dynamic_cast<MieInteractionType*>(nbiType);
335
336 mie_->addExplicitInteraction(atype1, atype2, mit->getSigma(),
337 mit->getEpsilon(), mit->getNrep(),
338 mit->getMatt());
339
340 vdwExplicit = true;
341 }
342
343 if (nbiType->isEAMTable() || nbiType->isEAMZhou()) {
344 // We found an explicit EAM interaction.
345 // override all other metallic entries for this pair of atom types:
346 for (it = interactions_[atid1][atid2].begin();
347 it != interactions_[atid1][atid2].end();) {
348 InteractionFamily ifam = (*it)->getFamily();
349 if (ifam == METALLIC_EMBEDDING_FAMILY) {
350 iHash_[atid1][atid2] ^= (*it)->getHash();
351 interactions_[atid1][atid2].erase(it++);
352 } else {
353 ++it;
354 }
355 }
356 interactions_[atid1][atid2].insert(eam_);
357 iHash_[atid1][atid2] |= EAM_INTERACTION;
358 metExplicit = true;
359 }
360
361 if (nbiType->isSC()) {
362 if (metExplicit) {
363 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
364 "InteractionManager::initialize found more than one "
365 "explicit\n"
366 "\tmetallic interaction for atom types %s - %s\n",
367 atype1->getName().c_str(), atype2->getName().c_str());
368 painCave.severity = OPENMD_ERROR;
369 painCave.isFatal = 1;
370 simError();
371 }
372 // We found an explicit Sutton-Chen interaction.
373 // override all other metallic entries for this pair of atom types:
374 for (it = interactions_[atid1][atid2].begin();
375 it != interactions_[atid1][atid2].end();) {
376 InteractionFamily ifam = (*it)->getFamily();
377 if (ifam == METALLIC_EMBEDDING_FAMILY) {
378 iHash_[atid1][atid2] ^= (*it)->getHash();
379 interactions_[atid1][atid2].erase(it++);
380 } else {
381 ++it;
382 }
383 }
384 interactions_[atid1][atid2].insert(sc_);
385 iHash_[atid1][atid2] |= SC_INTERACTION;
386 metExplicit = true;
387 }
388
389 if (nbiType->isMAW()) {
390 if (vdwExplicit) {
391 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
392 "InteractionManager::initialize found more than one "
393 "explicit\n"
394 "\tvan der Waals interaction for atom types %s - %s\n",
395 atype1->getName().c_str(), atype2->getName().c_str());
396 painCave.severity = OPENMD_ERROR;
397 painCave.isFatal = 1;
398 simError();
399 }
400 // We found an explicit MAW interaction.
401 // override all other vdw entries for this pair of atom types:
402 for (it = interactions_[atid1][atid2].begin();
403 it != interactions_[atid1][atid2].end();) {
404 InteractionFamily ifam = (*it)->getFamily();
405 if (ifam == VANDERWAALS_FAMILY) {
406 iHash_[atid1][atid2] ^= (*it)->getHash();
407 interactions_[atid1][atid2].erase(it++);
408 } else {
409 ++it;
410 }
411 }
412 interactions_[atid1][atid2].insert(maw_);
413 iHash_[atid1][atid2] |= MAW_INTERACTION;
414 MAWInteractionType* mit =
415 dynamic_cast<MAWInteractionType*>(nbiType);
416 maw_->addExplicitInteraction(atype1, atype2, mit->getD(),
417 mit->getBeta(), mit->getR(),
418 mit->getCA1(), mit->getCB1());
419 vdwExplicit = true;
420 }
421
422 if (nbiType->isInversePowerSeries()) {
423 if (vdwExplicit) {
424 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
425 "InteractionManager::initialize found more than one "
426 "explicit \n"
427 "\tvan der Waals interaction for atom types %s - %s\n",
428 atype1->getName().c_str(), atype2->getName().c_str());
429 painCave.severity = OPENMD_ERROR;
430 painCave.isFatal = 1;
431 simError();
432 }
433 // We found an explicit InversePowerSeries interaction.
434 // override all other vdw entries for this pair of atom types:
435 for (it = interactions_[atid1][atid2].begin();
436 it != interactions_[atid1][atid2].end();) {
437 InteractionFamily ifam = (*it)->getFamily();
438 if (ifam == VANDERWAALS_FAMILY) {
439 iHash_[atid1][atid2] ^= (*it)->getHash();
440 interactions_[atid1][atid2].erase(it++);
441 } else {
442 ++it;
443 }
444 }
445 interactions_[atid1][atid2].insert(inversePowerSeries_);
446 iHash_[atid1][atid2] |= INVERSEPOWERSERIES_INTERACTION;
448 dynamic_cast<InversePowerSeriesInteractionType*>(nbiType);
449
450 inversePowerSeries_->addExplicitInteraction(
451 atype1, atype2, ipsit->getPowers(), ipsit->getCoefficients());
452 vdwExplicit = true;
453 }
454 }
455 }
456 }
457
458 // Make sure every pair of atom types in this simulation has a
459 // non-bonded interaction. If not, just inform the user.
460
461 AtomTypeSet simTypes = info_->getSimulatedAtomTypes();
462 AtomTypeSet::iterator bt;
463
464 for (at = simTypes.begin(); at != simTypes.end(); ++at) {
465 atype1 = (*at);
466 atid1 = atype1->getIdent();
467 for (bt = at; bt != simTypes.end(); ++bt) {
468 atype2 = (*bt);
469 atid2 = atype2->getIdent();
470
471 if (interactions_[atid1][atid2].size() == 0) {
472 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
473 "InteractionManager could not find a matching non-bonded\n"
474 "\tinteraction for atom types %s - %s\n"
475 "\tProceeding without this interaction.\n",
476 atype1->getName().c_str(), atype2->getName().c_str());
477 painCave.severity = OPENMD_INFO;
478 painCave.isFatal = 0;
479 simError();
480 }
481 }
482 }
483
484 initialized_ = true;
485 }
486
487 void InteractionManager::setCutoffRadius(RealType rcut) {
488 electrostatic_->setCutoffRadius(rcut);
489 eam_->setCutoffRadius(rcut);
490 }
491
492 void InteractionManager::doPrePair(InteractionData& idat) {
493 if (!initialized_) initialize();
494
495 // excluded interaction, so just return
496 if (idat.excluded) return;
497
498 int& iHash = iHash_[idat.atid1][idat.atid2];
499
500 if ((iHash & EAM_INTERACTION) != 0) eam_->calcDensity(idat);
501 if ((iHash & SC_INTERACTION) != 0) sc_->calcDensity(idat);
502
503 // set<NonBondedInteraction*>::iterator it;
504 //
505 // for (it = interactions_[ idat.atypes ].begin();
506 // it != interactions_[ idat.atypes ].end(); ++it){
507 // if ((*it)->getFamily() == METALLIC_EMBEDDING_FAMILY) {
508 // dynamic_cast<MetallicInteraction*>(*it)->calcDensity(idat);
509 // }
510 // }
511
512 return;
513 }
514
515 void InteractionManager::doPreForce(SelfData& sdat) {
516 if (!initialized_) initialize();
517
518 int& sHash = sHash_[sdat.atid];
519
520 if ((sHash & EAM_INTERACTION) != 0) eam_->calcFunctional(sdat);
521 if ((sHash & SC_INTERACTION) != 0) sc_->calcFunctional(sdat);
522
523 // set<NonBondedInteraction*>::iterator it;
524 //
525 // for (it = interactions_[atid1][atid2].begin();
526 // it != interactions_[atid1][atid2].end(); ++it){
527 // if ((*it)->getFamily() == METALLIC_EMBEDDING_FAMILY) {
528 // dynamic_cast<MetallicInteraction*>(*it)->calcFunctional(sdat);
529 // }
530 // }
531
532 return;
533 }
534
535 void InteractionManager::doPair(InteractionData& idat) {
536 if (!initialized_) initialize();
537
538 int& iHash = iHash_[idat.atid1][idat.atid2];
539
540 if ((iHash & ELECTROSTATIC_INTERACTION) != 0)
541 electrostatic_->calcForce(idat);
542
543 // electrostatics still has to worry about indirect
544 // contributions from excluded pairs of atoms, but nothing else does:
545
546 if (idat.excluded) return;
547
548 if ((iHash & LJ_INTERACTION) != 0) lj_->calcForce(idat);
549 if ((iHash & GB_INTERACTION) != 0) gb_->calcForce(idat);
550 if ((iHash & STICKY_INTERACTION) != 0) sticky_->calcForce(idat);
551 if ((iHash & MORSE_INTERACTION) != 0) morse_->calcForce(idat);
552 if ((iHash & REPULSIVEPOWER_INTERACTION) != 0)
553 repulsivePower_->calcForce(idat);
554 if ((iHash & MIE_INTERACTION) != 0) mie_->calcForce(idat);
555 if ((iHash & EAM_INTERACTION) != 0) eam_->calcForce(idat);
556 if ((iHash & SC_INTERACTION) != 0) sc_->calcForce(idat);
557 if ((iHash & MAW_INTERACTION) != 0) maw_->calcForce(idat);
558 if ((iHash & INVERSEPOWERSERIES_INTERACTION) != 0)
559 inversePowerSeries_->calcForce(idat);
560
561 // set<NonBondedInteraction*>::iterator it;
562 //
563 // for (it = interactions_[ idat.atypes ].begin();
564 // it != interactions_[ idat.atypes ].end(); ++it) {
565 //
566 // if (!idat.excluded || (*it)->getFamily() == ELECTROSTATIC_FAMILY) {
567 // (*it)->calcForce(idat);
568 // }
569 // }
570
571 return;
572 }
573
574 void InteractionManager::doSelfCorrection(SelfData& sdat) {
575 if (!initialized_) initialize();
576
577 int& sHash = sHash_[sdat.atid];
578
579 if ((sHash & ELECTROSTATIC_INTERACTION) != 0) {
580 electrostatic_->calcSelfCorrection(sdat);
581 }
582
583 // set<NonBondedInteraction*>::iterator it;
584 //
585 // for (it = interactions_[atid1][atid2].begin();
586 // it != interactions_[atid1][atid2].end(); ++it){
587 // if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) {
588 // dynamic_cast<ElectrostaticInteraction*>(*it)->calcSelfCorrection(sdat);
589 // }
590 // }
591
592 return;
593 }
594
595 void InteractionManager::doSurfaceTerm(bool slabGeometry, int axis,
596 RealType& pot) {
597 if (!initialized_) initialize();
598 electrostatic_->calcSurfaceTerm(slabGeometry, axis, pot);
599 }
600
601 void InteractionManager::doReciprocalSpaceSum(RealType& pot) {
602 if (!initialized_) initialize();
603 electrostatic_->ReciprocalSpaceSum(pot);
604 }
605
606 RealType InteractionManager::getSuggestedCutoffRadius(int* atid) {
607 if (!initialized_) initialize();
608
609 AtomType* atype = typeMap_[*atid];
610
611 set<NonBondedInteractionPtr>::iterator it;
612 RealType cutoff = 0.0;
613
614 for (it = interactions_[*atid][*atid].begin();
615 it != interactions_[*atid][*atid].end(); ++it) {
616 cutoff =
617 max(cutoff, (*it)->getSuggestedCutoffRadius(make_pair(atype, atype)));
618 }
619 return cutoff;
620 }
621
622 RealType InteractionManager::getSuggestedCutoffRadius(AtomType* atype) {
623 if (!initialized_) initialize();
624
625 int atid = atype->getIdent();
626
627 set<NonBondedInteractionPtr>::iterator it;
628 RealType cutoff = 0.0;
629
630 for (it = interactions_[atid][atid].begin();
631 it != interactions_[atid][atid].end(); ++it) {
632 cutoff =
633 max(cutoff, (*it)->getSuggestedCutoffRadius(make_pair(atype, atype)));
634 }
635 return cutoff;
636 }
637} // namespace OpenMD
AtomType is what OpenMD looks to for unchanging data about an atom.
Definition AtomType.hpp:66
InversePowerSeriesInteractionType is a sum of powers in the inverse of r:
LennardJonesInteractionType is one of the basic interaction types.
MAWInteractionType (Metal-Angular-Water) is one of the basic Metal-to-NonMetal interaction types.
MieInteractionType is one of the basic interaction types.
MorseInteractionType is one of the basic non-bonded interactions.
NonBondedInteractionType class is responsible for keeping track of static (unchanging) parameters for...
RepulsivePowerInteractionType is one of the basic interaction types.
ForceField * getForceField()
Returns the force field.
Definition SimInfo.hpp:266
AtomTypeSet getSimulatedAtomTypes()
Returns the set of atom types present in this simulation.
Definition SimInfo.cpp:712
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
static const int ELECTROSTATIC_INTERACTION
Boolean flags for the iHash_ and sHash_ data structures.
InteractionFamily
The InteractionFamily enum.
@ VANDERWAALS_FAMILY
Long-range dispersion and short-range pauli repulsion.
@ METALLIC_EMBEDDING_FAMILY
Transition metal interactions involving electron density.
The InteractionData struct.
int atid1
atomType ident for atom 1
bool excluded
is this excluded from direct pairwise interactions? (some indirect interactions may still apply)
int atid2
atomType ident for atom 2
The SelfData struct.
int atid
atomType ident for the atom