ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/InteractionManager.cpp
Revision: 1656
Committed: Tue Oct 11 19:46:51 2011 UTC (13 years, 6 months ago) by jmichalk
File size: 15220 byte(s)
Log Message:
First few fixes for explicit interactions.

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

Properties

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