ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/InteractionManager.cpp
Revision: 1874
Committed: Wed May 15 15:09:35 2013 UTC (11 years, 11 months ago) by gezelter
File size: 15332 byte(s)
Log Message:
Fixed a bunch of cppcheck warnings.

File Contents

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

Properties

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