OpenMM
 All Classes Namespaces Functions Variables Pages
CustomAngleForce Class Reference

This class implements interactions between sets of three particles that depend on the angle between them. More...

+ Inheritance diagram for CustomAngleForce:

Public Member Functions

def getNumAngles
 getNumAngles(CustomAngleForce self) -> int More...
 
def getNumPerAngleParameters
 getNumPerAngleParameters(CustomAngleForce self) -> int More...
 
def getNumGlobalParameters
 getNumGlobalParameters(CustomAngleForce self) -> int More...
 
def getEnergyFunction
 getEnergyFunction(CustomAngleForce self) -> std::string const & More...
 
def setEnergyFunction
 setEnergyFunction(CustomAngleForce self, std::string const & energy) More...
 
def addPerAngleParameter
 addPerAngleParameter(CustomAngleForce self, std::string const & name) -> int More...
 
def getPerAngleParameterName
 getPerAngleParameterName(CustomAngleForce self, int index) -> std::string const & More...
 
def setPerAngleParameterName
 setPerAngleParameterName(CustomAngleForce self, int index, std::string const & name) More...
 
def addGlobalParameter
 addGlobalParameter(CustomAngleForce self, std::string const & name, double defaultValue) -> int More...
 
def getGlobalParameterName
 getGlobalParameterName(CustomAngleForce self, int index) -> std::string const & More...
 
def setGlobalParameterName
 setGlobalParameterName(CustomAngleForce self, int index, std::string const & name) More...
 
def getGlobalParameterDefaultValue
 getGlobalParameterDefaultValue(CustomAngleForce self, int index) -> double More...
 
def setGlobalParameterDefaultValue
 setGlobalParameterDefaultValue(CustomAngleForce self, int index, double defaultValue) More...
 
def addAngle
 addAngle(CustomAngleForce self, int particle1, int particle2, int particle3, vectord parameters) -> int More...
 
def getAngleParameters
 getAngleParameters(CustomAngleForce self, int index) More...
 
def setAngleParameters
 setAngleParameters(CustomAngleForce self, int index, int particle1, int particle2, int particle3, vectord parameters) More...
 
def updateParametersInContext
 updateParametersInContext(CustomAngleForce self, Context context) More...
 
def __init__
 init(OpenMM::CustomAngleForce self, std::string const & energy) -> CustomAngleForce init(OpenMM::CustomAngleForce self, CustomAngleForce other) -> CustomAngleForce More...
 
def __del__
 del(OpenMM::CustomAngleForce self) More...
 
- Public Member Functions inherited from Force
def __init__
 
def __del__
 del(OpenMM::Force self) More...
 
def getForceGroup
 getForceGroup(Force self) -> int More...
 
def setForceGroup
 setForceGroup(Force self, int group) More...
 
def __copy__
 
def __deepcopy__
 

Public Attributes

 this
 

Detailed Description

This class implements interactions between sets of three particles that depend on the angle between them.

Unlike HarmonicAngleForce, the functional form of the interaction is completely customizable, and may involve arbitrary algebraic expressions. In addition to the angle formed by the particles, it may depend on arbitrary global and per-angle parameters.

To use this class, create a CustomAngleForce object, passing an algebraic expression to the constructor that defines the interaction energy between each set of particles. The expression may depend on theta, the angle formed by the particles, as well as on any parameters you choose. Then call addPerAngleParameter() to define per-angle parameters, and addGlobalParameter() to define global parameters. The values of per-angle parameters are specified as part of the system definition, while values of global parameters may be modified during a simulation by calling Context::setParameter(). Finally, call addAngle() once for each angle. After an angle has been added, you can modify its parameters by calling setAngleParameters(). This will have no effect on Contexts that already exist unless you call updateParametersInContext().

As an example, the following code creates a CustomAngleForce that implements a harmonic potential:

CustomAngleForce* force = new CustomAngleForce("0.5*k*(theta-theta0)^2");

This force depends on two parameters: the spring constant k and equilibrium angle theta0. The following code defines these parameters:

force->addPerAngleParameter("k");
force->addPerAngleParameter("theta0");

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, sinh, cosh, tanh, erf, erfc, min, max, abs, step, delta. 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.

Constructor & Destructor Documentation

def __init__ (   self,
  args 
)

init(OpenMM::CustomAngleForce self, std::string const & energy) -> CustomAngleForce init(OpenMM::CustomAngleForce self, CustomAngleForce other) -> CustomAngleForce

Create a CustomAngleForce.

Parameters
energyan algebraic expression giving the interaction energy between three particles as a function of theta, the angle between them

References simtk.openmm.openmm.stripUnits().

def __del__ (   self)

del(OpenMM::CustomAngleForce self)

References simtk.openmm.openmm.stripUnits().

Member Function Documentation

def addAngle (   self,
  args 
)

addAngle(CustomAngleForce self, int particle1, int particle2, int particle3, vectord parameters) -> int

Add an angle term to the force field.

