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 1535 by gezelter, Fri Dec 31 18:31:56 2010 UTC vs.
Revision 1874 by gezelter, Wed May 15 15:09:35 2013 UTC

# Line 35 | Line 35
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).                        
38 > * [3]  Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).          
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] = {0.0, 0.0, 0.0, 0.0};
56 <  RealType InteractionManager::electrostaticScale_[4] = {0.0, 0.0, 0.0, 0.0};
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 <  map<int, AtomType*> InteractionManager::typeMap_;
63 <  map<pair<AtomType*, AtomType*>, set<NonBondedInteraction*> > InteractionManager::interactions_;
64 <
65 <  LJ* InteractionManager::lj_ = new LJ();
66 <  GB* InteractionManager::gb_ = new GB();
67 <  Sticky* InteractionManager::sticky_ = new Sticky();
68 <  Morse* InteractionManager::morse_ = new Morse();
69 <  EAM* InteractionManager::eam_ = new EAM();
70 <  SC* InteractionManager::sc_ = new SC();
71 <  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;
62 >  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    void InteractionManager::initialize() {
75 <    
75 >
76 >    if (initialized_) return;
77 >
78      ForceField* forceField_ = info_->getForceField();
79      
80      lj_->setForceField(forceField_);
# Line 87 | Line 83 | namespace OpenMD {
83      eam_->setForceField(forceField_);
84      sc_->setForceField(forceField_);
85      morse_->setForceField(forceField_);
86 +    electrostatic_->setSimInfo(info_);
87      electrostatic_->setForceField(forceField_);
88      maw_->setForceField(forceField_);
89 <
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
89 >    repulsivePower_->setForceField(forceField_);
90  
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
91      ForceField::AtomTypeContainer* atomTypes = forceField_->getAtomTypes();
92      ForceField::AtomTypeContainer::MapTypeIterator i1, i2;
93      AtomType* atype1;
94      AtomType* atype2;
95      pair<AtomType*, AtomType*> key;
118    pair<set<NonBondedInteraction*>::iterator, bool> ret;
96      
97      for (atype1 = atomTypes->beginType(i1); atype1 != NULL;
98           atype1 = atomTypes->nextType(i1)) {
99        
100        // add it to the map:
124      AtomTypeProperties atp = atype1->getATP();    
101        
102        pair<map<int,AtomType*>::iterator,bool> ret;    
103 <      ret = typeMap_.insert( pair<int, AtomType*>(atp.ident, atype1) );
103 >      ret = typeMap_.insert( pair<int, AtomType*>(atype1->getIdent(), atype1) );
104        if (ret.second == false) {
105          sprintf( painCave.errMsg,
106                   "InteractionManager already had a previous entry with ident %d\n",
107 <                 atp.ident);
107 >                 atype1->getIdent());
108          painCave.severity = OPENMD_INFO;
109          painCave.isFatal = 0;
110          simError();                
# Line 143 | Line 119 | namespace OpenMD {
119  
120        for( it2 = typeMap_.begin(); it2 != typeMap_.end(); ++it2) {        
121          atype2 = (*it2).second;
122 <        
147 <        bool vdwExplicit = false;
148 <        bool metExplicit = false;
149 <        bool hbExplicit = false;
150 <                      
122 >                              
123          key = make_pair(atype1, atype2);
124          
125          if (atype1->isLennardJones() && atype2->isLennardJones()) {
# Line 182 | Line 154 | namespace OpenMD {
154          
155          if (nbiType != NULL) {
156  
157 +          bool vdwExplicit = false;
158 +          bool metExplicit = false;
159 +          // bool hbExplicit = false;
160 +
161            if (nbiType->isLennardJones()) {
162              // We found an explicit Lennard-Jones interaction.  
163              // override all other vdw entries for this pair of atom types:
# Line 217 | Line 193 | namespace OpenMD {
193              interactions_[key].insert(morse_);
194              vdwExplicit = true;
195            }
196 +
197 +          if (nbiType->isRepulsivePower()) {
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 RepulsivePower 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(repulsivePower_);
217 +            vdwExplicit = true;
218 +          }
219            
220 +          
221            if (nbiType->isEAM()) {
222              // We found an explicit EAM interaction.  
223              // override all other metallic entries for this pair of atom types:
# Line 281 | Line 281 | namespace OpenMD {
281      }
282      
283      
284 <    // make sure every pair of atom types in this simulation has a
285 <    // non-bonded interaction:
284 >    // Make sure every pair of atom types in this simulation has a
285 >    // non-bonded interaction.  If not, just inform the user.
286  
287      set<AtomType*> simTypes = info_->getSimulatedAtomTypes();
288      set<AtomType*>::iterator it, jt;
289 +
290      for (it = simTypes.begin(); it != simTypes.end(); ++it) {
291        atype1 = (*it);
292 <      for (jt = simTypes.begin(); jt != simTypes.end(); ++jt) {
292 >      for (jt = it; jt != simTypes.end(); ++jt) {
293          atype2 = (*jt);
294          key = make_pair(atype1, atype2);
295          
296          if (interactions_[key].size() == 0) {
297            sprintf( painCave.errMsg,
298 <                   "InteractionManager unable to find an appropriate non-bonded\n"
299 <                   "\tinteraction for atom types %s - %s\n",
298 >                   "InteractionManager could not find a matching non-bonded\n"
299 >                   "\tinteraction for atom types %s - %s\n"
300 >                   "\tProceeding without this interaction.\n",
301                     atype1->getName().c_str(), atype2->getName().c_str());
302            painCave.severity = OPENMD_INFO;
303 <          painCave.isFatal = 1;
303 >          painCave.isFatal = 0;
304            simError();
305          }
306        }
307      }
308  
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
309      initialized_ = true;
310    }
321  
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    }
311  
312 <    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 <    }          
312 >  void InteractionManager::setCutoffRadius(RealType rcut) {
313      
314 <    if (simParams_->haveSwitchingFunctionType()) {
315 <      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 <    }
314 >    electrostatic_->setCutoffRadius(rcut);
315 >    eam_->setCutoffRadius(rcut);
316    }
317  
318 <  /**
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){
318 >  void InteractionManager::doPrePair(InteractionData idat){
319      
320      if (!initialized_) initialize();
321 <          
322 <    DensityData ddat;
321 >        
322 >    // excluded interaction, so just return
323 >    if (idat.excluded) return;
324  
496    ddat.atype1 = typeMap_[*atid1];
497    ddat.atype2 = typeMap_[*atid2];
498    ddat.rij = *rij;
499    ddat.rho_i_at_j = *rho_i_at_j;
500    ddat.rho_j_at_i = *rho_j_at_i;
501
502    pair<AtomType*, AtomType*> key = make_pair(ddat.atype1, ddat.atype2);
325      set<NonBondedInteraction*>::iterator it;
326  
327 <    for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
327 >    for (it = interactions_[ idat.atypes ].begin();
328 >         it != interactions_[ idat.atypes ].end(); ++it){
329        if ((*it)->getFamily() == METALLIC_FAMILY) {
330 <        dynamic_cast<MetallicInteraction*>(*it)->calcDensity(ddat);
330 >        dynamic_cast<MetallicInteraction*>(*it)->calcDensity(idat);
331        }
332      }
333      
334      return;    
335    }
336    
337 <  void InteractionManager::doPreForce(int *atid, RealType *rho, RealType *frho, RealType *dfrhodrho){
337 >  void InteractionManager::doPreForce(SelfData sdat){
338  
339      if (!initialized_) initialize();
340 <          
341 <    FunctionalData fdat;
519 <
520 <    fdat.atype = typeMap_[*atid];
521 <    fdat.rho = *rho;
522 <    fdat.frho = *frho;
523 <    fdat.dfrhodrho = *dfrhodrho;
524 <
525 <    pair<AtomType*, AtomType*> key = make_pair(fdat.atype, fdat.atype);
340 >    
341 >    pair<AtomType*, AtomType*> key = make_pair(sdat.atype, sdat.atype);
342      set<NonBondedInteraction*>::iterator it;
343      
344      for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
345        if ((*it)->getFamily() == METALLIC_FAMILY) {
346 <        dynamic_cast<MetallicInteraction*>(*it)->calcFunctional(fdat);
346 >        dynamic_cast<MetallicInteraction*>(*it)->calcFunctional(sdat);
347        }
348      }
349      
350      return;    
351    }
352  
353 <  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){
353 >  void InteractionManager::doPair(InteractionData idat){
354      
355      if (!initialized_) initialize();
540    
541    InteractionData idat;
542    
543    idat.atype1 = typeMap_[*atid1];
544    idat.atype2 = typeMap_[*atid2];
545    idat.d = Vector3d(d);
546    idat.rij = *r;
547    idat.r2 = *r2;
548    idat.rcut = *rcut;
549    idat.sw = *sw;
550    idat.vdwMult = vdwScale_[*topoDist];
551    idat.electroMult = electrostaticScale_[*topoDist];
552    idat.pot = *pot;
553    idat.vpair = *vpair;
554    idat.f1 = Vector3d(f1);
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;
356  
568    pair<AtomType*, AtomType*> key = make_pair(idat.atype1, idat.atype2);
357      set<NonBondedInteraction*>::iterator it;
358  
359 <    for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
360 <      (*it)->calcForce(idat);
573 <    
574 <    f1[0] = idat.f1.x();
575 <    f1[1] = idat.f1.y();
576 <    f1[2] = idat.f1.z();
577 <    
578 <    t1[0] = idat.t1.x();
579 <    t1[1] = idat.t1.y();
580 <    t1[2] = idat.t1.z();
581 <    
582 <    t2[0] = idat.t2.x();
583 <    t2[1] = idat.t2.y();
584 <    t2[2] = idat.t2.z();
359 >    for (it = interactions_[ idat.atypes ].begin();
360 >         it != interactions_[ idat.atypes ].end(); ++it) {
361  
362 <    return;    
363 <  }
362 >      // electrostatics still has to worry about indirect
363 >      // contributions from excluded pairs of atoms:
364  
365 <  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){
366 <
591 <    if (!initialized_) initialize();
592 <    
593 <    SkipCorrectionData skdat;
594 <    
595 <    skdat.atype1 = typeMap_[*atid1];
596 <    skdat.atype2 = typeMap_[*atid2];
597 <    skdat.d = Vector3d(d);
598 <    skdat.rij = *r;
599 <    skdat.skippedCharge1 = *skippedCharge1;
600 <    skdat.skippedCharge2 = *skippedCharge2;
601 <    skdat.sw = *sw;
602 <    skdat.electroMult = *electroMult;
603 <    skdat.pot = *pot;
604 <    skdat.vpair = *vpair;
605 <    skdat.f1 = Vector3d(f1);
606 <    skdat.eFrame1 = Mat3x3d(eFrame1);
607 <    skdat.eFrame2 = Mat3x3d(eFrame2);
608 <    skdat.t1 = Vector3d(t1);
609 <    skdat.t2 = Vector3d(t2);
610 <
611 <    pair<AtomType*, AtomType*> key = make_pair(skdat.atype1, skdat.atype2);
612 <    set<NonBondedInteraction*>::iterator it;
613 <
614 <    for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
615 <      if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) {
616 <        dynamic_cast<ElectrostaticInteraction*>(*it)->calcSkipCorrection(skdat);
365 >      if (!idat.excluded || (*it)->getFamily() == ELECTROSTATIC_FAMILY) {
366 >        (*it)->calcForce(idat);
367        }
368      }
369      
620    f1[0] = skdat.f1.x();
621    f1[1] = skdat.f1.y();
622    f1[2] = skdat.f1.z();
623    
624    t1[0] = skdat.t1.x();
625    t1[1] = skdat.t1.y();
626    t1[2] = skdat.t1.z();
627    
628    t2[0] = skdat.t2.x();
629    t2[1] = skdat.t2.y();
630    t2[2] = skdat.t2.z();
631
370      return;    
371    }
372  
373 <  void InteractionManager::doSelfCorrection(int *atid, RealType *eFrame, RealType *skippedCharge, RealType *pot, RealType *t){
373 >  void InteractionManager::doSelfCorrection(SelfData sdat){
374  
375      if (!initialized_) initialize();
376      
377 <    SelfCorrectionData scdat;
640 <    
641 <    scdat.atype = typeMap_[*atid];
642 <    scdat.eFrame = Mat3x3d(eFrame);
643 <    scdat.skippedCharge = *skippedCharge;
644 <    scdat.pot = *pot;
645 <    scdat.t = Vector3d(t);
646 <
647 <    pair<AtomType*, AtomType*> key = make_pair(scdat.atype, scdat.atype);
377 >    pair<AtomType*, AtomType*> key = make_pair(sdat.atype, sdat.atype);
378      set<NonBondedInteraction*>::iterator it;
379  
380      for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
381        if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) {
382 <        dynamic_cast<ElectrostaticInteraction*>(*it)->calcSelfCorrection(scdat);
382 >        dynamic_cast<ElectrostaticInteraction*>(*it)->calcSelfCorrection(sdat);
383        }
384      }
385 <        
656 <    t[0] = scdat.t.x();
657 <    t[1] = scdat.t.y();
658 <    t[2] = scdat.t.z();
659 <
385 >      
386      return;    
387    }
388  
663
389    RealType InteractionManager::getSuggestedCutoffRadius(int *atid) {
390      if (!initialized_) initialize();
391 <    
391 >
392      AtomType* atype = typeMap_[*atid];
393  
394      pair<AtomType*, AtomType*> key = make_pair(atype, atype);
# Line 671 | Line 396 | namespace OpenMD {
396      RealType cutoff = 0.0;
397      
398      for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
399 <      cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(atype, atype));  
399 >      cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(key));  
400      return cutoff;    
401    }
402  
# Line 683 | Line 408 | namespace OpenMD {
408      RealType cutoff = 0.0;
409      
410      for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
411 <      cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(atype, atype));  
411 >      cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(key));  
412      return cutoff;    
413    }
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
414   } //end namespace OpenMD
702
703 extern "C" {
704  
705 #define fortranDoPrePair FC_FUNC(do_prepair, DO_PREPAIR)
706 #define fortranDoPreForce FC_FUNC(do_preforce, DO_PREFORCE)
707 #define fortranDoPair FC_FUNC(do_pair, DO_PAIR)
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()->doPrePair(atid1, atid2, rij,
718                                                             rho_i_at_j,  
719                                                             rho_j_at_i);
720  }
721  void fortranDoPreForce(int *atid, RealType *rho, RealType *frho,
722                         RealType *dfrhodrho) {  
723    
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, 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()->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,
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()->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  
763  void fortranDoSelfCorrection(int *atid, RealType *eFrame, RealType *skippedCharge,
764                               RealType *pot, RealType *t) {
765    
766    return OpenMD::InteractionManager::Instance()->doSelfCorrection(atid,
767                                                                    eFrame,
768                                                                    skippedCharge,
769                                                                    pot, t);
770  }
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