OpenMM
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Pages
AmoebaVdwForce Class Reference

This class implements a buffered 14-7 potential used to model van der Waals forces. More...

#include <AmoebaVdwForce.h>

+ Inheritance diagram for AmoebaVdwForce:

Public Types

enum  NonbondedMethod { NoCutoff = 0, CutoffPeriodic = 1 }
 This is an enumeration of the different methods that may be used for handling long range nonbonded forces. More...
 

Public Member Functions

 AmoebaVdwForce ()
 Create an Amoeba VdwForce. More...
 
int getNumParticles () const
 Get the number of particles. More...
 
void setParticleParameters (int particleIndex, int parentIndex, double sigma, double epsilon, double reductionFactor)
 Set the force field parameters for a vdw particle. More...
 
void getParticleParameters (int particleIndex, int &parentIndex, double &sigma, double &epsilon, double &reductionFactor) const
 Get the force field parameters for a vdw particle. More...
 
int addParticle (int parentIndex, double sigma, double epsilon, double reductionFactor)
 Add the force field parameters for a vdw particle. More...
 
void setSigmaCombiningRule (const std::string &sigmaCombiningRule)
 Set sigma combining rule. More...
 
const std::string & getSigmaCombiningRule (void) const
 Get sigma combining rule. More...
 
void setEpsilonCombiningRule (const std::string &epsilonCombiningRule)
 Set epsilon combining rule. More...
 
const std::string & getEpsilonCombiningRule (void) const
 Get epsilon combining rule. More...
 
bool getUseDispersionCorrection () const
 Get whether to add a contribution to the energy that approximately represents the effect of VdW interactions beyond the cutoff distance. More...
 
void setUseDispersionCorrection (bool useCorrection)
 Set whether to add a contribution to the energy that approximately represents the effect of VdW interactions beyond the cutoff distance. More...
 
void setParticleExclusions (int particleIndex, const std::vector< int > &exclusions)
 Set exclusions for specified particle. More...
 
void getParticleExclusions (int particleIndex, std::vector< int > &exclusions) const
 Get exclusions for specified particle. More...
 
void setCutoff (double cutoff)
 Set the cutoff distance. More...
 
double getCutoff () const
 Get the cutoff distance. More...
 
NonbondedMethod getNonbondedMethod () const
 Get the method used for handling long range nonbonded interactions. More...
 
void setNonbondedMethod (NonbondedMethod method)
 Set the method used for handling long range nonbonded interactions. More...
 
void updateParametersInContext (Context &context)
 Update the per-particle parameters in a Context to match those stored in this Force object. More...
 
- Public Member Functions inherited from Force
 Force ()
 
virtual ~Force ()
 
int getForceGroup () const
 Get the force group this Force belongs to. More...
 
void setForceGroup (int group)
 Set the force group this Force belongs to. More...
 

Protected Member Functions

ForceImplcreateImpl () const
 When a Context is created, it invokes this method on each Force in the System. More...
 
- Protected Member Functions inherited from Force
ForceImplgetImplInContext (Context &context)
 Get the ForceImpl corresponding to this Force in a Context. More...
 
ContextImplgetContextImpl (Context &context)
 Get the ContextImpl corresponding to a Context. More...
 

Detailed Description

This class implements a buffered 14-7 potential used to model van der Waals forces.

To use it, create an AmoebaVdwForce object then call addParticle() once for each particle. After a particle has been added, you can modify its force field parameters by calling setParticleParameters(). This will have no effect on Contexts that already exist unless you call updateParametersInContext().

A unique feature of this class is that the interaction site for a particle does not need to be exactly at the particle's location. Instead, it can be placed a fraction of the distance from that particle to another one. This is typically done for hydrogens to place the interaction site slightly closer to the parent atom. The fraction is known as the "reduction factor", since it reduces the distance from the parent atom to the interaction site.

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.

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.

Constructor & Destructor Documentation

Create an Amoeba VdwForce.

Member Function Documentation

int addParticle ( int  parentIndex,
double  sigma,
double  epsilon,
double  reductionFactor 
)

Add the force field parameters for a vdw particle.