Parameters
particle1the index of the first particle connected by the angle
particle2the index of the second particle connected by the angle
particle3the index of the third particle connected by the angle
parametersthe list of parameters for the new angle

References simtk.openmm.openmm.stripUnits().

def addGlobalParameter (   self,
  args 
)

addGlobalParameter(CustomAngleForce self, std::string const & name, double defaultValue) -> int

Add a new global parameter that the interaction may depend on.

Parameters
namethe name of the parameter
defaultValuethe default value of the parameter

References simtk.openmm.openmm.stripUnits().

def addPerAngleParameter (   self,
  args 
)

addPerAngleParameter(CustomAngleForce self, std::string const & name) -> int

Add a new per-angle parameter that the interaction may depend on.

Parameters
namethe name of the parameter

References simtk.openmm.openmm.stripUnits().

def getAngleParameters (   self,
  args 
)

getAngleParameters(CustomAngleForce self, int index)

Get the force field parameters for an angle term.

Parameters
indexthe index of the angle for which to get parameters
particle1the index of the first particle connected by the angle
particle2the index of the second particle connected by the angle
particle3the index of the third particle connected by the angle
parametersthe list of parameters for the angle

References simtk.openmm.openmm.stripUnits().

def getEnergyFunction (   self)

getEnergyFunction(CustomAngleForce self) -> std::string const &

Get the algebraic expression that gives the interaction energy for each angle

References simtk.openmm.openmm.stripUnits().

def getGlobalParameterDefaultValue (   self,
  args 
)

getGlobalParameterDefaultValue(CustomAngleForce self, int index) -> double

Get the default value of a global parameter.

Parameters
indexthe index of the parameter for which to get the default value

References simtk.openmm.openmm.stripUnits().

def getGlobalParameterName (   self,
  args 
)

getGlobalParameterName(CustomAngleForce self, int index) -> std::string const &

Get the name of a global parameter.

Parameters
indexthe index of the parameter for which to get the name

References simtk.openmm.openmm.stripUnits().

def getNumAngles (   self)

getNumAngles(CustomAngleForce self) -> int

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

References simtk.openmm.openmm.stripUnits().

def getNumGlobalParameters (   self)

getNumGlobalParameters(CustomAngleForce self) -> int

Get the number of global parameters that the interaction depends on.

References simtk.openmm.openmm.stripUnits().

def getNumPerAngleParameters (   self)

getNumPerAngleParameters(CustomAngleForce self) -> int

Get the number of per-angle parameters that the interaction depends on.

References simtk.openmm.openmm.stripUnits().

def getPerAngleParameterName (   self,
  args 
)

getPerAngleParameterName(CustomAngleForce self, int index) -> std::string const &

Get the name of a per-angle parameter.

Parameters
indexthe index of the parameter for which to get the name

References simtk.openmm.openmm.stripUnits().

def setAngleParameters (   self,
  args 
)

setAngleParameters(CustomAngleForce self, int index, int particle1, int particle2, int particle3, vectord parameters)

Set the force field parameters for an angle term.

Parameters
indexthe index of the angle for which to set parameters
particle1the index of the first particle connected by the angle
particle2the index of the second particle connected by the angle
particle3the index of the third particle connected by the angle
parametersthe list of parameters for the angle

References simtk.openmm.openmm.stripUnits().

def setEnergyFunction (   self,
  args 
)

setEnergyFunction(CustomAngleForce self, std::string const & energy)

Set the algebraic expression that gives the interaction energy for each angle

References simtk.openmm.openmm.stripUnits().

def setGlobalParameterDefaultValue (   self,
  args 
)

setGlobalParameterDefaultValue(CustomAngleForce self, int index, double defaultValue)

Set the default value of a global parameter.

Parameters
indexthe index of the parameter for which to set the default value
namethe default value of the parameter

References simtk.openmm.openmm.stripUnits().

def setGlobalParameterName (   self,
  args 
)

setGlobalParameterName(CustomAngleForce self, int index, std::string const & name)

Set the name of a global parameter.

Parameters
indexthe index of the parameter for which to set the name
namethe name of the parameter

References simtk.openmm.openmm.stripUnits().

def setPerAngleParameterName (   self,
  args 
)

setPerAngleParameterName(CustomAngleForce self, int index, std::string const & name)

Set the name of a per-angle parameter.

Parameters
indexthe index of the parameter for which to set the name
namethe name of the parameter

References simtk.openmm.openmm.stripUnits().

def updateParametersInContext (   self,
  args 
)

updateParametersInContext(CustomAngleForce self, Context 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 updateParametersInState() to copy them over to the Context.

This method has several limitations. The only information it updates is the values of per-angle parameters. All other aspects of the Force (such as the energy function) are unaffected and can only be changed by reinitializing the Context. The set of particles involved in a angle cannot be changed, nor can new angles be added.

References simtk.openmm.openmm.stripUnits().

Member Data Documentation

this

The documentation for this class was generated from the following file: