OpenMM
 All Classes Namespaces Functions Variables Pages
Context Class Reference

A Context stores the complete state of a simulation. More...

Inherits _object.

Public Member Functions

def __del__
 del(OpenMM::Context self) More...
 
def getSystem
 getSystem(Context self) -> System More...
 
def getIntegrator
 getIntegrator(Context self) -> Integrator getIntegrator(Context self) -> Integrator More...
 
def getPlatform
 getPlatform(Context self) -> Platform getPlatform(Context self) -> Platform More...
 
def setTime
 setTime(Context self, double time) More...
 
def setPositions
 setPositions(self, positions) More...
 
def setVelocities
 setVelocities(self, velocities) More...
 
def setVelocitiesToTemperature
 setVelocitiesToTemperature(Context self, double temperature, int randomSeed=time(NULL)) setVelocitiesToTemperature(Context self, double temperature) More...
 
def getParameter
 getParameter(Context self, std::string const & name) -> double More...
 
def setParameter
 setParameter(Context self, std::string const & name, double value) More...
 
def setPeriodicBoxVectors
 setPeriodicBoxVectors(Context self, Vec3 const & a, Vec3 const & b, Vec3 const & c) More...
 
def applyConstraints
 applyConstraints(Context self, double tol) More...
 
def applyVelocityConstraints
 applyVelocityConstraints(Context self, double tol) More...
 
def computeVirtualSites
 computeVirtualSites(Context self) More...
 
def reinitialize
 reinitialize(Context self) More...
 
def getState
 getState(self, getPositions = False, getVelocities = False, getForces = False, getEnergy = False, getParameters = False, enforcePeriodicBox = False, groups = -1) -> State More...
 
def setState
 setState(Context self, State state) More...
 
def createCheckpoint
 createCheckpoint(Context self) -> std::string More...
 
def loadCheckpoint
 loadCheckpoint(Context self, std::string checkpoint) More...
 
def __init__
 init(OpenMM::Context self, System system, Integrator integrator) -> Context init(OpenMM::Context self, System system, Integrator integrator, Platform platform) -> Context init(OpenMM::Context self, System system, Integrator integrator, Platform platform, mapstringstring properties) -> Context init(OpenMM::Context self, Context other) -> Context More...
 

Public Attributes

 this
 

Detailed Description

A Context stores the complete state of a simulation.

More specifically, it includes:

  • The current time

  • The position of each particle

  • The velocity of each particle

  • The values of configurable parameters defined by Force objects in the System

You can retrieve a snapshot of the current state at any time by calling getState(). This allows you to record the state of the simulation at various points, either for analysis or for checkpointing. getState() can also be used to retrieve the current forces on each particle and the current energy of the System.

Constructor & Destructor Documentation

def __del__ (   self)

del(OpenMM::Context self)

References simtk.openmm.openmm.stripUnits().

def __init__ (   self,
  args 
)

init(OpenMM::Context self, System system, Integrator integrator) -> Context init(OpenMM::Context self, System system, Integrator integrator, Platform platform) -> Context init(OpenMM::Context self, System system, Integrator integrator, Platform platform, mapstringstring properties) -> Context init(OpenMM::Context self, Context other) -> Context

Construct a new Context in which to run a simulation, explicitly specifying what Platform should be used to perform calculations and the values of platform-specific properties.

Parameters
systemthe System which will be simulated
integratorthe Integrator which will be used to simulate the System
platformthe Platform to use for calculations
propertiesa set of values for platform-specific properties. Keys are the property names.

References simtk.openmm.openmm.stripUnits().

Member Function Documentation

def applyConstraints (   self,
  args 
)

applyConstraints(Context self, double tol)

Update the positions of particles so that all distance constraints are satisfied. This also recomputes the locations of all virtual sites.

Parameters
tolthe distance tolerance within which constraints must be satisfied.

References simtk.openmm.openmm.stripUnits().

def applyVelocityConstraints (   self,
  args 
)

applyVelocityConstraints(Context self, double tol)

Update the velocities of particles so the net velocity of each constrained distance is zero.

Parameters
tolthe velocity tolerance within which constraints must be satisfied.

References simtk.openmm.openmm.stripUnits().

def computeVirtualSites (   self)

computeVirtualSites(Context self)

Recompute the locations of all virtual sites. There is rarely a reason to call this, since virtual sites are also updated by applyConstraints(). This is only for the rare situations when you want to enforce virtual sites but not constraints.

References simtk.openmm.openmm.stripUnits().

def createCheckpoint (   self)

createCheckpoint(Context self) -> std::string

Create a checkpoint recording the current state of the Context. This should be treated as an opaque block of binary data. See loadCheckpoint() for more details.

Returns
a string containing the checkpoint data
def getIntegrator (   self,
  args 
)

getIntegrator(Context self) -> Integrator getIntegrator(Context self) -> Integrator

Get Integrator being used to by this context.

References simtk.openmm.openmm.stripUnits().

def getParameter (   self,
  args 
)

getParameter(Context self, std::string const & name) -> double

Get the value of an adjustable parameter defined by a Force object in the System.

Parameters
namethe name of the parameter to get

References simtk.openmm.openmm.stripUnits().

def getPlatform (   self,
  args 
)

getPlatform(Context self) -> Platform getPlatform(Context self) -> Platform

Get the Platform being used for calculations.

References simtk.openmm.openmm.stripUnits().

def getState (   self,
  getPositions = False,
  getVelocities = False,
  getForces = False,
  getEnergy = False,
  getParameters = False,
  enforcePeriodicBox = False,
  groups = -1 
)

getState(self, getPositions = False, getVelocities = False, getForces = False, getEnergy = False, getParameters = False, enforcePeriodicBox = False, groups = -1) -> State

Get a State object recording the current state information stored in this context.

Parameters
getPositions(bool=False) whether to store particle positions in the State
getVelocities(bool=False) whether to store particle velocities in the State
getForces(bool=False) whether to store the forces acting on particles in the State
getEnergy(bool=False) whether to store potential and kinetic energy in the State
getParameter(bool=False) whether to store context parameters in the State
enforcePeriodicBox(bool=False) if false, the position of each particle will be whatever position is stored in the Context, regardless of periodic boundary conditions. If true, particle positions will be translated so the center of every molecule lies in the same periodic box.
groups(int=-1) a set of bit flags for which force groups to include when computing forces and energies. Group i will be included if (groups&(1<<i)) != 0. The default value includes all groups.

References Context._getStateAsLists().

Referenced by Context.reinitialize().

def getSystem (   self)

getSystem(Context self) -> System

Get System being simulated in this context.

References simtk.openmm.openmm.stripUnits().

def loadCheckpoint (   self,
  args 
)

loadCheckpoint(Context self, std::string checkpoint)

Load a checkpoint that was written by createCheckpoint().

A checkpoint contains not only publicly visible data such as the particle positions and velocities, but also internal data such as the states of random number generators. Ideally, loading a checkpoint should restore the Context to an identical state to when it was written, such that continuing the simulation will produce an identical trajectory. This is not strictly guaranteed to be true, however, and should not be relied on. For most purposes, however, the internal state should be close enough to be reasonably considered equivalent.

A checkpoint contains data that is highly specific to the Context from which it was created. It depends on the details of the System, the Platform being used, and the hardware and software of the computer it was created on. If you try to load it on a computer with different hardware, or for a System that is different in any way, loading is likely to fail. Checkpoints created with different versions of OpenMM are also often incompatible. If a checkpoint cannot be loaded, that is signaled by throwing an exception.

Parameters
checkpoint(string) the checkpoint data to load
def reinitialize (   self)

reinitialize(Context self)

When a Context is created, it may cache information about the System being simulated and the Force objects contained in it. This means that, if the System or Forces are then modified, the Context might not see all of the changes. Call reinitialize() to force the Context to rebuild its internal representation of the System and pick up any changes that have been made.

This is an expensive operation, so you should try to avoid calling it too frequently.

References Context.getState(), and simtk.openmm.openmm.stripUnits().

def setParameter (   self,
  args 
)

setParameter(Context self, std::string const & name, double value)

Set the value of an adjustable parameter defined by a Force object in the System.

Parameters
namethe name of the parameter to set
valuethe value of the parameter

References simtk.openmm.openmm.stripUnits().

Referenced by Context.setState().

def setPeriodicBoxVectors (   self,
  args 
)

setPeriodicBoxVectors(Context self, Vec3 const & a, Vec3 const & b, Vec3 const & c)

Set the vectors defining the axes of the periodic box (measured in nm). They will affect any Force that uses periodic boundary conditions.

Currently, only rectangular boxes are supported. This means that a, b, and c must be aligned with the x, y, and z axes respectively. Future releases may support arbitrary triclinic boxes.

Parameters
athe vector defining the first edge of the periodic box
bthe vector defining the second edge of the periodic box
cthe vector defining the third edge of the periodic box

References simtk.openmm.openmm.stripUnits().

Referenced by Context.setState().

def setPositions (   self,
  args 
)

setPositions(self, positions)

Set the positions of all particles in the System (measured in nm). This method simply sets the positions without checking to see whether they satisfy distance constraints. If you want constraints to be enforced, call applyConstraints() after setting the positions.

Parameters
positionsa vector whose length equals the number of particles in the System. The i'th element contains the position of the i'th particle.

References simtk.openmm.openmm.stripUnits().

Referenced by Context.setState().

def setState (   self,
  state 
)

setState(Context self, State state)

Copy information from a State object into this Context. This restores the Context to approximately the same state it was in when the State was created. If the State does not include a piece of information (e.g. positions or velocities), that aspect of the Context is left unchanged.

Even when all possible information is included in the State, the effect of calling this method is still less complete than loadCheckpoint(). For example, it does not restore the internal states of random number generators. On the other hand, it has the advantage of not being hardware specific.

References Context.setParameter(), Context.setPeriodicBoxVectors(), Context.setPositions(), Context.setTime(), and Context.setVelocities().

def setTime (   self,
  args 
)

setTime(Context self, double time)

Set the current time of the simulation (in picoseconds).

References simtk.openmm.openmm.stripUnits().

Referenced by Context.setState().

def setVelocities (   self,
  args 
)

setVelocities(self, velocities)

Set the velocities of all particles in the System (measured in nm/picosecond).

Parameters
velocitiesa vector whose length equals the number of particles in the System. The i'th element contains the velocity of the i'th particle.

References simtk.openmm.openmm.stripUnits().

Referenced by Context.setState().

def setVelocitiesToTemperature (   self,
  args 
)

setVelocitiesToTemperature(Context self, double temperature, int randomSeed=time(NULL)) setVelocitiesToTemperature(Context self, double temperature)

Set the velocities of all particles in the System to random values chosen from a Boltzmann distribution at a given temperature.

Parameters
temperaturethe temperature for which to select the velocities (measured in Kelvin)
randomSeedthe random number seed to use when selecting velocities

References simtk.openmm.openmm.stripUnits().

Member Data Documentation

this

Referenced by System.__init__().


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