SlicedNonbondedForce

class SlicedNonbondedForce(*args)

Bases: NonbondedForce

This class implements sliced nonbonded interactions between particles, including a Coulomb force to represent electrostatics and a Lennard-Jones force to represent van der Waals interactions. It optionally supports periodic boundary conditions and cutoffs for long range interactions. Lennard-Jones interactions are calculated with the Lorentz-Berthelot combining rule: it uses the arithmetic mean of the sigmas and the geometric mean of the epsilons for the two interacting particles.

The particles are classified into \(n\) disjoint subsets, thus creating \(n(n+1)/2\) distinct types of particle pairs. Then, the total potential energy is sliced into this many parts by distinguishing the contributions of such pair types. Optionally, variables stored as an openmm.Context global parameters can multiply the Coulomb and/or Lennard-Jones parts of selected potential energy slices. Changes in the values of these global parameters via setParameter will automatically modify the total potential energy of a SlicedNonbondedForce. Derivatives with respect to scaling parameters can also be calculated upon request.

To use this class, create a SlicedNonbondedForce object, specifying the number \(n\) of particle subsets. Then, call addParticle() once for every single particle in the System to define its force field parameters. Alternatively, you can pass an existing openmm.NonbondedForce object when creating a SlicedNonbondedForce object, so that the latter inherits all properties of the former. After a particle has been added, you can modify its subset index from 0, the default value, to any positive integer lower than \(n\) by calling setParticleSubset(). You can also modify the force field parameters of an added particle by calling setParticleParameters(). These two methods will have no effect on Contexts that already exist unless you call updateParametersInContext().

SlicedNonbondedForce also lets you specify exceptions via addException() and modify their parameters via setExceptionParameters(). Exceptions are 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 exclude certain interactions from the force calculation.

By default, all energy slices of a SlicedNonbondedForce are added together, yielding the same total potential energy as a NonbondedForce with equal particle and exception parameters. In addition, you can select one or more slices and multiply their Coulomb and/or Lennard-Jones contributions by variables referred to as scaling parameters. For this, start by calling addGlobalParameter() to define a new Context parameter, then pass its name to addScalingParameter(). The same scaling parameter can affect multiple slices, and two scaling parameters can affect a single slice, but only if they multiply the Coulomb and Lennard-Jones parts of that slice separately.

With setScalingParameterDerivative(), you can request the derivative with respect to a scaling parameter. All requested derivatives can then be evaluated for a particular openmm.State by calling its getEnergyParameterDerivatives method. For this, such State must have been generated from an openmm.Context by passing True to the keyword getParameterDerivatives of the Context’s getState method.

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(). A context parameter cannot be used simultaneously to apply offsets and scale energy terms.

__init__(*args)

Overload 1:

Create a SlicedNonbondedForce.

Parameters:

numSubsets (int) – the number of particle subsets


Overload 2:

Create a SlicedNonbondedForce having the properties of an existing openmm.NonbondedForce.

Parameters:
  • force (openmm.NonbondedForce) – the NonbondedForce object from which to instantiate this SlicedNonbondedForce

  • numSubsets (int) – the number of particle subsets

Methods

addException(self, particle1, particle2, ...)

Add an interaction to the list of exceptions that should be calculated differently from other interactions.

addExceptionParameterOffset(self, parameter, ...)

Add an offset to the parameters of a particular exception, based on a global parameter.

addException_usingRMin(particle1, particle2, ...)

Add interaction exception using the product of the two atoms' elementary charges, rMin and epsilon, which is standard for AMBER force fields.

addGlobalParameter(self, name, defaultValue)

Add a new global parameter that parameter offsets may depend on.

addParticle(self, charge, sigma, epsilon)

Add the nonbonded force parameters for a particle.

addParticleParameterOffset(self, parameter, ...)

Add an offset to the per-particle parameters of a particular particle, based on a global parameter.

addParticle_usingRVdw(charge, rVDW, epsilon)

Add particle using elemetrary charge.

addScalingParameter(parameter, subset1, ...)

Add a scaling parameter to multiply a particular Coulomb slice.

addScalingParameterDerivative(parameter)

Request the derivative of this Force's energy with respect to a scaling parameter.

createExceptionsFromBonds(self, bonds, ...)

Identify exceptions based on the molecular topology.

getCutoffDistance(self)

Get the cutoff distance (in nm) being used for nonbonded interactions.

