NAMD
Functions
ComputeNonbondedAlch.h File Reference

Go to the source code of this file.

Functions

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)
 

Function Documentation

void vdw_energy ( const BigReal  A,
const BigReal  B,
const BigReal  r2,
BigReal U 
)
inline

Definition at line 25 of file ComputeNonbondedAlch.h.

References B.

Referenced by fep_vdw_forceandenergies().

26  {
27  const BigReal inv_r2 = 1. / r2;
28  const BigReal inv_r6 = inv_r2*inv_r2*inv_r2;
29  *U = inv_r6*(A*inv_r6 - B);
30 }
const BigReal A
const BigReal B
double BigReal
Definition: common.h:114
void vdw_forceandenergy ( const BigReal  A,
const BigReal  B,
const BigReal  r2,
BigReal U,
BigReal F 
)
inline

Copyright (c) 1995, 1996, 1997, 1998, 1999, 2000 by The Board of Trustees of the University of Illinois. All rights reserved.

Definition at line 16 of file ComputeNonbondedAlch.h.

References B.

Referenced by fep_vdw_forceandenergies(), and ti_vdw_force_energy_dUdl().

17  {
18  const BigReal inv_r2 = 1. / r2;
19  const BigReal inv_r6 = inv_r2*inv_r2*inv_r2;
20  *U = inv_r6*(A*inv_r6 - B);
21  *F = 6*inv_r2*(2*(*U) + B*inv_r6);
22 }
const BigReal A
const BigReal B
double BigReal
Definition: common.h:114
void vdw_fswitch_energy ( const BigReal  A,
const BigReal  B,
const BigReal  r2,
const BigReal  switchdist2,
const BigReal  cutoff2,
BigReal U 
)
inline

Definition at line 75 of file ComputeNonbondedAlch.h.

References cutoff2.

Referenced by fep_vdw_forceandenergies().

77  {
78  // TODO: Some rearrangement could be done to use rsqrt() instead?
79  const BigReal inv_r2 = 1. / r2;
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;
90 }
const BigReal A
__global__ void const int const TileList *__restrict__ TileExcl *__restrict__ const int *__restrict__ const int const float2 *__restrict__ cudaTextureObject_t const int *__restrict__ const float3 const float3 const float3 const float4 *__restrict__ const float cutoff2
const BigReal B
double BigReal
Definition: common.h:114
void vdw_fswitch_forceandenergy ( const BigReal  A,
const BigReal  B,
const BigReal  r2,
const BigReal  switchdist2,
const BigReal  cutoff2,
BigReal U,
BigReal F 
)
inline

Definition at line 56 of file ComputeNonbondedAlch.h.

References cutoff2.

Referenced by fep_vdw_forceandenergies(), and ti_vdw_force_energy_dUdl().

58  {
59  // TODO: Some rearrangement could be done to use rsqrt() instead?
60  const BigReal inv_r2 = 1. / r2;
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);
72 }
const BigReal A
__global__ void const int const TileList *__restrict__ TileExcl *__restrict__ const int *__restrict__ const int const float2 *__restrict__ cudaTextureObject_t const int *__restrict__ const float3 const float3 const float3 const float4 *__restrict__ const float cutoff2
const BigReal B
double BigReal
Definition: common.h:114
void vdw_fswitch_shift ( const BigReal  A,
const BigReal  B,
const BigReal  switchdist2,
const BigReal  cutoff2,
BigReal dU 
)
inline

Definition at line 97 of file ComputeNonbondedAlch.h.

References cutoff2.

Referenced by fep_vdw_forceandenergies(), and ti_vdw_force_energy_dUdl().

98  {
99  // TODO: Some rearrangement could be done to use rsqrt() instead.
100  const BigReal cutoff6 = cutoff2*cutoff2*cutoff2;
101  const BigReal switchdist6 = switchdist2*switchdist2*switchdist2;
102  const BigReal v_vdwa = -A / (cutoff6*switchdist6);
103  const BigReal v_vdwb = -B / sqrt(cutoff6*switchdist6);
104  *dU = v_vdwa - v_vdwb; //deltaV2 from Steinbach & Brooks
105 }
const BigReal A
__global__ void const int const TileList *__restrict__ TileExcl *__restrict__ const int *__restrict__ const int const float2 *__restrict__ cudaTextureObject_t const int *__restrict__ const float3 const float3 const float3 const float4 *__restrict__ const float cutoff2
const BigReal B
double BigReal
Definition: common.h:114
void vdw_switch ( const BigReal  r2,
const BigReal  switchdist2,
const BigReal  cutoff2,
const BigReal  switchfactor,
BigReal switchmul,
BigReal switchmul2 
)
inline

Definition at line 35 of file ComputeNonbondedAlch.h.

Referenced by fep_vdw_forceandenergies(), and ti_vdw_force_energy_dUdl().

37  {
38  if (r2 <= switchdist2) {
39  // "switching off" will always fall through here.
40  *switchmul = 1.0;
41  *switchmul2 = 0.0;
42  } else {
43  const BigReal cdiff2 = (cutoff2 - r2);
44  const BigReal sdiff2 = (r2 - switchdist2);
45  *switchmul = switchfactor*cdiff2*cdiff2*(cdiff2 + 3*sdiff2);
46  *switchmul2 = 12*switchfactor*cdiff2*sdiff2;
47  }
48 }
__global__ void const int const TileList *__restrict__ TileExcl *__restrict__ const int *__restrict__ const int const float2 *__restrict__ cudaTextureObject_t const int *__restrict__ const float3 const float3 const float3 const float4 *__restrict__ const float cutoff2
double BigReal
Definition: common.h:114