AmoebaVdwForce

class OpenMM::AmoebaVdwForce

This class models van der Waals forces in the AMOEBA force field. It can use either buffered 14-7 potential or a Lennard-Jones 12-6 potential.

This class can operate in two different modes. In one mode, force field parameters are defined for each particle. When two particles interact, a combining rule is used to calculate the interaction parameters based on the parameters for the two particles. To use the class in this mode, call the version of addParticle() that takes sigma and epsilon values. It should be called once for each particle in the System.

In the other mode, each particle has a type index, and parameters are specified for each type rather than each individual particle. By default this mode also uses a combining rule, but you can override it by defining alternate parameters to use for specific pairs of particle types. To use the class in this mode, call the version of addParticle() that takes a type index. It should be called once for each particle in the System. You also must call addParticleType() once for each type. If you wish to override the combining for particular pairs of types, do so by calling addTypePair().

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.

Support is also available for softcore interactions based on setting a per particle alchemical flag and setting the AmoebaVdwForce to use an “AlchemicalMethod” Context parameter “AmoebaVdwLambda” between 0.0 and 1.0.

The softcore functional form can be modified by setting the softcore power (default of 5) and the softcore alpha (default of 0,7). For more information on the softcore functional form see Eq. 2 from: Jiao, D.; Golubkov, P. A.; Darden, T. A.; Ren, P., Calculation of protein-ligand binding free energy by using a polarizable potential. Proc. Natl. Acad. Sci. U.S.A. 2008, 105 (17), 6290-6295.

Methods

AmoebaVdwForce

Create an Amoeba VdwForce.

getNumParticles

Get the number of particles

getNumParticleTypes

Get the number of particle types.

getNumTypePairs

Get the number of type pairs.

setParticleParameters

Set the force field parameters for a vdw particle.

getParticleParameters

Get the force field parameters for a vdw particle.

addParticle

Add the force field parameters for a vdw particle.

addParticle

Add the force field parameters for a vdw particle.

addParticleType

Add a particle type.

getParticleTypeParameters

Get the force field parameters for a particle type.

setParticleTypeParameters

Set the force field parameters for a particle type.

addTypePair

Add a type pair.

getTypePairParameters

Get the force field parameters for a type pair.

setTypePairParameters

Set the force field parameters for a type pair.

setSigmaCombiningRule

Set sigma combining rule

getSigmaCombiningRule

Get sigma combining rule

setEpsilonCombiningRule

Set epsilon combining rule

getEpsilonCombiningRule

Get epsilon combining rule

getUseDispersionCorrection

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

setUseDispersionCorrection

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

getUseParticleTypes

Get whether parameters were specified by particle or by particle type.

setParticleExclusions

Set exclusions for specified particle

getParticleExclusions

Get exclusions for specified particle

getCutoffDistance

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

setCutoffDistance

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

setCutoff

Set the cutoff distance.

getCutoff

Get the cutoff distance.

getNonbondedMethod

Get the method used for handling long range nonbonded interactions.

setNonbondedMethod

Set the method used for handling long range nonbonded interactions.

getPotentialFunction

Get the potential function to use.

setPotentialFunction

Set the potential function to use.

setSoftcorePower

Set the softcore power on lambda (default = 5).

getSoftcorePower

Get the softcore power on lambda.

setSoftcoreAlpha

Set the softcore alpha value (default = 0.7).

getSoftcoreAlpha

Get the softcore alpha value.

getAlchemicalMethod

Get the method used for alchemical interactions.

setAlchemicalMethod

Set the method used for handling long range nonbonded interactions.

updateParametersInContext

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

usesPeriodicBoundaryConditions

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

Enum: NonbondedMethod

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.

Enum: PotentialFunction

Buffered147

Use a buffered 14-7 potential. This is the default.

LennardJones

Use a Lennard-Jones 12-6 potential.

Enum: AlchemicalMethod

None

All vdW interactions are treated normally. This is the default.

Decouple

Maintain full strength vdW interactions between two alchemical particles.

Annihilate

Interactions between two alchemical particles are turned off at lambda=0.

const std::string &Lambda()

This is the name of the parameter which stores the current Amoeba vdW lambda value.

AmoebaVdwForce()

Create an Amoeba VdwForce.

int getNumParticles() const

Get the number of particles

int getNumParticleTypes() const

Get the number of particle types.

int getNumTypePairs() const

Get the number of type pairs.

void setParticleParameters(int particleIndex, int parentIndex, double sigma, double epsilon, double reductionFactor, bool isAlchemical = false, int typeIndex = -1)

Set the force field parameters for a vdw particle.

Parameters

  • particleIndex – the particle index

  • parentIndex – the index of the parent particle

  • sigma – vdw sigma

  • epsilon – vdw epsilon

  • reductionFactor – the fraction of the distance along the line from the parent particle to this particle at which the interaction site should be placed

  • isAlchemical – if true, this vdW particle is undergoing an alchemical change.

  • typeIndex – the index of the particle type for this particle

void getParticleParameters(int particleIndex, int &parentIndex, double &sigma, double &epsilon, double &reductionFactor, bool &isAlchemical, int &typeIndex) const

Get the force field parameters for a vdw particle.

Parameters

  • particleIndex – the particle index

  • parentIndex – [out] the index of the parent particle

  • sigma – [out] vdw sigma

  • epsilon – [out] vdw epsilon

  • reductionFactor – [out] the fraction of the distance along the line from the parent particle to this particle at which the interaction site should be placed

  • isAlchemical – [out] if true, this vdW particle is undergoing an alchemical change.

  • typeIndex – [out] the index of the particle type for this particle

int addParticle(int parentIndex, double sigma, double epsilon, double reductionFactor, bool isAlchemical = false)

Add the force field parameters for a vdw particle. This version is used when parameters are defined for each particle.

Parameters

  • parentIndex – the index of the parent particle

  • sigma – vdw sigma

  • epsilon – vdw epsilon

  • reductionFactor – the fraction of the distance along the line from the parent particle to this particle at which the interaction site should be placed

  • isAlchemical – if true, this vdW particle is undergoing an alchemical change.

Returns

index of added particle

int addParticle(int parentIndex, int typeIndex, double reductionFactor, bool isAlchemical = false)

Add the force field parameters for a vdw particle. This version is used when parameters are defined by particle type.

Parameters

  • parentIndex – the index of the parent particle

  • typeIndex – the index of the particle type for this particle

  • reductionFactor – the fraction of the distance along the line from the parent particle to this particle at which the interaction site should be placed

  • isAlchemical – if true, this vdW particle is undergoing an alchemical change.

Returns

index of added particle

int addParticleType(double sigma, double epsilon)

Add a particle type.

Parameters

  • sigma – the sigma value for particles of this type

  • epsilon – the epsilon value for particles of this type

Returns

the index of the particle type that was just added.

void getParticleTypeParameters(int typeIndex, double &sigma, double &epsilon) const

Get the force field parameters for a particle type.

Parameters

  • typeIndex – the index of the particle type

  • sigma – [out] the sigma value for particles of this type

  • epsilon – [out] the epsilon value for particles of this type

void setParticleTypeParameters(int typeIndex, double sigma, double epsilon)

Set the force field parameters for a particle type.

Parameters

  • typeIndex – the index of the particle type

  • sigma – the sigma value for particles of this type

  • epsilon – the epsilon value for particles of this type

int addTypePair(int type1, int type2, double sigma, double epsilon)

Add a type pair. This overrides the standard combining rule for interactions between particles of two particular types.

Parameters

  • type1 – the index of the first particle type

  • type2 – the index of the second particle type

  • sigma – the sigma value for interactions between particles of these two types

  • epsilon – the epsilon value for interactions between particles of these two types

Returns

the index of the type pair that was just added.

void getTypePairParameters(int pairIndex, int &type1, int &type2, double &sigma, double &epsilon) const

Get the force field parameters for a type pair. This overrides the standard combining rule for interactions between particles of two particular types.

Parameters

  • pairIndex – the index of the type pair

  • type1 – [out] the index of the first particle type

  • type2 – [out] the index of the second particle type

  • sigma – [out] the sigma value for interactions between particles of these two types

  • epsilon – [out] the epsilon value for interactions between particles of these two types

void setTypePairParameters(int pairIndex, int type1, int type2, double sigma, double epsilon)

Set the force field parameters for a type pair. This overrides the standard combining rule for interactions between particles of two particular types.

Parameters

  • pairIndex – the index of the type pair

  • type1 – the index of the first particle type

  • type2 – the index of the second particle type

  • sigma – the sigma value for interactions between particles of these two types

  • epsilon – the epsilon value for interactions between particles of these two types

void setSigmaCombiningRule(const std::string &sigmaCombiningRule)

Set sigma combining rule

Parameters

  • sigmaCombiningRule – sigma combining rule: ‘ARITHMETIC’, ‘GEOMETRIC’. ‘CUBIC-MEAN’

const std::string &getSigmaCombiningRule(void) const

Get sigma combining rule

Returns

sigmaCombiningRule sigma combining rule: ‘ARITHMETIC’, ‘GEOMETRIC’. ‘CUBIC-MEAN’

void setEpsilonCombiningRule(const std::string &epsilonCombiningRule)

Set epsilon combining rule

Parameters

  • epsilonCombiningRule – epsilon combining rule: ‘ARITHMETIC’, ‘GEOMETRIC’. ‘HARMONIC’, ‘W-H’, ‘HHG’

const std::string &getEpsilonCombiningRule(void) const

Get epsilon combining rule

Returns

epsilonCombiningRule epsilon combining rule: ‘ARITHMETIC’, ‘GEOMETRIC’. ‘HARMONIC’, ‘W-H’, ‘HHG’

bool getUseDispersionCorrection() const

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 setUseDispersionCorrection(bool 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.

bool getUseParticleTypes() const

Get whether parameters were specified by particle or by particle type.

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

Set exclusions for specified particle

Parameters

  • particleIndex – particle index

  • exclusions – vector of exclusions

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

Get exclusions for specified particle

Parameters

  • particleIndex – particle index

  • exclusions – [out] vector of exclusions

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

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 setCutoff(double cutoff)

Set the cutoff distance.

Deprecated

This method exists only for backward compatibility. Use setCutoffDistance() instead.

double getCutoff() const

Get the cutoff distance.

Deprecated

This method exists only for backward compatibility. Use getCutoffDistance() instead.

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.

PotentialFunction getPotentialFunction() const

Get the potential function to use.

void setPotentialFunction(PotentialFunction potential)

Set the potential function to use.

void setSoftcorePower(int n)

Set the softcore power on lambda (default = 5).

int getSoftcorePower() const

Get the softcore power on lambda.

void setSoftcoreAlpha(double alpha)

Set the softcore alpha value (default = 0.7).

double getSoftcoreAlpha() const

Get the softcore alpha value.

AlchemicalMethod getAlchemicalMethod() const

Get the method used for alchemical interactions.

void setAlchemicalMethod(AlchemicalMethod method)

Set the method used for handling long range nonbonded interactions.

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 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.

bool usesPeriodicBoundaryConditions() const

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

Returns

true if nonbondedMethod uses PBC and false otherwise