ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/nonbonded/InteractionManager.cpp
Revision: 2033
Committed: Sat Nov 1 14:12:16 2014 UTC (10 years, 6 months ago) by gezelter
File size: 20107 byte(s)
Log Message:
Fixed a selection issue in ParticleTimeCorrFunc, other whitespace stuff

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 1879 * [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 1879 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 gezelter 1895 int nTypes = atomTypes->size();
93 gezelter 1927 sHash_.resize(nTypes);
94 gezelter 1895 iHash_.resize(nTypes);
95     interactions_.resize(nTypes);
96 gezelter 1502 ForceField::AtomTypeContainer::MapTypeIterator i1, i2;
97     AtomType* atype1;
98     AtomType* atype2;
99 gezelter 1895 int atid1, atid2;
100    
101 gezelter 2033 // We only need to worry about the types that are actually in the
102     // simulation:
103    
104 gezelter 1895 set<AtomType*> atypes = info_->getSimulatedAtomTypes();
105    
106     lj_->setSimulatedAtomTypes(atypes);
107     gb_->setSimulatedAtomTypes(atypes);
108     sticky_->setSimulatedAtomTypes(atypes);
109     eam_->setSimulatedAtomTypes(atypes);
110     sc_->setSimulatedAtomTypes(atypes);
111     morse_->setSimulatedAtomTypes(atypes);
112     electrostatic_->setSimInfo(info_);
113     electrostatic_->setSimulatedAtomTypes(atypes);
114     maw_->setSimulatedAtomTypes(atypes);
115     repulsivePower_->setSimulatedAtomTypes(atypes);
116    
117     set<AtomType*>::iterator at;
118    
119     for (at = atypes.begin(); at != atypes.end(); ++at) {
120 gezelter 1935
121 gezelter 1895 atype1 = *at;
122     atid1 = atype1->getIdent();
123     iHash_[atid1].resize(nTypes);
124     interactions_[atid1].resize(nTypes);
125    
126     // add it to the map:
127 gezelter 1502 pair<map<int,AtomType*>::iterator,bool> ret;
128 gezelter 1895 ret = typeMap_.insert( pair<int, AtomType*>(atid1, atype1) );
129 gezelter 1502 if (ret.second == false) {
130     sprintf( painCave.errMsg,
131     "InteractionManager already had a previous entry with ident %d\n",
132 gezelter 1710 atype1->getIdent());
133 gezelter 1502 painCave.severity = OPENMD_INFO;
134     painCave.isFatal = 0;
135     simError();
136     }
137 gezelter 1935
138     if (atype1->isLennardJones()) {
139     sHash_[atid1] |= LJ_INTERACTION;
140     }
141     if (atype1->isElectrostatic()) {
142     sHash_[atid1] |= ELECTROSTATIC_INTERACTION;
143     }
144     if (atype1->isSticky()) {
145     sHash_[atid1] |= STICKY_INTERACTION;
146     }
147     if (atype1->isStickyPower()) {
148     sHash_[atid1] |= STICKY_INTERACTION;
149     }
150     if (atype1->isEAM()) {
151     sHash_[atid1] |= EAM_INTERACTION;
152     }
153     if (atype1->isSC()) {
154     sHash_[atid1] |= SC_INTERACTION;
155     }
156     if (atype1->isGayBerne()) {
157     sHash_[atid1] |= GB_INTERACTION;
158     }
159 gezelter 1502 }
160     // Now, iterate over all known types and add to the interaction map:
161    
162     map<int, AtomType*>::iterator it1, it2;
163     for (it1 = typeMap_.begin(); it1 != typeMap_.end(); ++it1) {
164     atype1 = (*it1).second;
165 gezelter 1895 atid1 = atype1->getIdent();
166 gezelter 1502
167     for( it2 = typeMap_.begin(); it2 != typeMap_.end(); ++it2) {
168     atype2 = (*it2).second;
169 gezelter 1895 atid2 = atype2->getIdent();
170 gezelter 1879
171 gezelter 1895 iHash_[atid1][atid2] = 0;
172 gezelter 1502
173     if (atype1->isLennardJones() && atype2->isLennardJones()) {
174 gezelter 1895 interactions_[atid1][atid2].insert(lj_);
175 gezelter 1927 iHash_[atid1][atid2] |= LJ_INTERACTION;
176 gezelter 1502 }
177     if (atype1->isElectrostatic() && atype2->isElectrostatic() ) {
178 gezelter 1929 interactions_[atid1][atid2].insert(electrostatic_);
179     iHash_[atid1][atid2] |= ELECTROSTATIC_INTERACTION;
180 gezelter 1502 }
181     if (atype1->isSticky() && atype2->isSticky() ) {
182 gezelter 1895 interactions_[atid1][atid2].insert(sticky_);
183 gezelter 1927 iHash_[atid1][atid2] |= STICKY_INTERACTION;
184 gezelter 1502 }
185     if (atype1->isStickyPower() && atype2->isStickyPower() ) {
186 gezelter 1895 interactions_[atid1][atid2].insert(sticky_);
187 gezelter 1927 iHash_[atid1][atid2] |= STICKY_INTERACTION;
188 gezelter 1502 }
189     if (atype1->isEAM() && atype2->isEAM() ) {
190 gezelter 1895 interactions_[atid1][atid2].insert(eam_);
191 gezelter 1927 iHash_[atid1][atid2] |= EAM_INTERACTION;
192 gezelter 1502 }
193     if (atype1->isSC() && atype2->isSC() ) {
194 gezelter 1895 interactions_[atid1][atid2].insert(sc_);
195 gezelter 1927 iHash_[atid1][atid2] |= SC_INTERACTION;
196 gezelter 1502 }
197     if (atype1->isGayBerne() && atype2->isGayBerne() ) {
198 gezelter 1895 interactions_[atid1][atid2].insert(gb_);
199 gezelter 1927 iHash_[atid1][atid2] |= GB_INTERACTION;
200 gezelter 1502 }
201     if ((atype1->isGayBerne() && atype2->isLennardJones())
202     || (atype1->isLennardJones() && atype2->isGayBerne())) {
203 gezelter 1895 interactions_[atid1][atid2].insert(gb_);
204 gezelter 1927 iHash_[atid1][atid2] |= GB_INTERACTION;
205 gezelter 1502 }
206    
207     // look for an explicitly-set non-bonded interaction type using the
208     // two atom types.
209     NonBondedInteractionType* nbiType = forceField_->getNonBondedInteractionType(atype1->getName(), atype2->getName());
210 gezelter 1535
211     if (nbiType != NULL) {
212 gezelter 1502
213 gezelter 1879 bool vdwExplicit = false;
214     bool metExplicit = false;
215     // bool hbExplicit = false;
216    
217 gezelter 1535 if (nbiType->isLennardJones()) {
218     // We found an explicit Lennard-Jones interaction.
219     // override all other vdw entries for this pair of atom types:
220     set<NonBondedInteraction*>::iterator it;
221 gezelter 1895 for (it = interactions_[atid1][atid2].begin();
222     it != interactions_[atid1][atid2].end(); ++it) {
223 gezelter 1535 InteractionFamily ifam = (*it)->getFamily();
224 gezelter 1879 if (ifam == VANDERWAALS_FAMILY) {
225 gezelter 1895 interactions_[atid1][atid2].erase(*it);
226     iHash_[atid1][atid2] ^= (*it)->getHash();
227 gezelter 1879 }
228 gezelter 1535 }
229 gezelter 1895 interactions_[atid1][atid2].insert(lj_);
230 gezelter 1927 iHash_[atid1][atid2] |= LJ_INTERACTION;
231 gezelter 1535 vdwExplicit = true;
232 gezelter 1502 }
233 gezelter 1535
234     if (nbiType->isMorse()) {
235     if (vdwExplicit) {
236     sprintf( painCave.errMsg,
237     "InteractionManager::initialize found more than one "
238     "explicit \n"
239     "\tvan der Waals 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 Morse interaction.
246     // override all other vdw entries for this pair of atom types:
247     set<NonBondedInteraction*>::iterator it;
248 gezelter 1895 for (it = interactions_[atid1][atid2].begin();
249     it != interactions_[atid1][atid2].end(); ++it) {
250 gezelter 1535 InteractionFamily ifam = (*it)->getFamily();
251 gezelter 1879 if (ifam == VANDERWAALS_FAMILY) {
252 gezelter 1895 interactions_[atid1][atid2].erase(*it);
253     iHash_[atid1][atid2] ^= (*it)->getHash();
254 gezelter 1879 }
255 gezelter 1535 }
256 gezelter 1895 interactions_[atid1][atid2].insert(morse_);
257 gezelter 1927 iHash_[atid1][atid2] |= MORSE_INTERACTION;
258 gezelter 1535 vdwExplicit = true;
259 gezelter 1502 }
260 jmichalk 1656
261     if (nbiType->isRepulsivePower()) {
262     if (vdwExplicit) {
263     sprintf( painCave.errMsg,
264     "InteractionManager::initialize found more than one "
265     "explicit \n"
266     "\tvan der Waals interaction for atom types %s - %s\n",
267     atype1->getName().c_str(), atype2->getName().c_str());
268     painCave.severity = OPENMD_ERROR;
269     painCave.isFatal = 1;
270     simError();
271     }
272     // We found an explicit RepulsivePower interaction.
273     // override all other vdw entries for this pair of atom types:
274     set<NonBondedInteraction*>::iterator it;
275 gezelter 1895 for (it = interactions_[atid1][atid2].begin();
276     it != interactions_[atid1][atid2].end(); ++it) {
277 jmichalk 1656 InteractionFamily ifam = (*it)->getFamily();
278 gezelter 1879 if (ifam == VANDERWAALS_FAMILY) {
279 gezelter 1895 interactions_[atid1][atid2].erase(*it);
280     iHash_[atid1][atid2] ^= (*it)->getHash();
281 gezelter 1879 }
282 jmichalk 1656 }
283 gezelter 1895 interactions_[atid1][atid2].insert(repulsivePower_);
284 gezelter 1927 iHash_[atid1][atid2] |= REPULSIVEPOWER_INTERACTION;
285 jmichalk 1656 vdwExplicit = true;
286     }
287 gezelter 1535
288 jmichalk 1656
289 gezelter 1535 if (nbiType->isEAM()) {
290     // We found an explicit EAM interaction.
291     // override all other metallic entries for this pair of atom types:
292     set<NonBondedInteraction*>::iterator it;
293 gezelter 1895 for (it = interactions_[atid1][atid2].begin();
294     it != interactions_[atid1][atid2].end(); ++it) {
295 gezelter 1535 InteractionFamily ifam = (*it)->getFamily();
296 gezelter 1879 if (ifam == METALLIC_FAMILY) {
297 gezelter 1895 interactions_[atid1][atid2].erase(*it);
298     iHash_[atid1][atid2] ^= (*it)->getHash();
299 gezelter 1879 }
300 gezelter 1535 }
301 gezelter 1895 interactions_[atid1][atid2].insert(eam_);
302 gezelter 1927 iHash_[atid1][atid2] |= EAM_INTERACTION;
303 gezelter 1535 metExplicit = true;
304 gezelter 1502 }
305 gezelter 1535
306     if (nbiType->isSC()) {
307     if (metExplicit) {
308     sprintf( painCave.errMsg,
309     "InteractionManager::initialize found more than one "
310     "explicit\n"
311     "\tmetallic interaction for atom types %s - %s\n",
312     atype1->getName().c_str(), atype2->getName().c_str());
313     painCave.severity = OPENMD_ERROR;
314     painCave.isFatal = 1;
315     simError();
316     }
317     // We found an explicit Sutton-Chen interaction.
318     // override all other metallic entries for this pair of atom types:
319     set<NonBondedInteraction*>::iterator it;
320 gezelter 1895 for (it = interactions_[atid1][atid2].begin();
321     it != interactions_[atid1][atid2].end(); ++it) {
322 gezelter 1535 InteractionFamily ifam = (*it)->getFamily();
323 gezelter 1879 if (ifam == METALLIC_FAMILY) {
324 gezelter 1895 interactions_[atid1][atid2].erase(*it);
325     iHash_[atid1][atid2] ^= (*it)->getHash();
326 gezelter 1879 }
327 gezelter 1535 }
328 gezelter 1895 interactions_[atid1][atid2].insert(sc_);
329 gezelter 1927 iHash_[atid1][atid2] |= SC_INTERACTION;
330 gezelter 1535 metExplicit = true;
331 gezelter 1502 }
332 gezelter 1535
333     if (nbiType->isMAW()) {
334     if (vdwExplicit) {
335     sprintf( painCave.errMsg,
336     "InteractionManager::initialize found more than one "
337     "explicit\n"
338     "\tvan der Waals interaction for atom types %s - %s\n",
339     atype1->getName().c_str(), atype2->getName().c_str());
340     painCave.severity = OPENMD_ERROR;
341     painCave.isFatal = 1;
342     simError();
343     }
344     // We found an explicit MAW interaction.
345     // override all other vdw entries for this pair of atom types:
346     set<NonBondedInteraction*>::iterator it;
347 gezelter 1895 for (it = interactions_[atid1][atid2].begin();
348     it != interactions_[atid1][atid2].end(); ++it) {
349 gezelter 1535 InteractionFamily ifam = (*it)->getFamily();
350 gezelter 1879 if (ifam == VANDERWAALS_FAMILY) {
351 gezelter 1895 interactions_[atid1][atid2].erase(*it);
352     iHash_[atid1][atid2] ^= (*it)->getHash();
353 gezelter 1879 }
354 gezelter 1535 }
355 gezelter 1895 interactions_[atid1][atid2].insert(maw_);
356 gezelter 1927 iHash_[atid1][atid2] |= MAW_INTERACTION;
357 gezelter 1535 vdwExplicit = true;
358     }
359 gezelter 1502 }
360     }
361     }
362    
363 gezelter 1535
364 gezelter 1583 // Make sure every pair of atom types in this simulation has a
365     // non-bonded interaction. If not, just inform the user.
366 gezelter 1535
367     set<AtomType*> simTypes = info_->getSimulatedAtomTypes();
368     set<AtomType*>::iterator it, jt;
369 gezelter 1583
370 gezelter 1535 for (it = simTypes.begin(); it != simTypes.end(); ++it) {
371     atype1 = (*it);
372 gezelter 1895 atid1 = atype1->getIdent();
373 gezelter 1583 for (jt = it; jt != simTypes.end(); ++jt) {
374 gezelter 1535 atype2 = (*jt);
375 gezelter 2005 atid2 = atype2->getIdent();
376 gezelter 1502
377 gezelter 1895 if (interactions_[atid1][atid2].size() == 0) {
378 gezelter 1502 sprintf( painCave.errMsg,
379 gezelter 1583 "InteractionManager could not find a matching non-bonded\n"
380     "\tinteraction for atom types %s - %s\n"
381     "\tProceeding without this interaction.\n",
382 gezelter 1502 atype1->getName().c_str(), atype2->getName().c_str());
383     painCave.severity = OPENMD_INFO;
384 gezelter 1583 painCave.isFatal = 0;
385 gezelter 1502 simError();
386     }
387     }
388     }
389 gezelter 1535
390     initialized_ = true;
391     }
392 gezelter 1584
393     void InteractionManager::setCutoffRadius(RealType rcut) {
394 gezelter 1587
395 gezelter 1584 electrostatic_->setCutoffRadius(rcut);
396 gezelter 1586 eam_->setCutoffRadius(rcut);
397 gezelter 1584 }
398    
399 gezelter 1895 void InteractionManager::doPrePair(InteractionData &idat){
400 gezelter 1502
401     if (!initialized_) initialize();
402 gezelter 1545
403 gezelter 1587 // excluded interaction, so just return
404     if (idat.excluded) return;
405    
406 gezelter 1895 int& iHash = iHash_[idat.atid1][idat.atid2];
407 gezelter 1504
408 gezelter 1927 if ((iHash & EAM_INTERACTION) != 0) eam_->calcDensity(idat);
409     if ((iHash & SC_INTERACTION) != 0) sc_->calcDensity(idat);
410 gezelter 1879
411     // set<NonBondedInteraction*>::iterator it;
412    
413     // for (it = interactions_[ idat.atypes ].begin();
414     // it != interactions_[ idat.atypes ].end(); ++it){
415     // if ((*it)->getFamily() == METALLIC_FAMILY) {
416     // dynamic_cast<MetallicInteraction*>(*it)->calcDensity(idat);
417     // }
418     // }
419 gezelter 1502
420     return;
421     }
422 gezelter 1504
423 gezelter 1895 void InteractionManager::doPreForce(SelfData &sdat){
424 gezelter 1502
425 gezelter 1504 if (!initialized_) initialize();
426 gezelter 1545
427 gezelter 1927 int& sHash = sHash_[sdat.atid];
428 gezelter 1879
429 jmichalk 2031 if ((sHash & EAM_INTERACTION) != 0) eam_->calcFunctional(sdat);
430 gezelter 1927 if ((sHash & SC_INTERACTION) != 0) sc_->calcFunctional(sdat);
431 gezelter 1879
432     // set<NonBondedInteraction*>::iterator it;
433 gezelter 1502
434 gezelter 1895 // for (it = interactions_[atid1][atid2].begin(); it != interactions_[atid1][atid2].end(); ++it){
435 gezelter 1879 // if ((*it)->getFamily() == METALLIC_FAMILY) {
436     // dynamic_cast<MetallicInteraction*>(*it)->calcFunctional(sdat);
437     // }
438     // }
439 gezelter 1504
440 gezelter 1502 return;
441     }
442    
443 gezelter 1895 void InteractionManager::doPair(InteractionData &idat){
444 gezelter 1504
445 gezelter 1502 if (!initialized_) initialize();
446 gezelter 1587
447 gezelter 1895 int& iHash = iHash_[idat.atid1][idat.atid2];
448 gezelter 1502
449 gezelter 1927 if ((iHash & ELECTROSTATIC_INTERACTION) != 0) electrostatic_->calcForce(idat);
450 gezelter 1879
451     // electrostatics still has to worry about indirect
452     // contributions from excluded pairs of atoms, but nothing else does:
453 gezelter 1502
454 gezelter 1879 if (idat.excluded) return;
455 gezelter 1502
456 gezelter 1927 if ((iHash & LJ_INTERACTION) != 0) lj_->calcForce(idat);
457     if ((iHash & GB_INTERACTION) != 0) gb_->calcForce(idat);
458     if ((iHash & STICKY_INTERACTION) != 0) sticky_->calcForce(idat);
459     if ((iHash & MORSE_INTERACTION) != 0) morse_->calcForce(idat);
460     if ((iHash & REPULSIVEPOWER_INTERACTION) != 0) repulsivePower_->calcForce(idat);
461     if ((iHash & EAM_INTERACTION) != 0) eam_->calcForce(idat);
462     if ((iHash & SC_INTERACTION) != 0) sc_->calcForce(idat);
463     if ((iHash & MAW_INTERACTION) != 0) maw_->calcForce(idat);
464 gezelter 1879
465     // set<NonBondedInteraction*>::iterator it;
466    
467     // for (it = interactions_[ idat.atypes ].begin();
468     // it != interactions_[ idat.atypes ].end(); ++it) {
469    
470     // // electrostatics still has to worry about indirect
471     // // contributions from excluded pairs of atoms:
472    
473     // if (!idat.excluded || (*it)->getFamily() == ELECTROSTATIC_FAMILY) {
474     // (*it)->calcForce(idat);
475     // }
476     // }
477 gezelter 1502
478 gezelter 1504 return;
479 gezelter 1502 }
480    
481 gezelter 1895 void InteractionManager::doSelfCorrection(SelfData &sdat){
482 gezelter 1502
483     if (!initialized_) initialize();
484 gezelter 1504
485 gezelter 1927 int& sHash = sHash_[sdat.atid];
486 gezelter 1502
487 gezelter 1927 if ((sHash & ELECTROSTATIC_INTERACTION) != 0) electrostatic_->calcSelfCorrection(sdat);
488 gezelter 1879
489     // set<NonBondedInteraction*>::iterator it;
490    
491 gezelter 1895 // for (it = interactions_[atid1][atid2].begin(); it != interactions_[atid1][atid2].end(); ++it){
492 gezelter 1879 // if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) {
493     // dynamic_cast<ElectrostaticInteraction*>(*it)->calcSelfCorrection(sdat);
494     // }
495     // }
496 gezelter 1536
497 gezelter 1504 return;
498 gezelter 1502 }
499    
500 gezelter 1925 void InteractionManager::doReciprocalSpaceSum(RealType &pot){
501 gezelter 1915 if (!initialized_) initialize();
502 gezelter 1921 electrostatic_->ReciprocalSpaceSum(pot);
503 gezelter 1915 }
504    
505 gezelter 1505 RealType InteractionManager::getSuggestedCutoffRadius(int *atid) {
506     if (!initialized_) initialize();
507 gezelter 1755
508 gezelter 1505 AtomType* atype = typeMap_[*atid];
509    
510     set<NonBondedInteraction*>::iterator it;
511     RealType cutoff = 0.0;
512    
513 gezelter 1921 for (it = interactions_[*atid][*atid].begin();
514     it != interactions_[*atid][*atid].end();
515 gezelter 1895 ++it)
516 gezelter 2033 cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(make_pair(atype,
517     atype)));
518 gezelter 1505 return cutoff;
519     }
520    
521 gezelter 1528 RealType InteractionManager::getSuggestedCutoffRadius(AtomType* atype) {
522     if (!initialized_) initialize();
523    
524 gezelter 1895 int atid = atype->getIdent();
525    
526 gezelter 1528 set<NonBondedInteraction*>::iterator it;
527     RealType cutoff = 0.0;
528    
529 gezelter 1921 for (it = interactions_[atid][atid].begin();
530     it != interactions_[atid][atid].end(); ++it)
531 gezelter 2033 cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(make_pair(atype,
532     atype)));
533 gezelter 1528 return cutoff;
534     }
535 gezelter 1502 } //end namespace OpenMD

Properties

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