ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/InteractionManager.cpp
Revision: 1868
Committed: Tue Apr 30 15:56:54 2013 UTC (12 years ago) by gezelter
File size: 15384 byte(s)
Log Message:
CLearing out some memory leaks

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

Properties

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