ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/InteractionManager.cpp
Revision: 1877
Committed: Thu Jun 6 15:43:35 2013 UTC (11 years, 10 months ago) by gezelter
File size: 17739 byte(s)
Log Message:
New electrostatic method, starting to do some performance tuning.

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 #include "nonbonded/InteractionManager.hpp"
44
45 namespace OpenMD {
46
47 InteractionManager::InteractionManager() {
48
49 initialized_ = false;
50
51 lj_ = new LJ();
52 gb_ = new GB();
53 sticky_ = new Sticky();
54 morse_ = new Morse();
55 repulsivePower_ = new RepulsivePower();
56 eam_ = new EAM();
57 sc_ = new SC();
58 electrostatic_ = new Electrostatic();
59 maw_ = new MAW();
60 }
61
62 InteractionManager::~InteractionManager() {
63 delete lj_;
64 delete gb_;
65 delete sticky_;
66 delete morse_;
67 delete repulsivePower_;
68 delete eam_;
69 delete sc_;
70 delete electrostatic_;
71 delete maw_;
72 }
73
74 void InteractionManager::initialize() {
75
76 if (initialized_) return;
77
78 ForceField* forceField_ = info_->getForceField();
79
80 lj_->setForceField(forceField_);
81 gb_->setForceField(forceField_);
82 sticky_->setForceField(forceField_);
83 eam_->setForceField(forceField_);
84 sc_->setForceField(forceField_);
85 morse_->setForceField(forceField_);
86 electrostatic_->setSimInfo(info_);
87 electrostatic_->setForceField(forceField_);
88 maw_->setForceField(forceField_);
89 repulsivePower_->setForceField(forceField_);
90
91 ForceField::AtomTypeContainer* atomTypes = forceField_->getAtomTypes();
92 ForceField::AtomTypeContainer::MapTypeIterator i1, i2;
93 AtomType* atype1;
94 AtomType* atype2;
95 pair<AtomType*, AtomType*> key;
96
97 for (atype1 = atomTypes->beginType(i1); atype1 != NULL;
98 atype1 = atomTypes->nextType(i1)) {
99
100 // add it to the map:
101
102 pair<map<int,AtomType*>::iterator,bool> ret;
103 ret = typeMap_.insert( pair<int, AtomType*>(atype1->getIdent(), atype1) );
104 if (ret.second == false) {
105 sprintf( painCave.errMsg,
106 "InteractionManager already had a previous entry with ident %d\n",
107 atype1->getIdent());
108 painCave.severity = OPENMD_INFO;
109 painCave.isFatal = 0;
110 simError();
111 }
112 }
113
114 // Now, iterate over all known types and add to the interaction map:
115
116 map<int, AtomType*>::iterator it1, it2;
117 for (it1 = typeMap_.begin(); it1 != typeMap_.end(); ++it1) {
118 atype1 = (*it1).second;
119
120 for( it2 = typeMap_.begin(); it2 != typeMap_.end(); ++it2) {
121 atype2 = (*it2).second;
122
123 key = make_pair(atype1, atype2);
124
125 iHash_[key] = 0;
126
127 if (atype1->isLennardJones() && atype2->isLennardJones()) {
128 interactions_[key].insert(lj_);
129 iHash_[key] |= LJ_PAIR;
130 }
131 if (atype1->isElectrostatic() && atype2->isElectrostatic() ) {
132 interactions_[key].insert(electrostatic_);
133 iHash_[key] |= ELECTROSTATIC_PAIR;
134 }
135 if (atype1->isSticky() && atype2->isSticky() ) {
136 interactions_[key].insert(sticky_);
137 iHash_[key] |= STICKY_PAIR;
138 }
139 if (atype1->isStickyPower() && atype2->isStickyPower() ) {
140 interactions_[key].insert(sticky_);
141 iHash_[key] |= STICKY_PAIR;
142 }
143 if (atype1->isEAM() && atype2->isEAM() ) {
144 interactions_[key].insert(eam_);
145 iHash_[key] |= EAM_PAIR;
146 }
147 if (atype1->isSC() && atype2->isSC() ) {
148 interactions_[key].insert(sc_);
149 iHash_[key] |= SC_PAIR;
150 }
151 if (atype1->isGayBerne() && atype2->isGayBerne() ) {
152 interactions_[key].insert(gb_);
153 iHash_[key] |= GB_PAIR;
154 }
155 if ((atype1->isGayBerne() && atype2->isLennardJones())
156 || (atype1->isLennardJones() && atype2->isGayBerne())) {
157 interactions_[key].insert(gb_);
158 iHash_[key] |= GB_PAIR;
159 }
160
161 // look for an explicitly-set non-bonded interaction type using the
162 // two atom types.
163 NonBondedInteractionType* nbiType = forceField_->getNonBondedInteractionType(atype1->getName(), atype2->getName());
164
165 if (nbiType != NULL) {
166
167 bool vdwExplicit = false;
168 bool metExplicit = false;
169 // bool hbExplicit = false;
170
171 if (nbiType->isLennardJones()) {
172 // We found an explicit Lennard-Jones interaction.
173 // override all other vdw entries for this pair of atom types:
174 set<NonBondedInteraction*>::iterator it;
175 for (it = interactions_[key].begin();
176 it != interactions_[key].end(); ++it) {
177 InteractionFamily ifam = (*it)->getFamily();
178 if (ifam == VANDERWAALS_FAMILY) {
179 interactions_[key].erase(*it);
180 // work on iHash here;
181 }
182 }
183 interactions_[key].insert(lj_);
184 iHash_[key] |= LJ_PAIR;
185 vdwExplicit = true;
186 }
187
188 if (nbiType->isMorse()) {
189 if (vdwExplicit) {
190 sprintf( painCave.errMsg,
191 "InteractionManager::initialize found more than one "
192 "explicit \n"
193 "\tvan der Waals interaction for atom types %s - %s\n",
194 atype1->getName().c_str(), atype2->getName().c_str());
195 painCave.severity = OPENMD_ERROR;
196 painCave.isFatal = 1;
197 simError();
198 }
199 // We found an explicit Morse interaction.
200 // override all other vdw entries for this pair of atom types:
201 set<NonBondedInteraction*>::iterator it;
202 for (it = interactions_[key].begin();
203 it != interactions_[key].end(); ++it) {
204 InteractionFamily ifam = (*it)->getFamily();
205 if (ifam == VANDERWAALS_FAMILY) {
206 interactions_[key].erase(*it);
207 // work on iHash here;
208 }
209 }
210 interactions_[key].insert(morse_);
211 iHash_[key] |= MORSE_PAIR;
212 vdwExplicit = true;
213 }
214
215 if (nbiType->isRepulsivePower()) {
216 if (vdwExplicit) {
217 sprintf( painCave.errMsg,
218 "InteractionManager::initialize found more than one "
219 "explicit \n"
220 "\tvan der Waals interaction for atom types %s - %s\n",
221 atype1->getName().c_str(), atype2->getName().c_str());
222 painCave.severity = OPENMD_ERROR;
223 painCave.isFatal = 1;
224 simError();
225 }
226 // We found an explicit RepulsivePower interaction.
227 // override all other vdw entries for this pair of atom types:
228 set<NonBondedInteraction*>::iterator it;
229 for (it = interactions_[key].begin();
230 it != interactions_[key].end(); ++it) {
231 InteractionFamily ifam = (*it)->getFamily();
232 if (ifam == VANDERWAALS_FAMILY) {
233 interactions_[key].erase(*it);
234 // work on iHash here;
235 }
236 }
237 interactions_[key].insert(repulsivePower_);
238 iHash_[key] |= REPULSIVEPOWER_PAIR;
239 vdwExplicit = true;
240 }
241
242
243 if (nbiType->isEAM()) {
244 // We found an explicit EAM interaction.
245 // override all other metallic entries for this pair of atom types:
246 set<NonBondedInteraction*>::iterator it;
247 for (it = interactions_[key].begin();
248 it != interactions_[key].end(); ++it) {
249 InteractionFamily ifam = (*it)->getFamily();
250 if (ifam == METALLIC_FAMILY) {
251 interactions_[key].erase(*it);
252 // work on iHash here;
253 }
254 }
255 interactions_[key].insert(eam_);
256 iHash_[key] |= EAM_PAIR;
257 metExplicit = true;
258 }
259
260 if (nbiType->isSC()) {
261 if (metExplicit) {
262 sprintf( painCave.errMsg,
263 "InteractionManager::initialize found more than one "
264 "explicit\n"
265 "\tmetallic interaction for atom types %s - %s\n",
266 atype1->getName().c_str(), atype2->getName().c_str());
267 painCave.severity = OPENMD_ERROR;
268 painCave.isFatal = 1;
269 simError();
270 }
271 // We found an explicit Sutton-Chen interaction.
272 // override all other metallic entries for this pair of atom types:
273 set<NonBondedInteraction*>::iterator it;
274 for (it = interactions_[key].begin();
275 it != interactions_[key].end(); ++it) {
276 InteractionFamily ifam = (*it)->getFamily();
277 if (ifam == METALLIC_FAMILY) {
278 interactions_[key].erase(*it);
279 // work on iHash here;
280 }
281 }
282 interactions_[key].insert(sc_);
283 iHash_[key] |= SC_PAIR;
284 metExplicit = true;
285 }
286
287 if (nbiType->isMAW()) {
288 if (vdwExplicit) {
289 sprintf( painCave.errMsg,
290 "InteractionManager::initialize found more than one "
291 "explicit\n"
292 "\tvan der Waals interaction for atom types %s - %s\n",
293 atype1->getName().c_str(), atype2->getName().c_str());
294 painCave.severity = OPENMD_ERROR;
295 painCave.isFatal = 1;
296 simError();
297 }
298 // We found an explicit MAW interaction.
299 // override all other vdw entries for this pair of atom types:
300 set<NonBondedInteraction*>::iterator it;
301 for (it = interactions_[key].begin();
302 it != interactions_[key].end(); ++it) {
303 InteractionFamily ifam = (*it)->getFamily();
304 if (ifam == VANDERWAALS_FAMILY) {
305 interactions_[key].erase(*it);
306 // work on iHash here;
307 }
308 }
309 interactions_[key].insert(maw_);
310 iHash_[key] |= MAW_PAIR;
311 vdwExplicit = true;
312 }
313 }
314 }
315 }
316
317
318 // Make sure every pair of atom types in this simulation has a
319 // non-bonded interaction. If not, just inform the user.
320
321 set<AtomType*> simTypes = info_->getSimulatedAtomTypes();
322 set<AtomType*>::iterator it, jt;
323
324 for (it = simTypes.begin(); it != simTypes.end(); ++it) {
325 atype1 = (*it);
326 for (jt = it; jt != simTypes.end(); ++jt) {
327 atype2 = (*jt);
328 key = make_pair(atype1, atype2);
329
330 if (interactions_[key].size() == 0) {
331 sprintf( painCave.errMsg,
332 "InteractionManager could not find a matching non-bonded\n"
333 "\tinteraction for atom types %s - %s\n"
334 "\tProceeding without this interaction.\n",
335 atype1->getName().c_str(), atype2->getName().c_str());
336 painCave.severity = OPENMD_INFO;
337 painCave.isFatal = 0;
338 simError();
339 }
340 }
341 }
342
343 initialized_ = true;
344 }
345
346 void InteractionManager::setCutoffRadius(RealType rcut) {
347
348 electrostatic_->setCutoffRadius(rcut);
349 eam_->setCutoffRadius(rcut);
350 }
351
352 void InteractionManager::doPrePair(InteractionData idat){
353
354 if (!initialized_) initialize();
355
356 // excluded interaction, so just return
357 if (idat.excluded) return;
358
359 int& iHash = iHash_[idat.atypes];
360
361 if ((iHash & EAM_PAIR) != 0) eam_->calcDensity(idat);
362 if ((iHash & SC_PAIR) != 0) sc_->calcDensity(idat);
363
364 // set<NonBondedInteraction*>::iterator it;
365
366 // for (it = interactions_[ idat.atypes ].begin();
367 // it != interactions_[ idat.atypes ].end(); ++it){
368 // if ((*it)->getFamily() == METALLIC_FAMILY) {
369 // dynamic_cast<MetallicInteraction*>(*it)->calcDensity(idat);
370 // }
371 // }
372
373 return;
374 }
375
376 void InteractionManager::doPreForce(SelfData sdat){
377
378 if (!initialized_) initialize();
379
380 // pair<AtomType*, AtomType*> key = make_pair(sdat.atype, sdat.atype);
381
382 int& iHash = iHash_[ make_pair(sdat.atype, sdat.atype) ];
383
384 if ((iHash & EAM_PAIR) != 0) eam_->calcFunctional(sdat);
385 if ((iHash & SC_PAIR) != 0) sc_->calcFunctional(sdat);
386
387 // set<NonBondedInteraction*>::iterator it;
388
389 // for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
390 // if ((*it)->getFamily() == METALLIC_FAMILY) {
391 // dynamic_cast<MetallicInteraction*>(*it)->calcFunctional(sdat);
392 // }
393 // }
394
395 return;
396 }
397
398 void InteractionManager::doPair(InteractionData idat){
399
400 if (!initialized_) initialize();
401
402 int& iHash = iHash_[idat.atypes];
403
404 if ((iHash & ELECTROSTATIC_PAIR) != 0) electrostatic_->calcForce(idat);
405
406 // electrostatics still has to worry about indirect
407 // contributions from excluded pairs of atoms, but nothing else does:
408
409 if (idat.excluded) return;
410
411 if ((iHash & LJ_PAIR) != 0) lj_->calcForce(idat);
412 if ((iHash & GB_PAIR) != 0) gb_->calcForce(idat);
413 if ((iHash & STICKY_PAIR) != 0) sticky_->calcForce(idat);
414 if ((iHash & MORSE_PAIR) != 0) morse_->calcForce(idat);
415 if ((iHash & REPULSIVEPOWER_PAIR) != 0) repulsivePower_->calcForce(idat);
416 if ((iHash & EAM_PAIR) != 0) eam_->calcForce(idat);
417 if ((iHash & SC_PAIR) != 0) sc_->calcForce(idat);
418 if ((iHash & MAW_PAIR) != 0) maw_->calcForce(idat);
419
420 // set<NonBondedInteraction*>::iterator it;
421
422 // for (it = interactions_[ idat.atypes ].begin();
423 // it != interactions_[ idat.atypes ].end(); ++it) {
424
425 // // electrostatics still has to worry about indirect
426 // // contributions from excluded pairs of atoms:
427
428 // if (!idat.excluded || (*it)->getFamily() == ELECTROSTATIC_FAMILY) {
429 // (*it)->calcForce(idat);
430 // }
431 // }
432
433 return;
434 }
435
436 void InteractionManager::doSelfCorrection(SelfData sdat){
437
438 if (!initialized_) initialize();
439
440 int& iHash = iHash_[ make_pair(sdat.atype, sdat.atype) ];
441
442 if ((iHash & ELECTROSTATIC_PAIR) != 0) electrostatic_->calcSelfCorrection(sdat);
443
444
445 // pair<AtomType*, AtomType*> key = make_pair(sdat.atype, sdat.atype);
446 // set<NonBondedInteraction*>::iterator it;
447
448 // for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
449 // if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) {
450 // dynamic_cast<ElectrostaticInteraction*>(*it)->calcSelfCorrection(sdat);
451 // }
452 // }
453
454 return;
455 }
456
457 RealType InteractionManager::getSuggestedCutoffRadius(int *atid) {
458 if (!initialized_) initialize();
459
460 AtomType* atype = typeMap_[*atid];
461
462 pair<AtomType*, AtomType*> key = make_pair(atype, atype);
463 set<NonBondedInteraction*>::iterator it;
464 RealType cutoff = 0.0;
465
466 for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
467 cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(key));
468 return cutoff;
469 }
470
471 RealType InteractionManager::getSuggestedCutoffRadius(AtomType* atype) {
472 if (!initialized_) initialize();
473
474 pair<AtomType*, AtomType*> key = make_pair(atype, atype);
475 set<NonBondedInteraction*>::iterator it;
476 RealType cutoff = 0.0;
477
478 for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
479 cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(key));
480 return cutoff;
481 }
482 } //end namespace OpenMD

Properties

Name Value
svn:eol-style native
svn:executable *