ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/nonbonded/InteractionManager.cpp
Revision: 1925
Committed: Wed Aug 7 15:24:16 2013 UTC (11 years, 8 months ago) by gezelter
File size: 19263 byte(s)
Log Message:
More ewald fixes, reporting reciprocal potential in stats.

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

Properties

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