ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/InteractionManager.cpp
Revision: 1586
Committed: Tue Jun 21 06:34:35 2011 UTC (13 years, 10 months ago) by gezelter
File size: 14210 byte(s)
Log Message:
bug fixes

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

Properties

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