DualAMDIntegrator

class simtk.openmm.amd.DualAMDIntegrator(dt, group, alphaTotal, ETotal, alphaGroup, EGroup)

DualAMDIntegrator implements a dual boost aMD integration algorithm.

This is similar to AMDIntegrator, but two different boosts are applied to the potential: one based on the total energy, and one based on the energy of a single force group (typically representing torsions).

For details, see Hamelberg et al., J. Chem. Phys. 127, 155102 (2007).

__init__(dt, group, alphaTotal, ETotal, alphaGroup, EGroup)

Create a DualAMDIntegrator.

Parameters:
  • dt (time) – The integration time step to use
  • group (int) – The force group to apply the second boost to
  • alphaTotal (energy) – The alpha parameter to use for the total energy
  • ETotal (energy) – The energy cutoff to use for the total energy
  • alphaGroup (energy) – The alpha parameter to use for the boosted force group
  • EGroup (energy) – The energy cutoff to use for the boosted force group

Methods

__init__(dt, group, alphaTotal, ETotal, ...) Create a DualAMDIntegrator.
addComputeGlobal((self, variable, ...) Add a step to the integration algorithm that computes a global value.
addComputePerDof((self, variable, ...) Add a step to the integration algorithm that computes a per-DOF value.
addComputeSum((self, variable, ...) Add a step to the integration algorithm that computes a sum over degrees of freedom.
addConstrainPositions((self) -> int) Add a step to the integration algorithm that updates particle positions so all constraints are satisfied.
addConstrainVelocities((self) -> int) Add a step to the integration algorithm that updates particle velocities so the net velocity along all constraints is 0.
addGlobalVariable((self, name, ...) Define a new global variable.
addPerDofVariable((self, name, ...) Define a new per-DOF variable.
addUpdateContextState((self) -> int) Add a step to the integration algorithm that allows Forces to update the context state.
beginIfBlock((self, condition) -> int) Add a step which begins a new “if” block.
beginWhileBlock((self, condition) -> int) Add a step which begins a new “while” block.
endBlock((self) -> int) Add a step which marks the end of the most recently begun “if” or “while” block.
getAlphaGroup() Get the value of alpha for the boosted force group.
getAlphaTotal() Get the value of alpha for the total energy.
getComputationStep(self, index) Get the details of a computation step that has been added to the integration algorithm.
getConstraintTolerance((self) -> double) Get the distance tolerance within which constraints are maintained, as a fraction of the constrained distance.
getEGroup() Get the energy threshold E for the boosted force group.
getETotal() Get the energy threshold E for the total energy.
getEffectiveEnergy(totalEnergy, groupEnergy) Given the actual potential energy of the system, return the value of the effective potential.
getGlobalVariable((self, index) -> double) Get the current value of a global variable.
getGlobalVariableByName((self, name) -> double) Get the current value of a global variable, specified by name.
getGlobalVariableName((self, ...) Get the name of a global variable.
getKineticEnergyExpression(...) Get the expression to use for computing the 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.
getPerDofVariable(self, index) getPerDofVariable(self, index) -> PyObject *
getPerDofVariableByName(self, name) Get the value of a per-DOF variable, specified by name.
getPerDofVariableName((self, ...) Get the name of a per-DOF variable.
getRandomNumberSeed((self) -> int) Get the random number seed.
getStepSize((self) -> double) Get the size of each time step, in picoseconds.
setAlphaGroup(alpha) Set the value of alpha for the boosted force group.
setAlphaTotal(alpha) Set the value of alpha for the total energy.
setConstraintTolerance(self, tol) Set the distance tolerance within which constraints are maintained, as a fraction of the constrained distance.
setEGroup(E) Set the energy threshold E for the boosted force group.
setETotal(E) Set the energy threshold E for the total energy.
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.
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
getAlphaTotal()

Get the value of alpha for the total energy.

setAlphaTotal(alpha)

Set the value of alpha for the total energy.

getETotal()

Get the energy threshold E for the total energy.

setETotal(E)

Set the energy threshold E for the total energy.

getAlphaGroup()

Get the value of alpha for the boosted force group.

setAlphaGroup(alpha)

Set the value of alpha for the boosted force group.

getEGroup()

Get the energy threshold E for the boosted force group.

setEGroup(E)

Set the energy threshold E for the boosted force group.

getEffectiveEnergy(totalEnergy, groupEnergy)

Given the actual potential energy of the system, return the value of the effective potential.

Parameters:
  • totalEnergy (energy) – the actual potential energy of the whole system
  • groupEnergy (energy) – the actual potential energy of the boosted force group
Returns:

the value of the effective potential

Return type:

energy

__delattr__

x.__delattr__(‘name’) <==> del x.name

__format__()

default object formatter

__getattribute__

x.__getattribute__(‘name’) <==> x.name

__hash__
__reduce__()

helper for pickle

__reduce_ex__()

helper for pickle

__sizeof__() → int

size of object in memory, in bytes

__str__
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

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

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
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
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