ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/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.
trunk/src/nonbonded/InteractionManager.cpp (file contents), Revision 1927 by gezelter, Wed Aug 14 20:19:19 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 +    repulsivePower_->setForceField(forceField_);
90  
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
91      ForceField::AtomTypeContainer* atomTypes = forceField_->getAtomTypes();
92 +    int nTypes = atomTypes->size();
93 +    sHash_.resize(nTypes);
94 +    iHash_.resize(nTypes);
95 +    interactions_.resize(nTypes);
96      ForceField::AtomTypeContainer::MapTypeIterator i1, i2;
97      AtomType* atype1;
98      AtomType* atype2;
99 <    pair<AtomType*, AtomType*> key;
100 <    pair<set<NonBondedInteraction*>::iterator, bool> ret;
99 >    int atid1, atid2;
100 >
101 >    // We only need to worry about the types that are actually in the simulation:
102 >
103 >    set<AtomType*> atypes = info_->getSimulatedAtomTypes();
104 >
105 >    lj_->setSimulatedAtomTypes(atypes);
106 >    gb_->setSimulatedAtomTypes(atypes);
107 >    sticky_->setSimulatedAtomTypes(atypes);
108 >    eam_->setSimulatedAtomTypes(atypes);
109 >    sc_->setSimulatedAtomTypes(atypes);
110 >    morse_->setSimulatedAtomTypes(atypes);
111 >    electrostatic_->setSimInfo(info_);
112 >    electrostatic_->setSimulatedAtomTypes(atypes);
113 >    maw_->setSimulatedAtomTypes(atypes);
114 >    repulsivePower_->setSimulatedAtomTypes(atypes);
115 >
116 >    set<AtomType*>::iterator at;
117 >
118 >    for (at = atypes.begin(); at != atypes.end(); ++at) {
119      
120 <    for (atype1 = atomTypes->beginType(i1); atype1 != NULL;
121 <         atype1 = atomTypes->nextType(i1)) {
120 >      //for (atype1 = atomTypes->beginType(i1); atype1 != NULL;
121 >      //   atype1 = atomTypes->nextType(i1)) {
122        
123 <      // add it to the map:
124 <      AtomTypeProperties atp = atype1->getATP();    
125 <      
123 >      atype1 = *at;
124 >      atid1 = atype1->getIdent();
125 >      iHash_[atid1].resize(nTypes);
126 >      interactions_[atid1].resize(nTypes);
127 >
128 >      // add it to the map:      
129        pair<map<int,AtomType*>::iterator,bool> ret;    
130 <      ret = typeMap_.insert( pair<int, AtomType*>(atp.ident, atype1) );
130 >      ret = typeMap_.insert( pair<int, AtomType*>(atid1, atype1) );
131        if (ret.second == false) {
132          sprintf( painCave.errMsg,
133                   "InteractionManager already had a previous entry with ident %d\n",
134 <                 atp.ident);
134 >                 atype1->getIdent());
135          painCave.severity = OPENMD_INFO;
136          painCave.isFatal = 0;
137          simError();                
138        }
139      }
140 <    
140 >
141 >    if (atype1->isLennardJones()) {
142 >      sHash_[atid1] |= LJ_INTERACTION;
143 >    }
144 >    if (atype1->isElectrostatic()) {
145 >      sHash_[atid1] |= ELECTROSTATIC_INTERACTION;
146 >    }
147 >    if (atype1->isSticky()) {
148 >      sHash_[atid1] |= STICKY_INTERACTION;
149 >    }
150 >    if (atype1->isStickyPower()) {
151 >      sHash_[atid1] |= STICKY_INTERACTION;
152 >    }
153 >    if (atype1->isEAM()) {
154 >      sHash_[atid1] |= EAM_INTERACTION;
155 >    }
156 >    if (atype1->isSC()) {
157 >      sHash_[atid1] |= SC_INTERACTION;
158 >    }
159 >    if (atype1->isGayBerne()) {
160 >      sHash_[atid1] |= GB_INTERACTION;
161 >    }
162 >  
163      // Now, iterate over all known types and add to the interaction map:
164      
165      map<int, AtomType*>::iterator it1, it2;
166      for (it1 = typeMap_.begin(); it1 != typeMap_.end(); ++it1) {
167        atype1 = (*it1).second;
168 +      atid1 = atype1->getIdent();
169  
170        for( it2 = typeMap_.begin(); it2 != typeMap_.end(); ++it2) {        
171          atype2 = (*it2).second;
172 +        atid2 = atype2->getIdent();
173 +                              
174 +        iHash_[atid1][atid2] = 0;
175          
147        bool vdwExplicit = false;
148        bool metExplicit = false;
149        bool hbExplicit = false;
150                      
151        key = make_pair(atype1, atype2);
152        
176          if (atype1->isLennardJones() && atype2->isLennardJones()) {
177 <          interactions_[key].insert(lj_);
177 >          interactions_[atid1][atid2].insert(lj_);
178 >          iHash_[atid1][atid2] |= LJ_INTERACTION;
179          }
180          if (atype1->isElectrostatic() && atype2->isElectrostatic() ) {
181 <          interactions_[key].insert(electrostatic_);
181 >          // Pairs of fluctuating density EAM atoms have their
182 >          // interactions handled via the EAM routines.  All other
183 >          // interactions with these atoms are handled via normal
184 >          // electrostatic channels:
185 >          if (!(atype1->isEAM() && atype2->isEAM())) {
186 >            interactions_[atid1][atid2].insert(electrostatic_);
187 >            iHash_[atid1][atid2] |= ELECTROSTATIC_INTERACTION;
188 >          }
189          }
190          if (atype1->isSticky() && atype2->isSticky() ) {
191 <          interactions_[key].insert(sticky_);
191 >          interactions_[atid1][atid2].insert(sticky_);
192 >          iHash_[atid1][atid2] |= STICKY_INTERACTION;
193          }
194          if (atype1->isStickyPower() && atype2->isStickyPower() ) {
195 <          interactions_[key].insert(sticky_);
195 >          interactions_[atid1][atid2].insert(sticky_);
196 >          iHash_[atid1][atid2] |= STICKY_INTERACTION;
197          }
198          if (atype1->isEAM() && atype2->isEAM() ) {
199 <          interactions_[key].insert(eam_);
199 >          interactions_[atid1][atid2].insert(eam_);
200 >          iHash_[atid1][atid2] |= EAM_INTERACTION;
201          }
202          if (atype1->isSC() && atype2->isSC() ) {
203 <          interactions_[key].insert(sc_);
203 >          interactions_[atid1][atid2].insert(sc_);
204 >          iHash_[atid1][atid2] |= SC_INTERACTION;
205          }
206          if (atype1->isGayBerne() && atype2->isGayBerne() ) {
207 <          interactions_[key].insert(gb_);
207 >          interactions_[atid1][atid2].insert(gb_);
208 >          iHash_[atid1][atid2] |= GB_INTERACTION;
209          }
210          if ((atype1->isGayBerne() && atype2->isLennardJones())
211              || (atype1->isLennardJones() && atype2->isGayBerne())) {
212 <          interactions_[key].insert(gb_);
212 >          interactions_[atid1][atid2].insert(gb_);
213 >          iHash_[atid1][atid2] |= GB_INTERACTION;
214          }
215          
216          // look for an explicitly-set non-bonded interaction type using the
# Line 182 | Line 219 | namespace OpenMD {
219          
220          if (nbiType != NULL) {
221  
222 +          bool vdwExplicit = false;
223 +          bool metExplicit = false;
224 +          // bool hbExplicit = false;
225 +
226            if (nbiType->isLennardJones()) {
227              // We found an explicit Lennard-Jones interaction.  
228              // override all other vdw entries for this pair of atom types:
229              set<NonBondedInteraction*>::iterator it;
230 <            for (it = interactions_[key].begin();
231 <                 it != interactions_[key].end(); ++it) {
230 >            for (it = interactions_[atid1][atid2].begin();
231 >                 it != interactions_[atid1][atid2].end(); ++it) {
232                InteractionFamily ifam = (*it)->getFamily();
233 <              if (ifam == VANDERWAALS_FAMILY) interactions_[key].erase(*it);
233 >              if (ifam == VANDERWAALS_FAMILY) {
234 >                interactions_[atid1][atid2].erase(*it);
235 >                iHash_[atid1][atid2] ^= (*it)->getHash();
236 >              }
237              }
238 <            interactions_[key].insert(lj_);
238 >            interactions_[atid1][atid2].insert(lj_);
239 >            iHash_[atid1][atid2] |= LJ_INTERACTION;
240              vdwExplicit = true;
241            }
242            
# Line 209 | Line 254 | namespace OpenMD {
254              // We found an explicit Morse interaction.  
255              // override all other vdw entries for this pair of atom types:
256              set<NonBondedInteraction*>::iterator it;
257 <            for (it = interactions_[key].begin();
258 <                 it != interactions_[key].end(); ++it) {
257 >            for (it = interactions_[atid1][atid2].begin();
258 >                 it != interactions_[atid1][atid2].end(); ++it) {
259                InteractionFamily ifam = (*it)->getFamily();
260 <              if (ifam == VANDERWAALS_FAMILY) interactions_[key].erase(*it);
260 >              if (ifam == VANDERWAALS_FAMILY) {
261 >                interactions_[atid1][atid2].erase(*it);
262 >                iHash_[atid1][atid2] ^= (*it)->getHash();
263 >              }
264 >            }
265 >            interactions_[atid1][atid2].insert(morse_);
266 >            iHash_[atid1][atid2] |= MORSE_INTERACTION;
267 >            vdwExplicit = true;
268 >          }
269 >
270 >          if (nbiType->isRepulsivePower()) {
271 >            if (vdwExplicit) {
272 >              sprintf( painCave.errMsg,
273 >                       "InteractionManager::initialize found more than one "
274 >                       "explicit \n"
275 >                       "\tvan der Waals interaction for atom types %s - %s\n",
276 >                       atype1->getName().c_str(), atype2->getName().c_str());
277 >              painCave.severity = OPENMD_ERROR;
278 >              painCave.isFatal = 1;
279 >              simError();
280              }
281 <            interactions_[key].insert(morse_);
281 >            // We found an explicit RepulsivePower interaction.  
282 >            // override all other vdw entries for this pair of atom types:
283 >            set<NonBondedInteraction*>::iterator it;
284 >            for (it = interactions_[atid1][atid2].begin();
285 >                 it != interactions_[atid1][atid2].end(); ++it) {
286 >              InteractionFamily ifam = (*it)->getFamily();
287 >              if (ifam == VANDERWAALS_FAMILY) {
288 >                interactions_[atid1][atid2].erase(*it);
289 >                iHash_[atid1][atid2] ^= (*it)->getHash();
290 >              }
291 >            }
292 >            interactions_[atid1][atid2].insert(repulsivePower_);
293 >            iHash_[atid1][atid2] |= REPULSIVEPOWER_INTERACTION;
294              vdwExplicit = true;
295            }
296            
297 +          
298            if (nbiType->isEAM()) {
299              // We found an explicit EAM interaction.  
300              // override all other metallic entries for this pair of atom types:
301              set<NonBondedInteraction*>::iterator it;
302 <            for (it = interactions_[key].begin();
303 <                 it != interactions_[key].end(); ++it) {
302 >            for (it = interactions_[atid1][atid2].begin();
303 >                 it != interactions_[atid1][atid2].end(); ++it) {
304                InteractionFamily ifam = (*it)->getFamily();
305 <              if (ifam == METALLIC_FAMILY) interactions_[key].erase(*it);
305 >              if (ifam == METALLIC_FAMILY) {
306 >                interactions_[atid1][atid2].erase(*it);
307 >                iHash_[atid1][atid2] ^= (*it)->getHash();
308 >              }
309              }
310 <            interactions_[key].insert(eam_);
310 >            interactions_[atid1][atid2].insert(eam_);
311 >            iHash_[atid1][atid2] |= EAM_INTERACTION;
312              metExplicit = true;
313            }
314            
# Line 245 | Line 326 | namespace OpenMD {
326              // We found an explicit Sutton-Chen interaction.  
327              // override all other metallic entries for this pair of atom types:
328              set<NonBondedInteraction*>::iterator it;
329 <            for (it = interactions_[key].begin();
330 <                 it != interactions_[key].end(); ++it) {
329 >            for (it = interactions_[atid1][atid2].begin();
330 >                 it != interactions_[atid1][atid2].end(); ++it) {
331                InteractionFamily ifam = (*it)->getFamily();
332 <              if (ifam == METALLIC_FAMILY) interactions_[key].erase(*it);
332 >              if (ifam == METALLIC_FAMILY) {
333 >                interactions_[atid1][atid2].erase(*it);
334 >                iHash_[atid1][atid2] ^= (*it)->getHash();
335 >              }
336              }
337 <            interactions_[key].insert(sc_);
337 >            interactions_[atid1][atid2].insert(sc_);
338 >            iHash_[atid1][atid2] |= SC_INTERACTION;
339              metExplicit = true;
340            }
341            
# Line 268 | Line 353 | namespace OpenMD {
353              // We found an explicit MAW interaction.  
354              // override all other vdw entries for this pair of atom types:
355              set<NonBondedInteraction*>::iterator it;
356 <            for (it = interactions_[key].begin();
357 <                 it != interactions_[key].end(); ++it) {
356 >            for (it = interactions_[atid1][atid2].begin();
357 >                 it != interactions_[atid1][atid2].end(); ++it) {
358                InteractionFamily ifam = (*it)->getFamily();
359 <              if (ifam == VANDERWAALS_FAMILY) interactions_[key].erase(*it);
359 >              if (ifam == VANDERWAALS_FAMILY) {
360 >                interactions_[atid1][atid2].erase(*it);
361 >                iHash_[atid1][atid2] ^= (*it)->getHash();
362 >              }
363              }
364 <            interactions_[key].insert(maw_);
364 >            interactions_[atid1][atid2].insert(maw_);
365 >            iHash_[atid1][atid2] |= MAW_INTERACTION;
366              vdwExplicit = true;
367            }        
368          }
# Line 281 | Line 370 | namespace OpenMD {
370      }
371      
372      
373 <    // make sure every pair of atom types in this simulation has a
374 <    // non-bonded interaction:
373 >    // Make sure every pair of atom types in this simulation has a
374 >    // non-bonded interaction.  If not, just inform the user.
375  
376      set<AtomType*> simTypes = info_->getSimulatedAtomTypes();
377      set<AtomType*>::iterator it, jt;
378 +
379      for (it = simTypes.begin(); it != simTypes.end(); ++it) {
380        atype1 = (*it);
381 <      for (jt = simTypes.begin(); jt != simTypes.end(); ++jt) {
381 >      atid1 = atype1->getIdent();
382 >      for (jt = it; jt != simTypes.end(); ++jt) {
383          atype2 = (*jt);
384 <        key = make_pair(atype1, atype2);
384 >        atid1 = atype1->getIdent();
385          
386 <        if (interactions_[key].size() == 0) {
386 >        if (interactions_[atid1][atid2].size() == 0) {
387            sprintf( painCave.errMsg,
388 <                   "InteractionManager unable to find an appropriate non-bonded\n"
389 <                   "\tinteraction for atom types %s - %s\n",
388 >                   "InteractionManager could not find a matching non-bonded\n"
389 >                   "\tinteraction for atom types %s - %s\n"
390 >                   "\tProceeding without this interaction.\n",
391                     atype1->getName().c_str(), atype2->getName().c_str());
392            painCave.severity = OPENMD_INFO;
393 <          painCave.isFatal = 1;
393 >          painCave.isFatal = 0;
394            simError();
395          }
396        }
397      }
398  
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
399      initialized_ = true;
400    }
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    }
401  
402 <    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 <    }          
402 >  void InteractionManager::setCutoffRadius(RealType rcut) {
403      
404 <    if (simParams_->haveSwitchingFunctionType()) {
405 <      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 <    }
404 >    electrostatic_->setCutoffRadius(rcut);
405 >    eam_->setCutoffRadius(rcut);
406    }
407  
408 <  /**
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){
408 >  void InteractionManager::doPrePair(InteractionData &idat){
409      
410      if (!initialized_) initialize();
411 <          
412 <    DensityData ddat;
411 >        
412 >    // excluded interaction, so just return
413 >    if (idat.excluded) return;
414  
415 <    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;
415 >    int& iHash = iHash_[idat.atid1][idat.atid2];
416  
417 <    pair<AtomType*, AtomType*> key = make_pair(ddat.atype1, ddat.atype2);
418 <    set<NonBondedInteraction*>::iterator it;
417 >    if ((iHash & EAM_INTERACTION) != 0) eam_->calcDensity(idat);
418 >    if ((iHash & SC_INTERACTION) != 0)  sc_->calcDensity(idat);
419  
420 <    for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
421 <      if ((*it)->getFamily() == METALLIC_FAMILY) {
422 <        dynamic_cast<MetallicInteraction*>(*it)->calcDensity(ddat);
423 <      }
424 <    }
420 >    // set<NonBondedInteraction*>::iterator it;
421 >
422 >    // for (it = interactions_[ idat.atypes ].begin();
423 >    //      it != interactions_[ idat.atypes ].end(); ++it){
424 >    //   if ((*it)->getFamily() == METALLIC_FAMILY) {
425 >    //     dynamic_cast<MetallicInteraction*>(*it)->calcDensity(idat);
426 >    //   }
427 >    // }
428      
429      return;    
430    }
431    
432 <  void InteractionManager::doPreForce(int *atid, RealType *rho, RealType *frho, RealType *dfrhodrho){
432 >  void InteractionManager::doPreForce(SelfData &sdat){
433  
434      if (!initialized_) initialize();
435 <          
436 <    FunctionalData fdat;
435 >    
436 >    int& sHash = sHash_[sdat.atid];
437  
438 <    fdat.atype = typeMap_[*atid];
439 <    fdat.rho = *rho;
522 <    fdat.frho = *frho;
523 <    fdat.dfrhodrho = *dfrhodrho;
438 >    if ((sHash & EAM_INTERACTION) != 0) eam_->calcFunctional(sdat);
439 >    if ((sHash & SC_INTERACTION) != 0)  sc_->calcFunctional(sdat);
440  
441 <    pair<AtomType*, AtomType*> key = make_pair(fdat.atype, fdat.atype);
526 <    set<NonBondedInteraction*>::iterator it;
441 >    // set<NonBondedInteraction*>::iterator it;
442      
443 <    for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
444 <      if ((*it)->getFamily() == METALLIC_FAMILY) {
445 <        dynamic_cast<MetallicInteraction*>(*it)->calcFunctional(fdat);
446 <      }
447 <    }
443 >    // for (it = interactions_[atid1][atid2].begin(); it != interactions_[atid1][atid2].end(); ++it){
444 >    //   if ((*it)->getFamily() == METALLIC_FAMILY) {
445 >    //     dynamic_cast<MetallicInteraction*>(*it)->calcFunctional(sdat);
446 >    //   }
447 >    // }
448      
449      return;    
450    }
451  
452 <  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){
452 >  void InteractionManager::doPair(InteractionData &idat){
453      
454      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;
455  
456 <    pair<AtomType*, AtomType*> key = make_pair(idat.atype1, idat.atype2);
569 <    set<NonBondedInteraction*>::iterator it;
456 >    int& iHash = iHash_[idat.atid1][idat.atid2];
457  
458 <    for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
459 <      (*it)->calcForce(idat);
460 <    
461 <    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();
458 >    if ((iHash & ELECTROSTATIC_INTERACTION) != 0) electrostatic_->calcForce(idat);
459 >      
460 >    // electrostatics still has to worry about indirect
461 >    // contributions from excluded pairs of atoms, but nothing else does:
462  
463 <    return;    
587 <  }
463 >    if (idat.excluded) return;
464  
465 <  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){
465 >    if ((iHash & LJ_INTERACTION) != 0)             lj_->calcForce(idat);
466 >    if ((iHash & GB_INTERACTION) != 0)             gb_->calcForce(idat);
467 >    if ((iHash & STICKY_INTERACTION) != 0)         sticky_->calcForce(idat);
468 >    if ((iHash & MORSE_INTERACTION) != 0)          morse_->calcForce(idat);
469 >    if ((iHash & REPULSIVEPOWER_INTERACTION) != 0) repulsivePower_->calcForce(idat);
470 >    if ((iHash & EAM_INTERACTION) != 0)            eam_->calcForce(idat);
471 >    if ((iHash & SC_INTERACTION) != 0)             sc_->calcForce(idat);
472 >    if ((iHash & MAW_INTERACTION) != 0)            maw_->calcForce(idat);
473  
474 <    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);
474 >    // set<NonBondedInteraction*>::iterator it;
475  
476 <    pair<AtomType*, AtomType*> key = make_pair(skdat.atype1, skdat.atype2);
477 <    set<NonBondedInteraction*>::iterator it;
476 >    // for (it = interactions_[ idat.atypes ].begin();
477 >    //      it != interactions_[ idat.atypes ].end(); ++it) {
478  
479 <    for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
480 <      if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) {
616 <        dynamic_cast<ElectrostaticInteraction*>(*it)->calcSkipCorrection(skdat);
617 <      }
618 <    }
619 <    
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();
479 >    //   // electrostatics still has to worry about indirect
480 >    //   // contributions from excluded pairs of atoms:
481  
482 +    //   if (!idat.excluded || (*it)->getFamily() == ELECTROSTATIC_FAMILY) {
483 +    //     (*it)->calcForce(idat);
484 +    //   }
485 +    // }
486 +    
487      return;    
488    }
489  
490 <  void InteractionManager::doSelfCorrection(int *atid, RealType *eFrame, RealType *skippedCharge, RealType *pot, RealType *t){
490 >  void InteractionManager::doSelfCorrection(SelfData &sdat){
491  
492      if (!initialized_) initialize();
493      
494 <    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);
494 >    int& sHash = sHash_[sdat.atid];
495  
496 <    pair<AtomType*, AtomType*> key = make_pair(scdat.atype, scdat.atype);
648 <    set<NonBondedInteraction*>::iterator it;
496 >    if ((sHash & ELECTROSTATIC_INTERACTION) != 0) electrostatic_->calcSelfCorrection(sdat);
497  
498 <    for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
651 <      if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) {
652 <        dynamic_cast<ElectrostaticInteraction*>(*it)->calcSelfCorrection(scdat);
653 <      }
654 <    }
655 <        
656 <    t[0] = scdat.t.x();
657 <    t[1] = scdat.t.y();
658 <    t[2] = scdat.t.z();
498 >    // set<NonBondedInteraction*>::iterator it;
499  
500 +    // for (it = interactions_[atid1][atid2].begin(); it != interactions_[atid1][atid2].end(); ++it){
501 +    //   if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) {
502 +    //     dynamic_cast<ElectrostaticInteraction*>(*it)->calcSelfCorrection(sdat);
503 +    //   }
504 +    // }
505 +      
506      return;    
507    }
508  
509 +  void InteractionManager::doReciprocalSpaceSum(RealType &pot){
510 +    if (!initialized_) initialize();
511 +    electrostatic_->ReciprocalSpaceSum(pot);
512 +  }
513  
514    RealType InteractionManager::getSuggestedCutoffRadius(int *atid) {
515      if (!initialized_) initialize();
516 <    
516 >
517      AtomType* atype = typeMap_[*atid];
518  
669    pair<AtomType*, AtomType*> key = make_pair(atype, atype);
519      set<NonBondedInteraction*>::iterator it;
520      RealType cutoff = 0.0;
521      
522 <    for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
523 <      cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(atype, atype));  
522 >    for (it = interactions_[*atid][*atid].begin();
523 >         it != interactions_[*atid][*atid].end();
524 >         ++it)
525 >      cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(make_pair(atype, atype)));  
526      return cutoff;    
527    }
528  
529    RealType InteractionManager::getSuggestedCutoffRadius(AtomType* atype) {
530      if (!initialized_) initialize();
531      
532 <    pair<AtomType*, AtomType*> key = make_pair(atype, atype);
532 >    int atid = atype->getIdent();
533 >
534      set<NonBondedInteraction*>::iterator it;
535      RealType cutoff = 0.0;
536      
537 <    for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
538 <      cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(atype, atype));  
537 >    for (it = interactions_[atid][atid].begin();
538 >         it != interactions_[atid][atid].end(); ++it)
539 >      cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(make_pair(atype, atype)));  
540      return cutoff;    
541    }
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
542   } //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