getEwaldErrorTolerance(self)

Get the error tolerance for Ewald summation.

getExceptionParameterOffset(self, index)

Get the offset added to the parameters of a particular exception, based on a global parameter.

getExceptionParameters(self, index)

Get the force field parameters for an interaction that should be calculated differently from others.

getExceptionsUsePeriodicBoundaryConditions(self)

Get whether periodic boundary conditions should be applied to exceptions.

getForceGroup(self)

Get the force group this Force belongs to.

getGlobalParameterDefaultValue(self, index)

Get the default value of a global parameter.

getGlobalParameterName(self, index)

Get the name of a global parameter.

getIncludeDirectSpace(self)

Get whether to include direct space interactions when calculating forces and energies.

getLJPMEParameters(self)

Get the parameters to use for dispersion term in LJ-PME calculations.

getLJPMEParametersInContext(context)

Get the PME parameters being used for the dispersion term for LJPME in a particular Context.

getName(self)

Get the name of this Force.

getNonbondedMethod(self)

Get the method used for handling long range nonbonded interactions.

getNonbondedMethodName()

Get the name of the method used for handling long range nonbonded interactions.

getNumExceptionParameterOffsets(self)

Get the number of exception parameter offsets that have been added.

getNumExceptions(self)

Get the number of special interactions that should be calculated differently from other interactions.

getNumGlobalParameters(self)

Get the number of global parameters that have been added.

getNumParticleParameterOffsets(self)

Get the number of particles parameter offsets that have been added.

getNumParticles(self)

Get the number of particles for which force field parameters have been defined.

getNumScalingParameterDerivatives()

Get the number of requested scaling parameter derivatives.

getNumScalingParameters()

Get the number of scaling parameters.

getNumSlices()

Get the number of slices determined by the specified number of particle subsets.

getNumSubsets()

Get the specified number of particle subsets.

getPMEParameters(self)

Get the parameters to use for PME calculations.

getPMEParametersInContext(context)

Get the parameters being used for PME in a particular Context.

getParticleParameterOffset(self, index)

Get the offset added to the per-particle parameters of a particular particle, based on a global parameter.

getParticleParameters(self, index)

Get the nonbonded force parameters for a particle.

getParticleSubset(index)

Get the subset to which a particle belongs.

getReactionFieldDielectric(self)

Get the dielectric constant to use for the solvent in the reaction field approximation.

getReciprocalSpaceForceGroup(self)

Get the force group that reciprocal space interactions for Ewald or PME are included in.

getScalingParameter(index, parameter, arg4)

Get the scaling parameter applied to a particular nonbonded slice.

getScalingParameterDerivativeName(index)

Get the name of the global parameter associated with a requested scaling parameter derivative.

getSwitchingDistance(self)

Get the distance at which the switching function begins to reduce the Lennard-Jones interaction.

getUseCudaFFT()

Get whether to use CUDA Toolkit's cuFFT library when executing in the CUDA platform.

getUseDispersionCorrection(self)

Get whether to add a contribution to the energy that approximately represents the effect of Lennard-Jones interactions beyond the cutoff distance.

getUseSwitchingFunction(self)

Get whether a switching function is applied to the Lennard-Jones interaction.

setCutoffDistance(self, distance)

Set the cutoff distance (in nm) being used for nonbonded interactions.

setEwaldErrorTolerance(self, tol)

Set the error tolerance for Ewald summation.

setExceptionParameterOffset(self, index, ...)

Set the offset added to the parameters of a particular exception, based on a global parameter.

setExceptionParameters(self, index, ...)

Set the force field parameters for an interaction that should be calculated differently from others.

setExceptionsUsePeriodicBoundaryConditions(...)

Set whether periodic boundary conditions should be applied to exceptions.

setForceGroup(self, group)

Set the force group this Force belongs to.

setGlobalParameterDefaultValue(self, index, ...)

Set the default value of a global parameter.

setGlobalParameterName(self, index, name)

Set the name of a global parameter.

setIncludeDirectSpace(self, include)

Set whether to include direct space interactions when calculating forces and energies.

setLJPMEParameters(self, alpha, nx, ny, nz)

Set the parameters to use for the dispersion term in LJPME calculations.

setName(self, name)

Set the name of this Force.

setNonbondedMethod(self, method)

Set the method used for handling long range nonbonded interactions.

setPMEParameters(self, alpha, nx, ny, nz)

