ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/InteractionManager.cpp
Revision: 1755
Committed: Thu Jun 14 01:58:35 2012 UTC (12 years, 10 months ago) by gezelter
File size: 15275 byte(s)
Log Message:
general bugfixes

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

Properties

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