MTSLangevinIntegrator

class openmm.mtsintegrator.MTSLangevinIntegrator(temperature, friction, dt, groups)

MTSLangevinIntegrator implements the BAOAB-RESPA multiple time step algorithm for constant temperature dynamics.

This integrator allows different forces to be evaluated at different frequencies, for example to evaluate the expensive, slowly changing forces less frequently than the inexpensive, quickly changing forces.

To use it, you must first divide your forces into two or more groups (by calling setForceGroup() on them) that should be evaluated at different frequencies. When you create the integrator, you provide a tuple for each group specifying the index of the force group and the frequency (as a fraction of the outermost time step) at which to evaluate it. For example:

integrator = MTSLangevinIntegrator(300*kelvin, 1/picosecond, 4*femtoseconds, [(0,1), (1,2), (2,8)])

This specifies that the outermost time step is 4 fs, so each step of the integrator will advance time by that much. It also says that force group 0 should be evaluated once per time step, force group 1 should be evaluated twice per time step (every 2 fs), and force group 2 should be evaluated eight times per time step (every 0.5 fs).

A common use of this algorithm is to evaluate reciprocal space nonbonded interactions less often than the bonded and direct space nonbonded interactions. The following example looks up the NonbondedForce, sets the reciprocal space interactions to their own force group, and then creates an integrator that evaluates them once every 4 fs, but all other interactions every 2 fs:

nonbonded = [f for f in system.getForces() if isinstance(f, NonbondedForce)][0]
nonbonded.setReciprocalSpaceForceGroup(1)
integrator = MTSLangevinIntegrator(300*kelvin, 1/picosecond, 4*femtoseconds, [(1,1), (0,2)])

For details, see Tuckerman et al., J. Chem. Phys. 97(3) pp. 1990-2001 (1992) and Lagardere et al., J. Phys. Chem. Lett. 10(10) pp. 2593-2599 (2019).

__init__(temperature, friction, dt, groups)

Create an MTSLangevinIntegrator.

Parameters
  • temperature (temperature) – the temperature of the heat bath

  • friction (1/temperature) – the friction coefficient which couples the system to the heat bath

  • dt (time) – The largest (outermost) integration time step to use

  • groups (list) – A list of tuples defining the force groups. The first element of each tuple is the force group index, and the second element is the number of times that force group should be evaluated in one time step.

Methods

__init__(temperature, friction, dt, groups)

Create an MTSLangevinIntegrator.

addComputeGlobal(self, variable, expression)

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

addComputePerDof(self, variable, expression)

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

addComputeSum(self, variable, expression)

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

addConstrainPositions(self)

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

addConstrainVelocities(self)

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

addGlobalVariable(self, name, initialValue)

Define a new global variable.

addPerDofVariable(self, name, initialValue)

Define a new per-DOF variable.

addTabulatedFunction(self, name, function)

Add a tabulated function that may appear in expressions.

addUpdateContextState(self)

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

beginIfBlock(self, condition)

Add a step which begins a new “if” block.

beginWhileBlock(self, condition)

Add a step which begins a new “while” block.

endBlock(self)

Add a step which marks the end of the most recently begun “if” or “while” block.

getComputationStep(self, index)

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

getConstraintTolerance(self)

Get the distance tolerance within which constraints are maintained, as a fraction of the constrained distance.

getGlobalVariable(self, index)

Get the current value of a global variable.

getGlobalVariableByName(self, name)

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

getGlobalVariableName(self, index)

Get the name of a global variable.

getIntegrationForceGroups(self)

Get which force groups to use for integration.

getKineticEnergyExpression(self)

Get the expression to use for computing the kinetic energy.

getNumComputations(self)

Get the number of computation steps that have been added.

getNumGlobalVariables(self)

Get the number of global variables that have been defined.

getNumPerDofVariables(self)

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

getNumTabulatedFunctions(self)

Get the number of tabulated functions that have been defined.

getPerDofVariable()

getPerDofVariableByName(self, name)

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

getPerDofVariableName(self, index)

Get the name of a per-DOF variable.

getRandomNumberSeed(self)

Get the random number seed.

getStepSize(self)

Get the size of each time step, in picoseconds.

getTabulatedFunction(-> TabulatedFunction)

Get a reference to a tabulated function that may appear in expressions.

getTabulatedFunctionName(self, index)

Get the name of a tabulated function that may appear in expressions.

setConstraintTolerance(self, tol)

Set the distance tolerance within which constraints are maintained, as a fraction of the constrained distance.

setGlobalVariable(self, index, value)

Set the value of a global variable.

setGlobalVariableByName(self, name, value)

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

setIntegrationForceGroups(groups)

Set which force groups to use for integration.

setKineticEnergyExpression(self, expression)

Set the expression to use for computing the kinetic energy.

setPerDofVariable(self, index, values)

Set the value of a per-DOF variable.

setPerDofVariableByName(self, name, values)

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

setRandomNumberSeed(self, seed)

Set the random number seed.

setStepSize(self, size)

Set the size of each time step, in picoseconds.

step(self, steps)

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

Attributes

BlockEnd

ComputeGlobal

ComputePerDof

ComputeSum

ConstrainPositions

ConstrainVelocities

IfBlockStart

UpdateContextState

WhileBlockStart

thisown

The membership flag

addComputeGlobal(self, variable, expression)int

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

Parameters
  • 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.

Returns

the index of the step that was added

Return type

int

addComputePerDof(self, variable, expression)int

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

Parameters
  • 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.

Returns

the index of the step that was added

Return type

int

addComputeSum(self, variable, expression)int

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

Parameters
  • 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.

Returns

the index of the step that was added

Return type

int

addConstrainPositions(self)int

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

Returns

the index of the step that was added

Return type

int

addConstrainVelocities(self)int

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

Returns

the index of the step that was added

Return type

int

addGlobalVariable(self, name, initialValue)int

Define a new global variable.

Parameters
  • name (string) – the name of the variable

  • initialValue (double) – the variable will initially be set to this value

Returns

the index of the variable that was added

Return type

int

addPerDofVariable(self, name, initialValue)int

Define a new per-DOF variable.

Parameters
  • name (string) – the name of the variable

  • initialValue (double) – the variable will initially be set to this value for all degrees of freedom

Returns

the index of the variable that was added

Return type

int

addTabulatedFunction(self, name, function)int

Add a tabulated function that may appear in expressions.

Parameters
  • name (string) – the name of the function as it appears in expressions

  • function (TabulatedFunction *) – a TabulatedFunction object defining the function. The TabulatedFunction should have been created on the heap with the “new” operator. The integrator takes over ownership of it, and deletes it when the integrator itself is deleted.

Returns

the index of the function that was added

Return type

int

addUpdateContextState(self)int

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

Returns

the index of the step that was added

Return type

int

beginIfBlock(self, condition)int

Add a step which begins a new “if” block.

Parameters

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.

Returns

the index of the step that was added

Return type

int

beginWhileBlock(self, condition)int

Add a step which begins a new “while” block.

Parameters

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.

Returns

the index of the step that was added

Return type

int

endBlock(self)int

Add a step which marks the end of the most recently begun “if” or “while” block.

Returns

the index of the step that was added

Return type

int

getComputationStep(self, index)

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

Parameters

index (int) – the index of the computation step to get

Returns

  • type (ComputationType) – the type of computation this step performs

  • variable (string) – 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 (string) – the expression this step evaluates. If this step does not evaluate an expression, this will be an empty string.

getConstraintTolerance(self)double

Get the distance tolerance within which constraints are maintained, as a fraction of the constrained distance.

getGlobalVariable(self, index)double

Get the current value of a global variable.

Parameters

index (int) – the index of the variable to get

Returns

the current value of the variable

Return type

double

getGlobalVariableByName(self, name)double

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

Parameters

name (string) – the name of the variable to get

Returns

the current value of the parameter

Return type

double

getGlobalVariableName(self, index)std::string const &

Get the name of a global variable.

Parameters

index (int) – the index of the variable to get

Returns

the name of the variable

Return type

string

getIntegrationForceGroups(self)int

Get which force groups to use for integration. By default, all force groups are included. This is interpreted as a set of bit flags: the forces from group i will be included if (groups&(1<<i)) != 0.

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.

getNumComputations(self)int

Get the number of computation steps that have been added.

getNumGlobalVariables(self)int

Get the number of global variables that have been defined.

getNumPerDofVariables(self)int

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

getNumTabulatedFunctions(self)int

Get the number of tabulated functions that have been defined.

getPerDofVariable(self, index)
getPerDofVariable(self, index)PyObject *
getPerDofVariableByName(self, name)

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

Parameters

name (string) – the name of the variable to get

Returns

values – the values of the variable for all degrees of freedom are stored into this

Return type

vector< Vec3 >

getPerDofVariableName(self, index)std::string const &

Get the name of a per-DOF variable.

Parameters

index (int) – the index of the variable to get

Returns

the name of the variable

Return type

string

getRandomNumberSeed(self)int

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

getStepSize(self)double

Get the size of each time step, in picoseconds. If this integrator uses variable time steps, the size of the most recent step is returned.

Returns

the step size, measured in ps

Return type

double

getTabulatedFunction(self, index)TabulatedFunction
getTabulatedFunction(self, index)TabulatedFunction

Get a reference to a tabulated function that may appear in expressions.

Parameters

index (int) – the index of the function to get

Returns

the TabulatedFunction object defining the function

Return type

TabulatedFunction

getTabulatedFunctionName(self, index)std::string const &

Get the name of a tabulated function that may appear in expressions.

Parameters

index (int) – the index of the function to get

Returns

the name of the function as it appears in expressions

Return type

string

setConstraintTolerance(self, tol)

Set the distance tolerance within which constraints are maintained, as a fraction of the constrained distance.

setGlobalVariable(self, index, value)

Set the value of a global variable.

Parameters
  • index (int) – the index of the variable to set

  • value (double) – the new value of the variable

setGlobalVariableByName(self, name, value)

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

Parameters
  • name (string) – the name of the variable to set

  • value (double) – the new value of the variable

setIntegrationForceGroups(groups)

Set which force groups to use for integration. By default, all force groups are included.

Parameters

groups (set or int) – a set of indices for which force groups to include when integrating the equations of motion. Alternatively, the groups can be passed as a single unsigned integer interpreted as a bitmask, in which case group i will be included if (groups&(1<<i)) != 0.

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.

setPerDofVariable(self, index, values)

Set the value of a per-DOF variable.

Parameters
  • index (int) – the index of the variable to set

  • values (vector< Vec3 >) – the new values of the variable for all degrees of freedom

setPerDofVariableByName(self, name, values)

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

Parameters
  • name (string) – the name of the variable to set

  • values (vector< Vec3 >) – the new values of the variable for all degrees of freedom

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.

setStepSize(self, size)

Set the size of each time step, in picoseconds. If this integrator uses variable time steps, the effect of calling this method is undefined, and it may simply be ignored.

Parameters

size (double) – the step size, measured in ps

step(self, steps)

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

Parameters

steps (int) – the number of time steps to take

property thisown

The membership flag