Set the parameters to use for PME calculations.

setParticleParameterOffset(self, index, ...)

Set the offset added to the per-particle parameters of a particular particle, based on a global parameter.

setParticleParameters(self, index, charge, ...)

Set the nonbonded force parameters for a particle.

setParticleSubset(index, subset)

Set the subset of a particle.

setReactionFieldDielectric(self, dielectric)

Set the dielectric constant to use for the solvent in the reaction field approximation.

setReciprocalSpaceForceGroup(self, group)

Set the force group that reciprocal space interactions for Ewald or PME are included in.

setScalingParameter(index, parameter, ...)

Modify an added scaling parameter.

setScalingParameterDerivative(index, parameter)

Set the name of the global parameter to associate with a requested scaling parameter derivative.

setSwitchingDistance(self, distance)

Set the distance at which the switching function begins to reduce the Lennard-Jones interaction.

setUseCuFFT(use)

Set whether whether to use CUDA Toolkit's cuFFT library when executing in the CUDA platform.

setUseDispersionCorrection(self, useCorrection)

Set whether to add a contribution to the energy that approximately represents the effect of Lennard-Jones interactions beyond the cutoff distance.

setUseSwitchingFunction(self, use)

Set whether a switching function is applied to the Lennard-Jones interaction.

updateParametersInContext(context)

Update the particle and exception parameters in a Context to match those stored in this Force object.

usesPeriodicBoundaryConditions(self)

Returns whether or not this force makes use of periodic boundary conditions.

Attributes

CutoffNonPeriodic

CutoffPeriodic

Ewald

LJPME

NoCutoff

PME

thisown

The membership flag

addException(self, particle1, particle2, chargeProd, sigma, epsilon, replace=False) int

Add an interaction to the list of exceptions that should be calculated differently from other interactions. If chargeProd and epsilon are both equal to 0, this will cause the interaction to be completely omitted from force and energy calculations.

Regardless of the NonbondedMethod used by this Force, cutoffs are never applied to exceptions. That is because they are primarily used for 1-4 interactions, which are really a type of bonded interaction and are parametrized together with the other bonded interactions.

In many cases, you can use createExceptionsFromBonds() rather than adding each exception explicitly.

Parameters:
  • particle1 (int) – the index of the first particle involved in the interaction

  • particle2 (int) – the index of the second particle involved in the interaction

  • chargeProd (double) – the scaled product of the atomic charges (i.e. the strength of the Coulomb interaction), measured in units of the proton charge squared

  • sigma (double) – the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm

  • epsilon (double) – the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol

  • replace (bool) – determines the behavior if there is already an exception for the same two particles. If true, the existing one is replaced. If false, an exception is thrown.

Returns:

int – the index of the exception that was added

addExceptionParameterOffset(self, parameter, exceptionIndex, chargeProdScale, sigmaScale, epsilonScale) int

Add an offset to the parameters of a particular exception, based on a global parameter.

Parameters:
  • parameter (string) – the name of the global parameter. It must have already been added with addGlobalParameter(). Its value can be modified at any time by calling Context::setParameter().

  • exceptionIndex (int) – the index of the exception whose parameters are affected

  • chargeProdScale (double) – this value multiplied by the parameter value is added to the exception’s charge product

  • sigmaScale (double) – this value multiplied by the parameter value is added to the exception’s sigma

  • epsilonScale (double) – this value multiplied by the parameter value is added to the exception’s epsilon

Returns:

int – the index of the offset that was added

addException_usingRMin(particle1, particle2, chargeProd, rMin, epsilon)

Add interaction exception using the product of the two atoms’ elementary charges, rMin and epsilon, which is standard for AMBER force fields. Note that rMin is the minimum energy point in the Lennard Jones potential. The conversion from sigma is: rMin = 2^1/6 * sigma.

addGlobalParameter(self, name, defaultValue) int

Add a new global parameter that parameter offsets may depend on. The default value provided to this method is the initial value of the parameter in newly created Contexts. You can change the value at any time by calling setParameter() on the Context.

Parameters:
  • name (string) – the name of the parameter

  • defaultValue (double) – the default value of the parameter

Returns:

int – the index of the parameter that was added

addParticle(self, charge, sigma, epsilon) int

Add the nonbonded force parameters for a particle. This should be called once for each particle in the System. When it is called for the i’th time, it specifies the parameters for the i’th particle. For calculating the Lennard-Jones interaction between two particles, the arithmetic mean of the sigmas and the geometric mean of the epsilons for the two interacting particles is used (the Lorentz-Berthelot combining rule).

Parameters:
  • charge (double) – the charge of the particle, measured in units of the proton charge

  • sigma (double) – the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm

  • epsilon (double) – the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol

Returns:

int – the index of the particle that was added

addParticleParameterOffset(self, parameter, particleIndex, chargeScale, sigmaScale, epsilonScale) int

Add an offset to the per-particle parameters of a particular particle, based on a global parameter.

Parameters:
  • parameter (string) – the name of the global parameter. It must have already been added with addGlobalParameter(). Its value can be modified at any time by calling Context::setParameter().

  • particleIndex (int) – the index of the particle whose parameters are affected

  • chargeScale (double) – this value multiplied by the parameter value is added to the particle’s charge

  • sigmaScale (double) – this value multiplied by the parameter value is added to the particle’s sigma

  • epsilonScale (double) – this value multiplied by the parameter value is added to the particle’s epsilon

Returns:

int – the index of the offset that was added

addParticle_usingRVdw(charge, rVDW, epsilon)

Add particle using elemetrary charge. Rvdw and epsilon, which is consistent with AMBER parameter file usage. Note that the sum of the radii of the two interacting atoms is the minimum energy point in the Lennard Jones potential and is often called rMin. The conversion from sigma follows: rVDW = 2^1/6 * sigma/2

addScalingParameter(parameter, subset1, subset2, includeCoulomb, includeLJ)

Add a scaling parameter to multiply a particular Coulomb slice. Its value will scale the Coulomb interactions between particles of a subset 1 with those of another (or the same) subset 2. The order of subset definition is irrelevant.

Parameters:
  • parameter (str) – the name of the global parameter. It must have already been added with addGlobalParameter(). Its value can be modified at any time by calling setParameter() on the openmm.Context

  • subset1 (int) – the index of a particle subset. Legal values are between 0 and the result of getNumSubsets()

  • subset2 (int) – the index of a particle subset. Legal values are between 0 and the result of getNumSubsets()

  • includeCoulomb (bool) – whether this scaling parameter applies to Coulomb interactions

  • includeLJ (bool) – whether this scaling parameter applies to Lennard-Jones interactions

Returns:

index (int) – the index of scaling parameter that was added

addScalingParameterDerivative(parameter)

Request the derivative of this Force’s energy with respect to a scaling parameter. This can be used to obtain the sum of particular energy slices. The parameter must have already been added with addGlobalParameter() and addScalingParameter().

Parameters:

parameter (str) – the name of the parameter

Returns:

index (int) – the index of scaling parameter derivative that was added

createExceptionsFromBonds(self, bonds, coulomb14Scale, lj14Scale)

Identify exceptions based on the molecular topology. Particles which are separated by one or two bonds are set to not interact at all, while pairs of particles separated by three bonds (known as “1-4 interactions”) have their Coulomb and Lennard-Jones interactions reduced by a fixed factor.

Parameters:
  • bonds (vector< std::pair< int, int > >) – the set of bonds based on which to construct exceptions. Each element specifies the indices of two particles that are bonded to each other.

  • coulomb14Scale (double) – pairs of particles separated by three bonds will have the strength of their Coulomb interaction multiplied by this factor

  • lj14Scale (double) – pairs of particles separated by three bonds will have the strength of their Lennard-Jones interaction multiplied by this factor

getCutoffDistance(self) double

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:

double – the cutoff distance, measured in nm

getEwaldErrorTolerance(self) double

Get the error tolerance for Ewald summation. This corresponds to the fractional error in the forces which is acceptable. This value is used to select the reciprocal space cutoff and separation parameter so that the average error level will be less than the tolerance. There is not a rigorous guarantee that all forces on all atoms will be less than the tolerance, however.

For PME calculations, if setPMEParameters() is used to set alpha to something other than 0, this value is ignored.

getExceptionParameterOffset(self, index)

Get the offset added to the parameters of a particular exception, based on a global parameter.

Parameters:
  • index (int) – the index of the offset to query, as returned by addExceptionParameterOffset()

  • parameter (string) – the name of the global parameter

  • exceptionIndex (int) – the index of the exception whose parameters are affected

  • chargeProdScale (double) – this value multiplied by the parameter value is added to the exception’s charge product

  • sigmaScale (double) – this value multiplied by the parameter value is added to the exception’s sigma

  • epsilonScale (double) – this value multiplied by the parameter value is added to the exception’s epsilon

