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

Comparing:
branches/development/src/nonbonded/Sticky.cpp (file contents), Revision 1665 by gezelter, Tue Nov 22 20:38:56 2011 UTC vs.
trunk/src/nonbonded/Sticky.cpp (file contents), Revision 2071 by gezelter, Sat Mar 7 21:41:51 2015 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).          
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   */
# Line 46 | Line 46
46   #include <cmath>
47   #include "nonbonded/Sticky.hpp"
48   #include "nonbonded/LJ.hpp"
49 + #include "types/StickyAdapter.hpp"
50   #include "utils/simError.h"
51  
52   using namespace std;
53   namespace OpenMD {
54    
55 <  Sticky::Sticky() : name_("Sticky"), initialized_(false), forceField_(NULL) {}
55 <  
56 <  StickyParam Sticky::getStickyParam(AtomType* atomType) {
55 >  Sticky::Sticky() : initialized_(false), forceField_(NULL), name_("Sticky") {}
56      
58    // Do sanity checking on the AtomType we were passed before
59    // building any data structures:
60    if (!atomType->isSticky() && !atomType->isStickyPower()) {
61      sprintf( painCave.errMsg,
62               "Sticky::getStickyParam was passed an atomType (%s) that does\n"
63               "\tnot appear to be a Sticky atom.\n",
64               atomType->getName().c_str());
65      painCave.severity = OPENMD_ERROR;
66      painCave.isFatal = 1;
67      simError();
68    }
69    
70    DirectionalAtomType* daType = dynamic_cast<DirectionalAtomType*>(atomType);
71    GenericData* data = daType->getPropertyByName("Sticky");
72    if (data == NULL) {
73      sprintf( painCave.errMsg, "Sticky::getStickyParam could not find\n"
74               "\tSticky parameters for atomType %s.\n",
75               daType->getName().c_str());
76      painCave.severity = OPENMD_ERROR;
77      painCave.isFatal = 1;
78      simError();
79    }
80    
81    StickyParamGenericData* stickyData = dynamic_cast<StickyParamGenericData*>(data);
82    if (stickyData == NULL) {
83      sprintf( painCave.errMsg,
84               "Sticky::getStickyParam could not convert GenericData to\n"
85               "\tStickyParamGenericData for atom type %s\n",
86               daType->getName().c_str());
87      painCave.severity = OPENMD_ERROR;
88      painCave.isFatal = 1;
89      simError();          
90    }
91    
92    return stickyData->getData();
93  }
94  
57    void Sticky::initialize() {    
58 <    
59 <    ForceFieldOptions& fopts = forceField_->getForceFieldOptions();
60 <    ForceField::AtomTypeContainer* atomTypes = forceField_->getAtomTypes();
61 <    ForceField::AtomTypeContainer::MapTypeIterator i;
100 <    AtomType* at;
58 >    Stypes.clear();
59 >    Stids.clear();
60 >    MixingMap.clear();
61 >    nSticky_=0;
62  
63 +    Stids.resize( forceField_->getNAtomType(), -1);
64 +
65      // Sticky handles all of the Sticky-Sticky interactions
66  
67 <    for (at = atomTypes->beginType(i); at != NULL;
68 <         at = atomTypes->nextType(i)) {
69 <      
107 <      if (at->isSticky() || at->isStickyPower())
108 <        addType(at);
67 >    set<AtomType*>::iterator at;
68 >    for (at = simTypes_.begin(); at != simTypes_.end(); ++at) {
69 >      if ((*at)->isSticky()) nSticky_++;
70      }
71 +
72 +    MixingMap.resize(nSticky_);
73 +
74 +    for (at = simTypes_.begin(); at != simTypes_.end(); ++at) {
75 +      if ((*at)->isSticky()) addType( *at );
76 +    }
77      
78      initialized_ = true;
79    }
80        
81    void Sticky::addType(AtomType* atomType){
82 +    StickyAdapter sticky1 = StickyAdapter(atomType);
83 +
84      // add it to the map:
116    AtomTypeProperties atp = atomType->getATP();    
85      
86 <    pair<map<int,AtomType*>::iterator,bool> ret;    
87 <    ret = StickyMap.insert( pair<int, AtomType*>(atp.ident, atomType) );
86 >    int atid = atomType->getIdent();
87 >    int stid = Stypes.size();
88 >  
89 >    pair<set<int>::iterator,bool> ret;    
90 >    ret = Stypes.insert( atid );
91      if (ret.second == false) {
92        sprintf( painCave.errMsg,
93                 "Sticky already had a previous entry with ident %d\n",
94 <               atp.ident);
94 >               atid) ;
95        painCave.severity = OPENMD_INFO;
96        painCave.isFatal = 0;
97        simError();        
98      }
128    
129    RealType w0i, v0i, v0pi, rli, rui, rlpi, rupi;
130    
131    StickyParam sticky1 = getStickyParam(atomType);
99  
100 +    Stids[atid] = stid;
101 +    MixingMap[stid].resize( nSticky_ );    
102 +
103 +
104      // Now, iterate over all known types and add to the mixing map:
105      
106 <    map<int, AtomType*>::iterator it;
107 <    for( it = StickyMap.begin(); it != StickyMap.end(); ++it) {
106 >    std::set<int>::iterator it;
107 >    for( it = Stypes.begin(); it != Stypes.end(); ++it) {
108        
138      AtomType* atype2 = (*it).second;
139    
140      StickyParam sticky2 = getStickyParam(atype2);
109  
110 +      int stid2 = Stids[ (*it) ];
111 +      AtomType* atype2 = forceField_->getAtomType( (*it) );
112 +      StickyAdapter sticky2 = StickyAdapter(atype2);
113 +
114        StickyInteractionData mixer;        
115            
116        // Mixing two different sticky types is silly, but if you want 2
# Line 146 | Line 118 | namespace OpenMD {
118        // Lorentz- Berthelot mixing rules (which happen to do the right thing
119        // when atomType and atype2 happen to be the same.
120        
121 <      mixer.rl   = 0.5 * ( sticky1.rl + sticky2.rl );
122 <      mixer.ru   = 0.5 * ( sticky1.ru + sticky2.ru );
123 <      mixer.rlp  = 0.5 * ( sticky1.rlp + sticky2.rlp );
124 <      mixer.rup  = 0.5 * ( sticky1.rup + sticky2.rup );
121 >      mixer.rl   = 0.5 * ( sticky1.getRl() + sticky2.getRl() );
122 >      mixer.ru   = 0.5 * ( sticky1.getRu() + sticky2.getRu() );
123 >      mixer.rlp  = 0.5 * ( sticky1.getRlp() + sticky2.getRlp() );
124 >      mixer.rup  = 0.5 * ( sticky1.getRup() + sticky2.getRup() );
125        mixer.rbig = max(mixer.ru, mixer.rup);
126 <      mixer.w0  = sqrt( sticky1.w0   * sticky2.w0  );
127 <      mixer.v0  = sqrt( sticky1.v0   * sticky2.v0  );
128 <      mixer.v0p = sqrt( sticky1.v0p  * sticky2.v0p );
129 <      mixer.isPower = atomType->isStickyPower() && atype2->isStickyPower();
126 >      mixer.w0  = sqrt( sticky1.getW0()   * sticky2.getW0()  );
127 >      mixer.v0  = sqrt( sticky1.getV0()   * sticky2.getV0()  );
128 >      mixer.v0p = sqrt( sticky1.getV0p()  * sticky2.getV0p() );
129 >      mixer.isPower = sticky1.isStickyPower() && sticky2.isStickyPower();
130  
131        CubicSpline* s = new CubicSpline();
132        s->addPoint(mixer.rl, 1.0);
# Line 166 | Line 138 | namespace OpenMD {
138        sp->addPoint(mixer.rup, 0.0);
139        mixer.sp = sp;
140  
141 <      
142 <      pair<AtomType*, AtomType*> key1, key2;
143 <      key1 = make_pair(atomType, atype2);
144 <      key2 = make_pair(atype2, atomType);
145 <      
174 <      MixingMap[key1] = mixer;
175 <      if (key2 != key1) {
176 <        MixingMap[key2] = mixer;
141 >      MixingMap[stid2].resize( nSticky_ );
142 >
143 >      MixingMap[stid][stid2] = mixer;
144 >      if (stid2 != stid) {
145 >        MixingMap[stid2][stid] = mixer;
146        }
147      }
148    }
# Line 191 | Line 160 | namespace OpenMD {
160    
161      if (!initialized_) initialize();
162      
163 <    map<pair<AtomType*, AtomType*>, StickyInteractionData>::iterator it;
164 <    it = MixingMap.find(idat.atypes);
165 <    if (it != MixingMap.end()) {
166 <
167 <      StickyInteractionData mixer = (*it).second;
163 >    StickyInteractionData &mixer = MixingMap[Stids[idat.atid1]][Stids[idat.atid2]];
164 >    
165 >    RealType w0  = mixer.w0;
166 >    RealType v0  = mixer.v0;
167 >    RealType v0p = mixer.v0p;
168 >    RealType rl  = mixer.rl;
169 >    RealType ru  = mixer.ru;
170 >    RealType rlp = mixer.rlp;
171 >    RealType rup = mixer.rup;
172 >    RealType rbig = mixer.rbig;
173 >    bool isPower = mixer.isPower;
174 >    
175 >    if ( *(idat.rij) <= rbig) {
176        
177 <      RealType w0  = mixer.w0;
178 <      RealType v0  = mixer.v0;
202 <      RealType v0p = mixer.v0p;
203 <      RealType rl  = mixer.rl;
204 <      RealType ru  = mixer.ru;
205 <      RealType rlp = mixer.rlp;
206 <      RealType rup = mixer.rup;
207 <      RealType rbig = mixer.rbig;
208 <      bool isPower = mixer.isPower;
177 >      RealType r3 = *(idat.r2) * *(idat.rij);
178 >      RealType r5 = r3 * *(idat.r2);
179        
180 <      if ( *(idat.rij) <= rbig) {
181 <        
182 <        RealType r3 = *(idat.r2) * *(idat.rij);
183 <        RealType r5 = r3 * *(idat.r2);
184 <        
185 <        RotMat3x3d A1trans = idat.A1->transpose();
186 <        RotMat3x3d A2trans = idat.A2->transpose();
187 <        
188 <        // rotate the inter-particle separation into the two different
189 <        // body-fixed coordinate systems:
190 <        
191 <        Vector3d ri = *(idat.A1) * *(idat.d);
192 <        
193 <        // negative sign because this is the vector from j to i:
194 <        
195 <        Vector3d rj = - *(idat.A2) * *(idat.d);
196 <        
197 <        RealType xi = ri.x();
198 <        RealType yi = ri.y();
199 <        RealType zi = ri.z();
200 <        
201 <        RealType xj = rj.x();
202 <        RealType yj = rj.y();
203 <        RealType zj = rj.z();
204 <        
205 <        RealType xi2 = xi * xi;
206 <        RealType yi2 = yi * yi;
207 <        RealType zi2 = zi * zi;
208 <        
209 <        RealType xj2 = xj * xj;
210 <        RealType yj2 = yj * yj;
211 <        RealType zj2 = zj * zj;    
212 <        
213 <        // calculate the switching info. from the splines
214 <        
215 <        RealType s = 0.0;
216 <        RealType dsdr = 0.0;
217 <        RealType sp = 0.0;
218 <        RealType dspdr = 0.0;
219 <        
220 <        if ( *(idat.rij) < ru) {
221 <          if ( *(idat.rij) < rl) {
252 <            s = 1.0;
253 <            dsdr = 0.0;
254 <          } else {          
255 <            // we are in the switching region
256 <            
257 <            pair<RealType, RealType> res = mixer.s->getValueAndDerivativeAt(*(idat.rij));
258 <            s = res.first;
259 <            dsdr = res.second;
260 <          }
261 <        }
262 <        
263 <        if (*(idat.rij) < rup) {
264 <          if ( *(idat.rij) < rlp) {
265 <            sp = 1.0;
266 <            dspdr = 0.0;
267 <          } else {
268 <            // we are in the switching region
269 <            
270 <            pair<RealType, RealType> res =mixer.sp->getValueAndDerivativeAt( *(idat.rij));
271 <            sp = res.first;
272 <            dspdr = res.second;
273 <          }
180 >      RotMat3x3d A1trans = idat.A1->transpose();
181 >      RotMat3x3d A2trans = idat.A2->transpose();
182 >      
183 >      // rotate the inter-particle separation into the two different
184 >      // body-fixed coordinate systems:
185 >      
186 >      Vector3d ri = *(idat.A1) * *(idat.d);
187 >      
188 >      // negative sign because this is the vector from j to i:
189 >      
190 >      Vector3d rj = - *(idat.A2) * *(idat.d);
191 >      
192 >      RealType xi = ri.x();
193 >      RealType yi = ri.y();
194 >      RealType zi = ri.z();
195 >      
196 >      RealType xj = rj.x();
197 >      RealType yj = rj.y();
198 >      RealType zj = rj.z();
199 >      
200 >      RealType xi2 = xi * xi;
201 >      RealType yi2 = yi * yi;
202 >      RealType zi2 = zi * zi;
203 >      
204 >      RealType xj2 = xj * xj;
205 >      RealType yj2 = yj * yj;
206 >      RealType zj2 = zj * zj;    
207 >      
208 >      // calculate the switching info. from the splines
209 >      
210 >      RealType s = 0.0;
211 >      RealType dsdr = 0.0;
212 >      RealType sp = 0.0;
213 >      RealType dspdr = 0.0;
214 >      
215 >      if ( *(idat.rij) < ru) {
216 >        if ( *(idat.rij) < rl) {
217 >          s = 1.0;
218 >          dsdr = 0.0;
219 >        } else {          
220 >          // we are in the switching region          
221 >          mixer.s->getValueAndDerivativeAt(*(idat.rij), s, dsdr);
222          }
223 <        
224 <        RealType wi = 2.0*(xi2-yi2)*zi / r3;
225 <        RealType wj = 2.0*(xj2-yj2)*zj / r3;
226 <        RealType w = wi+wj;
227 <        
280 <        
281 <        RealType zif = zi/ *(idat.rij)  - 0.6;
282 <        RealType zis = zi/ *(idat.rij)  + 0.8;
283 <        
284 <        RealType zjf = zj/ *(idat.rij)  - 0.6;
285 <        RealType zjs = zj/ *(idat.rij)  + 0.8;
286 <        
287 <        RealType wip = zif*zif*zis*zis - w0;
288 <        RealType wjp = zjf*zjf*zjs*zjs - w0;
289 <        RealType wp = wip + wjp;
290 <        
291 <        Vector3d dwi(4.0*xi*zi/r3  - 6.0*xi*zi*(xi2-yi2)/r5,
292 <                     - 4.0*yi*zi/r3  - 6.0*yi*zi*(xi2-yi2)/r5,
293 <                     2.0*(xi2-yi2)/r3  - 6.0*zi2*(xi2-yi2)/r5);
294 <        
295 <        Vector3d dwj(4.0*xj*zj/r3  - 6.0*xj*zj*(xj2-yj2)/r5,
296 <                     - 4.0*yj*zj/r3  - 6.0*yj*zj*(xj2-yj2)/r5,
297 <                     2.0*(xj2-yj2)/r3  - 6.0*zj2*(xj2-yj2)/r5);
298 <        
299 <        RealType uglyi = zif*zif*zis + zif*zis*zis;
300 <        RealType uglyj = zjf*zjf*zjs + zjf*zjs*zjs;
301 <
302 <        Vector3d dwip(-2.0*xi*zi*uglyi/r3,
303 <                      -2.0*yi*zi*uglyi/r3,
304 <                      2.0*(1.0/ *(idat.rij)  - zi2/r3)*uglyi);
305 <        
306 <        Vector3d dwjp(-2.0*xj*zj*uglyj/r3,
307 <                      -2.0*yj*zj*uglyj/r3,
308 <                      2.0*(1.0/ *(idat.rij)  - zj2/r3)*uglyj);
309 <        
310 <        Vector3d dwidu(4.0*(yi*zi2 + 0.5*yi*(xi2-yi2))/r3,
311 <                       4.0*(xi*zi2 - 0.5*xi*(xi2-yi2))/r3,
312 <                       - 8.0*xi*yi*zi/r3);
313 <        
314 <        Vector3d dwjdu(4.0*(yj*zj2 + 0.5*yj*(xj2-yj2))/r3,
315 <                       4.0*(xj*zj2 - 0.5*xj*(xj2-yj2))/r3,
316 <                       - 8.0*xj*yj*zj/r3);
317 <        
318 <        Vector3d dwipdu(2.0*yi*uglyi/ *(idat.rij) ,
319 <                        -2.0*xi*uglyi/ *(idat.rij) ,
320 <                        0.0);
321 <        
322 <        Vector3d dwjpdu(2.0*yj*uglyj/ *(idat.rij) ,
323 <                        -2.0*xj*uglyj/ *(idat.rij) ,
324 <                        0.0);
325 <        
326 <        if (isPower) {
327 <          RealType frac1 = 0.25;
328 <          RealType frac2 = 0.75;      
329 <          RealType wi2 = wi*wi;
330 <          RealType wj2 = wj*wj;
331 <          // sticky power has no w' function:
332 <          w = frac1 * wi * wi2 + frac2*wi + frac1*wj*wj2 + frac2*wj + v0p;
333 <          wp = 0.0;
334 <          dwi = frac1*3.0*wi2*dwi + frac2*dwi;
335 <          dwj = frac1*3.0*wj2*dwi + frac2*dwi;
336 <          dwip = V3Zero;
337 <          dwjp = V3Zero;
338 <          dwidu = frac1*3.0*wi2*dwidu + frac2*dwidu;
339 <          dwidu = frac1*3.0*wj2*dwjdu + frac2*dwjdu;
340 <          dwipdu = V3Zero;
341 <          dwjpdu = V3Zero;
342 <          sp = 0.0;
223 >      }
224 >      
225 >      if (*(idat.rij) < rup) {
226 >        if ( *(idat.rij) < rlp) {
227 >          sp = 1.0;
228            dspdr = 0.0;
229 +        } else {
230 +          // we are in the switching region          
231 +          mixer.sp->getValueAndDerivativeAt( *(idat.rij), sp, dspdr);
232          }
345        
346        *(idat.vpair) += 0.5*(v0*s*w + v0p*sp*wp);
347        (*(idat.pot))[HYDROGENBONDING_FAMILY] += 0.5*(v0*s*w + v0p*sp*wp)* *(idat.sw) ;
348        
349        // do the torques first since they are easy:
350        // remember that these are still in the body-fixed axes
351        
352        Vector3d ti = 0.5* *(idat.sw) *(v0*s*dwidu + v0p*sp*dwipdu);
353        Vector3d tj = 0.5* *(idat.sw) *(v0*s*dwjdu + v0p*sp*dwjpdu);
354        
355        // go back to lab frame using transpose of rotation matrix:
356        
357        *(idat.t1) += A1trans * ti;
358        *(idat.t2) += A2trans * tj;
359        
360        // Now, on to the forces:
361        
362        // first rotate the i terms back into the lab frame:
363        
364        Vector3d radcomi = (v0 * s * dwi + v0p * sp * dwip) *  *(idat.sw);
365        Vector3d radcomj = (v0 * s * dwj + v0p * sp * dwjp) *  *(idat.sw);
366        
367        Vector3d fii = A1trans * radcomi;
368        Vector3d fjj = A2trans * radcomj;
369        
370        // now assemble these with the radial-only terms:
371        
372        *(idat.f1) += 0.5 * ((v0*dsdr*w + v0p*dspdr*wp) * *(idat.d) /
373                             *(idat.rij)  + fii - fjj);
374        
233        }
234 +      
235 +      
236 +      RealType wi = 2.0*(xi2-yi2)*zi / r3;
237 +      RealType wj = 2.0*(xj2-yj2)*zj / r3;
238 +      RealType w = wi+wj;
239 +      
240 +      
241 +      RealType zif = zi/ *(idat.rij)  - 0.6;
242 +      RealType zis = zi/ *(idat.rij)  + 0.8;
243 +      
244 +      RealType zjf = zj/ *(idat.rij)  - 0.6;
245 +      RealType zjs = zj/ *(idat.rij)  + 0.8;
246 +      
247 +      RealType wip = zif*zif*zis*zis - w0;
248 +      RealType wjp = zjf*zjf*zjs*zjs - w0;
249 +      RealType wp = wip + wjp;
250 +      
251 +      Vector3d dwi(4.0*xi*zi/r3  - 6.0*xi*zi*(xi2-yi2)/r5,
252 +                   - 4.0*yi*zi/r3  - 6.0*yi*zi*(xi2-yi2)/r5,
253 +                   2.0*(xi2-yi2)/r3  - 6.0*zi2*(xi2-yi2)/r5);
254 +      
255 +      Vector3d dwj(4.0*xj*zj/r3  - 6.0*xj*zj*(xj2-yj2)/r5,
256 +                   - 4.0*yj*zj/r3  - 6.0*yj*zj*(xj2-yj2)/r5,
257 +                   2.0*(xj2-yj2)/r3  - 6.0*zj2*(xj2-yj2)/r5);
258 +      
259 +      RealType uglyi = zif*zif*zis + zif*zis*zis;
260 +      RealType uglyj = zjf*zjf*zjs + zjf*zjs*zjs;
261 +      
262 +      Vector3d dwip(-2.0*xi*zi*uglyi/r3,
263 +                    -2.0*yi*zi*uglyi/r3,
264 +                    2.0*(1.0/ *(idat.rij)  - zi2/r3)*uglyi);
265 +      
266 +      Vector3d dwjp(-2.0*xj*zj*uglyj/r3,
267 +                    -2.0*yj*zj*uglyj/r3,
268 +                    2.0*(1.0/ *(idat.rij)  - zj2/r3)*uglyj);
269 +      
270 +      Vector3d dwidu(4.0*(yi*zi2 + 0.5*yi*(xi2-yi2))/r3,
271 +                     4.0*(xi*zi2 - 0.5*xi*(xi2-yi2))/r3,
272 +                     - 8.0*xi*yi*zi/r3);
273 +      
274 +      Vector3d dwjdu(4.0*(yj*zj2 + 0.5*yj*(xj2-yj2))/r3,
275 +                     4.0*(xj*zj2 - 0.5*xj*(xj2-yj2))/r3,
276 +                     - 8.0*xj*yj*zj/r3);
277 +      
278 +      Vector3d dwipdu(2.0*yi*uglyi/ *(idat.rij) ,
279 +                      -2.0*xi*uglyi/ *(idat.rij) ,
280 +                      0.0);
281 +      
282 +      Vector3d dwjpdu(2.0*yj*uglyj/ *(idat.rij) ,
283 +                      -2.0*xj*uglyj/ *(idat.rij) ,
284 +                      0.0);
285 +      
286 +      if (isPower) {
287 +        cerr << "This is probably an error!\n";
288 +        RealType frac1 = 0.25;
289 +        RealType frac2 = 0.75;      
290 +        RealType wi2 = wi*wi;
291 +        RealType wj2 = wj*wj;
292 +        // sticky power has no w' function:
293 +        w = frac1 * wi * wi2 + frac2*wi + frac1*wj*wj2 + frac2*wj + v0p;
294 +        wp = 0.0;
295 +        dwi = frac1*RealType(3.0)*wi2*dwi + frac2*dwi;
296 +        dwj = frac1*RealType(3.0)*wj2*dwi + frac2*dwi;
297 +        dwip = V3Zero;
298 +        dwjp = V3Zero;
299 +        dwidu = frac1*RealType(3.0)*wi2*dwidu + frac2*dwidu;
300 +        dwidu = frac1*RealType(3.0)*wj2*dwjdu + frac2*dwjdu;
301 +        dwipdu = V3Zero;
302 +        dwjpdu = V3Zero;
303 +        sp = 0.0;
304 +        dspdr = 0.0;
305 +      }
306 +      
307 +      
308 +      
309 +      *(idat.vpair) += RealType(0.5)*(v0*s*w + v0p*sp*wp);
310 +      (*(idat.pot))[HYDROGENBONDING_FAMILY] += RealType(0.5)*(v0*s*w + v0p*sp*wp)* *(idat.sw) ;
311 +      
312 +      // do the torques first since they are easy:
313 +      // remember that these are still in the body-fixed axes
314 +      
315 +      Vector3d ti = RealType(0.5)* *(idat.sw) *(v0*s*dwidu + v0p*sp*dwipdu);
316 +      Vector3d tj = RealType(0.5)* *(idat.sw) *(v0*s*dwjdu + v0p*sp*dwjpdu);
317 +      
318 +      // go back to lab frame using transpose of rotation matrix:
319 +      
320 +      *(idat.t1) += A1trans * ti;
321 +      *(idat.t2) += A2trans * tj;
322 +      
323 +      // Now, on to the forces:
324 +      
325 +      // first rotate the i terms back into the lab frame:
326 +      
327 +      Vector3d radcomi = (v0 * s * dwi + v0p * sp * dwip) *  *(idat.sw);
328 +      Vector3d radcomj = (v0 * s * dwj + v0p * sp * dwjp) *  *(idat.sw);
329 +      
330 +      Vector3d fii = A1trans * radcomi;
331 +      Vector3d fjj = A2trans * radcomj;
332 +      
333 +      // now assemble these with the radial-only terms:
334 +      
335 +      *(idat.f1) += RealType(0.5) * ((v0*dsdr*w + v0p*dspdr*wp) * *(idat.d) /
336 +                                     *(idat.rij)  + fii - fjj);
337 +      
338      }
339      
340      return;      
341    }
342 <
342 >  
343    RealType Sticky::getSuggestedCutoffRadius(pair<AtomType*, AtomType*> atypes) {
344      if (!initialized_) initialize();  
345 <    map<pair<AtomType*, AtomType*>, StickyInteractionData>::iterator it;
346 <    it = MixingMap.find(atypes);
347 <    if (it == MixingMap.end())
348 <      return 0.0;
349 <    else  {
350 <      StickyInteractionData mixer = (*it).second;
351 <      return mixer.rbig;
345 >    int atid1 = atypes.first->getIdent();
346 >    int atid2 = atypes.second->getIdent();
347 >    int stid1 = Stids[atid1];
348 >    int stid2 = Stids[atid2];
349 >    
350 >    if (stid1 == -1 || stid2 == -1) return 0.0;
351 >    else {      
352 >      return MixingMap[stid1][stid2].rbig;
353      }
354    }
355   }

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines