CustomIntegrator Class Reference

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

Inheritance diagram for CustomIntegrator:
Integrator

List of all members.

Public Member Functions

def getNumGlobalVariables
 getNumGlobalVariables(self) -> int
def getNumPerDofVariables
 getNumPerDofVariables(self) -> int
def getNumComputations
 getNumComputations(self) -> int
def addGlobalVariable
 addGlobalVariable(self, string name, double initialValue) -> int
def getGlobalVariableName
 getGlobalVariableName(self, int index) -> string
def addPerDofVariable
 addPerDofVariable(self, string name, double initialValue) -> int
def getPerDofVariableName
 getPerDofVariableName(self, int index) -> string
def getGlobalVariable
 getGlobalVariable(self, int index) -> double
def setGlobalVariable
 setGlobalVariable(self, int index, double value)
def setGlobalVariableByName
 setGlobalVariableByName(self, string name, double value)
def setPerDofVariable
 setPerDofVariable(self, int index, std.vector<(Vec3,std.allocator<(Vec3)>)> values)
def setPerDofVariableByName
 setPerDofVariableByName(self, string name, std.vector<(Vec3,std.allocator<(Vec3)>)> values)
def addComputeGlobal
 addComputeGlobal(self, string variable, string expression) -> int
def addComputePerDof
 addComputePerDof(self, string variable, string expression) -> int
def addComputeSum
 addComputeSum(self, string variable, string expression) -> int
def addConstrainPositions
 addConstrainPositions(self) -> int
def addConstrainVelocities
 addConstrainVelocities(self) -> int
def addUpdateContextState
 addUpdateContextState(self) -> int
def getComputationStep
 getComputationStep(self, int index)
def getRandomNumberSeed
 getRandomNumberSeed(self) -> int
def setRandomNumberSeed
 setRandomNumberSeed(self, int seed)
def step
 step(self, int steps)
def getPerDofVariable
 getPerDofVariable(self, int index, std.vector<(Vec3,std.allocator<(Vec3)>)> values) getPerDofVariable(self, int index) -> PyObject
def __init__
 __init__(self, double stepSize) -> CustomIntegrator __init__(self, CustomIntegrator other) -> CustomIntegrator
def __del__
 __del__(self)

Public Attributes

 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:

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:

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

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. All trigonometric functions are defined in radians, and log is the natural logarithm. step(x) = 0 if x is less than 0, 1 otherwise. An expression may also involve intermediate quantities that are defined following the main expression, using ";" as a separator.


Member Function Documentation

def __del__ (   self  ) 

__del__(self)

Reimplemented from Integrator.

def __init__ (   self,
  args 
)

__init__(self, double stepSize) -> CustomIntegrator __init__(self, CustomIntegrator other) -> CustomIntegrator

Create a CustomIntegrator.

Parameters:
stepSize the step size with which to integrate the system (in picoseconds)
def addComputeGlobal (   self,
  args 
)

addComputeGlobal(self, string variable, string expression) -> int

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

Parameters:
variable the global variable to store the computed value into
expression a mathematical expression involving only global variables. In each integration step, its value is computed and stored into the specified variable.
def addComputePerDof (   self,
  args 
)

addComputePerDof(self, string variable, string expression) -> int

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

Parameters:
variable the per-DOF variable to store the computed value into
expression a 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.
def addComputeSum (   self,
  args 
)

addComputeSum(self, string variable, string expression) -> int

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

Parameters:
variable the global variable to store the computed value into
expression a 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.
def addConstrainPositions (   self  ) 

addConstrainPositions(self) -> int

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

def addConstrainVelocities (   self  ) 

addConstrainVelocities(self) -> int

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

def addGlobalVariable (   self,
  args 
)

addGlobalVariable(self, string name, double initialValue) -> int

Define a new global variable.

Parameters:
name the name of the variable
initialValue the variable will initially be set to this value
def addPerDofVariable (   self,
  args 
)

addPerDofVariable(self, string name, double initialValue) -> int

Define a new per-DOF variable.

Parameters:
name the name of the variable
initialValue the variable will initially be set to this value for all degrees of freedom
def addUpdateContextState (   self  ) 

addUpdateContextState(self) -> int

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

def getComputationStep (   self,
  args 
)

getComputationStep(self, int index)

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

Parameters:
index the index of the computation step to get
type on exit, the type of computation this step performs
variable on 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.
expression on exit, the expression this step evaluates. If this step does not evaluate an expression, this will be an empty string.
def getGlobalVariable (   self,
  args 
)

getGlobalVariable(self, int index) -> double

Get the current value of a global variable.

Parameters:
index the index of the variable to get
def getGlobalVariableName (   self,
  args 
)

getGlobalVariableName(self, int index) -> string

Get the name of a global variable.

Parameters:
index the index of the variable to get
def getNumComputations (   self  ) 

getNumComputations(self) -> int

Get the number of computation steps that have been added.

def getNumGlobalVariables (   self  ) 

getNumGlobalVariables(self) -> int

Get the number of global variables that have been defined.

def getNumPerDofVariables (   self  ) 

getNumPerDofVariables(self) -> int

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

def getPerDofVariable (   self,
  args 
)

getPerDofVariable(self, int index, std.vector<(Vec3,std.allocator<(Vec3)>)> values) getPerDofVariable(self, int index) -> PyObject

def getPerDofVariableName (   self,
  args 
)

getPerDofVariableName(self, int index) -> string

Get the name of a per-DOF variable.

Parameters:
index the index of the variable to get
def getRandomNumberSeed (   self  ) 

getRandomNumberSeed(self) -> int

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

def setGlobalVariable (   self,
  args 
)

setGlobalVariable(self, int index, double value)

Set the value of a global variable.

Parameters:
index the index of the variable to set
value the new value of the variable
def setGlobalVariableByName (   self,
  args 
)

setGlobalVariableByName(self, string name, double value)

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

Parameters:
name the name of the variable to set
value the new value of the variable
def setPerDofVariable (   self,
  args 
)

setPerDofVariable(self, int index, std.vector<(Vec3,std.allocator<(Vec3)>)> values)

Set the value of a per-DOF variable.

Parameters:
index the index of the variable to set
values the new values of the variable for all degrees of freedom
def setPerDofVariableByName (   self,
  args 
)

setPerDofVariableByName(self, string name, std.vector<(Vec3,std.allocator<(Vec3)>)> values)

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

Parameters:
name the name of the variable to set
values the new values of the variable for all degrees of freedom
def setRandomNumberSeed (   self,
  args 
)

setRandomNumberSeed(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.

def step (   self,
  args 
)

step(self, int steps)

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

Parameters:
steps the number of time steps to take

Reimplemented from Integrator.


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]
UpdateContextState = _openmm.CustomIntegrator_UpdateContextState [static]

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

Generated by  doxygen 1.6.2