getExceptionParameters(self, index)

Get the force field parameters for an interaction that should be calculated differently from others.

Parameters:

index (int) – the index of the interaction for which to get parameters

Returns:

  • particle1 (int) – the index of the first particle involved in the interaction

  • particle2 (int) – the index of the second particle involved in the interaction

  • chargeProd (double) – the scaled product of the atomic charges (i.e. the strength of the Coulomb interaction), measured in units of the proton charge squared

  • sigma (double) – the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm

  • epsilon (double) – the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol

getExceptionsUsePeriodicBoundaryConditions(self) bool

Get whether periodic boundary conditions should be applied to exceptions. Usually this is not appropriate, because exceptions are normally used to represent bonded interactions (1-2, 1-3, and 1-4 pairs), but there are situations when it does make sense. For example, you may want to simulate an infinite chain where one end of a molecule is bonded to the opposite end of the next periodic copy.

Regardless of this value, periodic boundary conditions are only applied to exceptions if they also are applied to other interactions. If the nonbonded method is NoCutoff or CutoffNonPeriodic, this value is ignored. Also note that cutoffs are never applied to exceptions, again because they are normally used to represent bonded interactions.

getForceGroup(self) int

Get the force group this Force belongs to.

getGlobalParameterDefaultValue(self, index) double

Get the default value of a global parameter.

Parameters:

index (int) – the index of the parameter for which to get the default value

Returns:

double – the parameter default value

getGlobalParameterName(self, index) std::string const &

Get the name of a global parameter.

Parameters:

index (int) – the index of the parameter for which to get the name

Returns:

string – the parameter name

getIncludeDirectSpace(self) bool

Get whether to include direct space interactions when calculating forces and energies. This is useful if you want to completely replace the direct space calculation, typically with a CustomNonbondedForce that computes it in a nonstandard way, while still using this object for the reciprocal space calculation.

getLJPMEParameters(self)

Get the parameters to use for dispersion term in LJ-PME calculations. If alpha is 0 (the default), these parameters are ignored and instead their values are chosen based on the Ewald error tolerance.

Returns:

  • alpha (double) – the separation parameter

  • nx (int) – the number of dispersion grid points along the X axis

  • ny (int) – the number of dispersion grid points along the Y axis

  • nz (int) – the number of dispersion grid points along the Z axis

getLJPMEParametersInContext(context)

Get the PME parameters being used for the dispersion term for LJPME in a particular Context. Because some platforms have restrictions on the allowed grid sizes, the values that are actually used may be slightly different from those specified with setPMEParameters(), or the standard values calculated based on the Ewald error tolerance. See the manual for details.

Parameters:

context (Context) – the Context for which to get the parameters

Returns:

  • alpha (double) – the separation parameter, measured in \(nm^{-1}\)

  • nx (int) – the number of grid points along the X axis

  • ny (int) – the number of grid points along the Y axis

  • nz (int) – the number of grid points along the Z axis

getName(self) std::string const &

Get the name of this Force. This is an arbitrary, user modifiable identifier. By default it equals the class name, but you can change it to anything useful.

getNonbondedMethod(self) OpenMM::NonbondedForce::NonbondedMethod

Get the method used for handling long range nonbonded interactions.

getNonbondedMethodName()

Get the name of the method used for handling long range nonbonded interactions.

getNumExceptionParameterOffsets(self) int

Get the number of exception parameter offsets that have been added.

getNumExceptions(self) int

Get the number of special interactions that should be calculated differently from other interactions.

getNumGlobalParameters(self) int

Get the number of global parameters that have been added.

getNumParticleParameterOffsets(self) int

Get the number of particles parameter offsets that have been added.

getNumParticles(self) int

Get the number of particles for which force field parameters have been defined.

getNumScalingParameterDerivatives()

Get the number of requested scaling parameter derivatives.

getNumScalingParameters()

Get the number of scaling parameters.

getNumSlices()

Get the number of slices determined by the specified number of particle subsets.

getNumSubsets()

Get the specified number of particle subsets.

getPMEParameters(self)

Get the parameters to use for PME calculations. If alpha is 0 (the default), these parameters are ignored and instead their values are chosen based on the Ewald error tolerance.

