OpenMM
|
This is an Integrator that can be used to implemented arbitrary, user defined integration algorithms. More...
#include <CustomIntegrator.h>
Public Types | |
enum | ComputationType { ComputeGlobal = 0, ComputePerDof = 1, ComputeSum = 2, ConstrainPositions = 3, ConstrainVelocities = 4, UpdateContextState = 5 } |
This is an enumeration of the different types of computations that may appear in an integration algorithm. More... | |
Public Member Functions | |
CustomIntegrator (double stepSize) | |
Create a CustomIntegrator. | |
int | getNumGlobalVariables () const |
Get the number of global variables that have been defined. | |
int | getNumPerDofVariables () const |
Get the number of per-DOF variables that have been defined. | |
int | getNumComputations () const |
Get the number of computation steps that have been added. | |
int | addGlobalVariable (const std::string &name, double initialValue) |
Define a new global variable. | |
const std::string & | getGlobalVariableName (int index) const |
Get the name of a global variable. | |
int | addPerDofVariable (const std::string &name, double initialValue) |
Define a new per-DOF variable. | |
const std::string & | getPerDofVariableName (int index) const |
Get the name of a per-DOF variable. | |
double | getGlobalVariable (int index) const |
Get the current value of a global variable. | |
void | setGlobalVariable (int index, double value) |
Set the value of a global variable. | |
void | setGlobalVariableByName (const std::string &name, double value) |
Set the value of a global variable, specified by name. | |
void | getPerDofVariable (int index, std::vector< Vec3 > &values) const |
Get the value of a per-DOF variable. | |
void | setPerDofVariable (int index, const std::vector< Vec3 > &values) |
Set the value of a per-DOF variable. | |
void | setPerDofVariableByName (const std::string &name, const std::vector< Vec3 > &values) |
Set the value of a per-DOF variable, specified by name. | |
int | addComputeGlobal (const std::string &variable, const std::string &expression) |
Add a step to the integration algorithm that computes a global value. | |
int | addComputePerDof (const std::string &variable, const std::string &expression) |
Add a step to the integration algorithm that computes a per-DOF value. | |
int | addComputeSum (const std::string &variable, const std::string &expression) |
Add a step to the integration algorithm that computes a sum over degrees of freedom. | |
int | addConstrainPositions () |
Add a step to the integration algorithm that updates particle positions so all constraints are satisfied. | |
int | addConstrainVelocities () |
Add a step to the integration algorithm that updates particle velocities so the net velocity along all constraints is 0. | |
int | addUpdateContextState () |
Add a step to the integration algorithm that allows Forces to update the context state. | |
void | getComputationStep (int index, ComputationType &type, std::string &variable, std::string &expression) const |
Get the details of a computation step that has been added to the integration algorithm. | |
const std::string & | getKineticEnergyExpression () const |
Get the expression to use for computing the kinetic energy. | |
void | setKineticEnergyExpression (const std::string &expression) |
Set the expression to use for computing the kinetic energy. | |
int | getRandomNumberSeed () const |
Get the random number seed. | |
void | setRandomNumberSeed (int seed) |
Set the random number seed. | |
void | step (int steps) |
Advance a simulation through time by taking a series of time steps. | |
Public Member Functions inherited from Integrator | |
Integrator () | |
virtual | ~Integrator () |
double | getStepSize () const |
Get the size of each time step, in picoseconds. | |
void | setStepSize (double size) |
Set the size of each time step, in picoseconds. | |
double | getConstraintTolerance () const |
Get the distance tolerance within which constraints are maintained, as a fraction of the constrained distance. | |
void | setConstraintTolerance (double tol) |
Set the distance tolerance within which constraints are maintained, as a fraction of the constrained distance. | |
Protected Member Functions | |
void | initialize (ContextImpl &context) |
This will be called by the Context when it is created. | |
void | cleanup () |
This will be called by the Context when it is destroyed to let the Integrator do any necessary cleanup. | |
void | stateChanged (State::DataType changed) |
When the user modifies the state, we need to mark that the forces need to be recalculated. | |
std::vector< std::string > | getKernelNames () |
Get the names of all Kernels used by this Integrator. | |
double | computeKineticEnergy () |
Compute the kinetic energy of the system at the current time. | |
Additional Inherited Members | |
Protected Attributes inherited from Integrator | |
ContextImpl * | context |
Context * | owner |
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");
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 mv2/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, 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. An expression may also involve intermediate quantities that are defined following the main expression, using ";" as a separator.
enum ComputationType |
This is an enumeration of the different types of computations that may appear in an integration algorithm.
CustomIntegrator | ( | double | stepSize | ) |
Create a CustomIntegrator.
stepSize | the step size with which to integrate the system (in picoseconds) |
int addComputeGlobal | ( | const std::string & | variable, |
const std::string & | expression | ||
) |
Add a step to the integration algorithm that computes a global value.
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. |
int addComputePerDof | ( | const std::string & | variable, |
const std::string & | expression | ||
) |
Add a step to the integration algorithm that computes a per-DOF value.
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. |
int addComputeSum | ( | const std::string & | variable, |
const std::string & | expression | ||
) |
Add a step to the integration algorithm that computes a sum over degrees of freedom.
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. |
int addConstrainPositions | ( | ) |
Add a step to the integration algorithm that updates particle positions so all constraints are satisfied.
int addConstrainVelocities | ( | ) |
Add a step to the integration algorithm that updates particle velocities so the net velocity along all constraints is 0.
int addGlobalVariable | ( | const std::string & | name, |
double | initialValue | ||
) |
Define a new global variable.
name | the name of the variable |
initialValue | the variable will initially be set to this value |
int addPerDofVariable | ( | const std::string & | name, |
double | initialValue | ||
) |
Define a new per-DOF variable.
name | the name of the variable |
initialValue | the variable will initially be set to this value for all degrees of freedom |
int addUpdateContextState | ( | ) |
Add a step to the integration algorithm that allows Forces to update the context state.
|
protectedvirtual |
This will be called by the Context when it is destroyed to let the Integrator do any necessary cleanup.
It will also get called again if the application calls reinitialize() on the Context.
Reimplemented from Integrator.
|
protectedvirtual |
Compute the kinetic energy of the system at the current time.
Implements Integrator.
void getComputationStep | ( | int | index, |
ComputationType & | type, | ||
std::string & | variable, | ||
std::string & | expression | ||
) | const |
Get the details of a computation step that has been added to the integration algorithm.
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. |
double getGlobalVariable | ( | int | index | ) | const |
Get the current value of a global variable.
index | the index of the variable to get |
const std::string& getGlobalVariableName | ( | int | index | ) | const |
Get the name of a global variable.
index | the index of the variable to get |
|
protectedvirtual |
Get the names of all Kernels used by this Integrator.
Implements Integrator.
const std::string& getKineticEnergyExpression | ( | ) | 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.
|
inline |
Get the number of computation steps that have been added.
|
inline |
Get the number of global variables that have been defined.
|
inline |
Get the number of per-DOF variables that have been defined.
void getPerDofVariable | ( | int | index, |
std::vector< Vec3 > & | values | ||
) | const |
Get the value of a per-DOF variable.
index | the index of the variable to get |
values | the values of the variable for all degrees of freedom are stored into this |
const std::string& getPerDofVariableName | ( | int | index | ) | const |
Get the name of a per-DOF variable.
index | the index of the variable to get |
|
inline |
Get the random number seed.
See setRandomNumberSeed() for details.
|
protectedvirtual |
This will be called by the Context when it is created.
It informs the Integrator of what context it will be integrating, and gives it a chance to do any necessary initialization. It will also get called again if the application calls reinitialize() on the Context.
Implements Integrator.
void setGlobalVariable | ( | int | index, |
double | value | ||
) |
Set the value of a global variable.
index | the index of the variable to set |
value | the new value of the variable |
void setGlobalVariableByName | ( | const std::string & | name, |
double | value | ||
) |
Set the value of a global variable, specified by name.
name | the name of the variable to set |
value | the new value of the variable |
void setKineticEnergyExpression | ( | const std::string & | 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.
void setPerDofVariable | ( | int | index, |
const std::vector< Vec3 > & | values | ||
) |
Set the value of a per-DOF variable.
index | the index of the variable to set |
values | the new values of the variable for all degrees of freedom |
void setPerDofVariableByName | ( | const std::string & | name, |
const std::vector< Vec3 > & | values | ||
) |
Set the value of a per-DOF variable, specified by name.
name | the name of the variable to set |
values | the new values of the variable for all degrees of freedom |
|
inline |
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.
|
protectedvirtual |
When the user modifies the state, we need to mark that the forces need to be recalculated.
Reimplemented from Integrator.
|
virtual |
Advance a simulation through time by taking a series of time steps.
steps | the number of time steps to take |
Implements Integrator.