ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/InteractionManager.cpp
Revision: 1550
Committed: Wed Apr 27 21:49:59 2011 UTC (14 years ago) by gezelter
File size: 22137 byte(s)
Log Message:
more fortran removal

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     InteractionManager* InteractionManager::_instance = NULL;
47     SimInfo* InteractionManager::info_ = NULL;
48 gezelter 1502 bool InteractionManager::initialized_ = false;
49 gezelter 1535
50     RealType InteractionManager::rCut_ = 0.0;
51     RealType InteractionManager::rSwitch_ = 0.0;
52     RealType InteractionManager::skinThickness_ = 0.0;
53     RealType InteractionManager::listRadius_ = 0.0;
54     CutoffMethod InteractionManager::cutoffMethod_ = SHIFTED_FORCE;
55     SwitchingFunctionType InteractionManager::sft_ = cubic;
56 gezelter 1536 RealType InteractionManager::vdwScale_[4] = {1.0, 0.0, 0.0, 0.0};
57     RealType InteractionManager::electrostaticScale_[4] = {1.0, 0.0, 0.0, 0.0};
58 gezelter 1535
59 gezelter 1502 map<int, AtomType*> InteractionManager::typeMap_;
60     map<pair<AtomType*, AtomType*>, set<NonBondedInteraction*> > InteractionManager::interactions_;
61 gezelter 1528
62     LJ* InteractionManager::lj_ = new LJ();
63     GB* InteractionManager::gb_ = new GB();
64     Sticky* InteractionManager::sticky_ = new Sticky();
65     Morse* InteractionManager::morse_ = new Morse();
66     EAM* InteractionManager::eam_ = new EAM();
67     SC* InteractionManager::sc_ = new SC();
68     Electrostatic* InteractionManager::electrostatic_ = new Electrostatic();
69 gezelter 1532 MAW* InteractionManager::maw_ = new MAW();
70 gezelter 1530 SwitchingFunction* InteractionManager::switcher_ = new SwitchingFunction();
71 gezelter 1502
72     InteractionManager* InteractionManager::Instance() {
73     if (!_instance) {
74     _instance = new InteractionManager();
75     }
76     return _instance;
77     }
78    
79     void InteractionManager::initialize() {
80    
81 gezelter 1535 ForceField* forceField_ = info_->getForceField();
82    
83 gezelter 1502 lj_->setForceField(forceField_);
84     gb_->setForceField(forceField_);
85     sticky_->setForceField(forceField_);
86     eam_->setForceField(forceField_);
87     sc_->setForceField(forceField_);
88     morse_->setForceField(forceField_);
89     electrostatic_->setForceField(forceField_);
90 gezelter 1532 maw_->setForceField(forceField_);
91 gezelter 1502
92 gezelter 1530 ForceFieldOptions& fopts = forceField_->getForceFieldOptions();
93 gezelter 1535
94 gezelter 1530 // Force fields can set options on how to scale van der Waals and electrostatic
95     // interactions for atoms connected via bonds, bends and torsions
96     // in this case the topological distance between atoms is:
97 gezelter 1536 // 0 = topologically unconnected
98 gezelter 1530 // 1 = bonded together
99     // 2 = connected via a bend
100     // 3 = connected via a torsion
101    
102 gezelter 1536 vdwScale_[0] = 1.0;
103 gezelter 1530 vdwScale_[1] = fopts.getvdw12scale();
104     vdwScale_[2] = fopts.getvdw13scale();
105     vdwScale_[3] = fopts.getvdw14scale();
106    
107 gezelter 1536 electrostaticScale_[0] = 1.0;
108 gezelter 1530 electrostaticScale_[1] = fopts.getelectrostatic12scale();
109     electrostaticScale_[2] = fopts.getelectrostatic13scale();
110 gezelter 1535 electrostaticScale_[3] = fopts.getelectrostatic14scale();
111 gezelter 1530
112 gezelter 1502 ForceField::AtomTypeContainer* atomTypes = forceField_->getAtomTypes();
113     ForceField::AtomTypeContainer::MapTypeIterator i1, i2;
114     AtomType* atype1;
115     AtomType* atype2;
116     pair<AtomType*, AtomType*> key;
117     pair<set<NonBondedInteraction*>::iterator, bool> ret;
118    
119     for (atype1 = atomTypes->beginType(i1); atype1 != NULL;
120     atype1 = atomTypes->nextType(i1)) {
121    
122     // add it to the map:
123     AtomTypeProperties atp = atype1->getATP();
124    
125     pair<map<int,AtomType*>::iterator,bool> ret;
126     ret = typeMap_.insert( pair<int, AtomType*>(atp.ident, atype1) );
127     if (ret.second == false) {
128     sprintf( painCave.errMsg,
129     "InteractionManager already had a previous entry with ident %d\n",
130     atp.ident);
131     painCave.severity = OPENMD_INFO;
132     painCave.isFatal = 0;
133     simError();
134     }
135     }
136    
137     // Now, iterate over all known types and add to the interaction map:
138    
139     map<int, AtomType*>::iterator it1, it2;
140     for (it1 = typeMap_.begin(); it1 != typeMap_.end(); ++it1) {
141     atype1 = (*it1).second;
142    
143     for( it2 = typeMap_.begin(); it2 != typeMap_.end(); ++it2) {
144     atype2 = (*it2).second;
145    
146     bool vdwExplicit = false;
147     bool metExplicit = false;
148     bool hbExplicit = false;
149    
150     key = make_pair(atype1, atype2);
151    
152     if (atype1->isLennardJones() && atype2->isLennardJones()) {
153     interactions_[key].insert(lj_);
154     }
155     if (atype1->isElectrostatic() && atype2->isElectrostatic() ) {
156     interactions_[key].insert(electrostatic_);
157     }
158     if (atype1->isSticky() && atype2->isSticky() ) {
159     interactions_[key].insert(sticky_);
160     }
161     if (atype1->isStickyPower() && atype2->isStickyPower() ) {
162     interactions_[key].insert(sticky_);
163     }
164     if (atype1->isEAM() && atype2->isEAM() ) {
165     interactions_[key].insert(eam_);
166     }
167     if (atype1->isSC() && atype2->isSC() ) {
168     interactions_[key].insert(sc_);
169     }
170     if (atype1->isGayBerne() && atype2->isGayBerne() ) {
171     interactions_[key].insert(gb_);
172     }
173     if ((atype1->isGayBerne() && atype2->isLennardJones())
174     || (atype1->isLennardJones() && atype2->isGayBerne())) {
175     interactions_[key].insert(gb_);
176     }
177    
178     // look for an explicitly-set non-bonded interaction type using the
179     // two atom types.
180     NonBondedInteractionType* nbiType = forceField_->getNonBondedInteractionType(atype1->getName(), atype2->getName());
181 gezelter 1535
182     if (nbiType != NULL) {
183 gezelter 1502
184 gezelter 1535 if (nbiType->isLennardJones()) {
185     // We found an explicit Lennard-Jones 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(lj_);
194     vdwExplicit = true;
195 gezelter 1502 }
196 gezelter 1535
197     if (nbiType->isMorse()) {
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 Morse 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(morse_);
217     vdwExplicit = true;
218 gezelter 1502 }
219 gezelter 1535
220     if (nbiType->isEAM()) {
221     // We found an explicit EAM interaction.
222     // override all other metallic entries for this pair of atom types:
223     set<NonBondedInteraction*>::iterator it;
224     for (it = interactions_[key].begin();
225     it != interactions_[key].end(); ++it) {
226     InteractionFamily ifam = (*it)->getFamily();
227     if (ifam == METALLIC_FAMILY) interactions_[key].erase(*it);
228     }
229     interactions_[key].insert(eam_);
230     metExplicit = true;
231 gezelter 1502 }
232 gezelter 1535
233     if (nbiType->isSC()) {
234     if (metExplicit) {
235     sprintf( painCave.errMsg,
236     "InteractionManager::initialize found more than one "
237     "explicit\n"
238     "\tmetallic interaction for atom types %s - %s\n",
239     atype1->getName().c_str(), atype2->getName().c_str());
240     painCave.severity = OPENMD_ERROR;
241     painCave.isFatal = 1;
242     simError();
243     }
244     // We found an explicit Sutton-Chen interaction.
245     // override all other metallic entries for this pair of atom types:
246     set<NonBondedInteraction*>::iterator it;
247     for (it = interactions_[key].begin();
248     it != interactions_[key].end(); ++it) {
249     InteractionFamily ifam = (*it)->getFamily();
250     if (ifam == METALLIC_FAMILY) interactions_[key].erase(*it);
251     }
252     interactions_[key].insert(sc_);
253     metExplicit = true;
254 gezelter 1502 }
255 gezelter 1535
256     if (nbiType->isMAW()) {
257     if (vdwExplicit) {
258     sprintf( painCave.errMsg,
259     "InteractionManager::initialize found more than one "
260     "explicit\n"
261     "\tvan der Waals interaction for atom types %s - %s\n",
262     atype1->getName().c_str(), atype2->getName().c_str());
263     painCave.severity = OPENMD_ERROR;
264     painCave.isFatal = 1;
265     simError();
266     }
267     // We found an explicit MAW interaction.
268     // override all other vdw entries for this pair of atom types:
269     set<NonBondedInteraction*>::iterator it;
270     for (it = interactions_[key].begin();
271     it != interactions_[key].end(); ++it) {
272     InteractionFamily ifam = (*it)->getFamily();
273     if (ifam == VANDERWAALS_FAMILY) interactions_[key].erase(*it);
274     }
275     interactions_[key].insert(maw_);
276     vdwExplicit = true;
277     }
278 gezelter 1502 }
279     }
280     }
281    
282 gezelter 1535
283     // make sure every pair of atom types in this simulation has a
284     // non-bonded interaction:
285    
286     set<AtomType*> simTypes = info_->getSimulatedAtomTypes();
287     set<AtomType*>::iterator it, jt;
288     for (it = simTypes.begin(); it != simTypes.end(); ++it) {
289     atype1 = (*it);
290     for (jt = simTypes.begin(); jt != simTypes.end(); ++jt) {
291     atype2 = (*jt);
292 gezelter 1502 key = make_pair(atype1, atype2);
293    
294     if (interactions_[key].size() == 0) {
295     sprintf( painCave.errMsg,
296     "InteractionManager unable to find an appropriate non-bonded\n"
297     "\tinteraction for atom types %s - %s\n",
298     atype1->getName().c_str(), atype2->getName().c_str());
299     painCave.severity = OPENMD_INFO;
300     painCave.isFatal = 1;
301     simError();
302     }
303     }
304     }
305 gezelter 1535
306     setupCutoffs();
307     setupSwitching();
308     setupNeighborlists();
309    
310 gezelter 1545 //int ljsp = cutoffMethod_ == SHIFTED_POTENTIAL ? 1 : 0;
311     //int ljsf = cutoffMethod_ == SHIFTED_FORCE ? 1 : 0;
312     //notifyFortranCutoffs(&rCut_, &rSwitch_, &ljsp, &ljsf);
313     //notifyFortranSkinThickness(&skinThickness_);
314 gezelter 1535
315     initialized_ = true;
316     }
317 gezelter 1502
318 gezelter 1535 /**
319     * setupCutoffs
320     *
321     * Sets the values of cutoffRadius and cutoffMethod
322     *
323     * cutoffRadius : realType
324     * If the cutoffRadius was explicitly set, use that value.
325     * If the cutoffRadius was not explicitly set:
326     * Are there electrostatic atoms? Use 12.0 Angstroms.
327     * No electrostatic atoms? Poll the atom types present in the
328     * simulation for suggested cutoff values (e.g. 2.5 * sigma).
329     * Use the maximum suggested value that was found.
330     *
331     * cutoffMethod : (one of HARD, SWITCHED, SHIFTED_FORCE, SHIFTED_POTENTIAL)
332     * If cutoffMethod was explicitly set, use that choice.
333     * If cutoffMethod was not explicitly set, use SHIFTED_FORCE
334     */
335     void InteractionManager::setupCutoffs() {
336    
337     Globals* simParams_ = info_->getSimParams();
338    
339     if (simParams_->haveCutoffRadius()) {
340     rCut_ = simParams_->getCutoffRadius();
341     } else {
342     if (info_->usesElectrostaticAtoms()) {
343     sprintf(painCave.errMsg,
344     "InteractionManager::setupCutoffs: No value was set for the cutoffRadius.\n"
345     "\tOpenMD will use a default value of 12.0 angstroms"
346     "\tfor the cutoffRadius.\n");
347     painCave.isFatal = 0;
348     painCave.severity = OPENMD_INFO;
349     simError();
350     rCut_ = 12.0;
351     } else {
352     RealType thisCut;
353     set<AtomType*>::iterator i;
354     set<AtomType*> atomTypes;
355     atomTypes = info_->getSimulatedAtomTypes();
356     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
357     thisCut = getSuggestedCutoffRadius((*i));
358     rCut_ = max(thisCut, rCut_);
359     }
360     sprintf(painCave.errMsg,
361     "InteractionManager::setupCutoffs: No value was set for the cutoffRadius.\n"
362     "\tOpenMD will use %lf angstroms.\n",
363     rCut_);
364     painCave.isFatal = 0;
365     painCave.severity = OPENMD_INFO;
366     simError();
367     }
368     }
369    
370     map<string, CutoffMethod> stringToCutoffMethod;
371     stringToCutoffMethod["HARD"] = HARD;
372     stringToCutoffMethod["SWITCHED"] = SWITCHED;
373     stringToCutoffMethod["SHIFTED_POTENTIAL"] = SHIFTED_POTENTIAL;
374     stringToCutoffMethod["SHIFTED_FORCE"] = SHIFTED_FORCE;
375    
376     if (simParams_->haveCutoffMethod()) {
377     string cutMeth = toUpperCopy(simParams_->getCutoffMethod());
378     map<string, CutoffMethod>::iterator i;
379     i = stringToCutoffMethod.find(cutMeth);
380     if (i == stringToCutoffMethod.end()) {
381     sprintf(painCave.errMsg,
382     "InteractionManager::setupCutoffs: Could not find chosen cutoffMethod %s\n"
383     "\tShould be one of: "
384     "HARD, SWITCHED, SHIFTED_POTENTIAL, or SHIFTED_FORCE\n",
385     cutMeth.c_str());
386     painCave.isFatal = 1;
387     painCave.severity = OPENMD_ERROR;
388     simError();
389     } else {
390     cutoffMethod_ = i->second;
391     }
392     } else {
393     sprintf(painCave.errMsg,
394     "InteractionManager::setupCutoffs: No value was set for the cutoffMethod.\n"
395     "\tOpenMD will use SHIFTED_FORCE.\n");
396     painCave.isFatal = 0;
397     painCave.severity = OPENMD_INFO;
398     simError();
399     cutoffMethod_ = SHIFTED_FORCE;
400     }
401     }
402    
403    
404     /**
405     * setupSwitching
406     *
407     * Sets the values of switchingRadius and
408     * If the switchingRadius was explicitly set, use that value (but check it)
409     * If the switchingRadius was not explicitly set: use 0.85 * cutoffRadius_
410     */
411     void InteractionManager::setupSwitching() {
412     Globals* simParams_ = info_->getSimParams();
413    
414     if (simParams_->haveSwitchingRadius()) {
415     rSwitch_ = simParams_->getSwitchingRadius();
416     if (rSwitch_ > rCut_) {
417     sprintf(painCave.errMsg,
418     "InteractionManager::setupSwitching: switchingRadius (%f) is larger than cutoffRadius(%f)\n",
419     rSwitch_, rCut_);
420     painCave.isFatal = 1;
421     painCave.severity = OPENMD_ERROR;
422     simError();
423     }
424     } else {
425     rSwitch_ = 0.85 * rCut_;
426     sprintf(painCave.errMsg,
427     "InteractionManager::setupSwitching: No value was set for the switchingRadius.\n"
428     "\tOpenMD will use a default value of 85 percent of the cutoffRadius.\n"
429     "\tswitchingRadius = %f. for this simulation\n", rSwitch_);
430     painCave.isFatal = 0;
431     painCave.severity = OPENMD_WARNING;
432     simError();
433     }
434    
435     if (simParams_->haveSwitchingFunctionType()) {
436     string funcType = simParams_->getSwitchingFunctionType();
437     toUpper(funcType);
438     if (funcType == "CUBIC") {
439     sft_ = cubic;
440     } else {
441     if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
442     sft_ = fifth_order_poly;
443     } else {
444     // throw error
445     sprintf( painCave.errMsg,
446     "InteractionManager::setupSwitching : Unknown switchingFunctionType. (Input file specified %s .)\n"
447     "\tswitchingFunctionType must be one of: "
448     "\"cubic\" or \"fifth_order_polynomial\".",
449     funcType.c_str() );
450     painCave.isFatal = 1;
451     painCave.severity = OPENMD_ERROR;
452     simError();
453     }
454     }
455     }
456 gezelter 1536
457     switcher_->setSwitchType(sft_);
458     switcher_->setSwitch(rSwitch_, rCut_);
459 gezelter 1535 }
460    
461     /**
462     * setupNeighborlists
463     *
464     * If the skinThickness was explicitly set, use that value (but check it)
465     * If the skinThickness was not explicitly set: use 1.0 angstroms
466     */
467     void InteractionManager::setupNeighborlists() {
468    
469     Globals* simParams_ = info_->getSimParams();
470    
471     if (simParams_->haveSkinThickness()) {
472     skinThickness_ = simParams_->getSkinThickness();
473     } else {
474     skinThickness_ = 1.0;
475     sprintf(painCave.errMsg,
476     "InteractionManager::setupNeighborlists: No value was set for the skinThickness.\n"
477     "\tOpenMD will use a default value of %f Angstroms\n"
478     "\tfor this simulation\n", skinThickness_);
479     painCave.severity = OPENMD_INFO;
480     painCave.isFatal = 0;
481     simError();
482     }
483    
484     listRadius_ = rCut_ + skinThickness_;
485     }
486    
487    
488 gezelter 1545 void InteractionManager::doPrePair(InteractionData idat){
489 gezelter 1502
490     if (!initialized_) initialize();
491 gezelter 1545
492 gezelter 1504 set<NonBondedInteraction*>::iterator it;
493    
494 gezelter 1545 for (it = interactions_[idat.atypes].begin();
495     it != interactions_[idat.atypes].end(); ++it){
496 gezelter 1504 if ((*it)->getFamily() == METALLIC_FAMILY) {
497 gezelter 1545 dynamic_cast<MetallicInteraction*>(*it)->calcDensity(idat);
498 gezelter 1504 }
499     }
500 gezelter 1502
501     return;
502     }
503 gezelter 1504
504 gezelter 1545 void InteractionManager::doPreForce(SelfData sdat){
505 gezelter 1502
506 gezelter 1504 if (!initialized_) initialize();
507 gezelter 1545
508     pair<AtomType*, AtomType*> key = make_pair(sdat.atype, sdat.atype);
509 gezelter 1504 set<NonBondedInteraction*>::iterator it;
510 gezelter 1502
511 gezelter 1504 for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
512     if ((*it)->getFamily() == METALLIC_FAMILY) {
513 gezelter 1545 dynamic_cast<MetallicInteraction*>(*it)->calcFunctional(sdat);
514 gezelter 1504 }
515     }
516    
517 gezelter 1502 return;
518     }
519    
520 gezelter 1545 void InteractionManager::doPair(InteractionData idat){
521 gezelter 1504
522 gezelter 1502 if (!initialized_) initialize();
523 gezelter 1545
524 gezelter 1502 set<NonBondedInteraction*>::iterator it;
525    
526 gezelter 1545 for (it = interactions_[idat.atypes].begin();
527     it != interactions_[idat.atypes].end(); ++it)
528 gezelter 1502 (*it)->calcForce(idat);
529    
530     return;
531     }
532    
533 gezelter 1545 void InteractionManager::doSkipCorrection(InteractionData idat){
534 gezelter 1502
535 gezelter 1545 if (!initialized_) initialize();
536 gezelter 1504
537     set<NonBondedInteraction*>::iterator it;
538    
539 gezelter 1545 for (it = interactions_[idat.atypes].begin();
540     it != interactions_[idat.atypes].end(); ++it){
541 gezelter 1504 if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) {
542 gezelter 1545 dynamic_cast<ElectrostaticInteraction*>(*it)->calcSkipCorrection(idat);
543 gezelter 1504 }
544     }
545 gezelter 1502
546 gezelter 1504 return;
547 gezelter 1502 }
548    
549 gezelter 1545 void InteractionManager::doSelfCorrection(SelfData sdat){
550 gezelter 1502
551     if (!initialized_) initialize();
552 gezelter 1504
553 gezelter 1545 pair<AtomType*, AtomType*> key = make_pair(sdat.atype, sdat.atype);
554 gezelter 1504 set<NonBondedInteraction*>::iterator it;
555 gezelter 1502
556 gezelter 1504 for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
557     if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) {
558 gezelter 1545 dynamic_cast<ElectrostaticInteraction*>(*it)->calcSelfCorrection(sdat);
559 gezelter 1504 }
560     }
561 gezelter 1536
562 gezelter 1504 return;
563 gezelter 1502 }
564    
565 gezelter 1505 RealType InteractionManager::getSuggestedCutoffRadius(int *atid) {
566     if (!initialized_) initialize();
567    
568     AtomType* atype = typeMap_[*atid];
569    
570     pair<AtomType*, AtomType*> key = make_pair(atype, atype);
571     set<NonBondedInteraction*>::iterator it;
572     RealType cutoff = 0.0;
573    
574     for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
575 gezelter 1545 cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(key));
576 gezelter 1505 return cutoff;
577     }
578    
579 gezelter 1528 RealType InteractionManager::getSuggestedCutoffRadius(AtomType* atype) {
580     if (!initialized_) initialize();
581    
582     pair<AtomType*, AtomType*> key = make_pair(atype, atype);
583     set<NonBondedInteraction*>::iterator it;
584     RealType cutoff = 0.0;
585    
586     for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
587 gezelter 1545 cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(key));
588 gezelter 1528 return cutoff;
589     }
590 gezelter 1502 } //end namespace OpenMD

Properties

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