Returns:

  • alpha (double) – the separation parameter

  • nx (int) – the number of grid points along the X axis

  • ny (int) – the number of grid points along the Y axis

  • nz (int) – the number of grid points along the Z axis

getPMEParametersInContext(context)

Get the parameters being used for PME in a particular Context. Because some platforms have restrictions on the allowed grid sizes, the values that are actually used may be slightly different from those specified with setPMEParameters(), or the standard values calculated based on the Ewald error tolerance. See the manual for details.

Parameters:

context (Context) – the Context for which to get the parameters

Returns:

  • alpha (double) – the separation parameter, measured in \(nm^{-1}\)

  • nx (int) – the number of grid points along the X axis

  • ny (int) – the number of grid points along the Y axis

  • nz (int) – the number of grid points along the Z axis

getParticleParameterOffset(self, index)

Get the offset added to the per-particle parameters of a particular particle, based on a global parameter.

Parameters:
  • index (int) – the index of the offset to query, as returned by addParticleParameterOffset()

  • parameter (string) – the name of the global parameter

  • particleIndex (int) – the index of the particle whose parameters are affected

  • chargeScale (double) – this value multiplied by the parameter value is added to the particle’s charge

  • sigmaScale (double) – this value multiplied by the parameter value is added to the particle’s sigma

  • epsilonScale (double) – this value multiplied by the parameter value is added to the particle’s epsilon

getParticleParameters(self, index)

Get the nonbonded force parameters for a particle.

Parameters:

index (int) – the index of the particle for which to get parameters

Returns:

  • charge (double) – the charge of the particle, measured in units of the proton charge

  • sigma (double) – the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm

  • epsilon (double) – the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol

getParticleSubset(index)

Get the subset to which a particle belongs.

Parameters:

index (int) – the index of the particle for which to get the subset

getReactionFieldDielectric(self) double

Get the dielectric constant to use for the solvent in the reaction field approximation.

getReciprocalSpaceForceGroup(self) int

Get the force group that reciprocal space interactions for Ewald or PME are included in. This allows multiple time step integrators to evaluate direct and reciprocal space interactions at different intervals: getForceGroup() specifies the group for direct space, and getReciprocalSpaceForceGroup() specifies the group for reciprocal space. If this is -1 (the default value), the same force group is used for reciprocal space as for direct space.

getScalingParameter(index, parameter, arg4)

Get the scaling parameter applied to a particular nonbonded slice.

Parameters:

index (int) – the index of the scaling parameter to query, as returned by addScalingParameter()

Returns:

  • parameter (str) – the name of the global parameter

  • subset1 (int) – the smallest index of the two particle subsets

  • subset2 (int) – the largest index of the two particle subsets

  • includeCoulomb (bool) – whether this scaling parameter applies to Coulomb interactions

  • includeLJ (bool) – whether this scaling parameter applies to Lennard-Jones interactions

getScalingParameterDerivativeName(index)

Get the name of the global parameter associated with a requested scaling parameter derivative.

Parameters:

index (int) – the index of the parameter derivative, between 0 and the result of getNumScalingParameterDerivatives()

getSwitchingDistance(self) double

Get the distance at which the switching function begins to reduce the Lennard-Jones interaction. This must be less than the cutoff distance.

getUseCudaFFT()

Get whether to use CUDA Toolkit’s cuFFT library when executing in the CUDA platform. The default value is False.

getUseDispersionCorrection(self) bool

Get whether to add a contribution to the energy that approximately represents the effect of Lennard-Jones 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.

getUseSwitchingFunction(self) bool

Get whether a switching function is applied to the Lennard-Jones interaction. If the nonbonded method is set to NoCutoff, this option is ignored.

