Class VdwForce

java.lang.Object
ffx.openmm.Force
ffx.openmm.amoeba.VdwForce
Direct Known Subclasses:
AmoebaVdwForce

public class VdwForce extends Force
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" -- either Decouple or Annihilate. For Decouple, two alchemical atoms interact normally. For Annihilate, all interactions involving an alchemical atom are influenced. The softcore state is specified by setting a single 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. https://www.pnas.org/content/105/17/6290.

  • Constructor Details

    • VdwForce

      public VdwForce()
      Create an Amoeba VdwForce.
  • Method Details

    • addParticle

      public int addParticle(int parentIndex, double sigma, double epsilon, double reductionFactor, int isAlchemical, double scaleFactor)
      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.
      scaleFactor - a scale factor to apply to all interactions involving this particle (used for CpHMD).
      Returns:
      index of added particle
    • addParticle

      public int addParticle(int parentIndex, int typeIndex, double reductionFactor, int isAlchemical, double scaleFactor)
      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.
      scaleFactor - a scale factor to apply to all interactions involving this particle (used for CpHMD).
      Returns:
      index of added particle
    • addParticleType

      public 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.
    • addTypePair

      public 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.
    • destroy

      public void destroy()
      Destroy the force.
      Specified by:
      destroy in class Force
    • getAlchemicalMethod

      public int getAlchemicalMethod()
      Get the alchemical method.
      Returns:
      The alchemical method.
    • getCutoff

      @Deprecated public double getCutoff()
      Deprecated.
      This method exists only for backward compatibility. Use getCutoffDistance() instead.
      Get the cutoff distance.
    • getCutoffDistance

      public double getCutoffDistance()
      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
    • getEpsilonCombiningRule

      public String getEpsilonCombiningRule()
      Get the epsilon combining rule.
      Returns:
      The epsilon combining rule.
    • getLambda

      public com.sun.jna.Pointer getLambda()
      Get the lambda parameter.
      Returns:
      The lambda parameter.
    • getNonbondedMethod

      public int getNonbondedMethod()
      Get the nonbonded method.
      Returns:
      The nonbonded method.
    • getNumParticles

      public int getNumParticles()
      Get the number of particles.
      Returns:
      The number of particles.
    • getNumParticleTypes

      public int getNumParticleTypes()
      Get the number of particle types.
      Returns:
      The number of particle types.
    • getNumTypePairs

      public int getNumTypePairs()
      Get the number of type pairs.
      Returns:
      The number of type pairs.
    • getParticleExclusions

      public IntArray getParticleExclusions(int i)
      Get the particle exclusions.
      Parameters:
      i - The index of the particle.
      Returns:
      An IntArray containing the exclusions.
    • getParticleParameters

      public void getParticleParameters(int index, com.sun.jna.ptr.IntByReference ired, com.sun.jna.ptr.DoubleByReference rad, com.sun.jna.ptr.DoubleByReference eps, com.sun.jna.ptr.DoubleByReference reductionFactor, com.sun.jna.ptr.IntByReference isAlchemical, com.sun.jna.ptr.IntByReference type, com.sun.jna.ptr.DoubleByReference scaleFactor)
      Get the particle parameters.
      Parameters:
      index - The index of the particle.
      ired - The index of the particle that this particle is reduced to (output).
      rad - The radius of the particle (output).
      eps - The epsilon of the particle (output).
      reductionFactor - The reduction factor (output).
      isAlchemical - Whether the particle is alchemical (output).
      type - The type of the particle (output).
      scaleFactor - The scale factor (output).
    • getParticleTypeParameters

      public void getParticleTypeParameters(int index, com.sun.jna.ptr.DoubleByReference rad, com.sun.jna.ptr.DoubleByReference eps)
      Get the particle type parameters.
      Parameters:
      index - The index of the particle type.
      rad - The radius of the particle type (output).
      eps - The epsilon of the particle type (output).
    • getPotentialFunction

      public int getPotentialFunction()
      Get the potential function.
      Returns:
      The potential function.
    • getSigmaCombiningRule

      public String getSigmaCombiningRule()
      Get the sigma combining rule.
      Returns:
      The sigma combining rule.
    • getSoftcoreAlpha

      public double getSoftcoreAlpha()
      Get the softcore alpha.
      Returns:
      The softcore alpha.
    • getSoftcorePower

      public int getSoftcorePower()
      Get the softcore power.
      Returns:
      The softcore power.
    • getTypePairParameters

      public void getTypePairParameters(int index, com.sun.jna.ptr.IntByReference type1, com.sun.jna.ptr.IntByReference type2, com.sun.jna.ptr.DoubleByReference rad, com.sun.jna.ptr.DoubleByReference eps)
      Get the type pair parameters.
      Parameters:
      index - The index of the type pair.
      type1 - The first type (output).
      type2 - The second type (output).
      rad - The radius (output).
      eps - The epsilon (output).
    • getUseDispersionCorrection

      public boolean getUseDispersionCorrection()
      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.
    • getUseParticleTypes

      public int getUseParticleTypes()
      Get whether to use particle types.
      Returns:
      1 if particle types are used, 0 otherwise.
    • setAlchemicalMethod

      public void setAlchemicalMethod(int method)
      Set the alchemical method.
      Parameters:
      method - The alchemical method.
    • setCutoff

      @Deprecated public void setCutoff(double cutoff)
      Deprecated.
      This method exists only for backward compatibility. Use setCutoffDistance() instead.
      Set the cutoff distance.
    • setCutoffDistance

      public 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
    • setEpsilonCombiningRule

      public void setEpsilonCombiningRule(String rule)
      Set the epsilon combining rule.
      Parameters:
      rule - The epsilon combining rule.
    • setLambdaName

      public void setLambdaName(String name)
      Set the lambda parameter name.
      Parameters:
      name - The name of the lambda parameter.
    • setNonbondedMethod

      public void setNonbondedMethod(int method)
      Set the nonbonded method.
      Parameters:
      method - The nonbonded method.
    • setParticleExclusions

      public void setParticleExclusions(int i, IntArray exclusions)
      Set the particle exclusions.
      Parameters:
      i - The index of the particle.
      exclusions - The exclusions.
    • setParticleParameters

      public void setParticleParameters(int index, int ired, double rad, double eps, double reductionFactor, int isAlchemical, int type, double scaleFactor)
      Set the particle parameters.
      Parameters:
      index - The index of the particle.
      ired - The index of the particle that this particle is reduced to.
      rad - The radius of the particle.
      eps - The epsilon of the particle.
      reductionFactor - The reduction factor.
      isAlchemical - Whether the particle is alchemical.
      type - The type of the particle.
      scaleFactor - The scale factor.
    • setParticleTypeParameters

      public void setParticleTypeParameters(int index, double rad, double eps)
      Set the particle type parameters.
      Parameters:
      index - The index of the particle type.
      rad - The radius of the particle type.
      eps - The epsilon of the particle type.
    • setPotentialFunction

      public void setPotentialFunction(int function)
      Set the potential function.
      Parameters:
      function - The potential function.
    • setSigmaCombiningRule

      public void setSigmaCombiningRule(String rule)
      Set the sigma combining rule.
      Parameters:
      rule - The sigma combining rule.
    • setSoftcoreAlpha

      public void setSoftcoreAlpha(double vdWSoftcoreAlpha)
      Set the softcore alpha.
      Parameters:
      vdWSoftcoreAlpha - The softcore alpha.
    • setSoftcorePower

      public void setSoftcorePower(int vdwSoftcorePower)
      Set the softcore power.
      Parameters:
      vdwSoftcorePower - The softcore power.
    • setTypePairParameters

      public void setTypePairParameters(int index, int type1, int type2, double rad, double eps)
      Set the type pair parameters.
      Parameters:
      index - The index of the type pair.
      type1 - The first type.
      type2 - The second type.
      rad - The radius.
      eps - The epsilon.
    • setUseDispersionCorrection

      public void setUseDispersionCorrection(boolean 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.
    • updateParametersInContext

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

    • usesPeriodicBoundaryConditions

      public boolean usesPeriodicBoundaryConditions()
      Returns whether or not this force makes use of periodic boundary conditions.
      Overrides:
      usesPeriodicBoundaryConditions in class Force
      Returns:
      true if nonbondedMethod uses PBC and false otherwise