openmmtools.cache.ContextCache

class openmmtools.cache.ContextCache(platform=None, **kwargs)[source]

LRU cache hosting the minimum amount of incompatible Contexts.

Two Contexts are compatible if they are in a compatible ThermodynamicState, and have compatible integrators. In general, two integrators are compatible if they have the same serialized state, but ContextCache can decide to store a single Context to optimize memory when two integrators differ by only few parameters that can be set after the Context is initialized. These parameters include all the global variables defined by a CustomIntegrator.

You can force ContextCache to consider an integrator global variable incompatible by adding it to the blacklist ContextCache.INCOMPATIBLE_INTEGRATOR_ATTRIBUTES. Similarly, you can add other attributes that should be considered compatible through the whitelist ContextCache.COMPATIBLE_INTEGRATOR_ATTRIBUTES. If an attribute in that dictionary is not found in the integrator, the cache will search for a corresponding getter and setter.

Parameters:
platform : simtk.openmm.Platform, optional

The OpenMM platform to use to create Contexts. If None, OpenMM tries to select the fastest one available (default is None).

**kwargs

Parameters to pass to the underlying LRUCache constructor such as capacity and time_to_live.

Warning

Python instance attributes are not copied when ContextCache.get_context() is called. You can force this by setting adding them to the whitelist ContextCache.COMPATIBLE_INTEGRATOR_ATTRIBUTES, but if modifying your Python attributes won’t modify the OpenMM serialization, this will likely cause problems so this is discouraged unless you know exactly what you are doing.

See also

LRUCache, states.ThermodynamicState.is_state_compatible

Examples

>>> from simtk import unit
>>> from openmmtools import testsystems
>>> from openmmtools.states import ThermodynamicState
>>> alanine = testsystems.AlanineDipeptideExplicit()
>>> thermodynamic_state = ThermodynamicState(alanine.system, 310*unit.kelvin)
>>> time_step = 1.0*unit.femtosecond

Two compatible thermodynamic states generate only a single cached Context. ContextCache can also (in few explicitly supported cases) recycle the same Context even if the integrators differ by some parameters.

>>> context_cache = ContextCache()
>>> context1, integrator1 = context_cache.get_context(thermodynamic_state,
...                                                   openmm.VerletIntegrator(time_step))
>>> thermodynamic_state.temperature = 300*unit.kelvin
>>> time_step2 = 2.0*unit.femtosecond
>>> context2, integrator2 = context_cache.get_context(thermodynamic_state,
...                                                   openmm.VerletIntegrator(time_step2))
>>> id(context1) == id(context2)
True
>>> len(context_cache)
1

When we switch to NPT the states are not compatible and so neither the Contexts are.

>>> integrator2 = openmm.VerletIntegrator(2.0*unit.femtosecond)
>>> thermodynamic_state_npt = copy.deepcopy(thermodynamic_state)
>>> thermodynamic_state_npt.pressure = 1.0*unit.atmosphere
>>> context3, integrator3 = context_cache.get_context(thermodynamic_state_npt,
...                                                   openmm.VerletIntegrator(time_step))
>>> id(context1) == id(context3)
False
>>> len(context_cache)
2

You can set a capacity and a time to live for contexts like in a normal LRUCache.

>>> context_cache = ContextCache(capacity=1, time_to_live=5)
>>> context2, integrator2 = context_cache.get_context(thermodynamic_state,
...                                                   openmm.VerletIntegrator(time_step))
>>> context3, integrator3 = context_cache.get_context(thermodynamic_state_npt,
...                                                   openmm.VerletIntegrator(time_step))
>>> len(context_cache)
1
Attributes:
platform

The OpenMM platform to use to create Contexts.

capacity

The maximum number of Context cached.

time_to_live

The Contexts expiration date in number of accesses to the LRUCache.

Methods

empty() Clear up cache and remove all Contexts.
get_context(thermodynamic_state[, integrator]) Return a context in the given thermodynamic state.
__init__(platform=None, **kwargs)[source]

x.__init__(…) initializes x; see help(type(x)) for signature

Methods

__init__([platform]) x.__init__(…) initializes x; see help(type(x)) for signature
empty() Clear up cache and remove all Contexts.
get_context(thermodynamic_state[, integrator]) Return a context in the given thermodynamic state.

Attributes

COMPATIBLE_INTEGRATOR_ATTRIBUTES
INCOMPATIBLE_INTEGRATOR_ATTRIBUTES
capacity The maximum number of Context cached.
platform The OpenMM platform to use to create Contexts.
time_to_live The Contexts expiration date in number of accesses to the LRUCache.
capacity

The maximum number of Context cached.

If None, the capacity is unlimited.

empty()[source]

Clear up cache and remove all Contexts.

get_context(thermodynamic_state, integrator=None)[source]

Return a context in the given thermodynamic state.

In general, the Context must be considered newly initialized. This means that positions and velocities must be set afterwards.

If the integrator is not provided, this will search the cache for any Context in the given ThermodynamicState, regardless of its integrator. In this case, the method guarantees that two consecutive calls with the same thermodynamic state will retrieve the same context.

This creates a new Context if no compatible one has been cached. If a compatible Context exists, the ThermodynamicState is applied to it, and the Context integrator state is changed to match the one passed as an argument. As a consequence, the returned integrator is guaranteed to be in the same state as the one provided, but it can be a different instance. This is to minimize the number of Contexts objects cached that use the same or very similar integrator.

Parameters:
thermodynamic_state : states.ThermodynamicState

The thermodynamic state of the system.

integrator : simtk.openmm.Integrator, optional

The integrator for the context (default is None).

Returns:
context : simtk.openmm.Context

The context in the given thermodynamic system.

context_integrator : simtk.openmm.Integrator

The integrator to be used to propagate the Context. Can be a difference instance from the one passed as an argument.

Warning

Python instance attributes are not copied when get_context() is called. You can force this by setting adding them to the whitelist ContextCache.COMPATIBLE_INTEGRATOR_ATTRIBUTES, but if modifying the attributes won’t modify the OpenMM serialization, this will likely cause problems so this is discouraged unless you know exactly what you’re doing.

platform

The OpenMM platform to use to create Contexts.

If None, OpenMM tries to select the fastest one available. This can be set only if the cache is empty.

time_to_live

The Contexts expiration date in number of accesses to the LRUCache.

If None, Contexts do not expire.