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

Comparing branches/development/src/nonbonded/GB.cpp (file contents):
Revision 1505 by gezelter, Sun Oct 3 22:18:59 2010 UTC vs.
Revision 1710 by gezelter, Fri May 18 21:44:02 2012 UTC

# Line 36 | Line 36
36   * [1]  Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).            
37   * [2]  Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).          
38   * [3]  Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).          
39 < * [4]  Vardeman & Gezelter, in progress (2009).                        
39 > * [4]  Kuang & Gezelter,  J. Chem. Phys. 133, 164101 (2010).
40 > * [5]  Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41   */
42  
43   #include <stdio.h>
# Line 45 | Line 46
46   #include <cmath>
47   #include "nonbonded/GB.hpp"
48   #include "utils/simError.h"
49 + #include "types/LennardJonesAdapter.hpp"
50 + #include "types/GayBerneAdapter.hpp"
51  
52   using namespace std;
53   namespace OpenMD {
54  
55 +  /* GB is the Gay-Berne interaction for ellipsoidal particles.  The original
56 +   * paper (for identical uniaxial particles) is:
57 +   *    J. G. Gay and B. J. Berne, J. Chem. Phys., 74, 3316-3319 (1981).
58 +   * A more-general GB potential for dissimilar uniaxial particles:
59 +   *    D. J. Cleaver, C. M. Care, M. P. Allen and M. P. Neal, Phys. Rev. E,
60 +   *    54, 559-567 (1996).
61 +   * Further parameterizations can be found in:
62 +   *    A. P. J. Emerson, G. R. Luckhurst and S. G. Whatling, Mol. Phys.,
63 +   *    82, 113-124 (1994).
64 +   * And a nice force expression:
65 +   *    G. R. Luckhurst and R. A. Stephens, Liq. Cryst. 8, 451-464 (1990).
66 +   * Even clearer force and torque expressions:
67 +   *    P. A. Golubkov and P. Y. Ren, J. Chem. Phys., 125, 64103 (2006).
68 +   * New expressions for cross interactions of strength parameters:
69 +   *    J. Wu, X. Zhen, H. Shen, G. Li, and P. Ren, J. Chem. Phys.,
70 +   *    135, 155104 (2011).
71 +   *
72 +   * In this version of the GB interaction, each uniaxial ellipsoidal type
73 +   * is described using a set of 6 parameters:
74 +   *  d:  range parameter for side-by-side (S) and cross (X) configurations
75 +   *  l:  range parameter for end-to-end (E) configuration
76 +   *  epsilon_X:  well-depth parameter for cross (X) configuration
77 +   *  epsilon_S:  well-depth parameter for side-by-side (S) configuration
78 +   *  epsilon_E:  well depth parameter for end-to-end (E) configuration
79 +   *  dw: "softness" of the potential
80 +   *
81 +   * Additionally, there are two "universal" paramters to govern the overall
82 +   * importance of the purely orientational (nu) and the mixed
83 +   * orientational / translational (mu) parts of strength of the interactions.
84 +   * These parameters have default or "canonical" values, but may be changed
85 +   * as a force field option:
86 +   * nu_: purely orientational part : defaults to 1
87 +   * mu_: mixed orientational / translational part : defaults to 2
88 +   */
89 +
90 +
91    GB::GB() : name_("GB"), initialized_(false), mu_(2.0), nu_(1.0), forceField_(NULL) {}
53  
54  GayBerneParam GB::getGayBerneParam(AtomType* atomType) {
92      
56    // Do sanity checking on the AtomType we were passed before
57    // building any data structures:
58    if (!atomType->isGayBerne()) {
59      sprintf( painCave.errMsg,
60               "GB::getGayBerneParam was passed an atomType (%s) that does\n"
61               "\tnot appear to be a Gay-Berne atom.\n",
62               atomType->getName().c_str());
63      painCave.severity = OPENMD_ERROR;
64      painCave.isFatal = 1;
65      simError();
66    }
67    
68    DirectionalAtomType* daType = dynamic_cast<DirectionalAtomType*>(atomType);
69    GenericData* data = daType->getPropertyByName("GayBerne");
70    if (data == NULL) {
71      sprintf( painCave.errMsg, "GB::getGayBerneParam could not find\n"
72               "\tGay-Berne parameters for atomType %s.\n",
73               daType->getName().c_str());
74      painCave.severity = OPENMD_ERROR;
75      painCave.isFatal = 1;
76      simError();
77    }
78    
79    GayBerneParamGenericData* gbData = dynamic_cast<GayBerneParamGenericData*>(data);
80    if (gbData == NULL) {
81      sprintf( painCave.errMsg,
82               "GB::getGayBerneParam could not convert GenericData to\n"
83               "\tGayBerneParamGenericData for atom type %s\n",
84               daType->getName().c_str());
85      painCave.severity = OPENMD_ERROR;
86      painCave.isFatal = 1;
87      simError();          
88    }
89    
90    return gbData->getData();
91  }
92
93  LJParam GB::getLJParam(AtomType* atomType) {
94    
95    // Do sanity checking on the AtomType we were passed before
96    // building any data structures:
97    if (!atomType->isLennardJones()) {
98      sprintf( painCave.errMsg,
99               "GB::getLJParam was passed an atomType (%s) that does not\n"
100               "\tappear to be a Lennard-Jones atom.\n",
101               atomType->getName().c_str());
102      painCave.severity = OPENMD_ERROR;
103      painCave.isFatal = 1;
104      simError();
105    }
106    
107    GenericData* data = atomType->getPropertyByName("LennardJones");
108    if (data == NULL) {
109      sprintf( painCave.errMsg, "GB::getLJParam could not find Lennard-Jones\n"
110               "\tparameters for atomType %s.\n", atomType->getName().c_str());
111      painCave.severity = OPENMD_ERROR;
112      painCave.isFatal = 1;
113      simError();
114    }
115    
116    LJParamGenericData* ljData = dynamic_cast<LJParamGenericData*>(data);
117    if (ljData == NULL) {
118      sprintf( painCave.errMsg,
119               "GB::getLJParam could not convert GenericData to LJParam for\n"
120               "\tatom type %s\n", atomType->getName().c_str());
121      painCave.severity = OPENMD_ERROR;
122      painCave.isFatal = 1;
123      simError();          
124    }
125    
126    return ljData->getData();
127  }
128  
129  RealType GB::getLJEpsilon(AtomType* atomType) {    
130    LJParam ljParam = getLJParam(atomType);
131    return ljParam.epsilon;
132  }
133  RealType GB::getLJSigma(AtomType* atomType) {    
134    LJParam ljParam = getLJParam(atomType);
135    return ljParam.sigma;
136  }
137  
93    void GB::initialize() {    
94      
95      ForceFieldOptions& fopts = forceField_->getForceFieldOptions();
# Line 150 | Line 105 | namespace OpenMD {
105      for (at = atomTypes->beginType(i); at != NULL;
106           at = atomTypes->nextType(i)) {
107        
108 <      if (at->isGayBerne() || at->isLennardJones())
108 >      LennardJonesAdapter lja = LennardJonesAdapter(at);
109 >      GayBerneAdapter gba = GayBerneAdapter(at);
110 >
111 >      if (gba.isGayBerne() || lja.isLennardJones())
112          addType(at);
113      }
114    
# Line 159 | Line 117 | namespace OpenMD {
117        
118    void GB::addType(AtomType* atomType){
119      // add it to the map:
162    AtomTypeProperties atp = atomType->getATP();    
120  
121      pair<map<int,AtomType*>::iterator,bool> ret;    
122 <    ret = GBMap.insert( pair<int, AtomType*>(atp.ident, atomType) );
122 >    ret = GBMap.insert( pair<int, AtomType*>(atomType->getIdent(), atomType) );
123      if (ret.second == false) {
124        sprintf( painCave.errMsg,
125                 "GB already had a previous entry with ident %d\n",
126 <               atp.ident);
126 >               atomType->getIdent() );
127        painCave.severity = OPENMD_INFO;
128        painCave.isFatal = 0;
129        simError();        
130      }
131      
132 <    RealType d1, l1, e1, er1, dw1;
133 <    
134 <    if (atomType->isGayBerne()) {
135 <      GayBerneParam gb1 = getGayBerneParam(atomType);
136 <      d1 = gb1.GB_d;
137 <      l1 = gb1.GB_l;
138 <      e1 = gb1.GB_eps;
139 <      er1 = gb1.GB_eps_ratio;
140 <      dw1 = gb1.GB_dw;
141 <    } else if (atomType->isLennardJones()) {
142 <      d1 = getLJSigma(atomType) / sqrt(2.0);
143 <      e1 = getLJEpsilon(atomType);
132 >    RealType d1, l1, eX1, eS1, eE1, dw1;
133 >        
134 >    LennardJonesAdapter lja1 = LennardJonesAdapter(atomType);
135 >    GayBerneAdapter gba1 = GayBerneAdapter(atomType);
136 >    if (gba1.isGayBerne()) {
137 >      d1 = gba1.getD();
138 >      l1 = gba1.getL();
139 >      eX1 = gba1.getEpsX();
140 >      eS1 = gba1.getEpsS();
141 >      eE1 = gba1.getEpsE();
142 >      dw1 = gba1.getDw();
143 >    } else if (lja1.isLennardJones()) {
144 >      d1 = lja1.getSigma() / sqrt(2.0);
145        l1 = d1;
146 <      er1 = 1.0;
146 >      eX1 = lja1.getEpsilon();
147 >      eS1 = eX1;
148 >      eE1 = eX1;
149        dw1 = 1.0;      
150      } else {
151        sprintf( painCave.errMsg,
# Line 204 | Line 164 | namespace OpenMD {
164      for( it = GBMap.begin(); it != GBMap.end(); ++it) {
165        
166        AtomType* atype2 = (*it).second;
167 +      LennardJonesAdapter lja2 = LennardJonesAdapter(atype2);
168 +      GayBerneAdapter gba2 = GayBerneAdapter(atype2);
169 +      RealType d2, l2, eX2, eS2, eE2, dw2;
170        
171 <      RealType d2, l2, e2, er2, dw2;
172 <      
173 <      if (atype2->isGayBerne()) {
174 <        GayBerneParam gb2 = getGayBerneParam(atype2);
175 <        d2 = gb2.GB_d;
176 <        l2 = gb2.GB_l;
177 <        e2 = gb2.GB_eps;
178 <        er2 = gb2.GB_eps_ratio;
179 <        dw2 = gb2.GB_dw;
217 <      } else if (atype2->isLennardJones()) {
218 <        d2 = getLJSigma(atype2) / sqrt(2.0);
219 <        e2 = getLJEpsilon(atype2);
171 >      if (gba2.isGayBerne()) {
172 >        d2 = gba2.getD();
173 >        l2 = gba2.getL();
174 >        eX2 = gba2.getEpsX();
175 >        eS2 = gba2.getEpsS();
176 >        eE2 = gba2.getEpsE();
177 >        dw2 = gba2.getDw();
178 >      } else if (lja2.isLennardJones()) {
179 >        d2 = lja2.getSigma() / sqrt(2.0);
180          l2 = d2;
181 <        er2 = 1.0;
181 >        eX2 = lja2.getEpsilon();
182 >        eS2 = eX2;
183 >        eE2 = eX2;
184          dw2 = 1.0;
185        }
186                        
187 <      GBInteractionData mixer;        
187 >      GBInteractionData mixer1, mixer2;    
188        
189        //  Cleaver paper uses sqrt of squares to get sigma0 for
190        //  mixed interactions.
191              
192 <      mixer.sigma0 = sqrt(d1*d1 + d2*d2);
193 <      mixer.xa2 = (l1*l1 - d1*d1)/(l1*l1 + d2*d2);
194 <      mixer.xai2 = (l2*l2 - d2*d2)/(l2*l2 + d1*d1);
195 <      mixer.x2 = (l1*l1 - d1*d1) * (l2*l2 - d2*d2) /
192 >      mixer1.sigma0 = sqrt(d1*d1 + d2*d2);
193 >      mixer1.xa2 = (l1*l1 - d1*d1)/(l1*l1 + d2*d2);
194 >      mixer1.xai2 = (l2*l2 - d2*d2)/(l2*l2 + d1*d1);
195 >      mixer1.x2 = (l1*l1 - d1*d1) * (l2*l2 - d2*d2) /
196          ((l2*l2 + d1*d1) * (l1*l1 + d2*d2));
197 +
198 +      mixer2.sigma0 = mixer1.sigma0;
199 +      // xa2 and xai2 for j-i pairs are reversed from the same i-j pairing.
200 +      // Swapping the particles reverses the anisotropy parameters:
201 +      mixer2.xa2 = mixer1.xai2;
202 +      mixer2.xai2 = mixer1.xa2;
203 +      mixer2.x2 = mixer1.x2;
204  
205        // assumed LB mixing rules for now:
206  
207 <      mixer.dw = 0.5 * (dw1 + dw2);
208 <      mixer.eps0 = sqrt(e1 * e2);
207 >      mixer1.dw = 0.5 * (dw1 + dw2);
208 >      mixer1.eps0 = sqrt(eX1 * eX2);
209 >
210 >      mixer2.dw = mixer1.dw;
211 >      mixer2.eps0 = mixer1.eps0;
212 >
213 >      RealType mi = RealType(1.0)/mu_;
214        
215 <      RealType er = sqrt(er1 * er2);
216 <      RealType ermu = pow(er,(1.0 / mu_));
217 <      RealType xp = (1.0 - ermu) / (1.0 + ermu);
218 <      RealType ap2 = 1.0 / (1.0 + ermu);
245 <      
246 <      mixer.xp2 = xp * xp;
247 <      mixer.xpap2 = xp * ap2;
248 <      mixer.xpapi2 = xp / ap2;
215 >      mixer1.xpap2  = (pow(eS1, mi) - pow(eE1, mi)) / (pow(eS1, mi) + pow(eE2, mi));
216 >      mixer1.xpapi2 = (pow(eS2, mi) - pow(eE2, mi)) / (pow(eS2, mi) + pow(eE1, mi));
217 >      mixer1.xp2    = (pow(eS1, mi) - pow(eE1, mi)) * (pow(eS2, mi) - pow(eE2, mi))  /
218 >        (pow(eS2, mi) + pow(eE1, mi)) / (pow(eS1, mi) + pow(eE2, mi)) ;
219  
220 +      // xpap2 and xpapi2 for j-i pairs are reversed from the same i-j pairing.
221 +      // Swapping the particles reverses the anisotropy parameters:
222 +      mixer2.xpap2 = mixer1.xpapi2;
223 +      mixer2.xpapi2 = mixer1.xpap2;
224 +      mixer2.xp2 = mixer1.xp2;
225 +
226        // only add this pairing if at least one of the atoms is a Gay-Berne atom
227  
228 <      if (atomType->isGayBerne() || atype2->isGayBerne()) {
228 >      if (gba1.isGayBerne() || gba2.isGayBerne()) {
229  
230          pair<AtomType*, AtomType*> key1, key2;
231          key1 = make_pair(atomType, atype2);
232          key2 = make_pair(atype2, atomType);
233          
234 <        MixingMap[key1] = mixer;
234 >        MixingMap[key1] = mixer1;
235          if (key2 != key1) {
236 <          MixingMap[key2] = mixer;
236 >          MixingMap[key2] = mixer2;
237          }
238        }
239      }      
240    }
241    
242 <  void GB::calcForce(InteractionData idat) {
242 >  void GB::calcForce(InteractionData &idat) {
243  
244      if (!initialized_) initialize();
245      
246 <    pair<AtomType*, AtomType*> key = make_pair(idat.atype1, idat.atype2);
271 <    GBInteractionData mixer = MixingMap[key];
246 >    GBInteractionData mixer = MixingMap[idat.atypes];
247  
248      RealType sigma0 = mixer.sigma0;
249      RealType dw     = mixer.dw;
# Line 280 | Line 255 | namespace OpenMD {
255      RealType xpap2  = mixer.xpap2;
256      RealType xpapi2 = mixer.xpapi2;
257  
258 <    Vector3d ul1 = idat.A1.getRow(2);
259 <    Vector3d ul2 = idat.A2.getRow(2);
258 >    // cerr << "atypes = " << idat.atypes.first->getName() << " " << idat.atypes.second->getName() << "\n";
259 >    // cerr << "sigma0 = " <<mixer.sigma0 <<"\n";
260 >    // cerr << "dw     = " <<mixer.dw <<"\n";
261 >    // cerr << "eps0   = " <<mixer.eps0 <<"\n";  
262 >    // cerr << "x2     = " <<mixer.x2 <<"\n";    
263 >    // cerr << "xa2    = " <<mixer.xa2 <<"\n";  
264 >    // cerr << "xai2   = " <<mixer.xai2 <<"\n";  
265 >    // cerr << "xp2    = " <<mixer.xp2 <<"\n";  
266 >    // cerr << "xpap2  = " <<mixer.xpap2 <<"\n";
267 >    // cerr << "xpapi2 = " <<mixer.xpapi2 <<"\n";
268  
269 <    RealType a, b, g;
269 >    Vector3d ul1 = idat.A1->getRow(2);
270 >    Vector3d ul2 = idat.A2->getRow(2);
271  
272 <    bool i_is_LJ = idat.atype1->isLennardJones();
273 <    bool j_is_LJ = idat.atype2->isLennardJones();
272 >    // cerr << "ul1 = " <<ul1<<"\n";
273 >    // cerr << "ul2 = " <<ul2<<"\n";
274  
275 +    RealType a, b, g;
276 +
277 +    // This is not good.  We should store this in the mixing map, and not
278 +    // query atom types in calc force.
279 +    bool i_is_LJ = idat.atypes.first->isLennardJones();
280 +    bool j_is_LJ = idat.atypes.second->isLennardJones();
281 +    
282      if (i_is_LJ) {
283        a = 0.0;
284        ul1 = V3Zero;
285      } else {
286 <      a = dot(idat.d, ul1);
286 >      a = dot(*(idat.d), ul1);
287      }
288  
289      if (j_is_LJ) {
290        b = 0.0;
291        ul2 = V3Zero;
292      } else {
293 <      b = dot(idat.d, ul2);
293 >      b = dot(*(idat.d), ul2);
294      }
295  
296      if (i_is_LJ || j_is_LJ)
# Line 307 | Line 298 | namespace OpenMD {
298      else
299        g = dot(ul1, ul2);
300  
301 <    RealType au = a / idat.rij;
302 <    RealType bu = b / idat.rij;
301 >    RealType au = a / *(idat.rij);
302 >    RealType bu = b / *(idat.rij);
303      
304      RealType au2 = au * au;
305      RealType bu2 = bu * bu;
306      RealType g2 = g * g;
307 <    
307 >
308      RealType H  = (xa2 * au2 + xai2 * bu2 - 2.0*x2*au*bu*g)  / (1.0 - x2*g2);
309      RealType Hp = (xpap2*au2 + xpapi2*bu2 - 2.0*xp2*au*bu*g) / (1.0 - xp2*g2);
310  
311 +    // cerr << "au2 = " << au2 << "\n";
312 +    // cerr << "bu2 = " << bu2 << "\n";
313 +    // cerr << "g2 = " << g2 << "\n";
314 +    // cerr << "H = " << H << "\n";
315 +    // cerr << "Hp = " << Hp << "\n";
316 +
317      RealType sigma = sigma0 / sqrt(1.0 - H);
318      RealType e1 = 1.0 / sqrt(1.0 - x2*g2);
319      RealType e2 = 1.0 - Hp;
320      RealType eps = eps0 * pow(e1,nu_) * pow(e2,mu_);
321 <    RealType BigR = dw*sigma0 / (idat.rij - sigma + dw*sigma0);
321 >    RealType BigR = dw*sigma0 / (*(idat.rij) - sigma + dw*sigma0);
322      
323      RealType R3 = BigR*BigR*BigR;
324      RealType R6 = R3*R3;
# Line 329 | Line 326 | namespace OpenMD {
326      RealType R12 = R6*R6;
327      RealType R13 = R6*R7;
328  
329 <    RealType U = idat.vdwMult * 4.0 * eps * (R12 - R6);
329 >    RealType U = *(idat.vdwMult) * 4.0 * eps * (R12 - R6);
330  
331      RealType s3 = sigma*sigma*sigma;
332      RealType s03 = sigma0*sigma0*sigma0;
333  
334 <    RealType pref1 = - idat.vdwMult * 8.0 * eps * mu_ * (R12 - R6) / (e2 * idat.rij);
334 >    // cerr << "vdwMult = " << *(idat.vdwMult) << "\n";
335 >    // cerr << "eps = " << eps <<"\n";
336 >    // cerr << "mu = " << mu_ << "\n";
337 >    // cerr << "R12 = " << R12 << "\n";
338 >    // cerr << "R6 = " << R6 << "\n";
339 >    // cerr << "R13 = " << R13 << "\n";
340 >    // cerr << "R7 = " << R7 << "\n";
341 >    // cerr << "e2 = " << e2 << "\n";
342 >    // cerr << "rij = " << *(idat.rij) << "\n";
343 >    // cerr << "s3 = " << s3 << "\n";
344 >    // cerr << "s03 = " << s03 << "\n";
345 >    // cerr << "dw = " << dw << "\n";
346  
347 <    RealType pref2 = idat.vdwMult * 8.0 * eps * s3 * (6.0*R13 - 3.0*R7) /(dw*idat.rij*s03);
347 >    RealType pref1 = - *(idat.vdwMult) * 8.0 * eps * mu_ * (R12 - R6) /
348 >      (e2 * *(idat.rij));
349  
350 <    RealType dUdr = - (pref1 * Hp + pref2 * (sigma0*sigma0*idat.rij/s3 + H));
350 >    RealType pref2 = *(idat.vdwMult) * 8.0 * eps * s3 * (6.0*R13 - 3.0*R7) /
351 >      (dw*  *(idat.rij) * s03);
352 >
353 >    RealType dUdr = - (pref1 * Hp + pref2 * (sigma0 * sigma0 *  
354 >                                             *(idat.rij) / s3 + H));
355      
356      RealType dUda = pref1 * (xpap2*au - xp2*bu*g) / (1.0 - xp2 * g2)
357        + pref2 * (xa2 * au - x2 *bu*g) / (1.0 - x2 * g2);
# Line 350 | Line 363 | namespace OpenMD {
363        + 8.0 * eps * mu_ * (R12 - R6) * (xp2*au*bu - Hp*xp2*g) /
364        (1.0 - xp2 * g2) / e2 + 8.0 * eps * s3 * (3.0 * R7 - 6.0 * R13) *
365        (x2 * au * bu - H * x2 * g) / (1.0 - x2 * g2) / (dw * s03);
353    
366  
367 <    Vector3d rhat = idat.d / idat.rij;  
368 <    Vector3d rxu1 = cross(idat.d, ul1);
369 <    Vector3d rxu2 = cross(idat.d, ul2);
367 >    // cerr << "pref = " << pref1 << " " << pref2 << "\n";
368 >    // cerr << "dU = " << dUdr << " " << dUda <<" " << dUdb << " " << dUdg << "\n";
369 >
370 >    Vector3d rhat = *(idat.d) / *(idat.rij);  
371 >    Vector3d rxu1 = cross(*(idat.d), ul1);
372 >    Vector3d rxu2 = cross(*(idat.d), ul2);
373      Vector3d uxu = cross(ul1, ul2);
359    
360    idat.pot += U*idat.sw;
361    idat.f1 += dUdr * rhat + dUda * ul1 + dUdb * ul2;    
362    idat.t1 += dUda * rxu1 - dUdg * uxu;
363    idat.t2 += dUdb * rxu2 - dUdg * uxu;
364    idat.vpair += U*idat.sw;
374  
375 +    (*(idat.pot))[VANDERWAALS_FAMILY] += U *  *(idat.sw);
376 +    *(idat.f1) += (dUdr * rhat + dUda * ul1 + dUdb * ul2) * *(idat.sw);
377 +    *(idat.t1) += (dUda * rxu1 - dUdg * uxu) * *(idat.sw);
378 +    *(idat.t2) += (dUdb * rxu2 + dUdg * uxu) * *(idat.sw);
379 +    *(idat.vpair) += U;
380 +
381 +    // cerr << "f1 term = " <<  (dUdr * rhat + dUda * ul1 + dUdb * ul2) * *(idat.sw) << "\n";
382 +    // cerr << "t1 term = " << (dUda * rxu1 - dUdg * uxu) * *(idat.sw) << "\n";
383 +    // cerr << "t2 term = " << (dUdb * rxu2 + dUdg * uxu) * *(idat.sw) << "\n";
384 +    // cerr << "vp term = " << U << "\n";
385 +
386      return;
387  
388    }
389  
390 <  RealType GB::getSuggestedCutoffRadius(AtomType* at1, AtomType* at2) {
390 >  RealType GB::getSuggestedCutoffRadius(pair<AtomType*, AtomType*> atypes) {
391      if (!initialized_) initialize();  
392  
393      RealType cut = 0.0;
394  
395 <    if (at1->isGayBerne()) {
396 <      GayBerneParam gb1 = getGayBerneParam(at1);
397 <      RealType d1 = gb1.GB_d;
398 <      RealType l1 = gb1.GB_l;
395 >    LennardJonesAdapter lja1 = LennardJonesAdapter(atypes.first);
396 >    GayBerneAdapter gba1 = GayBerneAdapter(atypes.first);
397 >    LennardJonesAdapter lja2 = LennardJonesAdapter(atypes.second);
398 >    GayBerneAdapter gba2 = GayBerneAdapter(atypes.second);
399 >
400 >    if (gba1.isGayBerne()) {
401 >      RealType d1 = gba1.getD();
402 >      RealType l1 = gba1.getL();
403        // sigma is actually sqrt(2)*l  for prolate ellipsoids
404 <      cut = max(cut, 2.5 * sqrt(2.0) * max(d1, l1));
405 <    } else if (at1->isLennardJones()) {
406 <      cut = max(cut, 2.5 * getLJSigma(at1));
404 >      cut = max(cut, RealType(2.5) * sqrt(RealType(2.0)) * max(d1, l1));
405 >    } else if (lja1.isLennardJones()) {
406 >      cut = max(cut, RealType(2.5) * lja1.getSigma());
407      }
408  
409 <    if (at2->isGayBerne()) {
410 <      GayBerneParam gb2 = getGayBerneParam(at2);
411 <      RealType d2 = gb2.GB_d;
412 <      RealType l2 = gb2.GB_l;
413 <      cut = max(cut, 2.5 * sqrt(2.0) * max(d2, l2));
414 <    } else if (at1->isLennardJones()) {
391 <      cut = max(cut, 2.5 * getLJSigma(at2));
409 >    if (gba2.isGayBerne()) {
410 >      RealType d2 = gba2.getD();
411 >      RealType l2 = gba2.getL();
412 >      cut = max(cut, RealType(2.5) * sqrt(RealType(2.0)) * max(d2, l2));
413 >    } else if (lja2.isLennardJones()) {
414 >      cut = max(cut, RealType(2.5) * lja2.getSigma());
415      }
416    
417      return cut;

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines