ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/nonbonded/InteractionManager.cpp
Revision: 1793
Committed: Fri Aug 31 21:16:10 2012 UTC (12 years, 8 months ago) by gezelter
File size: 15242 byte(s)
Log Message:
Cleaning up some warning messages on linux

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

Properties

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