OpenMM
 All Classes Namespaces Functions Variables Pages
CustomIntegrator Class Reference

This is an Integrator that can be used to implemented arbitrary, user defined integration algorithms. More...

+ Inheritance diagram for CustomIntegrator:

Public Member Functions

def getNumGlobalVariables
 getNumGlobalVariables(CustomIntegrator self) -> int More...
 
def getNumPerDofVariables
 getNumPerDofVariables(CustomIntegrator self) -> int More...
 
def getNumComputations
 getNumComputations(CustomIntegrator self) -> int More...
 
def addGlobalVariable
 addGlobalVariable(CustomIntegrator self, std::string const & name, double initialValue) -> int More...
 
def getGlobalVariableName
 getGlobalVariableName(CustomIntegrator self, int index) -> std::string const & More...
 
def addPerDofVariable
 addPerDofVariable(CustomIntegrator self, std::string const & name, double initialValue) -> int More...
 
def getPerDofVariableName
 getPerDofVariableName(CustomIntegrator self, int index) -> std::string const & More...
 
def getGlobalVariable
 getGlobalVariable(CustomIntegrator self, int index) -> double More...
 
def getGlobalVariableByName
 getGlobalVariableByName(CustomIntegrator self, std::string const & name) -> double More...
 
def setGlobalVariable
 setGlobalVariable(CustomIntegrator self, int index, double value) More...
 
def setGlobalVariableByName
 setGlobalVariableByName(CustomIntegrator self, std::string const & name, double value) More...
 
def getPerDofVariableByName
 getPerDofVariableByName(CustomIntegrator self, std::string const & name) More...
 
def setPerDofVariable
 setPerDofVariable(CustomIntegrator self, int index, std::vector< Vec3,std::allocator< Vec3 > > const & values) More...
 
def setPerDofVariableByName
 setPerDofVariableByName(CustomIntegrator self, std::string const & name, std::vector< Vec3,std::allocator< Vec3 > > const & values) More...
 
def addComputeGlobal
 addComputeGlobal(CustomIntegrator self, std::string const & variable, std::string const & expression) -> int More...
 
def addComputePerDof
 addComputePerDof(CustomIntegrator self, std::string const & variable, std::string const & expression) -> int More...
 
def addComputeSum
 addComputeSum(CustomIntegrator self, std::string const & variable, std::string const & expression) -> int More...
 
def addConstrainPositions
 addConstrainPositions(CustomIntegrator self) -> int More...
 
def addConstrainVelocities
 addConstrainVelocities(CustomIntegrator self) -> int More...
 
def addUpdateContextState
 addUpdateContextState(CustomIntegrator self) -> int More...
 
def getComputationStep
 getComputationStep(CustomIntegrator self, int index) More...
 
def getKineticEnergyExpression
 getKineticEnergyExpression(CustomIntegrator self) -> std::string const & More...
 
def setKineticEnergyExpression
 setKineticEnergyExpression(CustomIntegrator self, std::string const & expression) More...
 
def getRandomNumberSeed
 getRandomNumberSeed(CustomIntegrator self) -> int More...
 
def setRandomNumberSeed
 setRandomNumberSeed(CustomIntegrator self, int seed) More...
 
def step
 step(CustomIntegrator self, int steps) More...
 
def getPerDofVariable
 getPerDofVariable(CustomIntegrator self, int index) getPerDofVariable(CustomIntegrator self, int index) -> PyObject * More...
 
def __init__
 init(OpenMM::CustomIntegrator self, double stepSize) -> CustomIntegrator init(OpenMM::CustomIntegrator self, CustomIntegrator other) -> CustomIntegrator More...
 
def __del__
 del(OpenMM::CustomIntegrator self) More...
 
- Public Member Functions inherited from Integrator
def __init__
 
def __del__
 del(OpenMM::Integrator self) More...
 
def getStepSize
 getStepSize(Integrator self) -> double More...
 
def setStepSize
 setStepSize(Integrator self, double size) More...
 
def getConstraintTolerance
 getConstraintTolerance(Integrator self) -> double More...
 
def setConstraintTolerance
 setConstraintTolerance(Integrator self, double tol) More...
 
def step
 step(Integrator self, int steps) More...
 
def __getstate__
 
def __setstate__
 

Public Attributes

 this
 
- Public Attributes inherited from Integrator
 this
 

Static Public Attributes

 ComputeGlobal = _openmm.CustomIntegrator_ComputeGlobal
 
 ComputePerDof = _openmm.CustomIntegrator_ComputePerDof
 
 ComputeSum = _openmm.CustomIntegrator_ComputeSum
 
 ConstrainPositions = _openmm.CustomIntegrator_ConstrainPositions
 
 ConstrainVelocities = _openmm.CustomIntegrator_ConstrainVelocities
 
 UpdateContextState = _openmm.CustomIntegrator_UpdateContextState
 

Detailed Description

This is an Integrator that can be used to implemented arbitrary, user defined integration algorithms.

It is flexible enough to support a wide range of methods including both deterministic and stochastic integrators, Metropolized integrators, and integrators that must integrate additional quantities along with the particle positions and momenta.

To create an integration algorithm, you first define a set of variables the integrator will compute. Variables come in two types: global variables have a single value, while per-DOF variables have a value for every degree of freedom (x, y, or z coordinate of a particle). You can define as many variables as you want of each type. The value of any variable can be computed by the integration algorithm, or set directly by calling a method on the CustomIntegrator. All variables are persistent between integration steps; once a value is set, it keeps that value until it is changed by the user or recomputed in a later integration step.

Next, you define the algorithm as a series of computations. To execute a time step, the integrator performs the list of computations in order. Each computation updates the value of one global or per-DOF value. There are several types of computations that can be done:

  • Global: You provide a mathematical expression involving only global variables. It is evaluated and stored into a global variable.

  • Per-DOF: You provide a mathematical expression involving both global and per-DOF variables. It is evaluated once for every degree of freedom, and the values are stored into a per-DOF variable.

  • Sum: You provide a mathematical expression involving both global and per-DOF variables. It is evaluated once for every degree of freedom. All of those values are then added together, and the sum is stored into a global variable.

  • Constrain Positions: The particle positions are updated so that all distance constraints are satisfied.

  • Constrain Velocities: The particle velocities are updated so the net velocity along any constrained distance is 0.

Like all integrators, CustomIntegrator ignores any particle whose mass is 0. It is skipped when doing per-DOF computations, and is not included when computing sums over degrees of freedom.

In addition to the variables you define by calling addGlobalVariable() and addPerDofVariable(), the integrator provides the following pre-defined variables:

  • dt: (global) This is the step size being used by the integrator.

  • energy: (global, read-only) This is the current potential energy of the system.

  • energy0, energy1, energy2, ...: (global, read-only) This is similar to energy, but includes only the contribution from forces in one force group. A single computation step may only depend on a single energy variable (energy, energy0, energy1, etc.).

  • x: (per-DOF) This is the current value of the degree of freedom (the x, y, or z coordinate of a particle).

  • v: (per-DOF) This is the current velocity associated with the degree of freedom (the x, y, or z component of a particle's velocity).

  • f: (per-DOF, read-only) This is the current force acting on the degree of freedom (the x, y, or z component of the force on a particle).

  • f0, f1, f2, ...: (per-DOF, read-only) This is similar to f, but includes only the contribution from forces in one force group. A single computation step may only depend on a single force variable (f, f0, f1, etc.).

  • m: (per-DOF, read-only) This is the mass of the particle the degree of freedom is associated with.

  • uniform: (either global or per-DOF, read-only) This is a uniformly distributed random number between 0 and 1. Every time an expression is evaluated, a different value will be used. When used in a per-DOF expression, a different value will be used for every degree of freedom. Note, however, that if this variable appears multiple times in a single expression, the same value is used everywhere it appears in that expression.

  • gaussian: (either global or per-DOF, read-only) This is a Gaussian distributed random number with mean 0 and variance 1. Every time an expression is evaluated, a different value will be used. When used in a per-DOF expression, a different value will be used for every degree of freedom. Note, however, that if this variable appears multiple times in a single expression, the same value is used everywhere it appears in that expression.

  • A global variable is created for every adjustable parameter defined in the integrator's Context.

The following example uses a CustomIntegrator to implement a velocity Verlet integrator:

CustomIntegrator integrator(0.001);
integrator.addComputePerDof("v", "v+0.5*dt*f/m");
integrator.addComputePerDof("x", "x+dt*v");
integrator.addComputePerDof("v", "v+0.5*dt*f/m");

The first step updates the velocities based on the current forces. The second step updates the positions based on the new velocities, and the third step updates the velocities again. Although the first and third steps look identical, the forces used in them are different. You do not need to tell the integrator that; it will recognize that the positions have changed and know to recompute the forces automatically.

The above example has two problems. First, it does not respect distance constraints. To make the integrator work with constraints, you need to add extra steps to tell it when and how to apply them. Second, it never gives Forces an opportunity to update the context state. This should be done every time step so that, for example, an AndersenThermostat can randomize velocities or a MonteCarloBarostat can scale particle positions. You need to add a step to tell the integrator when to do this. The following example corrects both these problems, using the RATTLE algorithm to apply constraints:

CustomIntegrator integrator(0.001);
integrator.addPerDofVariable("x1", 0);
integrator.addUpdateContextState();
integrator.addComputePerDof("v", "v+0.5*dt*f/m");
integrator.addComputePerDof("x", "x+dt*v");
integrator.addComputePerDof("x1", "x");
integrator.addConstrainPositions();
integrator.addComputePerDof("v", "v+0.5*dt*f/m+(x-x1)/dt");
integrator.addConstrainVelocities();

CustomIntegrator can be used to implement multiple time step integrators. The following example shows an r-RESPA integrator. It assumes the quickly changing forces are in force group 0 and the slowly changing ones are in force group 1. It evaluates the "fast" forces four times as often as the "slow" forces.

CustomIntegrator integrator(0.004);
integrator.addComputePerDof("v", "v+0.5*dt*f1/m");
for (int i = 0; i < 4; i++) {
integrator.addComputePerDof("v", "v+0.5*(dt/4)*f0/m");
integrator.addComputePerDof("x", "x+(dt/4)*v");
integrator.addComputePerDof("v", "v+0.5*(dt/4)*f0/m");
}
integrator.addComputePerDof("v", "v+0.5*dt*f1/m");

An Integrator has one other job in addition to evolving the equations of motion: it defines how to compute the kinetic energy of the system. Depending on the integration method used, simply summing mv/2 over all degrees of freedom may not give the correct answer. For example, in a leapfrog integrator the velocities are "delayed" by half a time step, so the above formula would give the kinetic energy half a time step ago, not at the current time.

Call setKineticEnergyExpression() to set an expression for the kinetic energy. It is computed for every degree of freedom (excluding ones whose mass is 0) and the result is summed. The default expression is "m*v*v/2", which is correct for many integrators.

As example, the following line defines the correct way to compute kinetic energy when using a leapfrog algorithm:

integrator.setKineticEnergyExpression("m*v1*v1/2; v1=v+0.5*dt*f/m");

The kinetic energy expression may depend on the following pre-defined variables: x, v, f, m, dt. It also may depend on user-defined global and per-DOF variables, and on the values of adjustable parameters defined in the integrator's Context. It may not depend on any other variable, such as the potential energy, the force from a single force group, or a random number.

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, floor, ceil, step, delta, select. 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. select(x,y,z) = z if x = 0, y otherwise. An expression may also involve intermediate quantities that are defined following the main expression, using ";" as a separator.

Constructor & Destructor Documentation

def __init__ (   self,
  args 
)

init(OpenMM::CustomIntegrator self, double stepSize) -> CustomIntegrator init(OpenMM::CustomIntegrator self, CustomIntegrator other) -> CustomIntegrator

Create a CustomIntegrator.

Parameters
stepSizethe step size with which to integrate the system (in picoseconds)

References simtk.openmm.openmm.stripUnits().

def __del__ (   self)

del(OpenMM::CustomIntegrator self)

References simtk.openmm.openmm.stripUnits().

Member Function Documentation

def addComputeGlobal (   self,
  args 
)

addComputeGlobal(CustomIntegrator self, std::string const & variable, std::string const & expression) -> int

Add a step to the integration algorithm that computes a global value.

Parameters
variablethe global variable to store the computed value into
expressiona mathematical expression involving only global variables. In each integration step, its value is computed and stored into the specified variable.

References simtk.openmm.openmm.stripUnits().

Referenced by AMDForceGroupIntegrator.__init__(), and DualAMDIntegrator.__init__().

def addComputePerDof (   self,
  args 
)

addComputePerDof(CustomIntegrator self, std::string const & variable, std::string const & expression) -> int

Add a step to the integration algorithm that computes a per-DOF value.

Parameters
variablethe per-DOF variable to store the computed value into
expressiona mathematical expression involving both global and per-DOF variables. In each integration step, its value is computed for every degree of freedom and stored into the specified variable.

References simtk.openmm.openmm.stripUnits().

Referenced by AMDIntegrator.__init__(), MTSIntegrator.__init__(), AMDForceGroupIntegrator.__init__(), and DualAMDIntegrator.__init__().

def addComputeSum (   self,
  args 
)

addComputeSum(CustomIntegrator self, std::string const & variable, std::string const & expression) -> int

Add a step to the integration algorithm that computes a sum over degrees of freedom.

Parameters
variablethe global variable to store the computed value into
expressiona mathematical expression involving both global and per-DOF variables. In each integration step, its value is computed for every degree of freedom. Those values are then added together, and the sum is stored in the specified variable.

References simtk.openmm.openmm.stripUnits().

def addConstrainPositions (   self,
  args 
)

addConstrainPositions(CustomIntegrator self) -> int

Add a step to the integration algorithm that updates particle positions so all constraints are satisfied.

References simtk.openmm.openmm.stripUnits().

Referenced by AMDIntegrator.__init__(), MTSIntegrator.__init__(), AMDForceGroupIntegrator.__init__(), and DualAMDIntegrator.__init__().

def addConstrainVelocities (   self,
  args 
)

addConstrainVelocities(CustomIntegrator self) -> int

Add a step to the integration algorithm that updates particle velocities so the net velocity along all constraints is 0.

References simtk.openmm.openmm.stripUnits().

Referenced by MTSIntegrator.__init__().

def addGlobalVariable (   self,
  args 
)

addGlobalVariable(CustomIntegrator self, std::string const & name, double initialValue) -> int

Define a new global variable.

Parameters
namethe name of the variable
initialValuethe variable will initially be set to this value

References simtk.openmm.openmm.stripUnits().

Referenced by AMDIntegrator.__init__(), AMDForceGroupIntegrator.__init__(), and DualAMDIntegrator.__init__().

def addPerDofVariable (   self,
  args 
)

addPerDofVariable(CustomIntegrator self, std::string const & name, double initialValue) -> int

Define a new per-DOF variable.

Parameters
namethe name of the variable
initialValuethe variable will initially be set to this value for all degrees of freedom

References simtk.openmm.openmm.stripUnits().

Referenced by AMDIntegrator.__init__(), MTSIntegrator.__init__(), AMDForceGroupIntegrator.__init__(), and DualAMDIntegrator.__init__().

def addUpdateContextState (   self,
  args 
)

addUpdateContextState(CustomIntegrator self) -> int

Add a step to the integration algorithm that allows Forces to update the context state.

References simtk.openmm.openmm.stripUnits().

Referenced by AMDIntegrator.__init__(), MTSIntegrator.__init__(), AMDForceGroupIntegrator.__init__(), and DualAMDIntegrator.__init__().

def getComputationStep (   self,
  args 
)

getComputationStep(CustomIntegrator self, int index)

Get the details of a computation step that has been added to the integration algorithm.

Parameters
indexthe index of the computation step to get
typeon exit, the type of computation this step performs
variableon exit, the variable into which this step stores its result. If this step does not store a result in a variable, this will be an empty string.
expressionon exit, the expression this step evaluates. If this step does not evaluate an expression, this will be an empty string.

References simtk.openmm.openmm.stripUnits().

def getGlobalVariable (   self,
  args 
)

getGlobalVariable(CustomIntegrator self, int index) -> double

Get the current value of a global variable.

Parameters
indexthe index of the variable to get

References simtk.openmm.openmm.stripUnits().

Referenced by AMDIntegrator.getAlpha(), AMDForceGroupIntegrator.getAlphaGroup(), DualAMDIntegrator.getAlphaGroup(), DualAMDIntegrator.getAlphaTotal(), AMDIntegrator.getE(), AMDForceGroupIntegrator.getEGroup(), DualAMDIntegrator.getEGroup(), and DualAMDIntegrator.getETotal().

def getGlobalVariableByName (   self,
  args 
)

getGlobalVariableByName(CustomIntegrator self, std::string const & name) -> double

Get the current value of a global variable, specified by name.

Parameters
namethe name of the variable to get

References simtk.openmm.openmm.stripUnits().

def getGlobalVariableName (   self,
  args 
)

getGlobalVariableName(CustomIntegrator self, int index) -> std::string const &

Get the name of a global variable.

Parameters
indexthe index of the variable to get

References simtk.openmm.openmm.stripUnits().

def getKineticEnergyExpression (   self,
  args 
)

getKineticEnergyExpression(CustomIntegrator self) -> std::string const &

Get the expression to use for computing the kinetic energy. The expression is evaluated for every degree of freedom. Those values are then added together, and the sum is reported as the current kinetic energy.

References simtk.openmm.openmm.stripUnits().

def getNumComputations (   self,
  args 
)

getNumComputations(CustomIntegrator self) -> int

Get the number of computation steps that have been added.

References simtk.openmm.openmm.stripUnits().

def getNumGlobalVariables (   self,
  args 
)

getNumGlobalVariables(CustomIntegrator self) -> int

Get the number of global variables that have been defined.

References simtk.openmm.openmm.stripUnits().

def getNumPerDofVariables (   self,
  args 
)

getNumPerDofVariables(CustomIntegrator self) -> int

Get the number of per-DOF variables that have been defined.

References simtk.openmm.openmm.stripUnits().

def getPerDofVariable (   self,
  args 
)

getPerDofVariable(CustomIntegrator self, int index) getPerDofVariable(CustomIntegrator self, int index) -> PyObject *

def getPerDofVariableByName (   self,
  args 
)

getPerDofVariableByName(CustomIntegrator self, std::string const & name)

Get the value of a per-DOF variable, specified by name.

Parameters
namethe name of the variable to get
valuesthe values of the variable for all degrees of freedom are stored into this

References simtk.openmm.openmm.stripUnits().

def getPerDofVariableName (   self,
  args 
)

getPerDofVariableName(CustomIntegrator self, int index) -> std::string const &

Get the name of a per-DOF variable.

Parameters
indexthe index of the variable to get

References simtk.openmm.openmm.stripUnits().

def getRandomNumberSeed (   self,
  args 
)

getRandomNumberSeed(CustomIntegrator self) -> int

Get the random number seed. See setRandomNumberSeed() for details.

References simtk.openmm.openmm.stripUnits().

def setGlobalVariable (   self,
  args 
)

setGlobalVariable(CustomIntegrator self, int index, double value)

Set the value of a global variable.

Parameters
indexthe index of the variable to set
valuethe new value of the variable

References simtk.openmm.openmm.stripUnits().

Referenced by AMDIntegrator.setAlpha(), AMDForceGroupIntegrator.setAlphaGroup(), DualAMDIntegrator.setAlphaGroup(), DualAMDIntegrator.setAlphaTotal(), AMDIntegrator.setE(), AMDForceGroupIntegrator.setEGroup(), DualAMDIntegrator.setEGroup(), and DualAMDIntegrator.setETotal().

def setGlobalVariableByName (   self,
  args 
)

setGlobalVariableByName(CustomIntegrator self, std::string const & name, double value)

Set the value of a global variable, specified by name.

Parameters
namethe name of the variable to set
valuethe new value of the variable

References simtk.openmm.openmm.stripUnits().

def setKineticEnergyExpression (   self,
  args 
)

setKineticEnergyExpression(CustomIntegrator self, std::string const & expression)

Set the expression to use for computing the kinetic energy. The expression is evaluated for every degree of freedom. Those values are then added together, and the sum is reported as the current kinetic energy.

References simtk.openmm.openmm.stripUnits().

def setPerDofVariable (   self,
  args 
)

setPerDofVariable(CustomIntegrator self, int index, std::vector< Vec3,std::allocator< Vec3 > > const & values)

Set the value of a per-DOF variable.

Parameters
indexthe index of the variable to set
valuesthe new values of the variable for all degrees of freedom

References simtk.openmm.openmm.stripUnits().

def setPerDofVariableByName (   self,
  args 
)

setPerDofVariableByName(CustomIntegrator self, std::string const & name, std::vector< Vec3,std::allocator< Vec3 > > const & values)

Set the value of a per-DOF variable, specified by name.

Parameters
namethe name of the variable to set
valuesthe new values of the variable for all degrees of freedom

References simtk.openmm.openmm.stripUnits().

def setRandomNumberSeed (   self,
  args 
)

setRandomNumberSeed(CustomIntegrator self, int seed)

Set the random number seed. The precise meaning of this parameter is undefined, and is left up to each Platform to interpret in an appropriate way. It is guaranteed that if two simulations are run with different random number seeds, the sequence of random numbers will be different. On the other hand, no guarantees are made about the behavior of simulations that use the same seed. In particular, Platforms are permitted to use non-deterministic algorithms which produce different results on successive runs, even if those runs were initialized identically.

If seed is set to 0 (which is the default value assigned), a unique seed is chosen when a Context is created from this Force. This is done to ensure that each Context receives unique random seeds without you needing to set them explicitly.

References simtk.openmm.openmm.stripUnits().

def step (   self,
  args 
)

step(CustomIntegrator self, int steps)

Advance a simulation through time by taking a series of time steps.

Parameters
stepsthe number of time steps to take

References simtk.openmm.openmm.stripUnits().

Member Data Documentation

ComputeGlobal = _openmm.CustomIntegrator_ComputeGlobal
static
ComputePerDof = _openmm.CustomIntegrator_ComputePerDof
static
ComputeSum = _openmm.CustomIntegrator_ComputeSum
static
ConstrainPositions = _openmm.CustomIntegrator_ConstrainPositions
static
ConstrainVelocities = _openmm.CustomIntegrator_ConstrainVelocities
static
this

Referenced by System.__init__().

UpdateContextState = _openmm.CustomIntegrator_UpdateContextState
static

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