Package ffx.openmm

Class CustomTorsionForce

java.lang.Object
ffx.openmm.Force
ffx.openmm.CustomTorsionForce

public class CustomTorsionForce extends Force
This class implements interactions between sets of four particles that depend on the torsion angle between them. Unlike PeriodicTorsionForce, the functional form of the interaction is completely customizable, and may involve arbitrary algebraic expressions. In addition to the angle formed by the particles, it may depend on arbitrary global and per-torsion parameters.

To use this class, create a CustomTorsionForce object, passing an algebraic expression to the constructor that defines the interaction energy between each set of particles. The expression may depend on theta, the torsion angle formed by the particles, as well as on any parameters you choose. Then call addPerTorsionParameter() to define per-torsion parameters, and addGlobalParameter() to define global parameters. The values of per-torsion parameters are specified as part of the system definition, while values of global parameters may be modified during a simulation by calling Context::setParameter(). Finally, call addTorsion() once for each torsion. After an torsion has been added, you can modify its parameters by calling setTorsionParameters(). This will have no effect on Contexts that already exist unless you call updateParametersInContext(). Note that theta is guaranteed to be in the range [-pi,+pi], which may cause issues with force discontinuities if the energy function does not respect this domain.

As an example, the following code creates a CustomTorsionForce that implements a periodic potential:


 CustomTorsionForce force = new CustomTorsionForce("0.5*k*(1-cos(theta-theta0))");
 

This force depends on two parameters: the spring constant k and equilibrium angle theta0. The following code defines these parameters:


 force.addPerTorsionParameter("k");
 force.addPerTorsionParameter("theta0");
 

If a harmonic restraint is desired, it is important to be careful of the domain for theta, using an idiom like this:


 CustomTorsionForce force = new CustomTorsionForce("0.5*k*min(dtheta, 2*pi-dtheta)^2; dtheta = abs(theta-theta0); pi = 3.1415926535");
 

This class also has the ability to compute derivatives of the potential energy with respect to global parameters. Call addEnergyParameterDerivative() to request that the derivative with respect to a particular parameter be computed. You can then query its value in a Context by calling getState() on it.

Expressions may involve the operators + (add), - (subtract), * (multiply), / (divide), and ^ (power), and the following functions: sqrt, exp, log, sin, cos, sec, csc, tan, cot, asin, acos, atan, atan2, sinh, cosh, tanh, erf, erfc, min, max, abs, floor, ceil, step, delta, select. All trigonometric functions are defined in radians, and log is the natural logarithm. step(x) = 0 if x < 0, 1 otherwise. delta(x) = 1 if x = 0, 0 otherwise. select(x,y,z) = z if x = 0, y otherwise.

  • Constructor Details

    • CustomTorsionForce

      public CustomTorsionForce(String energy)
      Create a CustomTorsionForce.
      Parameters:
      energy - The algebraic expression that gives the interaction energy of each torsion as a function of theta, the torsion angle.
  • Method Details

    • addEnergyParameterDerivative

      public void addEnergyParameterDerivative(String name)
      Request that this Force compute the derivative of its energy with respect to a global parameter.
      Parameters:
      name - The name of the parameter.
    • addEnergyParameterDerivative

      public void addEnergyParameterDerivative(com.sun.jna.Pointer name)
      Request that this Force compute the derivative of its energy with respect to a global parameter.
      Parameters:
      name - The name of the parameter.
    • addGlobalParameter

      public int addGlobalParameter(String name, double defaultValue)
      Add a new global parameter that the interaction may depend on.
      Parameters:
      name - The name of the parameter.
      defaultValue - The default value of the parameter.
      Returns:
      The index of the parameter that was added.
    • addGlobalParameter

      public int addGlobalParameter(com.sun.jna.Pointer name, double defaultValue)
      Add a new global parameter that the interaction may depend on.
      Parameters:
      name - The name of the parameter.
      defaultValue - The default value of the parameter.
      Returns:
      The index of the parameter that was added.
    • addPerTorsionParameter

      public int addPerTorsionParameter(String name)
      Add a new per-torsion parameter that the interaction may depend on.
      Parameters:
      name - The name of the parameter.
      Returns:
      The index of the parameter that was added.
    • addPerTorsionParameter

      public int addPerTorsionParameter(com.sun.jna.Pointer name)
      Add a new per-torsion parameter that the interaction may depend on.
      Parameters:
      name - The name of the parameter.
      Returns:
      The index of the parameter that was added.
    • addTorsion

      public int addTorsion(int particle1, int particle2, int particle3, int particle4, com.sun.jna.ptr.PointerByReference parameters)
      Add a torsion to the Force.
      Parameters:
      particle1 - The index of the first particle forming the torsion.
      particle2 - The index of the second particle forming the torsion.
      particle3 - The index of the third particle forming the torsion.
      particle4 - The index of the fourth particle forming the torsion.
      parameters - The list of parameters for the new torsion.
      Returns:
      The index of the torsion that was added.
    • destroy

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

      public String getEnergyFunction()
      Get the algebraic expression that gives the interaction energy of each torsion.
      Returns:
      The energy expression.
    • getEnergyParameterDerivativeName

      public String getEnergyParameterDerivativeName(int index)
      Get the name of a parameter with respect to which the derivative of the energy should be computed.
      Parameters:
      index - The index of the parameter derivative, between 0 and getNumEnergyParameterDerivatives().
      Returns:
      The parameter name.
    • getGlobalParameterDefaultValue

      public double getGlobalParameterDefaultValue(int index)
      Get the default value of a global parameter.
      Parameters:
      index - The index of the parameter for which to get the default value.
      Returns:
      The parameter default value.
    • getGlobalParameterName

      public String getGlobalParameterName(int index)
      Get the name of a global parameter.
      Parameters:
      index - The index of the parameter for which to get the name.
      Returns:
      The parameter name.
    • getNumEnergyParameterDerivatives

      public int getNumEnergyParameterDerivatives()
      Get the number of parameters with respect to which the derivative of the energy should be computed.
      Returns:
      The number of parameters.
    • getNumGlobalParameters

      public int getNumGlobalParameters()
      Get the number of global parameters that the interaction depends on.
      Returns:
      The number of parameters.
    • getNumPerTorsionParameters

      public int getNumPerTorsionParameters()
      Get the number of per-torsion parameters that the interaction depends on.
      Returns:
      The number of parameters.
    • getNumTorsions

      public int getNumTorsions()
      Get the number of torsions for which force field parameters have been defined.
      Returns:
      The number of torsions.
    • getPerTorsionParameterName

      public String getPerTorsionParameterName(int index)
      Get the name of a per-torsion parameter.
      Parameters:
      index - The index of the parameter for which to get the name.
      Returns:
      The parameter name.
    • getTorsionParameters

      public void getTorsionParameters(int index, com.sun.jna.ptr.IntByReference particle1, com.sun.jna.ptr.IntByReference particle2, com.sun.jna.ptr.IntByReference particle3, com.sun.jna.ptr.IntByReference particle4, com.sun.jna.ptr.PointerByReference parameters)
      Get the force field parameters for a torsion.
      Parameters:
      index - The index of the torsion for which to get parameters.
      particle1 - The index of the first particle forming the torsion (output).
      particle2 - The index of the second particle forming the torsion (output).
      particle3 - The index of the third particle forming the torsion (output).
      particle4 - The index of the fourth particle forming the torsion (output).
      parameters - The list of parameters (output).
    • getTorsionParameters

      public void getTorsionParameters(int index, IntBuffer particle1, IntBuffer particle2, IntBuffer particle3, IntBuffer particle4, com.sun.jna.ptr.PointerByReference parameters)
      Get the force field parameters for a torsion.
      Parameters:
      index - The index of the torsion for which to get parameters.
      particle1 - The index of the first particle forming the torsion (output).
      particle2 - The index of the second particle forming the torsion (output).
      particle3 - The index of the third particle forming the torsion (output).
      particle4 - The index of the fourth particle forming the torsion (output).
      parameters - The list of parameters (output).
    • setEnergyFunction

      public void setEnergyFunction(String energy)
      Set the algebraic expression that gives the interaction energy of each torsion.
      Parameters:
      energy - The energy expression.
    • setEnergyFunction

      public void setEnergyFunction(com.sun.jna.Pointer energy)
      Set the algebraic expression that gives the interaction energy of each torsion.
      Parameters:
      energy - The energy expression.
    • setGlobalParameterDefaultValue

      public void setGlobalParameterDefaultValue(int index, double defaultValue)
      Set the default value of a global parameter.
      Parameters:
      index - The index of the parameter for which to set the default value.
      defaultValue - The default value of the parameter.
    • setGlobalParameterName

      public void setGlobalParameterName(int index, String name)
      Set the name of a global parameter.
      Parameters:
      index - The index of the parameter for which to set the name.
      name - The name of the parameter.
    • setGlobalParameterName

      public void setGlobalParameterName(int index, com.sun.jna.Pointer name)
      Set the name of a global parameter.
      Parameters:
      index - The index of the parameter for which to set the name.
      name - The name of the parameter.
    • setPerTorsionParameterName

      public void setPerTorsionParameterName(int index, String name)
      Set the name of a per-torsion parameter.
      Parameters:
      index - The index of the parameter for which to set the name.
      name - The name of the parameter.
    • setPerTorsionParameterName

      public void setPerTorsionParameterName(int index, com.sun.jna.Pointer name)
      Set the name of a per-torsion parameter.
      Parameters:
      index - The index of the parameter for which to set the name.
      name - The name of the parameter.
    • setTorsionParameters

      public void setTorsionParameters(int index, int particle1, int particle2, int particle3, int particle4, com.sun.jna.ptr.PointerByReference parameters)
      Set the force field parameters for a torsion.
      Parameters:
      index - The index of the torsion for which to set parameters.
      particle1 - The index of the first particle forming the torsion.
      particle2 - The index of the second particle forming the torsion.
      particle3 - The index of the third particle forming the torsion.
      particle4 - The index of the fourth particle forming the torsion.
      parameters - The list of parameters for the torsion.
    • setUsesPeriodicBoundaryConditions

      public void setUsesPeriodicBoundaryConditions(boolean periodic)
      Set whether this force should apply periodic boundary conditions when calculating displacements.
      Parameters:
      periodic - If true, periodic boundary conditions will be applied.
    • updateParametersInContext

      public void updateParametersInContext(Context context)
      Update the per-torsion parameters in a Context to match those stored in this Force object.
      Parameters:
      context - The Context in which to update the parameters.
    • usesPeriodicBoundaryConditions

      public boolean usesPeriodicBoundaryConditions()
      Check if the force uses periodic boundary conditions.
      Overrides:
      usesPeriodicBoundaryConditions in class Force
      Returns:
      True if the force uses periodic boundary conditions.