OpenMM
AmoebaVdwForce Class Reference

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

+ Inheritance diagram for AmoebaVdwForce:

List of all members.

Public Member Functions

def getNumParticles
 getNumParticles(self) -> int
def setParticleParameters
 Set the force field parameters for a vdw particle.
def getParticleParameters
 Get the force field parameters for a vdw particle.
def addParticle
 addParticle(self, parentIndex, sigma, epsilon, reductionFactor) -> int
def setSigmaCombiningRule
 Set sigma combining rule.
def getSigmaCombiningRule
 getSigmaCombiningRule(self) -> std::string const &
def setEpsilonCombiningRule
 Set epsilon combining rule.
def getEpsilonCombiningRule
 getEpsilonCombiningRule(self) -> std::string const &
def getUseDispersionCorrection
 getUseDispersionCorrection(self) -> bool
def setUseDispersionCorrection
 Set whether to add a contribution to the energy that approximately represents the effect of VdW interactions beyond the cutoff distance.
def setParticleExclusions
 Set exclusions for specified particle.
def getParticleExclusions
 Get exclusions for specified particle.
def setCutoff
 Set the cutoff distance.
def getCutoff
 getCutoff(self) -> double
def getNonbondedMethod
 getNonbondedMethod(self) -> OpenMM::AmoebaVdwForce::NonbondedMethod
def setNonbondedMethod
 Set the method used for handling long range nonbonded interactions.
def updateParametersInContext
 Update the per-particle parameters in a Context to match those stored in this Force object.
def usesPeriodicBoundaryConditions
 usesPeriodicBoundaryConditions(self) -> bool
def __init__
 __init__(self) -> AmoebaVdwForce __init__(self, other) -> AmoebaVdwForce

Public Attributes

 this

Static Public Attributes

 NoCutoff = _openmm.AmoebaVdwForce_NoCutoff
 CutoffPeriodic = _openmm.AmoebaVdwForce_CutoffPeriodic

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.


Constructor & Destructor Documentation

def __init__ (   self,
  args 
)

__init__(self) -> AmoebaVdwForce __init__(self, other) -> AmoebaVdwForce

Create an Amoeba VdwForce.


Member Function Documentation

def addParticle (   self,
  parentIndex,
  sigma,
  epsilon,
  reductionFactor 
)

addParticle(self, parentIndex, sigma, epsilon, reductionFactor) -> int

Add the force field parameters for a vdw particle.

Parameters:
parentIndex(int) the index of the parent particle
sigma(double) vdw sigma
epsilon(double) vdw epsilon
reductionFactor(double) the fraction of the distance along the line from the parent particle to this particle at which the interaction site should be placed
Returns:
(int) index of added particle
def getCutoff (   self)

getCutoff(self) -> double

Get the cutoff distance.

def getEpsilonCombiningRule (   self)

getEpsilonCombiningRule(self) -> std::string const &

Get epsilon combining rule

Returns:
(string) epsilonCombiningRule epsilon combining rule: 'ARITHMETIC', 'GEOMETRIC'. 'HARMONIC', 'HHG'
def getNonbondedMethod (   self)

getNonbondedMethod(self) -> OpenMM::AmoebaVdwForce::NonbondedMethod

Get the method used for handling long range nonbonded interactions.

def getNumParticles (   self)

getNumParticles(self) -> int

Get the number of particles

def getParticleExclusions (   self,
  particleIndex 
)

Get exclusions for specified particle.

Parameters:
particleIndex(int) particle index
Returns:
(vector< int >) vector of exclusions
def getParticleParameters (   self,
  particleIndex 
)

Get the force field parameters for a vdw particle.

Parameters:
particleIndex(int) the particle index
Returns:
(int) the index of the parent particle
(double) vdw sigma
(double) vdw epsilon
(double) the fraction of the distance along the line from the parent particle to this particle at which the interaction site should be placed
def getSigmaCombiningRule (   self)

getSigmaCombiningRule(self) -> std::string const &

Get sigma combining rule

Returns:
(string) sigmaCombiningRule sigma combining rule: 'ARITHMETIC', 'GEOMETRIC'. 'CUBIC-MEAN'

getUseDispersionCorrection(self) -> bool

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.

def setCutoff (   self,
  cutoff 
)

Set the cutoff distance.

def setEpsilonCombiningRule (   self,
  epsilonCombiningRule 
)

Set epsilon combining rule.

Parameters:
epsilonCombiningRule(string) epsilon combining rule: 'ARITHMETIC', 'GEOMETRIC'. 'HARMONIC', 'HHG'
def setNonbondedMethod (   self,
  method 
)

Set the method used for handling long range nonbonded interactions.

def setParticleExclusions (   self,
  particleIndex,
  exclusions 
)

Set exclusions for specified particle.

Parameters:
particleIndex(int) particle index
exclusions(vector< int >) vector of exclusions
def setParticleParameters (   self,
  particleIndex,
  parentIndex,
  sigma,
  epsilon,
  reductionFactor 
)

Set the force field parameters for a vdw particle.

Parameters:
particleIndex(int) the particle index
parentIndex(int) the index of the parent particle
sigma(double) vdw sigma
epsilon(double) vdw epsilon
reductionFactor(double) the fraction of the distance along the line from the parent particle to this particle at which the interaction site should be placed
def setSigmaCombiningRule (   self,
  sigmaCombiningRule 
)

Set sigma combining rule.

Parameters:
sigmaCombiningRule(string) sigma combining rule: 'ARITHMETIC', 'GEOMETRIC'. 'CUBIC-MEAN'
def setUseDispersionCorrection (   self,
  useCorrection 
)

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.

def updateParametersInContext (   self,
  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 updateParametersInContext() 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.

usesPeriodicBoundaryConditions(self) -> bool

Returns whether or not this force makes use of periodic boundary conditions.

Returns:
(bool) true if nonbondedMethod uses PBC and false otherwise

Reimplemented from Force.


Member Data Documentation

CutoffPeriodic = _openmm.AmoebaVdwForce_CutoffPeriodic [static]
NoCutoff = _openmm.AmoebaVdwForce_NoCutoff [static]

Reimplemented from Force.


The documentation for this class was generated from the following file:
 All Classes Functions Variables