HarmonicAngleForce

class simtk.openmm.openmm.HarmonicAngleForce(*args)

This class implements an interaction between groups of three particles that varies harmonically with the angle between them. To use it, create a HarmonicAngleForce object then call addAngle() once for each angle. After an angle has been added, you can modify its force field parameters by calling setAngleParameters(). This will have no effect on Contexts that already exist unless you call updateParametersInContext().

__init__(self) → HarmonicAngleForce

__init__(self, other) -> HarmonicAngleForce

Create a HarmonicAngleForce.

Methods

__init__((self) -> HarmonicAngleForce) __init__(self, other) -> HarmonicAngleForce
addAngle((self, particle1, particle2, ...) Add an angle term to the force field.
getAngleParameters(self, index) Get the force field parameters for an angle term.
getForceGroup((self) -> int) Get the force group this Force belongs to.
getNumAngles((self) -> int) Get the number of harmonic bond angle terms in the potential function
setAngleParameters(self, index, particle1, ...) Set the force field parameters for an angle term.
setForceGroup(self, group) Set the force group this Force belongs to.
setUsesPeriodicBoundaryConditions(self, periodic) Set whether this force should apply periodic boundary conditions when calculating displacements.
updateParametersInContext(self, context) Update the per-angle parameters in a Context to match those stored in this Force object.
usesPeriodicBoundaryConditions((self) -> bool) Returns whether or not this force makes use of periodic boundary conditions.
getNumAngles(self) → int

Get the number of harmonic bond angle terms in the potential function

addAngle(self, particle1, particle2, particle3, angle, k) → int

Add an angle term to the force field.

Parameters:
  • particle1 (int) – the index of the first particle forming the angle
  • particle2 (int) – the index of the second particle forming the angle
  • particle3 (int) – the index of the third particle forming the angle
  • angle (double) – the equilibrium angle, measured in radians
  • k (double) – the harmonic force constant for the angle, measured in kJ/mol/radian^2
Returns:

the index of the angle that was added

Return type:

int

getAngleParameters(self, index)

Get the force field parameters for an angle term.

Parameters:index (int) – the index of the angle for which to get parameters
Returns:
  • particle1 (int) – the index of the first particle forming the angle
  • particle2 (int) – the index of the second particle forming the angle
  • particle3 (int) – the index of the third particle forming the angle
  • angle (double) – the equilibrium angle, measured in radians
  • k (double) – the harmonic force constant for the angle, measured in kJ/mol/radian^2
setAngleParameters(self, index, particle1, particle2, particle3, angle, k)

Set the force field parameters for an angle term.

Parameters:
  • index (int) – the index of the angle for which to set parameters
  • particle1 (int) – the index of the first particle forming the angle
  • particle2 (int) – the index of the second particle forming the angle
  • particle3 (int) – the index of the third particle forming the angle
  • angle (double) – the equilibrium angle, measured in radians
  • k (double) – the harmonic force constant for the angle, measured in kJ/mol/radian^2
updateParametersInContext(self, context)

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

The only information this method updates is the values of per-angle parameters. The set of particles involved in a angle cannot be changed, nor can new angles be added.

setUsesPeriodicBoundaryConditions(self, periodic)

Set whether this force should apply periodic boundary conditions when calculating displacements. Usually this is not appropriate for bonded forces, but there are situations when it can be useful.

usesPeriodicBoundaryConditions(self) → bool

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

Returns:true if force uses PBC and false otherwise
Return type:bool
__copy__(self) → Force
getForceGroup(self) → int

Get the force group this Force belongs to.

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).