Parameters
parentIndexthe index of the parent particle
sigmavdw sigma
epsilonvdw epsilon
reductionFactorthe fraction of the distance along the line from the parent particle to this particle at which the interaction site should be placed
Returns
index of added particle
ForceImpl* createImpl ( ) const
protectedvirtual

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 getCutoff ( ) const

Get the cutoff distance.

const std::string& getEpsilonCombiningRule ( void  ) const

Get epsilon combining rule.

Returns
epsilonCombiningRule epsilon combining rule: 'ARITHMETIC', 'GEOMETRIC'. 'HARMONIC', 'HHG'
NonbondedMethod getNonbondedMethod ( ) const

Get the method used for handling long range nonbonded interactions.

int getNumParticles ( ) const
inline

Get the number of particles.

void getParticleExclusions ( int  particleIndex,
std::vector< int > &  exclusions 
) const

Get exclusions for specified particle.

Parameters
particleIndexparticle index
exclusionsvector of exclusions
void getParticleParameters ( int  particleIndex,
int &  parentIndex,
double &  sigma,
double &  epsilon,
double &  reductionFactor 
) const

Get the force field parameters for a vdw particle.

Parameters
particleIndexthe particle index
parentIndexthe index of the parent particle
sigmavdw sigma
epsilonvdw epsilon
reductionFactorthe fraction of the distance along the line from the parent particle to this particle at which the interaction site should be placed
const std::string& getSigmaCombiningRule ( void  ) const

Get sigma combining rule.

Returns
sigmaCombiningRule sigma combining rule: 'ARITHMETIC', 'GEOMETRIC'. 'CUBIC-MEAN'
bool getUseDispersionCorrection ( ) const
inline

Get whether to add a contribution to the energy that approximately represents the effect of VdW interactions beyond the cutoff distance.

The energy depends on the volume of the periodic box, and is only applicable when periodic boundary conditions are used. When running simulations at constant pressure, adding this contribution can improve the quality of results.

void setCutoff ( double  cutoff)

Set the cutoff distance.

void setEpsilonCombiningRule ( const std::string &  epsilonCombiningRule)

Set epsilon combining rule.

Parameters
epsilonCombiningRuleepsilon combining rule: 'ARITHMETIC', 'GEOMETRIC'. 'HARMONIC', 'HHG'
void setNonbondedMethod ( NonbondedMethod  method)

Set the method used for handling long range nonbonded interactions.

void setParticleExclusions ( int  particleIndex,
const std::vector< int > &  exclusions 
)

Set exclusions for specified particle.

Parameters
particleIndexparticle index
exclusionsvector of exclusions
void setParticleParameters ( int  particleIndex,
int  parentIndex,
double  sigma,
double  epsilon,
double  reductionFactor 
)

Set the force field parameters for a vdw particle.

Parameters
particleIndexthe particle index
parentIndexthe index of the parent particle
sigmavdw sigma
epsilonvdw epsilon
reductionFactorthe fraction of the distance along the line from the parent particle to this particle at which the interaction site should be placed
void setSigmaCombiningRule ( const std::string &  sigmaCombiningRule)

Set sigma combining rule.

Parameters
sigmaCombiningRulesigma combining rule: 'ARITHMETIC', 'GEOMETRIC'. 'CUBIC-MEAN'
void setUseDispersionCorrection ( bool  useCorrection)
inline

Set whether to add a contribution to the energy that approximately represents the effect of VdW interactions beyond the cutoff distance.

The energy depends on the volume of the periodic box, and is only applicable when periodic boundary conditions are used. When running simulations at constant pressure, adding this contribution can improve the quality of results.

void updateParametersInContext ( Context context)

Update the per-particle parameters in a Context to match those stored in this Force object.

This method provides an efficient method to update certain parameters in an existing Context without needing to reinitialize it. Simply call setParticleParameters() to modify this object's parameters, then call updateParametersInState() to copy them over to the Context.

The only information this method updates is the values of per-particle parameters. All other aspects of the Force (the nonbonded method, the cutoff distance, etc.) are unaffected and can only be changed by reinitializing the Context.


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