Go to the source code of this file.
 | 
| void  | vdw_forceandenergy (const BigReal A, const BigReal B, const BigReal r2, BigReal *U, BigReal *F) | 
|   | 
| void  | vdw_energy (const BigReal A, const BigReal B, const BigReal r2, BigReal *U) | 
|   | 
| void  | vdw_switch (const BigReal r2, const BigReal switchdist2, const BigReal cutoff2, const BigReal switchfactor, BigReal *switchmul, BigReal *switchmul2) | 
|   | 
| void  | vdw_fswitch_forceandenergy (const BigReal A, const BigReal B, const BigReal r2, const BigReal switchdist2, const BigReal cutoff2, BigReal *U, BigReal *F) | 
|   | 
| void  | vdw_fswitch_energy (const BigReal A, const BigReal B, const BigReal r2, const BigReal switchdist2, const BigReal cutoff2, BigReal *U) | 
|   | 
| void  | vdw_fswitch_shift (const BigReal A, const BigReal B, const BigReal switchdist2, const BigReal cutoff2, BigReal *dU) | 
|   | 
◆ vdw_energy()
◆ vdw_forceandenergy()
◆ vdw_fswitch_energy()
Definition at line 75 of file ComputeNonbondedAlch.h.
Referenced by fep_vdw_forceandenergies().
   80   const BigReal inv_r6 = inv_r2*inv_r2*inv_r2;
    81   const BigReal inv_r3 = sqrt(inv_r6);
    82   const BigReal inv_cutoff6 = 1. / (cutoff2*cutoff2*cutoff2);
    83   const BigReal inv_cutoff3 = sqrt(inv_cutoff6);
    84   const BigReal switchdist6 = switchdist2*switchdist2*switchdist2;
    85   const BigReal k_vdwa = A / (1 - switchdist6*inv_cutoff6);
    86   const BigReal k_vdwb = B / (1 - sqrt(switchdist6*inv_cutoff6));
    87   const BigReal tmpa = inv_r6 - inv_cutoff6;
    88   const BigReal tmpb = inv_r3 - inv_cutoff3;
    89   *U = k_vdwa*tmpa*tmpa - k_vdwb*tmpb*tmpb;
 
 
 
 
◆ vdw_fswitch_forceandenergy()
Definition at line 56 of file ComputeNonbondedAlch.h.
Referenced by fep_vdw_forceandenergies(), and ti_vdw_force_energy_dUdl().
   61   const BigReal inv_r6 = inv_r2*inv_r2*inv_r2;
    62   const BigReal inv_r3 = sqrt(inv_r6);
    63   const BigReal inv_cutoff6 = 1. / (cutoff2*cutoff2*cutoff2);
    64   const BigReal inv_cutoff3 = sqrt(inv_cutoff6);
    65   const BigReal switchdist6 = switchdist2*switchdist2*switchdist2;
    66   const BigReal k_vdwa = A / (1 - switchdist6*inv_cutoff6);
    67   const BigReal k_vdwb = B / (1 - sqrt(switchdist6*inv_cutoff6));
    68   const BigReal tmpa = inv_r6 - inv_cutoff6;
    69   const BigReal tmpb = inv_r3 - inv_cutoff3;
    70   *U = k_vdwa*tmpa*tmpa - k_vdwb*tmpb*tmpb;
    71   *F = 6*inv_r2*(2*k_vdwa*tmpa*inv_r6 - k_vdwb*tmpb*inv_r3);
 
 
 
 
◆ vdw_fswitch_shift()
◆ vdw_switch()