OpenMM
|
This is an Integrator that can be used to implemented arbitrary, user defined integration algorithms. More...
Public Member Functions | |
def | getNumGlobalVariables |
getNumGlobalVariables(self) -> int | |
def | getNumPerDofVariables |
getNumPerDofVariables(self) -> int | |
def | getNumComputations |
getNumComputations(self) -> int | |
def | addGlobalVariable |
addGlobalVariable(self, name, initialValue) -> int | |
def | getGlobalVariableName |
getGlobalVariableName(self, index) -> std::string const & | |
def | addPerDofVariable |
addPerDofVariable(self, name, initialValue) -> int | |
def | getPerDofVariableName |
getPerDofVariableName(self, index) -> std::string const & | |
def | getGlobalVariable |
getGlobalVariable(self, index) -> double | |
def | getGlobalVariableByName |
getGlobalVariableByName(self, name) -> double | |
def | setGlobalVariable |
Set the value of a global variable. | |
def | setGlobalVariableByName |
Set the value of a global variable, specified by name. | |
def | getPerDofVariableByName |
Get the value of a per-DOF variable, specified by name. | |
def | setPerDofVariable |
Set the value of a per-DOF variable. | |
def | setPerDofVariableByName |
Set the value of a per-DOF variable, specified by name. | |
def | addComputeGlobal |
addComputeGlobal(self, variable, expression) -> int | |
def | addComputePerDof |
addComputePerDof(self, variable, expression) -> int | |
def | addComputeSum |
addComputeSum(self, variable, expression) -> int | |
def | addConstrainPositions |
addConstrainPositions(self) -> int | |
def | addConstrainVelocities |
addConstrainVelocities(self) -> int | |
def | addUpdateContextState |
addUpdateContextState(self) -> int | |
def | beginIfBlock |
beginIfBlock(self, condition) -> int | |
def | beginWhileBlock |
beginWhileBlock(self, condition) -> int | |
def | endBlock |
endBlock(self) -> int | |
def | getComputationStep |
Get the details of a computation step that has been added to the integration algorithm. | |
def | getKineticEnergyExpression |
getKineticEnergyExpression(self) -> std::string const & | |
def | setKineticEnergyExpression |
Set the expression to use for computing the kinetic energy. | |
def | getRandomNumberSeed |
getRandomNumberSeed(self) -> int | |
def | setRandomNumberSeed |
Set the random number seed. | |
def | step |
Advance a simulation through time by taking a series of time steps. | |
def | getPerDofVariable |
getPerDofVariable(self, index) getPerDofVariable(self, index) -> PyObject * | |
def | __init__ |
__init__(self, stepSize) -> CustomIntegrator __init__(self, other) -> CustomIntegrator | |
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 | |
BeginIfBlock = _openmm.CustomIntegrator_BeginIfBlock | |
BeginWhileBlock = _openmm.CustomIntegrator_BeginWhileBlock | |
EndBlock = _openmm.CustomIntegrator_EndBlock |
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.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");
The sequence of computations in a CustomIntegrator can include flow control in the form of "if" and "while" blocks. The computations inside an "if" block are executed either zero or one times, depending on whether a condition is true. The computations inside a "while" block are executed repeatedly for as long as the condition remains true. Be very careful when writing "while" blocks; there is nothing to stop you from creating an infinite loop!
For example, suppose you are writing a Monte Carlo algorithm. Assume you have already computed a new set of particle coordinates "xnew" and a step acceptance probability "acceptanceProbability". The following lines use an "if" block to decide whether to accept the step, and if it is accepted, store the new positions into "x".
integrator.beginIfBlock("uniform < acceptanceProbability"); integrator.computePerDof("x", "xnew"); integrator.endBlock();
The condition in an "if" or "while" block is evaluated globally, so it may only involve global variables, not per-DOF ones. It may use any of the following comparison operators: =, <. >, !=, <=, >=. Blocks may be nested inside each other.
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.
def __init__ | ( | self, | |
args | |||
) |
__init__(self, stepSize) -> CustomIntegrator __init__(self, other) -> CustomIntegrator
Create a CustomIntegrator.
stepSize | (double) the step size with which to integrate the system (in picoseconds) |
def addComputeGlobal | ( | self, | |
variable, | |||
expression | |||
) |
addComputeGlobal(self, variable, expression) -> int
Add a step to the integration algorithm that computes a global value.
variable | (string) the global variable to store the computed value into |
expression | (string) a mathematical expression involving only global variables. In each integration step, its value is computed and stored into the specified variable. |
def addComputePerDof | ( | self, | |
variable, | |||
expression | |||
) |
addComputePerDof(self, variable, expression) -> int
Add a step to the integration algorithm that computes a per-DOF value.
variable | (string) the per-DOF variable to store the computed value into |
expression | (string) 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, | |
variable, | |||
expression | |||
) |
addComputeSum(self, variable, expression) -> int
Add a step to the integration algorithm that computes a sum over degrees of freedom.
variable | (string) the global variable to store the computed value into |
expression | (string) 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, | |
name, | |||
initialValue | |||
) |
addGlobalVariable(self, name, initialValue) -> int
Define a new global variable.
name | (string) the name of the variable |
initialValue | (double) the variable will initially be set to this value |
def addPerDofVariable | ( | self, | |
name, | |||
initialValue | |||
) |
addPerDofVariable(self, name, initialValue) -> int
Define a new per-DOF variable.
name | (string) the name of the variable |
initialValue | (double) 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 beginIfBlock | ( | self, | |
condition | |||
) |
beginIfBlock(self, condition) -> int
Add a step which begins a new "if" block.
condition | (string) a mathematical expression involving a comparison operator and global variables. All steps between this one and the end of the block are executed only if the condition is true. |
def beginWhileBlock | ( | self, | |
condition | |||
) |
beginWhileBlock(self, condition) -> int
Add a step which begins a new "while" block.
condition | (string) a mathematical expression involving a comparison operator and global variables. All steps between this one and the end of the block are executed repeatedly as long as the condition remains true. |
def endBlock | ( | self | ) |
endBlock(self) -> int
Add a step which marks the end of the most recently begun "if" or "while" block.
def getComputationStep | ( | self, | |
index | |||
) |
Get the details of a computation step that has been added to the integration algorithm.
index | (int) the index of the computation step to get |
def getGlobalVariable | ( | self, | |
index | |||
) |
getGlobalVariable(self, index) -> double
Get the current value of a global variable.
index | (int) the index of the variable to get |
def getGlobalVariableByName | ( | self, | |
name | |||
) |
getGlobalVariableByName(self, name) -> double
Get the current value of a global variable, specified by name.
name | (string) the name of the variable to get |
def getGlobalVariableName | ( | self, | |
index | |||
) |
getGlobalVariableName(self, index) -> std::string const &
Get the name of a global variable.
index | (int) the index of the variable to get |
def getKineticEnergyExpression | ( | self | ) |
getKineticEnergyExpression(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.
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, index) getPerDofVariable(self, index) -> PyObject *
def getPerDofVariableByName | ( | self, | |
name | |||
) |
Get the value of a per-DOF variable, specified by name.
name | (string) the name of the variable to get |
def getPerDofVariableName | ( | self, | |
index | |||
) |
getPerDofVariableName(self, index) -> std::string const &
Get the name of a per-DOF variable.
index | (int) 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, | |
index, | |||
value | |||
) |
Set the value of a global variable.
index | (int) the index of the variable to set |
value | (double) the new value of the variable |
def setGlobalVariableByName | ( | self, | |
name, | |||
value | |||
) |
Set the value of a global variable, specified by name.
name | (string) the name of the variable to set |
value | (double) the new value of the variable |
def setKineticEnergyExpression | ( | self, | |
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.
def setPerDofVariable | ( | self, | |
index, | |||
values | |||
) |
Set the value of a per-DOF variable.
index | (int) the index of the variable to set |
values | (vector< Vec3 >) the new values of the variable for all degrees of freedom |
def setPerDofVariableByName | ( | self, | |
name, | |||
values | |||
) |
Set the value of a per-DOF variable, specified by name.
name | (string) the name of the variable to set |
values | (vector< Vec3 >) the new values of the variable for all degrees of freedom |
def setRandomNumberSeed | ( | self, | |
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.
def step | ( | self, | |
steps | |||
) |
Advance a simulation through time by taking a series of time steps.
steps | (int) the number of time steps to take |
Reimplemented from Integrator.
BeginIfBlock = _openmm.CustomIntegrator_BeginIfBlock [static] |
BeginWhileBlock = _openmm.CustomIntegrator_BeginWhileBlock [static] |
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] |
EndBlock = _openmm.CustomIntegrator_EndBlock [static] |
Reimplemented from Integrator.
UpdateContextState = _openmm.CustomIntegrator_UpdateContextState [static] |