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 1688 by gezelter, Wed Mar 14 17:56:01 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 49 | Line 50 | namespace OpenMD {
50   using namespace std;
51   namespace OpenMD {
52  
53 +  /* GB is the Gay-Berne interaction for ellipsoidal particles.  The original
54 +   * paper (for identical uniaxial particles) is:
55 +   *    J. G. Gay and B. J. Berne, J. Chem. Phys., 74, 3316-3319 (1981).
56 +   * A more-general GB potential for dissimilar uniaxial particles:
57 +   *    D. J. Cleaver, C. M. Care, M. P. Allen and M. P. Neal, Phys. Rev. E,
58 +   *    54, 559-567 (1996).
59 +   * Further parameterizations can be found in:
60 +   *    A. P. J. Emerson, G. R. Luckhurst and S. G. Whatling, Mol. Phys.,
61 +   *    82, 113-124 (1994).
62 +   * And a nice force expression:
63 +   *    G. R. Luckhurst and R. A. Stephens, Liq. Cryst. 8, 451-464 (1990).
64 +   * Even clearer force and torque expressions:
65 +   *    P. A. Golubkov and P. Y. Ren, J. Chem. Phys., 125, 64103 (2006).
66 +   * New expressions for cross interactions of strength parameters:
67 +   *    J. Wu, X. Zhen, H. Shen, G. Li, and P. Ren, J. Chem. Phys.,
68 +   *    135, 155104 (2011).
69 +   *
70 +   * In this version of the GB interaction, each uniaxial ellipsoidal type
71 +   * is described using a set of 6 parameters:
72 +   *  d:  range parameter for side-by-side (S) and cross (X) configurations
73 +   *  l:  range parameter for end-to-end (E) configuration
74 +   *  epsilon_X:  well-depth parameter for cross (X) configuration
75 +   *  epsilon_S:  well-depth parameter for side-by-side (S) configuration
76 +   *  epsilon_E:  well depth parameter for end-to-end (E) configuration
77 +   *  dw: "softness" of the potential
78 +   *
79 +   * Additionally, there are two "universal" paramters to govern the overall
80 +   * importance of the purely orientational (nu) and the mixed
81 +   * orientational / translational (mu) parts of strength of the interactions.
82 +   * These parameters have default or "canonical" values, but may be changed
83 +   * as a force field option:
84 +   * nu_: purely orientational part : defaults to 1
85 +   * mu_: mixed orientational / translational part : defaults to 2
86 +   */
87 +
88 +
89    GB::GB() : name_("GB"), initialized_(false), mu_(2.0), nu_(1.0), forceField_(NULL) {}
90    
91    GayBerneParam GB::getGayBerneParam(AtomType* atomType) {
# Line 172 | Line 209 | namespace OpenMD {
209        simError();        
210      }
211      
212 <    RealType d1, l1, e1, er1, dw1;
212 >    RealType d1, l1, eX1, eS1, eE1, dw1;
213      
214      if (atomType->isGayBerne()) {
215        GayBerneParam gb1 = getGayBerneParam(atomType);
216        d1 = gb1.GB_d;
217        l1 = gb1.GB_l;
218 <      e1 = gb1.GB_eps;
219 <      er1 = gb1.GB_eps_ratio;
218 >      eX1 = gb1.GB_eps_X;
219 >      eS1 = gb1.GB_eps_S;
220 >      eE1 = gb1.GB_eps_E;
221        dw1 = gb1.GB_dw;
222      } else if (atomType->isLennardJones()) {
223        d1 = getLJSigma(atomType) / sqrt(2.0);
186      e1 = getLJEpsilon(atomType);
224        l1 = d1;
225 <      er1 = 1.0;
225 >      eX1 = getLJEpsilon(atomType);
226 >      eS1 = eX1;
227 >      eE1 = eX1;
228        dw1 = 1.0;      
229      } else {
230        sprintf( painCave.errMsg,
# Line 205 | Line 244 | namespace OpenMD {
244        
245        AtomType* atype2 = (*it).second;
246        
247 <      RealType d2, l2, e2, er2, dw2;
247 >      RealType d2, l2, eX2, eS2, eE2, dw2;
248        
249        if (atype2->isGayBerne()) {
250          GayBerneParam gb2 = getGayBerneParam(atype2);
251          d2 = gb2.GB_d;
252          l2 = gb2.GB_l;
253 <        e2 = gb2.GB_eps;
254 <        er2 = gb2.GB_eps_ratio;
253 >        eX2 = gb2.GB_eps_X;
254 >        eS2 = gb2.GB_eps_S;
255 >        eE2 = gb2.GB_eps_E;
256          dw2 = gb2.GB_dw;
257        } else if (atype2->isLennardJones()) {
258          d2 = getLJSigma(atype2) / sqrt(2.0);
219        e2 = getLJEpsilon(atype2);
259          l2 = d2;
260 <        er2 = 1.0;
260 >        eX2 = getLJEpsilon(atype2);
261 >        eS2 = eX2;
262 >        eE2 = eX2;
263          dw2 = 1.0;
264        }
265                        
266 <      GBInteractionData mixer;        
266 >      GBInteractionData mixer1, mixer2;    
267        
268        //  Cleaver paper uses sqrt of squares to get sigma0 for
269        //  mixed interactions.
270              
271 <      mixer.sigma0 = sqrt(d1*d1 + d2*d2);
272 <      mixer.xa2 = (l1*l1 - d1*d1)/(l1*l1 + d2*d2);
273 <      mixer.xai2 = (l2*l2 - d2*d2)/(l2*l2 + d1*d1);
274 <      mixer.x2 = (l1*l1 - d1*d1) * (l2*l2 - d2*d2) /
271 >      mixer1.sigma0 = sqrt(d1*d1 + d2*d2);
272 >      mixer1.xa2 = (l1*l1 - d1*d1)/(l1*l1 + d2*d2);
273 >      mixer1.xai2 = (l2*l2 - d2*d2)/(l2*l2 + d1*d1);
274 >      mixer1.x2 = (l1*l1 - d1*d1) * (l2*l2 - d2*d2) /
275          ((l2*l2 + d1*d1) * (l1*l1 + d2*d2));
276 +
277 +      mixer2.sigma0 = mixer1.sigma0;
278 +      // xa2 and xai2 for j-i pairs are reversed from the same i-j pairing.
279 +      // Swapping the particles reverses the anisotropy parameters:
280 +      mixer2.xa2 = mixer1.xai2;
281 +      mixer2.xai2 = mixer1.xa2;
282 +      mixer2.x2 = mixer1.x2;
283  
284        // assumed LB mixing rules for now:
285  
286 <      mixer.dw = 0.5 * (dw1 + dw2);
287 <      mixer.eps0 = sqrt(e1 * e2);
286 >      mixer1.dw = 0.5 * (dw1 + dw2);
287 >      mixer1.eps0 = sqrt(eX1 * eX2);
288 >
289 >      mixer2.dw = mixer1.dw;
290 >      mixer2.eps0 = mixer1.eps0;
291 >
292 >      RealType mi = RealType(1.0)/mu_;
293        
294 <      RealType er = sqrt(er1 * er2);
295 <      RealType ermu = pow(er,(1.0 / mu_));
296 <      RealType xp = (1.0 - ermu) / (1.0 + ermu);
297 <      RealType ap2 = 1.0 / (1.0 + ermu);
245 <      
246 <      mixer.xp2 = xp * xp;
247 <      mixer.xpap2 = xp * ap2;
248 <      mixer.xpapi2 = xp / ap2;
294 >      mixer1.xpap2  = (pow(eS1, mi) - pow(eE1, mi)) / (pow(eS1, mi) + pow(eE2, mi));
295 >      mixer1.xpapi2 = (pow(eS2, mi) - pow(eE2, mi)) / (pow(eS2, mi) + pow(eE1, mi));
296 >      mixer1.xp2    = (pow(eS1, mi) - pow(eE1, mi)) * (pow(eS2, mi) - pow(eE2, mi))  /
297 >        (pow(eS2, mi) + pow(eE1, mi)) / (pow(eS1, mi) + pow(eE2, mi)) ;
298  
299 +      // xpap2 and xpapi2 for j-i pairs are reversed from the same i-j pairing.
300 +      // Swapping the particles reverses the anisotropy parameters:
301 +      mixer2.xpap2 = mixer1.xpapi2;
302 +      mixer2.xpapi2 = mixer1.xpap2;
303 +      mixer2.xp2 = mixer1.xp2;
304 +
305        // only add this pairing if at least one of the atoms is a Gay-Berne atom
306  
307        if (atomType->isGayBerne() || atype2->isGayBerne()) {
# Line 255 | Line 310 | namespace OpenMD {
310          key1 = make_pair(atomType, atype2);
311          key2 = make_pair(atype2, atomType);
312          
313 <        MixingMap[key1] = mixer;
313 >        MixingMap[key1] = mixer1;
314          if (key2 != key1) {
315 <          MixingMap[key2] = mixer;
315 >          MixingMap[key2] = mixer2;
316          }
317        }
318      }      
319    }
320    
321 <  void GB::calcForce(InteractionData idat) {
321 >  void GB::calcForce(InteractionData &idat) {
322  
323      if (!initialized_) initialize();
324      
325 <    pair<AtomType*, AtomType*> key = make_pair(idat.atype1, idat.atype2);
271 <    GBInteractionData mixer = MixingMap[key];
325 >    GBInteractionData mixer = MixingMap[idat.atypes];
326  
327      RealType sigma0 = mixer.sigma0;
328      RealType dw     = mixer.dw;
# Line 280 | Line 334 | namespace OpenMD {
334      RealType xpap2  = mixer.xpap2;
335      RealType xpapi2 = mixer.xpapi2;
336  
337 <    Vector3d ul1 = idat.A1.getRow(2);
338 <    Vector3d ul2 = idat.A2.getRow(2);
337 >    // cerr << "atypes = " << idat.atypes.first->getName() << " " << idat.atypes.second->getName() << "\n";
338 >    // cerr << "sigma0 = " <<mixer.sigma0 <<"\n";
339 >    // cerr << "dw     = " <<mixer.dw <<"\n";
340 >    // cerr << "eps0   = " <<mixer.eps0 <<"\n";  
341 >    // cerr << "x2     = " <<mixer.x2 <<"\n";    
342 >    // cerr << "xa2    = " <<mixer.xa2 <<"\n";  
343 >    // cerr << "xai2   = " <<mixer.xai2 <<"\n";  
344 >    // cerr << "xp2    = " <<mixer.xp2 <<"\n";  
345 >    // cerr << "xpap2  = " <<mixer.xpap2 <<"\n";
346 >    // cerr << "xpapi2 = " <<mixer.xpapi2 <<"\n";
347  
348 <    RealType a, b, g;
348 >    Vector3d ul1 = idat.A1->getRow(2);
349 >    Vector3d ul2 = idat.A2->getRow(2);
350  
351 <    bool i_is_LJ = idat.atype1->isLennardJones();
352 <    bool j_is_LJ = idat.atype2->isLennardJones();
351 >    // cerr << "ul1 = " <<ul1<<"\n";
352 >    // cerr << "ul2 = " <<ul2<<"\n";
353  
354 +    RealType a, b, g;
355 +
356 +    bool i_is_LJ = idat.atypes.first->isLennardJones();
357 +    bool j_is_LJ = idat.atypes.second->isLennardJones();
358 +    
359      if (i_is_LJ) {
360        a = 0.0;
361        ul1 = V3Zero;
362      } else {
363 <      a = dot(idat.d, ul1);
363 >      a = dot(*(idat.d), ul1);
364      }
365  
366      if (j_is_LJ) {
367        b = 0.0;
368        ul2 = V3Zero;
369      } else {
370 <      b = dot(idat.d, ul2);
370 >      b = dot(*(idat.d), ul2);
371      }
372  
373      if (i_is_LJ || j_is_LJ)
# Line 307 | Line 375 | namespace OpenMD {
375      else
376        g = dot(ul1, ul2);
377  
378 <    RealType au = a / idat.rij;
379 <    RealType bu = b / idat.rij;
378 >    RealType au = a / *(idat.rij);
379 >    RealType bu = b / *(idat.rij);
380      
381      RealType au2 = au * au;
382      RealType bu2 = bu * bu;
383      RealType g2 = g * g;
384 <    
384 >
385      RealType H  = (xa2 * au2 + xai2 * bu2 - 2.0*x2*au*bu*g)  / (1.0 - x2*g2);
386      RealType Hp = (xpap2*au2 + xpapi2*bu2 - 2.0*xp2*au*bu*g) / (1.0 - xp2*g2);
387  
388 +    // cerr << "au2 = " << au2 << "\n";
389 +    // cerr << "bu2 = " << bu2 << "\n";
390 +    // cerr << "g2 = " << g2 << "\n";
391 +    // cerr << "H = " << H << "\n";
392 +    // cerr << "Hp = " << Hp << "\n";
393 +
394      RealType sigma = sigma0 / sqrt(1.0 - H);
395      RealType e1 = 1.0 / sqrt(1.0 - x2*g2);
396      RealType e2 = 1.0 - Hp;
397      RealType eps = eps0 * pow(e1,nu_) * pow(e2,mu_);
398 <    RealType BigR = dw*sigma0 / (idat.rij - sigma + dw*sigma0);
398 >    RealType BigR = dw*sigma0 / (*(idat.rij) - sigma + dw*sigma0);
399      
400      RealType R3 = BigR*BigR*BigR;
401      RealType R6 = R3*R3;
# Line 329 | Line 403 | namespace OpenMD {
403      RealType R12 = R6*R6;
404      RealType R13 = R6*R7;
405  
406 <    RealType U = idat.vdwMult * 4.0 * eps * (R12 - R6);
406 >    RealType U = *(idat.vdwMult) * 4.0 * eps * (R12 - R6);
407  
408      RealType s3 = sigma*sigma*sigma;
409      RealType s03 = sigma0*sigma0*sigma0;
410  
411 <    RealType pref1 = - idat.vdwMult * 8.0 * eps * mu_ * (R12 - R6) / (e2 * idat.rij);
411 >    // cerr << "vdwMult = " << *(idat.vdwMult) << "\n";
412 >    // cerr << "eps = " << eps <<"\n";
413 >    // cerr << "mu = " << mu_ << "\n";
414 >    // cerr << "R12 = " << R12 << "\n";
415 >    // cerr << "R6 = " << R6 << "\n";
416 >    // cerr << "R13 = " << R13 << "\n";
417 >    // cerr << "R7 = " << R7 << "\n";
418 >    // cerr << "e2 = " << e2 << "\n";
419 >    // cerr << "rij = " << *(idat.rij) << "\n";
420 >    // cerr << "s3 = " << s3 << "\n";
421 >    // cerr << "s03 = " << s03 << "\n";
422 >    // cerr << "dw = " << dw << "\n";
423  
424 <    RealType pref2 = idat.vdwMult * 8.0 * eps * s3 * (6.0*R13 - 3.0*R7) /(dw*idat.rij*s03);
424 >    RealType pref1 = - *(idat.vdwMult) * 8.0 * eps * mu_ * (R12 - R6) /
425 >      (e2 * *(idat.rij));
426  
427 <    RealType dUdr = - (pref1 * Hp + pref2 * (sigma0*sigma0*idat.rij/s3 + H));
427 >    RealType pref2 = *(idat.vdwMult) * 8.0 * eps * s3 * (6.0*R13 - 3.0*R7) /
428 >      (dw*  *(idat.rij) * s03);
429 >
430 >    RealType dUdr = - (pref1 * Hp + pref2 * (sigma0 * sigma0 *  
431 >                                             *(idat.rij) / s3 + H));
432      
433      RealType dUda = pref1 * (xpap2*au - xp2*bu*g) / (1.0 - xp2 * g2)
434        + pref2 * (xa2 * au - x2 *bu*g) / (1.0 - x2 * g2);
# Line 350 | Line 440 | namespace OpenMD {
440        + 8.0 * eps * mu_ * (R12 - R6) * (xp2*au*bu - Hp*xp2*g) /
441        (1.0 - xp2 * g2) / e2 + 8.0 * eps * s3 * (3.0 * R7 - 6.0 * R13) *
442        (x2 * au * bu - H * x2 * g) / (1.0 - x2 * g2) / (dw * s03);
353    
443  
444 <    Vector3d rhat = idat.d / idat.rij;  
445 <    Vector3d rxu1 = cross(idat.d, ul1);
446 <    Vector3d rxu2 = cross(idat.d, ul2);
444 >    // cerr << "pref = " << pref1 << " " << pref2 << "\n";
445 >    // cerr << "dU = " << dUdr << " " << dUda <<" " << dUdb << " " << dUdg << "\n";
446 >
447 >    Vector3d rhat = *(idat.d) / *(idat.rij);  
448 >    Vector3d rxu1 = cross(*(idat.d), ul1);
449 >    Vector3d rxu2 = cross(*(idat.d), ul2);
450      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;
451  
452 +    (*(idat.pot))[VANDERWAALS_FAMILY] += U *  *(idat.sw);
453 +    *(idat.f1) += (dUdr * rhat + dUda * ul1 + dUdb * ul2) * *(idat.sw);
454 +    *(idat.t1) += (dUda * rxu1 - dUdg * uxu) * *(idat.sw);
455 +    *(idat.t2) += (dUdb * rxu2 + dUdg * uxu) * *(idat.sw);
456 +    *(idat.vpair) += U;
457 +
458 +    // cerr << "f1 term = " <<  (dUdr * rhat + dUda * ul1 + dUdb * ul2) * *(idat.sw) << "\n";
459 +    // cerr << "t1 term = " << (dUda * rxu1 - dUdg * uxu) * *(idat.sw) << "\n";
460 +    // cerr << "t2 term = " << (dUdb * rxu2 + dUdg * uxu) * *(idat.sw) << "\n";
461 +    // cerr << "vp term = " << U << "\n";
462 +
463      return;
464  
465    }
466  
467 <  RealType GB::getSuggestedCutoffRadius(AtomType* at1, AtomType* at2) {
467 >  RealType GB::getSuggestedCutoffRadius(pair<AtomType*, AtomType*> atypes) {
468      if (!initialized_) initialize();  
469  
470      RealType cut = 0.0;
471  
472 <    if (at1->isGayBerne()) {
473 <      GayBerneParam gb1 = getGayBerneParam(at1);
472 >    if (atypes.first->isGayBerne()) {
473 >      GayBerneParam gb1 = getGayBerneParam(atypes.first);
474        RealType d1 = gb1.GB_d;
475        RealType l1 = gb1.GB_l;
476        // sigma is actually sqrt(2)*l  for prolate ellipsoids
477 <      cut = max(cut, 2.5 * sqrt(2.0) * max(d1, l1));
478 <    } else if (at1->isLennardJones()) {
479 <      cut = max(cut, 2.5 * getLJSigma(at1));
477 >      cut = max(cut, RealType(2.5) * sqrt(RealType(2.0)) * max(d1, l1));
478 >    } else if (atypes.first->isLennardJones()) {
479 >      cut = max(cut, RealType(2.5) * getLJSigma(atypes.first));
480      }
481  
482 <    if (at2->isGayBerne()) {
483 <      GayBerneParam gb2 = getGayBerneParam(at2);
482 >    if (atypes.second->isGayBerne()) {
483 >      GayBerneParam gb2 = getGayBerneParam(atypes.second);
484        RealType d2 = gb2.GB_d;
485        RealType l2 = gb2.GB_l;
486 <      cut = max(cut, 2.5 * sqrt(2.0) * max(d2, l2));
487 <    } else if (at1->isLennardJones()) {
488 <      cut = max(cut, 2.5 * getLJSigma(at2));
486 >      cut = max(cut, RealType(2.5) * sqrt(RealType(2.0)) * max(d2, l2));
487 >    } else if (atypes.second->isLennardJones()) {
488 >      cut = max(cut, RealType(2.5) * getLJSigma(atypes.second));
489      }
490    
491      return cut;

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines