Class CustomNonbondedForce
To use this class, create a CustomNonbondedForce object, passing an algebraic expression to the constructor that defines the interaction energy between each pair of particles. The expression may depend on r, the distance between the particles, as well as on any parameters you choose. Then call addPerParticleParameter() to define per-particle parameters, and addGlobalParameter() to define global parameters. The values of per-particle parameters are specified as part of the system definition, while values of global parameters may be modified during a simulation by calling Context::setParameter().
Next, call addParticle() once for each particle in the System to set the values of its per-particle parameters. 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().
CustomNonbondedForce also lets you specify "exclusions", particular pairs of particles whose interactions should be omitted from force and energy calculations. This is most often used for particles that are bonded to each other.
As an example, the following code creates a CustomNonbondedForce that implements a 12-6 Lennard-Jones potential:
CustomNonbondedForce force = new CustomNonbondedForce("4*epsilon*((sigma/r)^12-(sigma/r)^6); sigma=0.5*(sigma1+sigma2); epsilon=sqrt(epsilon1*epsilon2)");
This force depends on two parameters: sigma and epsilon. The following code defines these as per-particle parameters:
force.addPerParticleParameter("sigma");
force.addPerParticleParameter("epsilon");
The expression must be symmetric with respect to the two particles. It typically will only be evaluated once for each pair of particles, and no guarantee is made about which particle will be identified as "particle 1". In the above example, the energy only depends on the products sigma1*sigma2 and epsilon1*epsilon2, both of which are unchanged if the labels 1 and 2 are reversed. In contrast, if it depended on the difference sigma1-sigma2, the results would be undefined, because reversing the labels 1 and 2 would change the energy.
The energy also may depend on "computed values". These are similar to per-particle parameters, but instead of being specified in advance, their values are computed based on global and per-particle parameters. For example, the following code uses a global parameter (lambda) to interpolate between two different sigma values for each particle (sigmaA and sigmaB).
CustomNonbondedForce force = new CustomNonbondedForce("4*epsilon*((sigma/r)^12-(sigma/r)^6); sigma=0.5*(sigma1+sigma2); epsilon=sqrt(epsilon1*epsilon2)");
force.addComputedValue("sigma", "(1-lambda)*sigmaA + lambda*sigmaB");
force.addGlobalParameter("lambda", 0);
force.addPerParticleParameter("sigmaA");
force.addPerParticleParameter("sigmaB");
force.addPerParticleParameter("epsilon");
You could, of course, embed the computation of sigma directly into the energy expression, but then it would need to be repeated for every interaction. By separating it out as a computed value, it only needs to be computed once for each particle instead of once for each interaction, thus saving computation time.
CustomNonbondedForce can operate in two modes. By default, it computes the interaction of every particle in the System with every other particle. Alternatively, you can restrict it to only a subset of particle pairs. To do this, specify one or more "interaction groups". An interaction group consists of two sets of particles that should interact with each other. Every particle in the first set interacts with every particle in the second set. For example, you might use this feature to compute a solute-solvent interaction energy, while omitting all interactions between two solute atoms or two solvent atoms.
-
Field Summary
-
Constructor Summary
Constructors -
Method Summary
Modifier and TypeMethodDescriptionint
addComputedValue
(String name, String expression) Add a computed value to the force.void
Add an energy parameter derivative to the force.int
addExclusion
(int particle1, int particle2) Add an exclusion to the force.int
addFunction
(String name, com.sun.jna.ptr.PointerByReference function, double min, double max) Add a tabulated function that may appear in the energy expression.int
addGlobalParameter
(String name, double value) Add a global parameter to the force.int
addInteractionGroup
(IntSet group1, IntSet group2) Add an interaction group to the force.int
addParticle
(DoubleArray parameters) Add a particle to the force.int
Add a per-particle parameter to the force.int
addTabulatedFunction
(String name, com.sun.jna.ptr.PointerByReference function) Add a tabulated function that may appear in the energy expression.void
createExclusionsFromBonds
(BondArray bonds, int bondCutoff) Create exclusions based on the molecular topology.void
destroy()
Destroy the force.void
getComputedValueParameters
(int index, com.sun.jna.ptr.PointerByReference name, com.sun.jna.ptr.PointerByReference expression) Get the parameters for a computed value.double
Get the cutoff distance.Get the energy expression for the force.getEnergyParameterDerivativeName
(int index) Get the name of a parameter with respect to which the derivative of the energy should be computed.void
getExclusionParticles
(int index, com.sun.jna.ptr.IntByReference particle1, com.sun.jna.ptr.IntByReference particle2) Get the particles in an exclusion.void
getExclusionParticles
(int index, IntBuffer particle1, IntBuffer particle2) Get the particles in an exclusion.void
getFunctionParameters
(int index, com.sun.jna.ptr.PointerByReference name, com.sun.jna.ptr.PointerByReference function, com.sun.jna.ptr.DoubleByReference min, com.sun.jna.ptr.DoubleByReference max) Get the parameters for a tabulated function.void
getFunctionParameters
(int index, com.sun.jna.ptr.PointerByReference name, com.sun.jna.ptr.PointerByReference function, DoubleBuffer min, DoubleBuffer max) Get the parameters for a tabulated function.double
getGlobalParameterDefaultValue
(int index) Get the default value of a global parameter.getGlobalParameterName
(int index) Get the name of a global parameter.void
getInteractionGroupParameters
(int index, com.sun.jna.ptr.PointerByReference group1, com.sun.jna.ptr.PointerByReference group2) Get the parameters for an interaction group.int
Get the nonbonded method.int
Get the number of computed values.int
Get the number of parameters with respect to which the derivative of the energy should be computed.int
Get the number of exclusions.int
Deprecated.This method exists only for backward compatibility.int
Get the number of global parameters.int
Get the number of interaction groups.int
Get the number of particles.int
Get the number of per-particle parameters.int
Get the number of tabulated functions.void
getParticleParameters
(int index, com.sun.jna.ptr.PointerByReference parameters) Get the parameters for a particle.getPerParticleParameterName
(int index) Get the name of a per-particle parameter.double
Get the switching distance.com.sun.jna.ptr.PointerByReference
getTabulatedFunction
(int index) Get a reference to a tabulated function that may appear in the energy expression.getTabulatedFunctionName
(int index) Get the name of a tabulated function that may appear in the energy expression.int
Get whether to use the long range correction.int
Get whether to use a switching function.void
setComputedValueParameters
(int index, String name, String expression) Set the parameters for a computed value.void
setCutoffDistance
(double distance) Set the cutoff distance.void
setEnergyFunction
(String expression) Set the energy expression for the force.void
setExclusionParticles
(int index, int particle1, int particle2) Set the particles in an exclusion.void
setFunctionParameters
(int index, String name, com.sun.jna.ptr.PointerByReference function, double min, double max) Set the parameters for a tabulated function.void
setGlobalParameterDefaultValue
(int index, double value) Set the default value of a global parameter.void
setGlobalParameterName
(int index, String name) Set the name of a global parameter.void
setInteractionGroupParameters
(int index, com.sun.jna.ptr.PointerByReference group1, com.sun.jna.ptr.PointerByReference group2) Set the parameters for an interaction group.void
setNonbondedMethod
(int method) Set the nonbonded method.void
setParticleParameters
(int index, com.sun.jna.ptr.PointerByReference parameters) Set the parameters for a particle.void
setPerParticleParameterName
(int index, String name) Set the name of a per-particle parameter.void
setSwitchingDistance
(double distance) Set the switching distance.void
setUseLongRangeCorrection
(int use) Set whether to use the long range correction.void
setUseSwitchingFunction
(int use) Set whether to use a switching function.void
updateParametersInContext
(Context context) Update the per-particle parameters and tabulated functions in a Context to match those stored in this Force object.boolean
Check if the force uses periodic boundary conditions.Methods inherited from class ffx.openmm.Force
getForceGroup, getForceIndex, getName, getPointer, setForceGroup, setForceIndex, setName
-
Constructor Details
-
CustomNonbondedForce
Constructor.- Parameters:
energy
- The energy expression for the force.
-
-
Method Details
-
addComputedValue
Add a computed value to the force.- Parameters:
name
- The name of the computed value.expression
- The expression for computing the value.- Returns:
- The index of the computed value that was added.
-
addEnergyParameterDerivative
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.
-
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
public int addFunction(String name, com.sun.jna.ptr.PointerByReference function, double min, double max) 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.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
Add a global parameter to the force.- Parameters:
name
- The name of the parameter.value
- The default value of the parameter.- Returns:
- The index of the parameter that was added.
-
addInteractionGroup
Add an interaction group to the force.- Parameters:
group1
- The set of particles in the first group.group2
- The set of particles in the second group.- Returns:
- The index of the interaction group that was added.
-
addParticle
Add a particle to the force.- Parameters:
parameters
- The parameters for the new particle.- Returns:
- The index of the particle that was added.
-
addPerParticleParameter
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
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.
-
createExclusionsFromBonds
Create exclusions based on the molecular topology.- Parameters:
bonds
- The bonds to use for determining exclusions.bondCutoff
- The number of bonds within which particles should be excluded.
-
destroy
public void destroy()Destroy the force. -
getComputedValueParameters
public void getComputedValueParameters(int index, com.sun.jna.ptr.PointerByReference name, com.sun.jna.ptr.PointerByReference expression) Get the parameters for a computed value.- Parameters:
index
- The index of the computed value to get.name
- The name of the computed value (output).expression
- The expression for computing the value (output).
-
getCutoffDistance
public double getCutoffDistance()Get the cutoff distance.- Returns:
- The cutoff distance, measured in nm.
-
getEnergyFunction
Get the energy expression for the force.- Returns:
- The energy expression for the force.
-
getEnergyParameterDerivativeName
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.
-
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 to get.particle1
- The index of the first particle in the exclusion (output).particle2
- The index of the second particle in the exclusion (output).
-
getExclusionParticles
Get the particles in an exclusion.- Parameters:
index
- The index of the exclusion to get.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 function, 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 to get.name
- The name of the function as it appears in expressions (output).function
- A TabulatedFunction object defining 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 function, DoubleBuffer min, DoubleBuffer max) Get the parameters for a tabulated function.- Parameters:
index
- The index of the function to get.name
- The name of the function as it appears in expressions (output).function
- A TabulatedFunction object defining 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 to get.- Returns:
- The default value of the parameter.
-
getGlobalParameterName
Get the name of a global parameter.- Parameters:
index
- The index of the parameter to get.- Returns:
- The name of the parameter.
-
getInteractionGroupParameters
public void getInteractionGroupParameters(int index, com.sun.jna.ptr.PointerByReference group1, com.sun.jna.ptr.PointerByReference group2) Get the parameters for an interaction group.- Parameters:
index
- The index of the interaction group to get.group1
- The set of particles in the first group (output).group2
- The set of particles in the second group (output).
-
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.
-
getNumExclusions
public int getNumExclusions()Get the number of exclusions.- Returns:
- The number of exclusions.
-
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.
-
getNumInteractionGroups
public int getNumInteractionGroups()Get the number of interaction groups.- Returns:
- The number of interaction groups.
-
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, com.sun.jna.ptr.PointerByReference parameters) Get the parameters for a particle.- Parameters:
index
- The index of the particle to get.parameters
- The parameters for the particle (output).
-
getPerParticleParameterName
Get the name of a per-particle parameter.- Parameters:
index
- The index of the parameter to get.- Returns:
- The name of the parameter.
-
getSwitchingDistance
public double getSwitchingDistance()Get the switching distance.- Returns:
- The switching distance, measured in nm.
-
getTabulatedFunction
public com.sun.jna.ptr.PointerByReference getTabulatedFunction(int index) Get a reference to a tabulated function that may appear in the energy expression.- Parameters:
index
- The index of the function to get.- Returns:
- The TabulatedFunction object defining the function.
-
getTabulatedFunctionName
Get the name of a tabulated function that may appear in the energy expression.- Parameters:
index
- The index of the function to get.- Returns:
- The name of the function as it appears in expressions.
-
getUseLongRangeCorrection
public int getUseLongRangeCorrection()Get whether to use the long range correction.- Returns:
- 1 if the long range correction is used, 0 otherwise.
-
getUseSwitchingFunction
public int getUseSwitchingFunction()Get whether to use a switching function.- Returns:
- 1 if a switching function is used, 0 otherwise.
-
setComputedValueParameters
Set the parameters for a computed value.- Parameters:
index
- The index of the computed value to set.name
- The name of the computed value.expression
- The expression for computing the value.
-
setCutoffDistance
public void setCutoffDistance(double distance) Set the cutoff distance.- Parameters:
distance
- The cutoff distance, measured in nm.
-
setEnergyFunction
Set the energy expression for the force.- Parameters:
expression
- The energy expression for the force.
-
setExclusionParticles
public void setExclusionParticles(int index, int particle1, int particle2) Set the particles in an exclusion.- Parameters:
index
- The index of the exclusion to set.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 function, double min, double max) Set the parameters for a tabulated function.- Parameters:
index
- The index of the function to set.name
- The name of the function as it appears in expressions.function
- A TabulatedFunction object defining 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 value) Set the default value of a global parameter.- Parameters:
index
- The index of the parameter to set.value
- The default value of the parameter.
-
setGlobalParameterName
Set the name of a global parameter.- Parameters:
index
- The index of the parameter to set.name
- The name of the parameter.
-
setInteractionGroupParameters
public void setInteractionGroupParameters(int index, com.sun.jna.ptr.PointerByReference group1, com.sun.jna.ptr.PointerByReference group2) Set the parameters for an interaction group.- Parameters:
index
- The index of the interaction group to set.group1
- The set of particles in the first group.group2
- The set of particles in the second group.
-
setNonbondedMethod
public void setNonbondedMethod(int method) Set the nonbonded method.- Parameters:
method
- The nonbonded method.
-
setParticleParameters
public void setParticleParameters(int index, com.sun.jna.ptr.PointerByReference parameters) Set the parameters for a particle.- Parameters:
index
- The index of the particle to set.parameters
- The parameters for the particle.
-
setPerParticleParameterName
Set the name of a per-particle parameter.- Parameters:
index
- The index of the parameter to set.name
- The name of the parameter.
-
setSwitchingDistance
public void setSwitchingDistance(double distance) Set the switching distance.- Parameters:
distance
- The switching distance, measured in nm.
-
setUseLongRangeCorrection
public void setUseLongRangeCorrection(int use) Set whether to use the long range correction.- Parameters:
use
- 1 to use the long range correction, 0 otherwise.
-
setUseSwitchingFunction
public void setUseSwitchingFunction(int use) Set whether to use a switching function.- Parameters:
use
- 1 to use a switching function, 0 otherwise.
-
updateParametersInContext
Update the per-particle parameters and tabulated functions 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 classForce
- Returns:
- True if the force uses periodic boundary conditions.
-