ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/InteractionManager.cpp
(Generate patch)

Comparing branches/development/src/nonbonded/InteractionManager.cpp (file contents):
Revision 1504 by gezelter, Sat Oct 2 20:41:53 2010 UTC vs.
Revision 1535 by gezelter, Fri Dec 31 18:31:56 2010 UTC

# Line 40 | Line 40
40   */
41  
42   #include "nonbonded/InteractionManager.hpp"
43 + #include "UseTheForce/doForces_interface.h"
44  
45   namespace OpenMD {
46 <
46 >
47 >  InteractionManager* InteractionManager::_instance = NULL;
48 >  SimInfo* InteractionManager::info_ = NULL;
49    bool InteractionManager::initialized_ = false;
50 <  ForceField* InteractionManager::forceField_ = NULL;  
51 <  InteractionManager* InteractionManager::_instance = NULL;
50 >
51 >  RealType InteractionManager::rCut_ = 0.0;
52 >  RealType InteractionManager::rSwitch_ = 0.0;
53 >  RealType InteractionManager::skinThickness_ = 0.0;
54 >  RealType InteractionManager::listRadius_ = 0.0;
55 >  CutoffMethod InteractionManager::cutoffMethod_ = SHIFTED_FORCE;
56 >  SwitchingFunctionType InteractionManager::sft_ = cubic;
57 >  RealType InteractionManager::vdwScale_[4] = {0.0, 0.0, 0.0, 0.0};
58 >  RealType InteractionManager::electrostaticScale_[4] = {0.0, 0.0, 0.0, 0.0};
59 >
60    map<int, AtomType*> InteractionManager::typeMap_;
61    map<pair<AtomType*, AtomType*>, set<NonBondedInteraction*> > InteractionManager::interactions_;
62 +
63 +  LJ* InteractionManager::lj_ = new LJ();
64 +  GB* InteractionManager::gb_ = new GB();
65 +  Sticky* InteractionManager::sticky_ = new Sticky();
66 +  Morse* InteractionManager::morse_ = new Morse();
67 +  EAM* InteractionManager::eam_ = new EAM();
68 +  SC* InteractionManager::sc_ = new SC();
69 +  Electrostatic* InteractionManager::electrostatic_ = new Electrostatic();
70 +  MAW* InteractionManager::maw_ = new MAW();
71 +  SwitchingFunction* InteractionManager::switcher_ = new SwitchingFunction();
72  
73    InteractionManager* InteractionManager::Instance() {
74      if (!_instance) {
# Line 58 | Line 79 | namespace OpenMD {
79  
80    void InteractionManager::initialize() {
81      
82 <    lj_ = new LJ();
83 <    gb_ = new GB();
63 <    sticky_ = new Sticky();
64 <    eam_ = new EAM();
65 <    sc_ = new SC();
66 <    morse_ = new Morse();
67 <    electrostatic_ = new Electrostatic();
68 <
82 >    ForceField* forceField_ = info_->getForceField();
83 >    
84      lj_->setForceField(forceField_);
85      gb_->setForceField(forceField_);
86      sticky_->setForceField(forceField_);
# Line 73 | Line 88 | namespace OpenMD {
88      sc_->setForceField(forceField_);
89      morse_->setForceField(forceField_);
90      electrostatic_->setForceField(forceField_);
91 +    maw_->setForceField(forceField_);
92  
93 +    ForceFieldOptions& fopts = forceField_->getForceFieldOptions();
94 +
95 +    // Force fields can set options on how to scale van der Waals and electrostatic
96 +    // interactions for atoms connected via bonds, bends and torsions
97 +    // in this case the topological distance between atoms is:
98 +    // 0 = the atom itself
99 +    // 1 = bonded together
100 +    // 2 = connected via a bend
101 +    // 3 = connected via a torsion
102 +
103 +    vdwScale_[0] = 0.0;
104 +    vdwScale_[1] = fopts.getvdw12scale();
105 +    vdwScale_[2] = fopts.getvdw13scale();
106 +    vdwScale_[3] = fopts.getvdw14scale();
107 +
108 +    electrostaticScale_[0] = 0.0;
109 +    electrostaticScale_[1] = fopts.getelectrostatic12scale();
110 +    electrostaticScale_[2] = fopts.getelectrostatic13scale();
111 +    electrostaticScale_[3] = fopts.getelectrostatic14scale();    
112 +
113      ForceField::AtomTypeContainer* atomTypes = forceField_->getAtomTypes();
114      ForceField::AtomTypeContainer::MapTypeIterator i1, i2;
115      AtomType* atype1;
# Line 143 | Line 179 | namespace OpenMD {
179          // look for an explicitly-set non-bonded interaction type using the
180          // two atom types.
181          NonBondedInteractionType* nbiType = forceField_->getNonBondedInteractionType(atype1->getName(), atype2->getName());
182 +        
183 +        if (nbiType != NULL) {
184  
185 <        if (nbiType->isLennardJones()) {
186 <          // We found an explicit Lennard-Jones interaction.  
187 <          // override all other vdw entries for this pair of atom types:
188 <          set<NonBondedInteraction*>::iterator it;
189 <          for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it) {
190 <            InteractionFamily ifam = (*it)->getFamily();
191 <            if (ifam == VANDERWAALS_FAMILY) interactions_[key].erase(*it);
185 >          if (nbiType->isLennardJones()) {
186 >            // We found an explicit Lennard-Jones 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(lj_);
195 >            vdwExplicit = true;
196            }
197 <          interactions_[key].insert(lj_);
198 <          vdwExplicit = true;
199 <        }
200 <
201 <        if (nbiType->isMorse()) {
202 <          if (vdwExplicit) {
203 <            sprintf( painCave.errMsg,
204 <                     "InteractionManager::initialize found more than one explicit\n"
205 <                     "\tvan der Waals interaction for atom types %s - %s\n",
206 <                     atype1->getName().c_str(), atype2->getName().c_str());
207 <            painCave.severity = OPENMD_ERROR;
208 <            painCave.isFatal = 1;
209 <            simError();
197 >          
198 >          if (nbiType->isMorse()) {
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 Morse 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(morse_);
218 >            vdwExplicit = true;
219            }
220 <          // We found an explicit Morse interaction.  
221 <          // override all other vdw entries for this pair of atom types:
222 <          set<NonBondedInteraction*>::iterator it;
223 <          for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it) {
224 <            InteractionFamily ifam = (*it)->getFamily();
225 <            if (ifam == VANDERWAALS_FAMILY) interactions_[key].erase(*it);
226 <          }
227 <          interactions_[key].insert(morse_);
228 <          vdwExplicit = true;
229 <        }
230 <
231 <        if (nbiType->isEAM()) {
181 <          // We found an explicit EAM interaction.  
182 <          // override all other metallic entries for this pair of atom types:
183 <          set<NonBondedInteraction*>::iterator it;
184 <          for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it) {
185 <            InteractionFamily ifam = (*it)->getFamily();
186 <            if (ifam == METALLIC_FAMILY) interactions_[key].erase(*it);
220 >          
221 >          if (nbiType->isEAM()) {
222 >            // We found an explicit EAM interaction.  
223 >            // override all other metallic entries for this pair of atom types:
224 >            set<NonBondedInteraction*>::iterator it;
225 >            for (it = interactions_[key].begin();
226 >                 it != interactions_[key].end(); ++it) {
227 >              InteractionFamily ifam = (*it)->getFamily();
228 >              if (ifam == METALLIC_FAMILY) interactions_[key].erase(*it);
229 >            }
230 >            interactions_[key].insert(eam_);
231 >            metExplicit = true;
232            }
233 <          interactions_[key].insert(eam_);
234 <          metExplicit = true;
235 <        }
236 <
237 <        if (nbiType->isSC()) {
238 <          if (metExplicit) {
239 <            sprintf( painCave.errMsg,
240 <                     "InteractionManager::initialize found more than one explicit\n"
241 <                     "\tmetallic interaction for atom types %s - %s\n",
242 <                     atype1->getName().c_str(), atype2->getName().c_str());
243 <            painCave.severity = OPENMD_ERROR;
244 <            painCave.isFatal = 1;
245 <            simError();
233 >          
234 >          if (nbiType->isSC()) {
235 >            if (metExplicit) {
236 >              sprintf( painCave.errMsg,
237 >                       "InteractionManager::initialize found more than one "
238 >                       "explicit\n"
239 >                       "\tmetallic 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 Sutton-Chen interaction.  
246 >            // override all other metallic entries for this pair of atom types:
247 >            set<NonBondedInteraction*>::iterator it;
248 >            for (it = interactions_[key].begin();
249 >                 it != interactions_[key].end(); ++it) {
250 >              InteractionFamily ifam = (*it)->getFamily();
251 >              if (ifam == METALLIC_FAMILY) interactions_[key].erase(*it);
252 >            }
253 >            interactions_[key].insert(sc_);
254 >            metExplicit = true;
255            }
256 <          // We found an explicit Sutton-Chen interaction.  
257 <          // override all other metallic entries for this pair of atom types:
258 <          set<NonBondedInteraction*>::iterator it;
259 <          for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it) {
260 <            InteractionFamily ifam = (*it)->getFamily();
261 <            if (ifam == METALLIC_FAMILY) interactions_[key].erase(*it);
262 <          }
263 <          interactions_[key].insert(sc_);
264 <          metExplicit = true;
256 >          
257 >          if (nbiType->isMAW()) {
258 >            if (vdwExplicit) {
259 >              sprintf( painCave.errMsg,
260 >                       "InteractionManager::initialize found more than one "
261 >                       "explicit\n"
262 >                       "\tvan der Waals interaction for atom types %s - %s\n",
263 >                       atype1->getName().c_str(), atype2->getName().c_str());
264 >              painCave.severity = OPENMD_ERROR;
265 >              painCave.isFatal = 1;
266 >              simError();
267 >            }
268 >            // We found an explicit MAW interaction.  
269 >            // override all other vdw entries for this pair of atom types:
270 >            set<NonBondedInteraction*>::iterator it;
271 >            for (it = interactions_[key].begin();
272 >                 it != interactions_[key].end(); ++it) {
273 >              InteractionFamily ifam = (*it)->getFamily();
274 >              if (ifam == VANDERWAALS_FAMILY) interactions_[key].erase(*it);
275 >            }
276 >            interactions_[key].insert(maw_);
277 >            vdwExplicit = true;
278 >          }        
279          }
280        }
281      }
282      
283 <    // make sure every pair of atom types has a non-bonded interaction:
284 <    for (atype1 = atomTypes->beginType(i1); atype1 != NULL;
285 <         atype1 = atomTypes->nextType(i1)) {
286 <      for (atype2 = atomTypes->beginType(i2); atype2 != NULL;
287 <           atype2 = atomTypes->nextType(i2)) {
283 >    
284 >    // make sure every pair of atom types in this simulation has a
285 >    // non-bonded interaction:
286 >
287 >    set<AtomType*> simTypes = info_->getSimulatedAtomTypes();
288 >    set<AtomType*>::iterator it, jt;
289 >    for (it = simTypes.begin(); it != simTypes.end(); ++it) {
290 >      atype1 = (*it);
291 >      for (jt = simTypes.begin(); jt != simTypes.end(); ++jt) {
292 >        atype2 = (*jt);
293          key = make_pair(atype1, atype2);
294          
295          if (interactions_[key].size() == 0) {
# Line 230 | Line 303 | namespace OpenMD {
303          }
304        }
305      }
306 <  }
306 >
307 >    setupCutoffs();
308 >    setupSwitching();
309 >    setupNeighborlists();
310 >    notifyFortranSkinThickness(&skinThickness_);
311 >
312 >    int ljsp = cutoffMethod_ == SHIFTED_POTENTIAL ? 1 : 0;
313 >    int ljsf = cutoffMethod_ == SHIFTED_FORCE ? 1 : 0;
314 >    notifyFortranCutoffs(&rCut_, &rSwitch_, &ljsp, &ljsf);
315 >
316 >    int isError;
317 >    initFortranFF(&isError);
318 >
319 >    initialized_ = true;
320 >  }
321    
322 <  void InteractionManager::do_prepair(int *atid1, int *atid2, RealType *rij, RealType *rho_i_at_j, RealType *rho_j_at_i){
322 >  /**
323 >   * setupCutoffs
324 >   *
325 >   * Sets the values of cutoffRadius and cutoffMethod
326 >   *
327 >   * cutoffRadius : realType
328 >   *  If the cutoffRadius was explicitly set, use that value.
329 >   *  If the cutoffRadius was not explicitly set:
330 >   *      Are there electrostatic atoms?  Use 12.0 Angstroms.
331 >   *      No electrostatic atoms?  Poll the atom types present in the
332 >   *      simulation for suggested cutoff values (e.g. 2.5 * sigma).
333 >   *      Use the maximum suggested value that was found.
334 >   *
335 >   * cutoffMethod : (one of HARD, SWITCHED, SHIFTED_FORCE, SHIFTED_POTENTIAL)
336 >   *      If cutoffMethod was explicitly set, use that choice.
337 >   *      If cutoffMethod was not explicitly set, use SHIFTED_FORCE
338 >   */
339 >  void InteractionManager::setupCutoffs() {
340 >    
341 >    Globals* simParams_ = info_->getSimParams();
342 >    
343 >    if (simParams_->haveCutoffRadius()) {
344 >      rCut_ = simParams_->getCutoffRadius();
345 >    } else {      
346 >      if (info_->usesElectrostaticAtoms()) {
347 >        sprintf(painCave.errMsg,
348 >                "InteractionManager::setupCutoffs: No value was set for the cutoffRadius.\n"
349 >                "\tOpenMD will use a default value of 12.0 angstroms"
350 >                "\tfor the cutoffRadius.\n");
351 >        painCave.isFatal = 0;
352 >        painCave.severity = OPENMD_INFO;
353 >        simError();
354 >        rCut_ = 12.0;
355 >      } else {
356 >        RealType thisCut;
357 >        set<AtomType*>::iterator i;
358 >        set<AtomType*> atomTypes;
359 >        atomTypes = info_->getSimulatedAtomTypes();        
360 >        for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
361 >          thisCut = getSuggestedCutoffRadius((*i));
362 >          rCut_ = max(thisCut, rCut_);
363 >        }
364 >        sprintf(painCave.errMsg,
365 >                "InteractionManager::setupCutoffs: No value was set for the cutoffRadius.\n"
366 >                "\tOpenMD will use %lf angstroms.\n",
367 >                rCut_);
368 >        painCave.isFatal = 0;
369 >        painCave.severity = OPENMD_INFO;
370 >        simError();
371 >      }            
372 >    }
373 >
374 >    map<string, CutoffMethod> stringToCutoffMethod;
375 >    stringToCutoffMethod["HARD"] = HARD;
376 >    stringToCutoffMethod["SWITCHED"] = SWITCHED;
377 >    stringToCutoffMethod["SHIFTED_POTENTIAL"] = SHIFTED_POTENTIAL;    
378 >    stringToCutoffMethod["SHIFTED_FORCE"] = SHIFTED_FORCE;
379 >  
380 >    if (simParams_->haveCutoffMethod()) {
381 >      string cutMeth = toUpperCopy(simParams_->getCutoffMethod());
382 >      map<string, CutoffMethod>::iterator i;
383 >      i = stringToCutoffMethod.find(cutMeth);
384 >      if (i == stringToCutoffMethod.end()) {
385 >        sprintf(painCave.errMsg,
386 >                "InteractionManager::setupCutoffs: Could not find chosen cutoffMethod %s\n"
387 >                "\tShould be one of: "
388 >                "HARD, SWITCHED, SHIFTED_POTENTIAL, or SHIFTED_FORCE\n",
389 >                cutMeth.c_str());
390 >        painCave.isFatal = 1;
391 >        painCave.severity = OPENMD_ERROR;
392 >        simError();
393 >      } else {
394 >        cutoffMethod_ = i->second;
395 >      }
396 >    } else {
397 >      sprintf(painCave.errMsg,
398 >              "InteractionManager::setupCutoffs: No value was set for the cutoffMethod.\n"
399 >              "\tOpenMD will use SHIFTED_FORCE.\n");
400 >        painCave.isFatal = 0;
401 >        painCave.severity = OPENMD_INFO;
402 >        simError();
403 >        cutoffMethod_ = SHIFTED_FORCE;        
404 >    }
405 >  }
406 >
407 >
408 >  /**
409 >   * setupSwitching
410 >   *
411 >   * Sets the values of switchingRadius and
412 >   *  If the switchingRadius was explicitly set, use that value (but check it)
413 >   *  If the switchingRadius was not explicitly set: use 0.85 * cutoffRadius_
414 >   */
415 >  void InteractionManager::setupSwitching() {
416 >    Globals* simParams_ = info_->getSimParams();
417 >
418 >    if (simParams_->haveSwitchingRadius()) {
419 >      rSwitch_ = simParams_->getSwitchingRadius();
420 >      if (rSwitch_ > rCut_) {        
421 >        sprintf(painCave.errMsg,
422 >                "InteractionManager::setupSwitching: switchingRadius (%f) is larger than cutoffRadius(%f)\n",
423 >                rSwitch_, rCut_);
424 >        painCave.isFatal = 1;
425 >        painCave.severity = OPENMD_ERROR;
426 >        simError();
427 >      }
428 >    } else {      
429 >      rSwitch_ = 0.85 * rCut_;
430 >      sprintf(painCave.errMsg,
431 >              "InteractionManager::setupSwitching: No value was set for the switchingRadius.\n"
432 >              "\tOpenMD will use a default value of 85 percent of the cutoffRadius.\n"
433 >              "\tswitchingRadius = %f. for this simulation\n", rSwitch_);
434 >      painCave.isFatal = 0;
435 >      painCave.severity = OPENMD_WARNING;
436 >      simError();
437 >    }          
438      
439 +    if (simParams_->haveSwitchingFunctionType()) {
440 +      string funcType = simParams_->getSwitchingFunctionType();
441 +      toUpper(funcType);
442 +      if (funcType == "CUBIC") {
443 +        sft_ = cubic;
444 +      } else {
445 +        if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
446 +          sft_ = fifth_order_poly;
447 +        } else {
448 +          // throw error        
449 +          sprintf( painCave.errMsg,
450 +                   "InteractionManager::setupSwitching : Unknown switchingFunctionType. (Input file specified %s .)\n"
451 +                   "\tswitchingFunctionType must be one of: "
452 +                   "\"cubic\" or \"fifth_order_polynomial\".",
453 +                   funcType.c_str() );
454 +          painCave.isFatal = 1;
455 +          painCave.severity = OPENMD_ERROR;
456 +          simError();
457 +        }          
458 +      }
459 +    }
460 +  }
461 +
462 +  /**
463 +   * setupNeighborlists
464 +   *
465 +   *  If the skinThickness was explicitly set, use that value (but check it)
466 +   *  If the skinThickness was not explicitly set: use 1.0 angstroms
467 +   */
468 +  void InteractionManager::setupNeighborlists() {  
469 +
470 +    Globals* simParams_ = info_->getSimParams();    
471 +  
472 +    if (simParams_->haveSkinThickness()) {
473 +      skinThickness_ = simParams_->getSkinThickness();
474 +    } else {      
475 +      skinThickness_ = 1.0;
476 +      sprintf(painCave.errMsg,
477 +              "InteractionManager::setupNeighborlists: No value was set for the skinThickness.\n"
478 +              "\tOpenMD will use a default value of %f Angstroms\n"
479 +              "\tfor this simulation\n", skinThickness_);
480 +      painCave.severity = OPENMD_INFO;
481 +      painCave.isFatal = 0;
482 +      simError();
483 +    }            
484 +
485 +    listRadius_ = rCut_ + skinThickness_;
486 +
487 +  }
488 +
489 +
490 +  void InteractionManager::doPrePair(int *atid1, int *atid2, RealType *rij, RealType *rho_i_at_j, RealType *rho_j_at_i){
491 +    
492      if (!initialized_) initialize();
493            
494      DensityData ddat;
# Line 256 | Line 511 | namespace OpenMD {
511      return;    
512    }
513    
514 <  void InteractionManager::do_preforce(int *atid, RealType *rho, RealType *frho, RealType *dfrhodrho){
514 >  void InteractionManager::doPreForce(int *atid, RealType *rho, RealType *frho, RealType *dfrhodrho){
515  
516      if (!initialized_) initialize();
517            
# Line 279 | Line 534 | namespace OpenMD {
534      return;    
535    }
536  
537 <  void InteractionManager::do_pair(int *atid1, int *atid2, RealType *d, RealType *r, RealType *r2, RealType *rcut, RealType *sw, RealType *vdwMult,RealType *electroMult, RealType *pot, RealType *vpair, RealType *f1, RealType *eFrame1, RealType *eFrame2, RealType *A1, RealType *A2, RealType *t1, RealType *t2, RealType *rho1, RealType *rho2, RealType *dfrho1, RealType *dfrho2, RealType *fshift1, RealType *fshift2){
537 >  void InteractionManager::doPair(int *atid1, int *atid2, RealType *d, RealType *r, RealType *r2, RealType *rcut, RealType *sw, int *topoDist, RealType *pot, RealType *vpair, RealType *f1, RealType *eFrame1, RealType *eFrame2, RealType *A1, RealType *A2, RealType *t1, RealType *t2, RealType *rho1, RealType *rho2, RealType *dfrho1, RealType *dfrho2, RealType *fshift1, RealType *fshift2){
538      
539      if (!initialized_) initialize();
540      
# Line 292 | Line 547 | namespace OpenMD {
547      idat.r2 = *r2;
548      idat.rcut = *rcut;
549      idat.sw = *sw;
550 <    idat.vdwMult = *vdwMult;
551 <    idat.electroMult = *electroMult;
550 >    idat.vdwMult = vdwScale_[*topoDist];
551 >    idat.electroMult = electrostaticScale_[*topoDist];
552      idat.pot = *pot;
553      idat.vpair = *vpair;
554      idat.f1 = Vector3d(f1);
# Line 331 | Line 586 | namespace OpenMD {
586      return;    
587    }
588  
589 <  void InteractionManager::do_skip_correction(int *atid1, int *atid2, RealType *d, RealType *r, RealType *skippedCharge1, RealType *skippedCharge2, RealType *sw, RealType *electroMult, RealType *pot, RealType *vpair, RealType *f1, RealType *eFrame1, RealType *eFrame2, RealType *t1, RealType *t2){
589 >  void InteractionManager::doSkipCorrection(int *atid1, int *atid2, RealType *d, RealType *r, RealType *skippedCharge1, RealType *skippedCharge2, RealType *sw, RealType *electroMult, RealType *pot, RealType *vpair, RealType *f1, RealType *eFrame1, RealType *eFrame2, RealType *t1, RealType *t2){
590  
591      if (!initialized_) initialize();
592      
# Line 377 | Line 632 | namespace OpenMD {
632      return;    
633    }
634  
635 <  void InteractionManager::do_self_correction(int *atid, RealType *eFrame, RealType *skippedCharge, RealType *pot, RealType *t){
635 >  void InteractionManager::doSelfCorrection(int *atid, RealType *eFrame, RealType *skippedCharge, RealType *pot, RealType *t){
636  
637      if (!initialized_) initialize();
638      
# Line 405 | Line 660 | namespace OpenMD {
660      return;    
661    }
662  
663 +
664 +  RealType InteractionManager::getSuggestedCutoffRadius(int *atid) {
665 +    if (!initialized_) initialize();
666 +    
667 +    AtomType* atype = typeMap_[*atid];
668 +
669 +    pair<AtomType*, AtomType*> key = make_pair(atype, atype);
670 +    set<NonBondedInteraction*>::iterator it;
671 +    RealType cutoff = 0.0;
672 +    
673 +    for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
674 +      cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(atype, atype));  
675 +    return cutoff;    
676 +  }
677 +
678 +  RealType InteractionManager::getSuggestedCutoffRadius(AtomType* atype) {
679 +    if (!initialized_) initialize();
680 +    
681 +    pair<AtomType*, AtomType*> key = make_pair(atype, atype);
682 +    set<NonBondedInteraction*>::iterator it;
683 +    RealType cutoff = 0.0;
684 +    
685 +    for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
686 +      cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(atype, atype));  
687 +    return cutoff;    
688 +  }
689 +
690 +
691 +  void InteractionManager::setSwitch(RealType *rIn, RealType *rOut) {
692 +    switcher_->setSwitch(*rIn, *rOut);    
693 +  }
694 +
695 +  void InteractionManager::getSwitch(RealType *r2, RealType *sw, RealType *dswdr, RealType *r,
696 +                                     int *in_switching_region) {
697 +    bool isr = switcher_->getSwitch(*r2, *sw, *dswdr, *r);    
698 +    *in_switching_region = (int)isr;
699 +  }
700 +
701   } //end namespace OpenMD
702  
703   extern "C" {
# Line 415 | Line 708 | extern "C" {
708   #define fortranDoSkipCorrection FC_FUNC(do_skip_correction, DO_SKIP_CORRECTION)
709   #define fortranDoSelfCorrection FC_FUNC(do_self_correction, DO_SELF_CORRECTION)
710   #define fortranGetCutoff FC_FUNC(get_cutoff, GET_CUTOFF)
711 + #define fortranSetSwitch FC_FUNC(set_switch, SET_SWITCH)
712 + #define fortranGetSwitch FC_FUNC(get_switch, GET_SWITCH)
713  
714    void fortranDoPrePair(int *atid1, int *atid2, RealType *rij,
715                          RealType *rho_i_at_j, RealType *rho_j_at_i) {
716              
717 <    return OpenMD::InteractionManager::Instance()->do_prepair(atid1, atid2, rij,
718 <                                                              rho_i_at_j,  
719 <                                                              rho_j_at_i);
717 >    return OpenMD::InteractionManager::Instance()->doPrePair(atid1, atid2, rij,
718 >                                                             rho_i_at_j,  
719 >                                                             rho_j_at_i);
720    }
721 <  void fortranDoPreforce(int *atid, RealType *rho, RealType *frho,
721 >  void fortranDoPreForce(int *atid, RealType *rho, RealType *frho,
722                           RealType *dfrhodrho) {  
723      
724 <    return OpenMD::InteractionManager::Instance()->do_preforce(atid, rho, frho,
725 <                                                               dfrhodrho);    
724 >    return OpenMD::InteractionManager::Instance()->doPreForce(atid, rho, frho,
725 >                                                              dfrhodrho);    
726    }
727    
728    void fortranDoPair(int *atid1, int *atid2, RealType *d, RealType *r,
729 <                     RealType *r2, RealType *rcut, RealType *sw,
730 <                     RealType *vdwMult, RealType *electroMult, RealType *pot,
436 <                     RealType *vpair, RealType *f1, RealType *eFrame1,
729 >                     RealType *r2, RealType *rcut, RealType *sw, int *topoDist,
730 >                     RealType *pot, RealType *vpair, RealType *f1, RealType *eFrame1,
731                       RealType *eFrame2, RealType *A1, RealType *A2,
732                       RealType *t1, RealType *t2, RealType *rho1, RealType *rho2,
733                       RealType *dfrho1, RealType *dfrho2, RealType *fshift1,
734                       RealType *fshift2){
735      
736 <    return OpenMD::InteractionManager::Instance()->do_pair(atid1, atid2, d, r, r2, rcut,
737 <                                                           sw, vdwMult, electroMult, pot,
738 <                                                           vpair, f1, eFrame1, eFrame2,
739 <                                                           A1, A2, t1, t2, rho1, rho2,
740 <                                                           dfrho1, dfrho2, fshift1, fshift2);
736 >    return OpenMD::InteractionManager::Instance()->doPair(atid1, atid2, d, r,
737 >                                                          r2, rcut, sw, topoDist,
738 >                                                          pot, vpair, f1,
739 >                                                          eFrame1, eFrame2,
740 >                                                          A1, A2, t1, t2, rho1,
741 >                                                          rho2, dfrho1, dfrho2,
742 >                                                          fshift1, fshift2);
743    }
744 <
745 <  void fortranDoSkipCorrection(int *atid1, int *atid2, RealType *d, RealType *r,
746 <                               RealType *skippedCharge1, RealType *skippedCharge2,
747 <                               RealType *sw, RealType *electroMult, RealType *pot,
744 >  
745 >  void fortranDoSkipCorrection(int *atid1, int *atid2, RealType *d,
746 >                               RealType *r, RealType *skippedCharge1,
747 >                               RealType *skippedCharge2, RealType *sw,
748 >                               RealType *electroMult, RealType *pot,
749                                 RealType *vpair, RealType *f1,
750                                 RealType *eFrame1, RealType *eFrame2,
751                                 RealType *t1, RealType *t2){
752      
753 <    return OpenMD::InteractionManager::Instance()->do_skip_correction(atid1, atid2, d, r,
754 <                                                                      skippedCharge1,
755 <                                                                      skippedCharge2,
756 <                                                                      sw, electroMult, pot,
757 <                                                                      vpair, f1, eFrame1,
758 <                                                                      eFrame2, t1, t2);
753 >    return OpenMD::InteractionManager::Instance()->doSkipCorrection(atid1,
754 >                                                                    atid2, d,
755 >                                                                    r,
756 >                                                                    skippedCharge1,
757 >                                                                    skippedCharge2,
758 >                                                                    sw, electroMult, pot,
759 >                                                                    vpair, f1, eFrame1,
760 >                                                                    eFrame2, t1, t2);
761    }
762 <
762 >  
763    void fortranDoSelfCorrection(int *atid, RealType *eFrame, RealType *skippedCharge,
764                                 RealType *pot, RealType *t) {
765      
766 <    return OpenMD::InteractionManager::Instance()->do_self_correction(atid, eFrame,
767 <                                                                      skippedCharge,
768 <                                                                      pot, t);
766 >    return OpenMD::InteractionManager::Instance()->doSelfCorrection(atid,
767 >                                                                    eFrame,
768 >                                                                    skippedCharge,
769 >                                                                    pot, t);
770    }
771 <  RealType fortranGetCutoff() {
772 <    return OpenMD::InteractionManager::Instance()->getCutoff();
771 >  RealType fortranGetCutoff(int *atid) {    
772 >    return OpenMD::InteractionManager::Instance()->getSuggestedCutoffRadius(atid);
773    }
774 +
775 +  void fortranGetSwitch(RealType *r2, RealType *sw, RealType *dswdr, RealType *r,
776 +                        int *in_switching_region) {
777 +    
778 +    return OpenMD::InteractionManager::Instance()->getSwitch(r2, sw, dswdr, r,
779 +                                                             in_switching_region);
780 +  }
781 +
782 +  void fortranSetSwitch(RealType *rIn, RealType *rOut) {    
783 +    return OpenMD::InteractionManager::Instance()->setSwitch(rIn, rOut);
784 +  }
785 +  
786   }

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines