NonbondedForce Class Reference

This class implements nonbonded interactions between particles, including a Coulomb force to represent electrostatics and a Lennard-Jones force to represent van der Waals interactions. More...

#include <NonbondedForce.h>

Inheritance diagram for NonbondedForce:
Force

List of all members.

Classes

class  ExceptionInfo
 This is an internal class used to record information about an exception.
class  ParticleInfo
 This is an internal class used to record information about a particle.

Public Types

enum  NonbondedMethod {
  NoCutoff = 0, CutoffNonPeriodic = 1, CutoffPeriodic = 2, Ewald = 3,
  PME = 4
}
 

This is an enumeration of the different methods that may be used for handling long range nonbonded forces.

More...

Public Member Functions

 NonbondedForce ()
 Create a NonbondedForce.
int getNumParticles () const
 Get the number of particles for which force field parameters have been defined.
int getNumExceptions () const
 Get the number of special interactions that should be calculated differently from other interactions.
NonbondedMethod getNonbondedMethod () const
 Get the method used for handling long range nonbonded interactions.
void setNonbondedMethod (NonbondedMethod method)
 Set the method used for handling long range nonbonded interactions.
double getCutoffDistance () const
 Get the cutoff distance (in nm) being used for nonbonded interactions.
void setCutoffDistance (double distance)
 Set the cutoff distance (in nm) being used for nonbonded interactions.
double getReactionFieldDielectric () const
 Get the dielectric constant to use for the solvent in the reaction field approximation.
void setReactionFieldDielectric (double dielectric)
 Set the dielectric constant to use for the solvent in the reaction field approximation.
double getEwaldErrorTolerance () const
 Get the error tolerance for Ewald summation.
void setEwaldErrorTolerance (double tol)
 Get the error tolerance for Ewald summation.
int addParticle (double charge, double sigma, double epsilon)
 Add the nonbonded force parameters for a particle.
void getParticleParameters (int index, double &charge, double &sigma, double &epsilon) const
 Get the nonbonded force parameters for a particle.
void setParticleParameters (int index, double charge, double sigma, double epsilon)
 Set the nonbonded force parameters for a particle.
int addException (int particle1, int particle2, double chargeProd, double sigma, double epsilon, bool replace=false)
 Add an interaction to the list of exceptions that should be calculated differently from other interactions.
void getExceptionParameters (int index, int &particle1, int &particle2, double &chargeProd, double &sigma, double &epsilon) const
 Get the force field parameters for an interaction that should be calculated differently from others.
void setExceptionParameters (int index, int particle1, int particle2, double chargeProd, double sigma, double epsilon)
 Set the force field parameters for an interaction that should be calculated differently from others.
void createExceptionsFromBonds (const std::vector< std::pair< int, int > > &bonds, double coulomb14Scale, double lj14Scale)
 Identify exceptions based on the molecular topology.

Protected Member Functions

ForceImplcreateImpl ()
 When a Context is created, it invokes this method on each Force in the System.

Detailed Description

This class implements nonbonded interactions between particles, including a Coulomb force to represent electrostatics and a Lennard-Jones force to represent van der Waals interactions.

It optionally supports periodic boundary conditions and cutoffs for long range interactions. Lennard-Jones interactions are calculated with the Lorentz-Bertelot combining rule: it uses the arithmetic mean of the sigmas and the geometric mean of the epsilons for the two interacting particles.

To use this class, create a NonbondedForce object, then call addParticle() once for each particle in the System to define its parameters. The number of particles for which you define nonbonded parameters must be exactly equal to the number of particles in the System, or else an exception will be thrown when you try to create a Context. After a particle has been added, you can modify its force field parameters by calling setParticleParameters().

NonbondedForce also lets you specify "exceptions", particular pairs of particles whose interactions should be computed based on different parameters than those defined for the individual particles. This can be used to completely exclude certain interactions from the force calculation, or to alter how they interact with each other.

Many molecular force fields omit Coulomb and Lennard-Jones interactions between particles separated by one or two bonds, while using modified parameters for those separated by three bonds (known as "1-4 interactions"). This class provides a convenience method for this case called createExceptionsFromBonds(). You pass to it a list of bonds and the scale factors to use for 1-4 interactions. It identifies all pairs of particles which are separated by 1, 2, or 3 bonds, then automatically creates exceptions for them.


Member Enumeration Documentation

This is an enumeration of the different methods that may be used for handling long range nonbonded forces.

Enumerator:
NoCutoff 

No cutoff is applied to nonbonded interactions.

The full set of N^2 interactions is computed exactly. This necessarily means that periodic boundary conditions cannot be used. This is the default.

CutoffNonPeriodic 

Interactions beyond the cutoff distance are ignored.

Coulomb interactions closer than the cutoff distance are modified using the reaction field method.

CutoffPeriodic 

Periodic boundary conditions are used, so that each particle interacts only with the nearest periodic copy of each other particle.

Interactions beyond the cutoff distance are ignored. Coulomb interactions closer than the cutoff distance are modified using the reaction field method.

Ewald 

Periodic boundary conditions are used, and Ewald summation is used to compute the interaction of each particle with all periodic copies of every other particle.

PME 

Periodic boundary conditions are used, and Particle-Mesh Ewald (PME) summation is used to compute the interaction of each particle with all periodic copies of every other particle.


Constructor & Destructor Documentation

NonbondedForce (  ) 

Create a NonbondedForce.


Member Function Documentation

int addException ( int  particle1,
int  particle2,
double  chargeProd,
double  sigma,
double  epsilon,
bool  replace = false 
)

Add an interaction to the list of exceptions that should be calculated differently from other interactions.

If chargeProd and epsilon are both equal to 0, this will cause the interaction to be completely omitted from force and energy calculations.

In many cases, you can use createExceptionsFromBonds() rather than adding each exception explicitly.

Parameters:
particle1 the index of the first particle involved in the interaction
particle2 the index of the second particle involved in the interaction
chargeProd the scaled product of the atomic charges (i.e. the strength of the Coulomb interaction), measured in units of the proton charge squared
sigma the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm
epsilon the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol
replace determines the behavior if there is already an exception for the same two particles. If true, the existing one is replaced. If false, an exception is thrown.
Returns:
the index of the exception that was added
int addParticle ( double  charge,
double  sigma,
double  epsilon 
)

Add the nonbonded force parameters for a particle.

This should be called once for each particle in the System. When it is called for the i'th time, it specifies the parameters for the i'th particle. For calculating the Lennard-Jones interaction between two particles, the arithmetic mean of the sigmas and the geometric mean of the epsilons for the two interacting particles is used (the Lorentz-Bertelot combining rule).

Parameters:
charge the charge of the particle, measured in units of the proton charge
sigma the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm
epsilon the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol
Returns:
the index of the particle that was added
void createExceptionsFromBonds ( const std::vector< std::pair< int, int > > &  bonds,
double  coulomb14Scale,
double  lj14Scale 
)

Identify exceptions based on the molecular topology.

Particles which are separated by one or two bonds are set to not interact at all, while pairs of particles separated by three bonds (known as "1-4 interactions") have their Coulomb and Lennard-Jones interactions reduced by a fixed factor.

Parameters:
bonds the set of bonds based on which to construct exceptions. Each element specifies the indices of two particles that are bonded to each other.
coulomb14Scale pairs of particles separated by three bonds will have the strength of their Coulomb interaction multiplied by this factor
lj14Scale pairs of particles separated by three bonds will have the strength of their Lennard-Jones interaction multiplied by this factor
ForceImpl* createImpl (  )  [protected, virtual]

When a Context is created, it invokes this method on each Force in the System.

It should create a new ForceImpl object which can be used by the context for calculating forces. The ForceImpl will be deleted automatically when the Context is deleted.

Implements Force.

double getCutoffDistance (  )  const

Get the cutoff distance (in nm) being used for nonbonded interactions.

If the NonbondedMethod in use is NoCutoff, this value will have no effect.

Returns:
the cutoff distance, measured in nm
double getEwaldErrorTolerance (  )  const

Get the error tolerance for Ewald summation.

This corresponds to the fractional error in the forces which is acceptable. This value is used to select the reciprocal space cutoff and separation parameter so that the average error level will be less than the tolerance. There is not a rigorous guarantee that all forces on all atoms will be less than the tolerance, however.

void getExceptionParameters ( int  index,
int &  particle1,
int &  particle2,
double &  chargeProd,
double &  sigma,
double &  epsilon 
) const

Get the force field parameters for an interaction that should be calculated differently from others.

Parameters:
index the index of the interaction for which to get parameters
particle1 the index of the first particle involved in the interaction
particle2 the index of the second particle involved in the interaction
chargeProd the scaled product of the atomic charges (i.e. the strength of the Coulomb interaction), measured in units of the proton charge squared
sigma the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm
epsilon the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol
NonbondedMethod getNonbondedMethod (  )  const

Get the method used for handling long range nonbonded interactions.

int getNumExceptions (  )  const [inline]

Get the number of special interactions that should be calculated differently from other interactions.

int getNumParticles (  )  const [inline]

Get the number of particles for which force field parameters have been defined.

void getParticleParameters ( int  index,
double &  charge,
double &  sigma,
double &  epsilon 
) const

Get the nonbonded force parameters for a particle.

Parameters:
index the index of the particle for which to get parameters
charge the charge of the particle, measured in units of the proton charge
sigma the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm
epsilon the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol
double getReactionFieldDielectric (  )  const

Get the dielectric constant to use for the solvent in the reaction field approximation.

void setCutoffDistance ( double  distance  ) 

Set the cutoff distance (in nm) being used for nonbonded interactions.

If the NonbondedMethod in use is NoCutoff, this value will have no effect.

Parameters:
distance the cutoff distance, measured in nm
void setEwaldErrorTolerance ( double  tol  ) 

Get the error tolerance for Ewald summation.

This corresponds to the fractional error in the forces which is acceptable. This value is used to select the reciprocal space cutoff and separation parameter so that the average error level will be less than the tolerance. There is not a rigorous guarantee that all forces on all atoms will be less than the tolerance, however.

void setExceptionParameters ( int  index,
int  particle1,
int  particle2,
double  chargeProd,
double  sigma,
double  epsilon 
)

Set the force field parameters for an interaction that should be calculated differently from others.

If chargeProd and epsilon are both equal to 0, this will cause the interaction to be completely omitted from force and energy calculations.

Parameters:
index the index of the interaction for which to get parameters
particle1 the index of the first particle involved in the interaction
particle2 the index of the second particle involved in the interaction
chargeProd the scaled product of the atomic charges (i.e. the strength of the Coulomb interaction), measured in units of the proton charge squared
sigma the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm
epsilon the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol
void setNonbondedMethod ( NonbondedMethod  method  ) 

Set the method used for handling long range nonbonded interactions.

void setParticleParameters ( int  index,
double  charge,
double  sigma,
double  epsilon 
)

Set the nonbonded force parameters for a particle.

When calculating the Lennard-Jones interaction between two particles, it uses the arithmetic mean of the sigmas and the geometric mean of the epsilons for the two interacting particles (the Lorentz-Bertelot combining rule).

Parameters:
index the index of the particle for which to set parameters
charge the charge of the particle, measured in units of the proton charge
sigma the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm
epsilon the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol
void setReactionFieldDielectric ( double  dielectric  ) 

Set the dielectric constant to use for the solvent in the reaction field approximation.


The documentation for this class was generated from the following file:

Generated by  doxygen 1.6.2