ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/InteractionManager.cpp
Revision: 1825
Committed: Wed Jan 9 19:27:52 2013 UTC (12 years, 3 months ago) by gezelter
File size: 15156 byte(s)
Log Message:
Deleting unused variables

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     * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (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     void InteractionManager::initialize() {
63 gezelter 1755
64     if (initialized_) return;
65    
66 gezelter 1535 ForceField* forceField_ = info_->getForceField();
67    
68 gezelter 1502 lj_->setForceField(forceField_);
69     gb_->setForceField(forceField_);
70     sticky_->setForceField(forceField_);
71     eam_->setForceField(forceField_);
72     sc_->setForceField(forceField_);
73     morse_->setForceField(forceField_);
74 gezelter 1584 electrostatic_->setSimInfo(info_);
75 gezelter 1502 electrostatic_->setForceField(forceField_);
76 gezelter 1532 maw_->setForceField(forceField_);
77 jmichalk 1656 repulsivePower_->setForceField(forceField_);
78 gezelter 1502
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 gezelter 1710 ret = typeMap_.insert( pair<int, AtomType*>(atype1->getIdent(), atype1) );
93 gezelter 1502 if (ret.second == false) {
94     sprintf( painCave.errMsg,
95     "InteractionManager already had a previous entry with ident %d\n",
96 gezelter 1710 atype1->getIdent());
97 gezelter 1502 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 gezelter 1825 // bool hbExplicit = false;
115 gezelter 1502
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 gezelter 1535
148     if (nbiType != NULL) {
149 gezelter 1502
150 gezelter 1535 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 gezelter 1502 }
162 gezelter 1535
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 gezelter 1502 }
185 jmichalk 1656
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 gezelter 1535
209 jmichalk 1656
210 gezelter 1535 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 gezelter 1502 }
222 gezelter 1535
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 gezelter 1502 }
245 gezelter 1535
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 gezelter 1502 }
269     }
270     }
271    
272 gezelter 1535
273 gezelter 1583 // Make sure every pair of atom types in this simulation has a
274     // non-bonded interaction. If not, just inform the user.
275 gezelter 1535
276     set<AtomType*> simTypes = info_->getSimulatedAtomTypes();
277     set<AtomType*>::iterator it, jt;
278 gezelter 1583
279 gezelter 1535 for (it = simTypes.begin(); it != simTypes.end(); ++it) {
280     atype1 = (*it);
281 gezelter 1583 for (jt = it; jt != simTypes.end(); ++jt) {
282 gezelter 1535 atype2 = (*jt);
283 gezelter 1502 key = make_pair(atype1, atype2);
284    
285     if (interactions_[key].size() == 0) {
286     sprintf( painCave.errMsg,
287 gezelter 1583 "InteractionManager could not find a matching non-bonded\n"
288     "\tinteraction for atom types %s - %s\n"
289     "\tProceeding without this interaction.\n",
290 gezelter 1502 atype1->getName().c_str(), atype2->getName().c_str());
291     painCave.severity = OPENMD_INFO;
292 gezelter 1583 painCave.isFatal = 0;
293 gezelter 1502 simError();
294     }
295     }
296     }
297 gezelter 1535
298     initialized_ = true;
299     }
300 gezelter 1584
301     void InteractionManager::setCutoffRadius(RealType rcut) {
302 gezelter 1587
303 gezelter 1584 electrostatic_->setCutoffRadius(rcut);
304 gezelter 1586 eam_->setCutoffRadius(rcut);
305 gezelter 1584 }
306    
307 gezelter 1545 void InteractionManager::doPrePair(InteractionData idat){
308 gezelter 1502
309     if (!initialized_) initialize();
310 gezelter 1545
311 gezelter 1587 // excluded interaction, so just return
312     if (idat.excluded) return;
313    
314 gezelter 1504 set<NonBondedInteraction*>::iterator it;
315    
316 gezelter 1571 for (it = interactions_[ idat.atypes ].begin();
317     it != interactions_[ idat.atypes ].end(); ++it){
318 gezelter 1504 if ((*it)->getFamily() == METALLIC_FAMILY) {
319 gezelter 1545 dynamic_cast<MetallicInteraction*>(*it)->calcDensity(idat);
320 gezelter 1504 }
321     }
322 gezelter 1502
323     return;
324     }
325 gezelter 1504
326 gezelter 1545 void InteractionManager::doPreForce(SelfData sdat){
327 gezelter 1502
328 gezelter 1504 if (!initialized_) initialize();
329 gezelter 1545
330     pair<AtomType*, AtomType*> key = make_pair(sdat.atype, sdat.atype);
331 gezelter 1504 set<NonBondedInteraction*>::iterator it;
332 gezelter 1502
333 gezelter 1504 for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
334     if ((*it)->getFamily() == METALLIC_FAMILY) {
335 gezelter 1545 dynamic_cast<MetallicInteraction*>(*it)->calcFunctional(sdat);
336 gezelter 1504 }
337     }
338    
339 gezelter 1502 return;
340     }
341    
342 gezelter 1545 void InteractionManager::doPair(InteractionData idat){
343 gezelter 1504
344 gezelter 1502 if (!initialized_) initialize();
345 gezelter 1587
346 gezelter 1502 set<NonBondedInteraction*>::iterator it;
347    
348 gezelter 1571 for (it = interactions_[ idat.atypes ].begin();
349 gezelter 1587 it != interactions_[ idat.atypes ].end(); ++it) {
350 gezelter 1502
351 gezelter 1587 // electrostatics still has to worry about indirect
352     // contributions from excluded pairs of atoms:
353 gezelter 1502
354 gezelter 1587 if (!idat.excluded || (*it)->getFamily() == ELECTROSTATIC_FAMILY) {
355     (*it)->calcForce(idat);
356 gezelter 1504 }
357     }
358 gezelter 1502
359 gezelter 1504 return;
360 gezelter 1502 }
361    
362 gezelter 1545 void InteractionManager::doSelfCorrection(SelfData sdat){
363 gezelter 1502
364     if (!initialized_) initialize();
365 gezelter 1504
366 gezelter 1545 pair<AtomType*, AtomType*> key = make_pair(sdat.atype, sdat.atype);
367 gezelter 1504 set<NonBondedInteraction*>::iterator it;
368 gezelter 1502
369 gezelter 1504 for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
370     if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) {
371 gezelter 1545 dynamic_cast<ElectrostaticInteraction*>(*it)->calcSelfCorrection(sdat);
372 gezelter 1504 }
373     }
374 gezelter 1536
375 gezelter 1504 return;
376 gezelter 1502 }
377    
378 gezelter 1505 RealType InteractionManager::getSuggestedCutoffRadius(int *atid) {
379     if (!initialized_) initialize();
380 gezelter 1755
381 gezelter 1505 AtomType* atype = typeMap_[*atid];
382    
383     pair<AtomType*, AtomType*> key = make_pair(atype, atype);
384     set<NonBondedInteraction*>::iterator it;
385     RealType cutoff = 0.0;
386    
387     for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
388 gezelter 1545 cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(key));
389 gezelter 1505 return cutoff;
390     }
391    
392 gezelter 1528 RealType InteractionManager::getSuggestedCutoffRadius(AtomType* atype) {
393     if (!initialized_) initialize();
394    
395     pair<AtomType*, AtomType*> key = make_pair(atype, atype);
396     set<NonBondedInteraction*>::iterator it;
397     RealType cutoff = 0.0;
398    
399     for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
400 gezelter 1545 cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(key));
401 gezelter 1528 return cutoff;
402     }
403 gezelter 1502 } //end namespace OpenMD

Properties

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