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 1536 by gezelter, Wed Jan 5 14:49:05 2011 UTC vs.
Revision 1787 by gezelter, Wed Aug 29 18:13:11 2012 UTC

# Line 36 | Line 36
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).                        
39 > * [4]  Kuang & Gezelter,  J. Chem. Phys. 133, 164101 (2010).
40 > * [5]  Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41   */
42  
43   #include "nonbonded/InteractionManager.hpp"
43 #include "UseTheForce/doForces_interface.h"
44  
45   namespace OpenMD {
46  
47 <  InteractionManager* InteractionManager::_instance = NULL;
48 <  SimInfo* InteractionManager::info_ = NULL;
49 <  bool InteractionManager::initialized_ = false;
47 >  InteractionManager::InteractionManager() {
48  
49 <  RealType InteractionManager::rCut_ = 0.0;
50 <  RealType InteractionManager::rSwitch_ = 0.0;
51 <  RealType InteractionManager::skinThickness_ = 0.0;
52 <  RealType InteractionManager::listRadius_ = 0.0;
53 <  CutoffMethod InteractionManager::cutoffMethod_ = SHIFTED_FORCE;
54 <  SwitchingFunctionType InteractionManager::sft_ = cubic;
55 <  RealType InteractionManager::vdwScale_[4] = {1.0, 0.0, 0.0, 0.0};
56 <  RealType InteractionManager::electrostaticScale_[4] = {1.0, 0.0, 0.0, 0.0};
57 <
58 <  map<int, AtomType*> InteractionManager::typeMap_;
59 <  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) {
75 <      _instance = new InteractionManager();
76 <    }
77 <    return _instance;
49 >    initialized_ = false;
50 >        
51 >    lj_ = new LJ();
52 >    gb_ = new GB();
53 >    sticky_ = new Sticky();
54 >    morse_ = new Morse();
55 >    repulsivePower_ = new RepulsivePower();
56 >    eam_ = new EAM();
57 >    sc_ = new SC();
58 >    electrostatic_ = new Electrostatic();
59 >    maw_ = new MAW();
60    }
61  
62    void InteractionManager::initialize() {
63 <    
63 >
64 >    if (initialized_) return;
65 >
66      ForceField* forceField_ = info_->getForceField();
67      
68      lj_->setForceField(forceField_);
# Line 87 | Line 71 | namespace OpenMD {
71      eam_->setForceField(forceField_);
72      sc_->setForceField(forceField_);
73      morse_->setForceField(forceField_);
74 +    electrostatic_->setSimInfo(info_);
75      electrostatic_->setForceField(forceField_);
76      maw_->setForceField(forceField_);
77 <
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 = topologically unconnected
99 <    // 1 = bonded together
100 <    // 2 = connected via a bend
101 <    // 3 = connected via a torsion
102 <
103 <    vdwScale_[0] = 1.0;
104 <    vdwScale_[1] = fopts.getvdw12scale();
105 <    vdwScale_[2] = fopts.getvdw13scale();
106 <    vdwScale_[3] = fopts.getvdw14scale();
107 <
108 <    electrostaticScale_[0] = 1.0;
109 <    electrostaticScale_[1] = fopts.getelectrostatic12scale();
110 <    electrostaticScale_[2] = fopts.getelectrostatic13scale();
111 <    electrostaticScale_[3] = fopts.getelectrostatic14scale();    
77 >    repulsivePower_->setForceField(forceField_);
78  
79      ForceField::AtomTypeContainer* atomTypes = forceField_->getAtomTypes();
80      ForceField::AtomTypeContainer::MapTypeIterator i1, i2;
# Line 121 | Line 87 | namespace OpenMD {
87           atype1 = atomTypes->nextType(i1)) {
88        
89        // add it to the map:
124      AtomTypeProperties atp = atype1->getATP();    
90        
91        pair<map<int,AtomType*>::iterator,bool> ret;    
92 <      ret = typeMap_.insert( pair<int, AtomType*>(atp.ident, atype1) );
92 >      ret = typeMap_.insert( pair<int, AtomType*>(atype1->getIdent(), atype1) );
93        if (ret.second == false) {
94          sprintf( painCave.errMsg,
95                   "InteractionManager already had a previous entry with ident %d\n",
96 <                 atp.ident);
96 >                 atype1->getIdent());
97          painCave.severity = OPENMD_INFO;
98          painCave.isFatal = 0;
99          simError();                
# Line 217 | Line 182 | namespace OpenMD {
182              interactions_[key].insert(morse_);
183              vdwExplicit = true;
184            }
185 +
186 +          if (nbiType->isRepulsivePower()) {
187 +            if (vdwExplicit) {
188 +              sprintf( painCave.errMsg,
189 +                       "InteractionManager::initialize found more than one "
190 +                       "explicit \n"
191 +                       "\tvan der Waals interaction for atom types %s - %s\n",
192 +                       atype1->getName().c_str(), atype2->getName().c_str());
193 +              painCave.severity = OPENMD_ERROR;
194 +              painCave.isFatal = 1;
195 +              simError();
196 +            }
197 +            // We found an explicit RepulsivePower interaction.  
198 +            // override all other vdw entries for this pair of atom types:
199 +            set<NonBondedInteraction*>::iterator it;
200 +            for (it = interactions_[key].begin();
201 +                 it != interactions_[key].end(); ++it) {
202 +              InteractionFamily ifam = (*it)->getFamily();
203 +              if (ifam == VANDERWAALS_FAMILY) interactions_[key].erase(*it);
204 +            }
205 +            interactions_[key].insert(repulsivePower_);
206 +            vdwExplicit = true;
207 +          }
208            
209 +          
210            if (nbiType->isEAM()) {
211              // We found an explicit EAM interaction.  
212              // override all other metallic entries for this pair of atom types:
# Line 281 | Line 270 | namespace OpenMD {
270      }
271      
272      
273 <    // make sure every pair of atom types in this simulation has a
274 <    // non-bonded interaction:
273 >    // Make sure every pair of atom types in this simulation has a
274 >    // non-bonded interaction.  If not, just inform the user.
275  
276      set<AtomType*> simTypes = info_->getSimulatedAtomTypes();
277      set<AtomType*>::iterator it, jt;
278 +
279      for (it = simTypes.begin(); it != simTypes.end(); ++it) {
280        atype1 = (*it);
281 <      for (jt = simTypes.begin(); jt != simTypes.end(); ++jt) {
281 >      for (jt = it; jt != simTypes.end(); ++jt) {
282          atype2 = (*jt);
283          key = make_pair(atype1, atype2);
284          
285          if (interactions_[key].size() == 0) {
286            sprintf( painCave.errMsg,
287 <                   "InteractionManager unable to find an appropriate non-bonded\n"
288 <                   "\tinteraction for atom types %s - %s\n",
287 >                   "InteractionManager could not find a matching non-bonded\n"
288 >                   "\tinteraction for atom types %s - %s\n"
289 >                   "\tProceeding without this interaction.\n",
290                     atype1->getName().c_str(), atype2->getName().c_str());
291            painCave.severity = OPENMD_INFO;
292 <          painCave.isFatal = 1;
292 >          painCave.isFatal = 0;
293            simError();
294          }
295        }
296      }
297  
307    setupCutoffs();
308    setupSwitching();
309    setupNeighborlists();
310
311    int ljsp = cutoffMethod_ == SHIFTED_POTENTIAL ? 1 : 0;
312    int ljsf = cutoffMethod_ == SHIFTED_FORCE ? 1 : 0;
313    notifyFortranCutoffs(&rCut_, &rSwitch_, &ljsp, &ljsf);
314    notifyFortranSkinThickness(&skinThickness_);
315
298      initialized_ = true;
317  }
318  
319  /**
320   * setupCutoffs
321   *
322   * Sets the values of cutoffRadius and cutoffMethod
323   *
324   * cutoffRadius : realType
325   *  If the cutoffRadius was explicitly set, use that value.
326   *  If the cutoffRadius was not explicitly set:
327   *      Are there electrostatic atoms?  Use 12.0 Angstroms.
328   *      No electrostatic atoms?  Poll the atom types present in the
329   *      simulation for suggested cutoff values (e.g. 2.5 * sigma).
330   *      Use the maximum suggested value that was found.
331   *
332   * cutoffMethod : (one of HARD, SWITCHED, SHIFTED_FORCE, SHIFTED_POTENTIAL)
333   *      If cutoffMethod was explicitly set, use that choice.
334   *      If cutoffMethod was not explicitly set, use SHIFTED_FORCE
335   */
336  void InteractionManager::setupCutoffs() {
337    
338    Globals* simParams_ = info_->getSimParams();
339    
340    if (simParams_->haveCutoffRadius()) {
341      rCut_ = simParams_->getCutoffRadius();
342    } else {      
343      if (info_->usesElectrostaticAtoms()) {
344        sprintf(painCave.errMsg,
345                "InteractionManager::setupCutoffs: No value was set for the cutoffRadius.\n"
346                "\tOpenMD will use a default value of 12.0 angstroms"
347                "\tfor the cutoffRadius.\n");
348        painCave.isFatal = 0;
349        painCave.severity = OPENMD_INFO;
350        simError();
351        rCut_ = 12.0;
352      } else {
353        RealType thisCut;
354        set<AtomType*>::iterator i;
355        set<AtomType*> atomTypes;
356        atomTypes = info_->getSimulatedAtomTypes();        
357        for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
358          thisCut = getSuggestedCutoffRadius((*i));
359          rCut_ = max(thisCut, rCut_);
360        }
361        sprintf(painCave.errMsg,
362                "InteractionManager::setupCutoffs: No value was set for the cutoffRadius.\n"
363                "\tOpenMD will use %lf angstroms.\n",
364                rCut_);
365        painCave.isFatal = 0;
366        painCave.severity = OPENMD_INFO;
367        simError();
368      }            
369    }
370
371    map<string, CutoffMethod> stringToCutoffMethod;
372    stringToCutoffMethod["HARD"] = HARD;
373    stringToCutoffMethod["SWITCHED"] = SWITCHED;
374    stringToCutoffMethod["SHIFTED_POTENTIAL"] = SHIFTED_POTENTIAL;    
375    stringToCutoffMethod["SHIFTED_FORCE"] = SHIFTED_FORCE;
376  
377    if (simParams_->haveCutoffMethod()) {
378      string cutMeth = toUpperCopy(simParams_->getCutoffMethod());
379      map<string, CutoffMethod>::iterator i;
380      i = stringToCutoffMethod.find(cutMeth);
381      if (i == stringToCutoffMethod.end()) {
382        sprintf(painCave.errMsg,
383                "InteractionManager::setupCutoffs: Could not find chosen cutoffMethod %s\n"
384                "\tShould be one of: "
385                "HARD, SWITCHED, SHIFTED_POTENTIAL, or SHIFTED_FORCE\n",
386                cutMeth.c_str());
387        painCave.isFatal = 1;
388        painCave.severity = OPENMD_ERROR;
389        simError();
390      } else {
391        cutoffMethod_ = i->second;
392      }
393    } else {
394      sprintf(painCave.errMsg,
395              "InteractionManager::setupCutoffs: No value was set for the cutoffMethod.\n"
396              "\tOpenMD will use SHIFTED_FORCE.\n");
397        painCave.isFatal = 0;
398        painCave.severity = OPENMD_INFO;
399        simError();
400        cutoffMethod_ = SHIFTED_FORCE;        
401    }
299    }
300  
301 <
405 <  /**
406 <   * setupSwitching
407 <   *
408 <   * Sets the values of switchingRadius and
409 <   *  If the switchingRadius was explicitly set, use that value (but check it)
410 <   *  If the switchingRadius was not explicitly set: use 0.85 * cutoffRadius_
411 <   */
412 <  void InteractionManager::setupSwitching() {
413 <    Globals* simParams_ = info_->getSimParams();
414 <
415 <    if (simParams_->haveSwitchingRadius()) {
416 <      rSwitch_ = simParams_->getSwitchingRadius();
417 <      if (rSwitch_ > rCut_) {        
418 <        sprintf(painCave.errMsg,
419 <                "InteractionManager::setupSwitching: switchingRadius (%f) is larger than cutoffRadius(%f)\n",
420 <                rSwitch_, rCut_);
421 <        painCave.isFatal = 1;
422 <        painCave.severity = OPENMD_ERROR;
423 <        simError();
424 <      }
425 <    } else {      
426 <      rSwitch_ = 0.85 * rCut_;
427 <      sprintf(painCave.errMsg,
428 <              "InteractionManager::setupSwitching: No value was set for the switchingRadius.\n"
429 <              "\tOpenMD will use a default value of 85 percent of the cutoffRadius.\n"
430 <              "\tswitchingRadius = %f. for this simulation\n", rSwitch_);
431 <      painCave.isFatal = 0;
432 <      painCave.severity = OPENMD_WARNING;
433 <      simError();
434 <    }          
301 >  void InteractionManager::setCutoffRadius(RealType rcut) {
302      
303 <    if (simParams_->haveSwitchingFunctionType()) {
304 <      string funcType = simParams_->getSwitchingFunctionType();
438 <      toUpper(funcType);
439 <      if (funcType == "CUBIC") {
440 <        sft_ = cubic;
441 <      } else {
442 <        if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
443 <          sft_ = fifth_order_poly;
444 <        } else {
445 <          // throw error        
446 <          sprintf( painCave.errMsg,
447 <                   "InteractionManager::setupSwitching : Unknown switchingFunctionType. (Input file specified %s .)\n"
448 <                   "\tswitchingFunctionType must be one of: "
449 <                   "\"cubic\" or \"fifth_order_polynomial\".",
450 <                   funcType.c_str() );
451 <          painCave.isFatal = 1;
452 <          painCave.severity = OPENMD_ERROR;
453 <          simError();
454 <        }          
455 <      }
456 <    }
457 <
458 <    switcher_->setSwitchType(sft_);
459 <    switcher_->setSwitch(rSwitch_, rCut_);
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_;
303 >    electrostatic_->setCutoffRadius(rcut);
304 >    eam_->setCutoffRadius(rcut);
305    }
306  
307 <
489 <  void InteractionManager::doPrePair(int *atid1, int *atid2, RealType *rij, RealType *rho_i_at_j, RealType *rho_j_at_i){
307 >  void InteractionManager::doPrePair(InteractionData idat){
308      
309      if (!initialized_) initialize();
310 <          
311 <    DensityData ddat;
310 >        
311 >    // excluded interaction, so just return
312 >    if (idat.excluded) return;
313  
495    ddat.atype1 = typeMap_[*atid1];
496    ddat.atype2 = typeMap_[*atid2];
497    ddat.rij = *rij;
498    ddat.rho_i_at_j = *rho_i_at_j;
499    ddat.rho_j_at_i = *rho_j_at_i;
500
501    pair<AtomType*, AtomType*> key = make_pair(ddat.atype1, ddat.atype2);
314      set<NonBondedInteraction*>::iterator it;
315  
316 <    for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
316 >    for (it = interactions_[ idat.atypes ].begin();
317 >         it != interactions_[ idat.atypes ].end(); ++it){
318        if ((*it)->getFamily() == METALLIC_FAMILY) {
319 <        dynamic_cast<MetallicInteraction*>(*it)->calcDensity(ddat);
319 >        dynamic_cast<MetallicInteraction*>(*it)->calcDensity(idat);
320        }
321      }
322      
323      return;    
324    }
325    
326 <  void InteractionManager::doPreForce(int *atid, RealType *rho, RealType *frho, RealType *dfrhodrho){
326 >  void InteractionManager::doPreForce(SelfData sdat){
327  
328      if (!initialized_) initialize();
329 <          
330 <    FunctionalData fdat;
518 <
519 <    fdat.atype = typeMap_[*atid];
520 <    fdat.rho = *rho;
521 <    fdat.frho = *frho;
522 <    fdat.dfrhodrho = *dfrhodrho;
523 <
524 <    pair<AtomType*, AtomType*> key = make_pair(fdat.atype, fdat.atype);
329 >    
330 >    pair<AtomType*, AtomType*> key = make_pair(sdat.atype, sdat.atype);
331      set<NonBondedInteraction*>::iterator it;
332      
333      for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
334        if ((*it)->getFamily() == METALLIC_FAMILY) {
335 <        dynamic_cast<MetallicInteraction*>(*it)->calcFunctional(fdat);
335 >        dynamic_cast<MetallicInteraction*>(*it)->calcFunctional(sdat);
336        }
337      }
338      
339      return;    
340    }
341  
342 <  void InteractionManager::doPair(int *atid1, int *atid2, RealType *d, RealType *r, RealType *r2, RealType *sw, int *topoDist, RealType *A1, RealType *A2, RealType *eFrame1, RealType *eFrame2, RealType *vpair, RealType *pot, RealType  *f1, RealType *t1, RealType *t2, RealType *rho1, RealType *rho2, RealType *dfrho1, RealType *dfrho2, RealType *fshift1, RealType *fshift2){
342 >  void InteractionManager::doPair(InteractionData idat){
343      
344      if (!initialized_) initialize();
539    
540    InteractionData idat;
541    
542    idat.atype1 = typeMap_[*atid1];
543    idat.atype2 = typeMap_[*atid2];
544    idat.d = Vector3d(d);
545    idat.rij = *r;
546    idat.r2 = *r2;
547    idat.sw = *sw;
548    idat.vdwMult = vdwScale_[*topoDist];
549    idat.electroMult = electrostaticScale_[*topoDist];
550    for (int i = 0; i < 4; i++) {
551      idat.vpair[i] = vpair[i];
552      idat.pot[i] = pot[i];
553    }
554    idat.f1 = Vector3d(f1[0], f1[1], f1[2]);
555    idat.eFrame1 = Mat3x3d(eFrame1);
556    idat.eFrame2 = Mat3x3d(eFrame2);
557    idat.A1 = RotMat3x3d(A1);
558    idat.A2 = RotMat3x3d(A2);
559    idat.t1 = Vector3d(t1);
560    idat.t2 = Vector3d(t2);
561    idat.rho1 = *rho1;
562    idat.rho2 = *rho2;
563    idat.dfrho1 = *dfrho1;
564    idat.dfrho2 = *dfrho2;
565    idat.fshift1 = *fshift1;
566    idat.fshift2 = *fshift2;
345  
568    pair<AtomType*, AtomType*> key = make_pair(idat.atype1, idat.atype2);
346      set<NonBondedInteraction*>::iterator it;
347  
348 <    for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
349 <      (*it)->calcForce(idat);
573 <    
574 <    for (int i = 0; i < 4; i++) {
575 <      vpair[i] = idat.vpair[i];
576 <      pot[i] = idat.pot[i];
577 <    }
348 >    for (it = interactions_[ idat.atypes ].begin();
349 >         it != interactions_[ idat.atypes ].end(); ++it) {
350  
351 <    for (int i = 0; i < 3; i++) {
352 <      f1[i] = idat.f1[i];
581 <      t1[i] = idat.t1[i];
582 <      t2[i] = idat.t2[i];
583 <    }
351 >      // electrostatics still has to worry about indirect
352 >      // contributions from excluded pairs of atoms:
353  
354 <    return;    
355 <  }
587 <
588 <  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){
589 <
590 <    if (!initialized_) initialize();
591 <    
592 <    SkipCorrectionData skdat;
593 <    
594 <    skdat.atype1 = typeMap_[*atid1];
595 <    skdat.atype2 = typeMap_[*atid2];
596 <    skdat.d = Vector3d(d);
597 <    skdat.rij = *r;
598 <    skdat.skippedCharge1 = *skippedCharge1;
599 <    skdat.skippedCharge2 = *skippedCharge2;
600 <    skdat.sw = *sw;
601 <    skdat.electroMult = *electroMult;
602 <    for (int i = 0; i < 4; i++) {
603 <      skdat.vpair[i] = vpair[i];
604 <      skdat.pot[i] = pot[i];
605 <    }
606 <    skdat.f1 = Vector3d(f1);
607 <    skdat.eFrame1 = Mat3x3d(eFrame1);
608 <    skdat.eFrame2 = Mat3x3d(eFrame2);
609 <    skdat.t1 = Vector3d(t1);
610 <    skdat.t2 = Vector3d(t2);
611 <
612 <    pair<AtomType*, AtomType*> key = make_pair(skdat.atype1, skdat.atype2);
613 <    set<NonBondedInteraction*>::iterator it;
614 <
615 <    for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
616 <      if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) {
617 <        dynamic_cast<ElectrostaticInteraction*>(*it)->calcSkipCorrection(skdat);
354 >      if (!idat.excluded || (*it)->getFamily() == ELECTROSTATIC_FAMILY) {
355 >        (*it)->calcForce(idat);
356        }
357      }
358      
621    for (int i = 0; i < 4; i++) {
622      vpair[i] = skdat.vpair[i];
623      pot[i] = skdat.pot[i];
624    }
625
626    for (int i = 0; i < 3; i++) {
627      f1[i] = skdat.f1[i];
628      t1[i] = skdat.t1[i];
629      t2[i] = skdat.t2[i];
630    }
631
359      return;    
360    }
361  
362 <  void InteractionManager::doSelfCorrection(int *atid, RealType *eFrame, RealType *skippedCharge, RealType *pot, RealType *t){
362 >  void InteractionManager::doSelfCorrection(SelfData sdat){
363  
364      if (!initialized_) initialize();
365      
366 <    SelfCorrectionData scdat;
640 <    
641 <    scdat.atype = typeMap_[*atid];
642 <    scdat.eFrame = Mat3x3d(eFrame);
643 <    scdat.skippedCharge = *skippedCharge;
644 <    for (int i = 0; i < 4; i++){
645 <      scdat.pot[i] = pot[i];
646 <    }    
647 <    scdat.t = Vector3d(t);
648 <
649 <    pair<AtomType*, AtomType*> key = make_pair(scdat.atype, scdat.atype);
366 >    pair<AtomType*, AtomType*> key = make_pair(sdat.atype, sdat.atype);
367      set<NonBondedInteraction*>::iterator it;
368  
369      for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
370        if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) {
371 <        dynamic_cast<ElectrostaticInteraction*>(*it)->calcSelfCorrection(scdat);
371 >        dynamic_cast<ElectrostaticInteraction*>(*it)->calcSelfCorrection(sdat);
372        }
373      }
374        
658    t[0] = scdat.t.x();
659    t[1] = scdat.t.y();
660    t[2] = scdat.t.z();
661
375      return;    
376    }
377  
665
378    RealType InteractionManager::getSuggestedCutoffRadius(int *atid) {
379      if (!initialized_) initialize();
380 <    
380 >
381      AtomType* atype = typeMap_[*atid];
382  
383      pair<AtomType*, AtomType*> key = make_pair(atype, atype);
# Line 673 | Line 385 | namespace OpenMD {
385      RealType cutoff = 0.0;
386      
387      for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
388 <      cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(atype, atype));  
388 >      cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(key));  
389      return cutoff;    
390    }
391  
# Line 685 | Line 397 | namespace OpenMD {
397      RealType cutoff = 0.0;
398      
399      for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
400 <      cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(atype, atype));  
400 >      cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(key));  
401      return cutoff;    
402    }
691
692  void InteractionManager::getSwitch(RealType *r2, RealType *sw, RealType *dswdr, RealType *r,
693                                     int *in_switching_region) {
694    bool isr = switcher_->getSwitch( *r2 , *sw, *dswdr, *r);    
695    *in_switching_region = (int)isr;
696  }
697
403   } //end namespace OpenMD
699
700 extern "C" {
701  
702 #define fortranDoPrePair FC_FUNC(do_prepair, DO_PREPAIR)
703 #define fortranDoPreForce FC_FUNC(do_preforce, DO_PREFORCE)
704 #define fortranDoPair FC_FUNC(do_pair, DO_PAIR)
705 #define fortranDoSkipCorrection FC_FUNC(do_skip_correction, DO_SKIP_CORRECTION)
706 #define fortranDoSelfCorrection FC_FUNC(do_self_correction, DO_SELF_CORRECTION)
707 #define fortranGetCutoff FC_FUNC(get_cutoff, GET_CUTOFF)
708 #define fortranGetSwitch FC_FUNC(get_switch, GET_SWITCH)
709
710  void fortranDoPrePair(int *atid1, int *atid2, RealType *rij,
711                        RealType *rho_i_at_j, RealType *rho_j_at_i) {
712            
713    return OpenMD::InteractionManager::Instance()->doPrePair(atid1, atid2, rij,
714                                                             rho_i_at_j,  
715                                                             rho_j_at_i);
716  }
717  void fortranDoPreForce(int *atid, RealType *rho, RealType *frho,
718                         RealType *dfrhodrho) {  
719    
720    return OpenMD::InteractionManager::Instance()->doPreForce(atid, rho, frho,
721                                                              dfrhodrho);    
722  }
723  
724  void fortranDoPair(int *atid1, int *atid2, RealType *d, RealType *r, RealType *r2,
725                     RealType *sw, int *topoDist, RealType *A1, RealType *A2,
726                     RealType *eFrame1, RealType *eFrame2,
727                     RealType *vpair, RealType *pot,
728                     RealType *f1, RealType *t1, RealType *t2,
729                     RealType *rho1, RealType *rho2, RealType *dfrho1, RealType *dfrho2,
730                     RealType *fshift1, RealType *fshift2){
731    
732    return OpenMD::InteractionManager::Instance()->doPair(atid1, atid2, d, r, r2,
733                                                          sw, topoDist, A1, A2,
734                                                          eFrame1, eFrame2,
735                                                          vpair, pot,
736                                                          f1, t1, t2,
737                                                          rho1, rho2,
738                                                          dfrho1, dfrho2,
739                                                          fshift1, fshift2);
740  }
741  
742  void fortranDoSkipCorrection(int *atid1, int *atid2, RealType *d,
743                               RealType *r, RealType *skippedCharge1,
744                               RealType *skippedCharge2, RealType *sw,
745                               RealType *electroMult, RealType *pot,
746                               RealType *vpair, RealType *f1,
747                               RealType *eFrame1, RealType *eFrame2,
748                               RealType *t1, RealType *t2){
749    
750    return OpenMD::InteractionManager::Instance()->doSkipCorrection(atid1,
751                                                                    atid2, d,
752                                                                    r,
753                                                                    skippedCharge1,
754                                                                    skippedCharge2,
755                                                                    sw, electroMult, pot,
756                                                                    vpair, f1, eFrame1,
757                                                                    eFrame2, t1, t2);
758  }
759  
760  void fortranDoSelfCorrection(int *atid, RealType *eFrame, RealType *skippedCharge,
761                               RealType *pot, RealType *t) {
762    
763    return OpenMD::InteractionManager::Instance()->doSelfCorrection(atid,
764                                                                    eFrame,
765                                                                    skippedCharge,
766                                                                    pot, t);
767  }
768  RealType fortranGetCutoff(int *atid) {    
769    return OpenMD::InteractionManager::Instance()->getSuggestedCutoffRadius(atid);
770  }
771
772  void fortranGetSwitch(RealType *r2, RealType *sw, RealType *dswdr, RealType *r,
773                        int *in_switching_region) {    
774    return OpenMD::InteractionManager::Instance()->getSwitch(r2, sw, dswdr, r,
775                                                             in_switching_region);
776  }  
777 }

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines