Package ffx.openmm

Class CustomGBForce

java.lang.Object
ffx.openmm.Force
ffx.openmm.CustomGBForce
Direct Known Subclasses:
FixedChargeGBForce

public class CustomGBForce extends Force
This class implements complex, multiple stage nonbonded interactions between particles. It is designed primarily for implementing Generalized Born implicit solvation models, although it is not strictly limited to that purpose. The interaction is specified as a series of computations, each defined by an arbitrary algebraic expression. It also allows tabulated functions to be defined and used with the computations. It optionally supports periodic boundary conditions and cutoffs for long range interactions.

The computation consists of calculating some number of per-particle computed values, followed by one or more energy terms. A computed value is a scalar value that is computed for each particle in the system. It may depend on an arbitrary set of global and per-particle parameters, and well as on other computed values that have been calculated before it. Once all computed values have been calculated, the energy terms and their derivatives are evaluated to determine the system energy and particle forces. The energy terms may depend on global parameters, per-particle parameters, and per-particle computed values.

When specifying a computed value or energy term, you provide an algebraic expression to evaluate and a computation type describing how the expression is to be evaluated. There are two main types of computations:

  • Single Particle: The expression is evaluated once for each particle in the System. In the case of a computed value, this means the value for a particle depends only on other properties of that particle (its position, parameters, and other computed values). In the case of an energy term, it means each particle makes an independent contribution to the System energy.
  • Particle Pairs: The expression is evaluated for every pair of particles in the system. In the case of a computed value, the value for a particular particle is calculated by pairing it with every other particle in the system, evaluating the expression for each pair, and summing them. For an energy term, each particle pair makes an independent contribution to the System energy. (Note that energy terms are assumed to be symmetric with respect to the two interacting particles, and therefore are evaluated only once per pair. In contrast, expressions for computed values need not be symmetric and therefore are calculated twice for each pair: once when calculating the value for the first particle, and again when calculating the value for the second particle.)

Be aware that, although this class is extremely general in the computations it can define, particular Platforms may only support more restricted types of computations. In particular, all currently existing Platforms require that the first computed value must be a particle pair computation, and all computed values after the first must be single particle computations. This is sufficient for most Generalized Born models, but might not permit some other types of calculations to be implemented.

This is a complicated class to use, and an example may help to clarify it. The following code implements the OBC variant of the GB/SA solvation model, using the ACE approximation to estimate surface area:

   
    CustomGBForce* custom = new CustomGBForce();
    custom->addPerParticleParameter("q");
    custom->addPerParticleParameter("radius");
    custom->addPerParticleParameter("scale");
    custom->addGlobalParameter("solventDielectric", obc->getSolventDielectric());
    custom->addGlobalParameter("soluteDielectric", obc->getSoluteDielectric());
    custom->addComputedValue("I", "step(r+sr2-or1)*0.5*(1/L-1/U+0.25*(1/U^2-1/L^2)*(r-sr2*sr2/r)+0.5*log(L/U)/r+C);"
                                  "U=r+sr2;"
                                  "C=2*(1/or1-1/L)*step(sr2-r-or1);"
                                  "L=max(or1, D);"
                                  "D=abs(r-sr2);"
                                  "sr2 = scale2*or2;"
                                  "or1 = radius1-0.009; or2 = radius2-0.009", CustomGBForce::ParticlePairNoExclusions);
    custom->addComputedValue("B", "1/(1/or-tanh(1*psi-0.8*psi^2+4.85*psi^3)/radius);"
                                  "psi=I*or; or=radius-0.009", CustomGBForce::SingleParticle);
    custom->addEnergyTerm("28.3919551*(radius+0.14)^2*(radius/B)^6-0.5*138.935456*(1/soluteDielectric-1/solventDielectric)*q^2/B",
                          CustomGBForce::SingleParticle);
    custom->addEnergyTerm("-138.935456*(1/soluteDielectric-1/solventDielectric)*q1*q2/f;"
                          "f=sqrt(r^2+B1*B2*exp(-r^2/(4*B1*B2)))", CustomGBForce::ParticlePair);
   
 

It begins by defining three per-particle parameters (charge, atomic radius, and scale factor) and two global parameters (the dielectric constants for the solute and solvent). It then defines a computed value "I" of type ParticlePair. The expression for evaluating it is a complicated function of the distance between each pair of particles (r), their atomic radii (radius1 and radius2), and their scale factors (scale1 and scale2). Very roughly speaking, it is a measure of the distance between each particle and other nearby particles.

Next a computation is defined for the Born Radius (B). It is computed independently for each particle, and is a function of that particle's atomic radius and the intermediate value I defined above.

Finally, two energy terms are defined. The first one is computed for each particle and represents the surface area term, as well as the self interaction part of the polarization energy. The second term is calculated for each pair of particles, and represents the screening of electrostatic interactions by the solvent.

After defining the force as shown above, you should then call addParticle() once for each particle in the System to set the values of its per-particle parameters (q, radius, and scale). The number of particles for which you set parameters must be exactly equal to the number of particles in the System, or else an exception will be thrown when you try to create a Context. After a particle has been added, you can modify its parameters by calling setParticleParameters(). This will have no effect on Contexts that already exist unless you call updateParametersInContext().

CustomGBForce also lets you specify "exclusions", particular pairs of particles whose interactions should be omitted from calculations. This is most often used for particles that are bonded to each other. Even if you specify exclusions, however, you can use the computation type ParticlePairNoExclusions to indicate that exclusions should not be applied to a particular piece of the computation.

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 is less than 0, 1 otherwise. delta(x) = 1 if x is 0, 0 otherwise. select(x,y,z) = z if x = 0, y otherwise. In expressions for particle pair calculations, the names of per-particle parameters and computed values have the suffix "1" or "2" appended to them to indicate the values for the two interacting particles. As seen in the above example, an expression may also involve intermediate quantities that are defined following the main expression, using ";" as a separator.

In addition, you can call addTabulatedFunction() to define a new function based on tabulated values. You specify the function by creating a TabulatedFunction object. That function can then appear in expressions.

  • Constructor Details

    • CustomGBForce

      public CustomGBForce()
      Create a CustomGBForce.
  • Method Details

    • addComputedValue

      public int addComputedValue(String name, String expression, int type)
      Add a computed value to the force.
      Parameters:
      name - The name of the computed value.
      expression - The expression for computing the value.
      type - The type of computation to perform.
      Returns:
      The index of the computed value that was added.
    • addEnergyParameterDerivative

      public void addEnergyParameterDerivative(String name)
      Add an energy parameter derivative to the force.
      Parameters:
      name - The name of the parameter to compute the derivative of the energy with respect to.
    • addEnergyTerm

      public int addEnergyTerm(String expression, int type)
      Add an energy term to the force.
      Parameters:
      expression - The expression for computing the energy term.
      type - The type of computation to perform.
      Returns:
      The index of the energy term that was added.
    • addExclusion

      public int addExclusion(int particle1, int particle2)
      Add an exclusion to the force.
      Parameters:
      particle1 - The index of the first particle in the exclusion.
      particle2 - The index of the second particle in the exclusion.
      Returns:
      The index of the exclusion that was added.
    • addFunction

      @Deprecated public int addFunction(String name, com.sun.jna.ptr.PointerByReference values, double min, double max)
      Deprecated.
      This method exists only for backward compatibility. Use addTabulatedFunction() instead.
      Add a tabulated function that may appear in the energy expression.
      Parameters:
      name - The name of the function as it appears in expressions.
      values - The tabulated values of the function.
      min - The minimum value of the independent variable for which the function is defined.
      max - The maximum value of the independent variable for which the function is defined.
      Returns:
      The index of the function that was added.
    • addGlobalParameter

      public int addGlobalParameter(String name, double defaultValue)
      Add a global parameter to the force.
      Parameters:
      name - The name of the parameter.
      defaultValue - The default value of the parameter.
      Returns:
      The index of the parameter that was added.
    • addParticle

      public int addParticle(DoubleArray particleParameters)
      Add a particle to the force.
      Parameters:
      particleParameters - The parameters for the particle.
      Returns:
      The index of the particle that was added.
    • addPerParticleParameter

      public int addPerParticleParameter(String name)
      Add a per-particle parameter to the force.
      Parameters:
      name - The name of the parameter.
      Returns:
      The index of the parameter that was added.
    • addTabulatedFunction

      public int addTabulatedFunction(String name, com.sun.jna.ptr.PointerByReference function)
      Add a tabulated function that may appear in the energy expression.
      Parameters:
      name - The name of the function as it appears in expressions.
      function - A TabulatedFunction object defining the function.
      Returns:
      The index of the function that was added.
    • destroy

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

      public void getComputedValueParameters(int index, com.sun.jna.ptr.PointerByReference name, com.sun.jna.ptr.PointerByReference expression, IntBuffer type)
      Get the parameters for a computed value.
      Parameters:
      index - The index of the computed value.
      name - The name of the computed value (output).
      expression - The expression for computing the value (output).
      type - The type of computation to perform (output).
    • getComputedValueParameters

      public void getComputedValueParameters(int index, com.sun.jna.ptr.PointerByReference name, com.sun.jna.ptr.PointerByReference expression, com.sun.jna.ptr.IntByReference type)
      Get the parameters for a computed value.
      Parameters:
      index - The index of the computed value.
      name - The name of the computed value (output).
      expression - The expression for computing the value (output).
      type - The type of computation to perform (output).
    • getCutoffDistance

      public double getCutoffDistance()
      Get the cutoff distance.
      Returns:
      The cutoff distance.
    • 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.
      Returns:
      The name of the parameter.
    • getEnergyTermParameters

      public void getEnergyTermParameters(int index, com.sun.jna.ptr.PointerByReference expression, IntBuffer type)
      Get the parameters for an energy term.
      Parameters:
      index - The index of the energy term.
      expression - The expression for computing the energy term (output).
      type - The type of computation to perform (output).
    • getEnergyTermParameters

      public void getEnergyTermParameters(int index, com.sun.jna.ptr.PointerByReference expression, com.sun.jna.ptr.IntByReference type)
      Get the parameters for an energy term.
      Parameters:
      index - The index of the energy term.
      expression - The expression for computing the energy term (output).
      type - The type of computation to perform (output).
    • getExclusionParticles

      public void getExclusionParticles(int index, IntBuffer particle1, IntBuffer particle2)
      Get the particles in an exclusion.
      Parameters:
      index - The index of the exclusion.
      particle1 - The index of the first particle in the exclusion (output).
      particle2 - The index of the second particle in the exclusion (output).
    • getExclusionParticles

      public void getExclusionParticles(int index, com.sun.jna.ptr.IntByReference particle1, com.sun.jna.ptr.IntByReference particle2)
      Get the particles in an exclusion.
      Parameters:
      index - The index of the exclusion.
      particle1 - The index of the first particle in the exclusion (output).
      particle2 - The index of the second particle in the exclusion (output).
    • getFunctionParameters

      public void getFunctionParameters(int index, com.sun.jna.ptr.PointerByReference name, com.sun.jna.ptr.PointerByReference values, DoubleBuffer min, DoubleBuffer max)
      Get the parameters for a tabulated function.
      Parameters:
      index - The index of the function.
      name - The name of the function (output).
      values - The tabulated values of the function (output).
      min - The minimum value of the independent variable for which the function is defined (output).
      max - The maximum value of the independent variable for which the function is defined (output).
    • getFunctionParameters

      public void getFunctionParameters(int index, com.sun.jna.ptr.PointerByReference name, com.sun.jna.ptr.PointerByReference values, com.sun.jna.ptr.DoubleByReference min, com.sun.jna.ptr.DoubleByReference max)
      Get the parameters for a tabulated function.
      Parameters:
      index - The index of the function.
      name - The name of the function (output).
      values - The tabulated values of the function (output).
      min - The minimum value of the independent variable for which the function is defined (output).
      max - The maximum value of the independent variable for which the function is defined (output).
    • getGlobalParameterDefaultValue

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

      public String getGlobalParameterName(int index)
      Get the name of a global parameter.
      Parameters:
      index - The index of the parameter.
      Returns:
      The name of the parameter.
    • getNonbondedMethod

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

      public int getNumComputedValues()
      Get the number of computed values.
      Returns:
      The number of computed values.
    • 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.
    • getNumEnergyTerms

      public int getNumEnergyTerms()
      Get the number of energy terms.
      Returns:
      The number of energy terms.
    • getNumExclusions

      public int getNumExclusions()
      Get the number of exclusions.
      Returns:
      The number of exclusions.
    • getNumFunctions

      @Deprecated public int getNumFunctions()
      Deprecated.
      This method exists only for backward compatibility. Use getNumTabulatedFunctions() instead.
      Get the number of tabulated functions.
      Returns:
      The number of tabulated functions.
    • getNumGlobalParameters

      public int getNumGlobalParameters()
      Get the number of global parameters.
      Returns:
      The number of global parameters.
    • getNumParticles

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

      public int getNumPerParticleParameters()
      Get the number of per-particle parameters.
      Returns:
      The number of per-particle parameters.
    • getNumTabulatedFunctions

      public int getNumTabulatedFunctions()
      Get the number of tabulated functions.
      Returns:
      The number of tabulated functions.
    • getParticleParameters

      public void getParticleParameters(int index, DoubleArray particleParameters)
      Get the parameters for a particle.
      Parameters:
      index - The index of the particle.
      particleParameters - The parameters for the particle (output).
    • getPerParticleParameterName

      public String getPerParticleParameterName(int index)
      Get the name of a per-particle parameter.
      Parameters:
      index - The index of the parameter.
      Returns:
      The name of the parameter.
    • getTabulatedFunction

      public com.sun.jna.ptr.PointerByReference getTabulatedFunction(int index)
      Get a reference to a tabulated function.
      Parameters:
      index - The index of the function.
      Returns:
      A reference to the function.
    • getTabulatedFunctionName

      public String getTabulatedFunctionName(int index)
      Get the name of a tabulated function.
      Parameters:
      index - The index of the function.
      Returns:
      The name of the function.
    • setComputedValueParameters

      public void setComputedValueParameters(int index, String name, String expression, int type)
      Set the parameters for a computed value.
      Parameters:
      index - The index of the computed value.
      name - The name of the computed value.
      expression - The expression for computing the value.
      type - The type of computation to perform.
    • setCutoffDistance

      public void setCutoffDistance(double distance)
      Set the cutoff distance.
      Parameters:
      distance - The cutoff distance.
    • setEnergyTermParameters

      public void setEnergyTermParameters(int index, String expression, int type)
      Set the parameters for an energy term.
      Parameters:
      index - The index of the energy term.
      expression - The expression for computing the energy term.
      type - The type of computation to perform.
    • setExclusionParticles

      public void setExclusionParticles(int index, int particle1, int particle2)
      Set the particles in an exclusion.
      Parameters:
      index - The index of the exclusion.
      particle1 - The index of the first particle in the exclusion.
      particle2 - The index of the second particle in the exclusion.
    • setFunctionParameters

      public void setFunctionParameters(int index, String name, com.sun.jna.ptr.PointerByReference values, double min, double max)
      Set the parameters for a tabulated function.
      Parameters:
      index - The index of the function.
      name - The name of the function.
      values - The tabulated values of the function.
      min - The minimum value of the independent variable for which the function is defined.
      max - The maximum value of the independent variable for which the function is defined.
    • setGlobalParameterDefaultValue

      public void setGlobalParameterDefaultValue(int index, double defaultValue)
      Set the default value of a global parameter.
      Parameters:
      index - The index of the parameter.
      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.
      name - The name of the parameter.
    • setNonbondedMethod

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

      public void setParticleParameters(int index, DoubleArray particleParameters)
      Set the parameters for a particle.
      Parameters:
      index - The index of the particle.
      particleParameters - The parameters for the particle.
    • setPerParticleParameterName

      public void setPerParticleParameterName(int index, String name)
      Set the name of a per-particle parameter.
      Parameters:
      index - The index of the parameter.
      name - The name of the parameter.
    • updateParametersInContext

      public void updateParametersInContext(Context context)
      Update the parameters in the context.
      Parameters:
      context - The context.
    • 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.