Class NonbondedForce
- Direct Known Subclasses:
FixedChargeNonbondedForce
To use this class, create a NonbondedForce object, then call addParticle() once for each particle in the System to define its parameters. The number of particles for which you define nonbonded 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 force field parameters by calling setParticleParameters(). This will have no effect on Contexts that already exist unless you call updateParametersInContext().
NonbondedForce also lets you specify "exceptions", particular pairs of particles whose interactions should be computed based on different parameters than those defined for the individual particles. This can be used to completely exclude certain interactions from the force calculation, or to alter how they interact with each other.
Many molecular force fields omit Coulomb and Lennard-Jones interactions between particles separated by one or two bonds, while using modified parameters for those separated by three bonds (known as "1-4 interactions"). This class provides a convenience method for this case called createExceptionsFromBonds(). You pass to it a list of bonds and the scale factors to use for 1-4 interactions. It identifies all pairs of particles which are separated by 1, 2, or 3 bonds, then automatically creates exceptions for them.
When using a cutoff, by default Lennard-Jones interactions are sharply truncated at the cutoff distance. Optionally you can instead use a switching function to make the interaction smoothly go to zero over a finite distance range. To enable this, call setUseSwitchingFunction(). You must also call setSwitchingDistance() to specify the distance at which the interaction should begin to decrease. The switching distance must be less than the cutoff distance.
Another optional feature of this class (enabled by default) is to add a contribution to the energy which approximates the effect of all Lennard-Jones interactions beyond the cutoff in a periodic system. When running a simulation at constant pressure, this can improve the quality of the result. Call setUseDispersionCorrection() to set whether this should be used.
In some applications, it is useful to be able to inexpensively change the parameters of small groups of particles. Usually this is done to interpolate between two sets of parameters. For example, a titratable group might have two states it can exist in, each described by a different set of parameters for the atoms that make up the group. You might then want to smoothly interpolate between the two states. This is done by first calling addGlobalParameter() to define a Context parameter, then addParticleParameterOffset() to create a "parameter offset" that depends on the Context parameter. Each offset defines the following:
- A Context parameter used to interpolate between the states.
- A single particle whose parameters are influenced by the Context parameter.
- Three scale factors (chargeScale, sigmaScale, and epsilonScale) that specify how the Context parameter affects the particle.
The "effective" parameters for a particle (those used to compute forces) are given by
charge = baseCharge + param*chargeScale
sigma = baseSigma + param*sigmaScale
epsilon = baseEpsilon + param*epsilonScale
where the "base" values are the ones specified by addParticle() and "param" is the current value
of the Context parameter. A single Context parameter can apply offsets to multiple particles,
and multiple parameters can be used to apply offsets to the same particle. Parameters can also be used
to modify exceptions in exactly the same way by calling addExceptionParameterOffset().-
Field Summary
-
Constructor Summary
Constructors -
Method Summary
Modifier and TypeMethodDescriptionint
addException
(int particle1, int particle2, double chargeProd, double sigma, double epsilon, boolean replace) Add an exception to the force.int
addParticle
(double charge, double sigma, double eps) Add a particle.void
createExceptionsFromBonds
(BondArray bondArray, double coulomb14Scale, double lj14Scale) Create exceptions from bonds.void
destroy()
Destroy the force.double
Get the cutoff distance.void
getExceptionParameters
(int index, com.sun.jna.ptr.IntByReference particle1, com.sun.jna.ptr.IntByReference particle2, com.sun.jna.ptr.DoubleByReference chargeProd, com.sun.jna.ptr.DoubleByReference sigma, com.sun.jna.ptr.DoubleByReference eps) Get the exception parameters.int
Get the number of exceptions.int
Get the number of particles.void
getParticleParameters
(int index, com.sun.jna.ptr.DoubleByReference charge, com.sun.jna.ptr.DoubleByReference sigma, com.sun.jna.ptr.DoubleByReference eps) Get the particle parameters.void
setCutoffDistance
(double cutoffDistance) Set the cutoff distance.void
setExceptionParameters
(int index, int particle1, int particle2, double chargeProd, double sigma, double eps) Set the exception parameters.void
setNonbondedMethod
(int method) Set the nonbonded method.void
setParticleParameters
(int index, double charge, double sigma, double eps) Set the particle parameters.void
setPMEParameters
(double aEwald, int nx, int ny, int nz) Set the PME parameters.void
setSwitchingDistance
(double switchingDistance) Set the switching distance.void
setUseDispersionCorrection
(int useDispersionCorrection) Set if a dispersion correction will be used.void
setUseSwitchingFunction
(int useSwitchingFunction) Set if a switching function will be used.void
updateParametersInContext
(Context context) Update the parameters in the context.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
-
NonbondedForce
public NonbondedForce()Create a new NonbondedForce.
-
-
Method Details
-
addException
public int addException(int particle1, int particle2, double chargeProd, double sigma, double epsilon, boolean replace) Add an exception to the force.- Parameters:
particle1
- The index of the first particle.particle2
- The index of the second particle.chargeProd
- The charge product.sigma
- The sigma vdW parameter.epsilon
- The epsilon vdW parameter.replace
- If true, replace any existing exception.- Returns:
- The index of the exception that was added.
-
addParticle
public int addParticle(double charge, double sigma, double eps) Add a particle.- Parameters:
charge
- The atomic charge.sigma
- The vdW sigma.eps
- The vdW eps.- Returns:
- The index of the particle that was added.
-
createExceptionsFromBonds
Create exceptions from bonds.- Parameters:
bondArray
- The bond array.coulomb14Scale
- The coulomb 1-4 scale.lj14Scale
- The LJ 1-4 scale.
-
destroy
public void destroy()Destroy the force. -
getCutoffDistance
public double getCutoffDistance()Get the cutoff distance.- Returns:
- The cutoff distance.
-
getExceptionParameters
public void getExceptionParameters(int index, com.sun.jna.ptr.IntByReference particle1, com.sun.jna.ptr.IntByReference particle2, com.sun.jna.ptr.DoubleByReference chargeProd, com.sun.jna.ptr.DoubleByReference sigma, com.sun.jna.ptr.DoubleByReference eps) Get the exception parameters.- Parameters:
index
- The exception index.particle1
- The first particle.particle2
- The second particle.chargeProd
- The charge product.sigma
- The sigma vdW parameter.eps
- The eps vdW parameter.
-
getNumExceptions
public int getNumExceptions()Get the number of exceptions.- Returns:
- The number of exceptions.
-
getNumParticles
public int getNumParticles()Get the number of particles.- Returns:
- The number of particles.
-
getParticleParameters
public void getParticleParameters(int index, com.sun.jna.ptr.DoubleByReference charge, com.sun.jna.ptr.DoubleByReference sigma, com.sun.jna.ptr.DoubleByReference eps) Get the particle parameters.- Parameters:
index
- The particle index.charge
- The atomic charge.sigma
- The vdW sigma.eps
- The vdW eps.
-
setCutoffDistance
public void setCutoffDistance(double cutoffDistance) Set the cutoff distance.- Parameters:
cutoffDistance
- The cutoff distance.
-
setExceptionParameters
public void setExceptionParameters(int index, int particle1, int particle2, double chargeProd, double sigma, double eps) Set the exception parameters.- Parameters:
index
- The exception index.particle1
- The first particle.particle2
- The second particle.chargeProd
- The charge product.sigma
- The sigma vdW parameter.eps
- The eps vdW parameter.
-
setNonbondedMethod
public void setNonbondedMethod(int method) Set the nonbonded method.- Parameters:
method
- The nonbonded method.
-
setPMEParameters
public void setPMEParameters(double aEwald, int nx, int ny, int nz) Set the PME parameters.- Parameters:
aEwald
- The Ewald alpha.nx
- The PME grid size in x.ny
- The PME grid size in y.nz
- The PME grid size in z.
-
setParticleParameters
public void setParticleParameters(int index, double charge, double sigma, double eps) Set the particle parameters.- Parameters:
index
- The particle index.charge
- The atomic charge.sigma
- The vdW sigma.eps
- The vdW eps.
-
setSwitchingDistance
public void setSwitchingDistance(double switchingDistance) Set the switching distance.- Parameters:
switchingDistance
- The switching distance.
-
setUseDispersionCorrection
public void setUseDispersionCorrection(int useDispersionCorrection) Set if a dispersion correction will be used.- Parameters:
useDispersionCorrection
- The dispersion correction flag.
-
setUseSwitchingFunction
public void setUseSwitchingFunction(int useSwitchingFunction) Set if a switching function will be used.- Parameters:
useSwitchingFunction
- The switching distance flag.
-
updateParametersInContext
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 classForce
- Returns:
- True if the force uses periodic boundary conditions.
-