setCutoffDistance(self, 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 (double) – the cutoff distance, measured in nm

setEwaldErrorTolerance(self, tol)

Set the error tolerance for Ewald summation. This corresponds to the fractional error in the forces which is acceptable. This value is used to select the reciprocal space cutoff and separation parameter so that the average error level will be less than the tolerance. There is not a rigorous guarantee that all forces on all atoms will be less than the tolerance, however.

For PME calculations, if setPMEParameters() is used to set alpha to something other than 0, this value is ignored.

setExceptionParameterOffset(self, index, parameter, exceptionIndex, chargeProdScale, sigmaScale, epsilonScale)

Set the offset added to the parameters of a particular exception, based on a global parameter.

Parameters:
  • index (int) – the index of the offset to modify, as returned by addExceptionParameterOffset()

  • parameter (string) – the name of the global parameter. It must have already been added with addGlobalParameter(). Its value can be modified at any time by calling Context::setParameter().

  • exceptionIndex (int) – the index of the exception whose parameters are affected

  • chargeProdScale (double) – this value multiplied by the parameter value is added to the exception’s charge product

  • sigmaScale (double) – this value multiplied by the parameter value is added to the exception’s sigma

  • epsilonScale (double) – this value multiplied by the parameter value is added to the exception’s epsilon

setExceptionParameters(self, index, particle1, particle2, chargeProd, sigma, epsilon)

Set the force field parameters for an interaction that should be calculated differently from others. If chargeProd and epsilon are both equal to 0, this will cause the interaction to be completely omitted from force and energy calculations.

Regardless of the NonbondedMethod used by this Force, cutoffs are never applied to exceptions. That is because they are primarily used for 1-4 interactions, which are really a type of bonded interaction and are parametrized together with the other bonded interactions.

Parameters:
  • index (int) – the index of the interaction for which to get parameters

  • particle1 (int) – the index of the first particle involved in the interaction

  • particle2 (int) – the index of the second particle involved in the interaction

  • chargeProd (double) – the scaled product of the atomic charges (i.e. the strength of the Coulomb interaction), measured in units of the proton charge squared

  • sigma (double) – the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm

  • epsilon (double) – the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol

setExceptionsUsePeriodicBoundaryConditions(self, periodic)

Set whether periodic boundary conditions should be applied to exceptions. Usually this is not appropriate, because exceptions are normally used to represent bonded interactions (1-2, 1-3, and 1-4 pairs), but there are situations when it does make sense. For example, you may want to simulate an infinite chain where one end of a molecule is bonded to the opposite end of the next periodic copy.

Regardless of this value, periodic boundary conditions are only applied to exceptions if they also get applied to other interactions. If the nonbonded method is NoCutoff or CutoffNonPeriodic, this value is ignored. Also note that cutoffs are never applied to exceptions, again because they are normally used to represent bonded interactions.

setForceGroup(self, group)

Set the force group this Force belongs to.

Parameters:

group (int) – the group index. Legal values are between 0 and 31 (inclusive).

setGlobalParameterDefaultValue(self, index, defaultValue)

Set the default value of a global parameter.

Parameters:
  • index (int) – the index of the parameter for which to set the default value

  • defaultValue (double) – the default value of the parameter

setGlobalParameterName(self, index, name)

Set the name of a global parameter.

Parameters:
  • index (int) – the index of the parameter for which to set the name

  • name (string) – the name of the parameter

setIncludeDirectSpace(self, include)

Set whether to include direct space interactions when calculating forces and energies. This is useful if you want to completely replace the direct space calculation, typically with a CustomNonbondedForce that computes it in a nonstandard way, while still using this object for the reciprocal space calculation.

setLJPMEParameters(self, alpha, nx, ny, nz)

Set the parameters to use for the dispersion term in LJPME calculations. If alpha is 0 (the default), these parameters are ignored and instead their values are chosen based on the Ewald error tolerance.

Parameters:
  • alpha (double) – the separation parameter

  • nx (int) – the number of grid points along the X axis

  • ny (int) – the number of grid points along the Y axis

  • nz (int) – the number of grid points along the Z axis

setName(self, name)

Set the name of this Force. This is an arbitrary, user modifiable identifier. By default it equals the class name, but you can change it to anything useful.

setNonbondedMethod(self, method)

Set the method used for handling long range nonbonded interactions.

setPMEParameters(self, alpha, nx, ny, nz)

Set the parameters to use for PME calculations. If alpha is 0 (the default), these parameters are ignored and instead their values are chosen based on the Ewald error tolerance.

Parameters:
  • alpha (double) – the separation parameter

  • nx (int) – the number of grid points along the X axis

  • ny (int) – the number of grid points along the Y axis

  • nz (int) – the number of grid points along the Z axis

setParticleParameterOffset(self, index, parameter, particleIndex, chargeScale, sigmaScale, epsilonScale)

Set the offset added to the per-particle parameters of a particular particle, based on a global parameter.

Parameters:
  • index (int) – the index of the offset to modify, as returned by addParticleParameterOffset()

  • parameter (string) – the name of the global parameter. It must have already been added with addGlobalParameter(). Its value can be modified at any time by calling Context::setParameter().

  • particleIndex (int) – the index of the particle whose parameters are affected

  • chargeScale (double) – this value multiplied by the parameter value is added to the particle’s charge

  • sigmaScale (double) – this value multiplied by the parameter value is added to the particle’s sigma

  • epsilonScale (double) – this value multiplied by the parameter value is added to the particle’s epsilon

setParticleParameters(self, index, charge, sigma, epsilon)

Set the nonbonded force parameters for a particle. When calculating the Lennard-Jones interaction between two particles, it uses the arithmetic mean of the sigmas and the geometric mean of the epsilons for the two interacting particles (the Lorentz-Berthelot combining rule).

Parameters:
  • index (int) – the index of the particle for which to set parameters

  • charge (double) – the charge of the particle, measured in units of the proton charge

  • sigma (double) – the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm

  • epsilon (double) – the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol

setParticleSubset(index, subset)

Set the subset of a particle.

Parameters:
  • index (int) – the index of the particle for which to set the subset

  • subset (int) – the subset to which this particle belongs

setReactionFieldDielectric(self, dielectric)

Set the dielectric constant to use for the solvent in the reaction field approximation.

setReciprocalSpaceForceGroup(self, group)

Set the force group that reciprocal space interactions for Ewald or PME are included in. This allows multiple time step integrators to evaluate direct and reciprocal space interactions at different intervals: setForceGroup() specifies the group for direct space, and setReciprocalSpaceForceGroup() specifies the group for reciprocal space. If this is -1 (the default value), the same force group is used for reciprocal space as for direct space.

Parameters:

group (int) – the group index. Legal values are between 0 and 31 (inclusive), or -1 to use the same force group that is specified for direct space.

setScalingParameter(index, parameter, subset1, subset2, includeCoulomb, includeLJ)

Modify an added scaling parameter.

Parameters:
  • index (int) – the index of the scaling parameter to modify, as returned by addExceptionChargeOffset()

  • parameter (str) – the name of the global parameter. It must have already been added with addGlobalParameter(). Its value can be modified at any time by calling setParameter() on the openmm.Context

  • subset1 (int) – the index of a particle subset. Legal values are between 0 and the result of getNumSubsets()

  • subset2 (int) – the index of a particle subset. Legal values are between 0 and the result of getNumSubsets()

  • includeCoulomb (bool) – whether this scaling parameter applies to Coulomb interactions

  • includeLJ (bool) – whether this scaling parameter applies to Lennard-Jones interactions

setScalingParameterDerivative(index, parameter)

Set the name of the global parameter to associate with a requested scaling parameter derivative.

Parameters:
  • index (int) – the index of the parameter derivative, between 0 and getNumScalingParameterDerivatives`

  • parameter (str) – the name of the parameter

setSwitchingDistance(self, distance)

Set the distance at which the switching function begins to reduce the Lennard-Jones interaction. This must be less than the cutoff distance.

setUseCuFFT(use)

Set whether whether to use CUDA Toolkit’s cuFFT library when executing in the CUDA platform. This choice has no effect when using other platforms or when the CUDA Toolkit is version 7.0 or older.

Parameters:

use (bool) – whether to use the cuFFT library

setUseDispersionCorrection(self, useCorrection)

Set whether to add a contribution to the energy that approximately represents the effect of Lennard-Jones 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.

setUseSwitchingFunction(self, use)

Set whether a switching function is applied to the Lennard-Jones interaction. If the nonbonded method is set to NoCutoff, this option is ignored.

updateParametersInContext(context)

Update the particle and exception 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() and setExceptionParameters() to modify this object’s parameters, then call updateParametersInContext() to copy them over to the Context.

This method has several limitations. The only information it updates is the parameters of particles and exceptions. All other aspects of the Force (the nonbonded method, the cutoff distance, etc.) are unaffected and can only be changed by reinitializing the Context. Furthermore, only the chargeProd, sigma, and epsilon values of an exception can be changed; the pair of particles involved in the exception cannot change. Finally, this method cannot be used to add new particles or exceptions, only to change the parameters of existing ones.

Parameters:

context (Context) – the Context in which to update the parameters

usesPeriodicBoundaryConditions(self) bool

Returns whether or not this force makes use of periodic boundary conditions.

Returns:

bool – true if force uses PBC and false otherwise