openmmtools.mcmc.HMCMove

class openmmtools.mcmc.HMCMove(timestep=Quantity(value=1.0, unit=femtosecond), n_steps=1000, **kwargs)[source]

Hybrid Monte Carlo dynamics.

This move assigns a velocity from the Maxwell-Boltzmann distribution and executes a number of velocity Verlet steps to propagate dynamics.

Parameters:
timestep : simtk.unit.Quantity, optional

The timestep to use for HMC dynamics, which uses velocity Verlet following velocity randomization (time units, default is 1*simtk.unit.femtosecond)

n_steps : int, optional

The number of dynamics steps to take before Metropolis acceptance/rejection (default is 1000).

context_cache : openmmtools.cache.ContextCache, optional

The ContextCache to use for Context creation. If None, the global cache openmmtools.cache.global_context_cache is used (default is None).

Examples

First we need to create the thermodynamic state and the sampler state to propagate. Here we create an alanine dipeptide system in vacuum.

>>> from simtk import unit
>>> from openmmtools import testsystems
>>> from openmmtools.states import ThermodynamicState, SamplerState
>>> test = testsystems.AlanineDipeptideVacuum()
>>> sampler_state = SamplerState(positions=test.positions)
>>> thermodynamic_state = ThermodynamicState(system=test.system, temperature=298*unit.kelvin)

Create a GHMC move with default parameters.

>>> move = HMCMove()

or create a GHMC move with specified parameters.

>>> move = HMCMove(timestep=0.5*unit.femtoseconds, n_steps=10)

Perform one update of the sampler state. The sampler state is updated with the new state.

>>> move.apply(thermodynamic_state, sampler_state)
>>> np.allclose(sampler_state.positions, test.positions)
False

The same move can be applied to a different state, here an ideal gas.

>>> test = testsystems.IdealGas()
>>> sampler_state = SamplerState(positions=test.positions)
>>> thermodynamic_state = ThermodynamicState(system=test.system,
...                                          temperature=298*unit.kelvin)
>>> move.apply(thermodynamic_state, sampler_state)
>>> np.allclose(sampler_state.positions, test.positions)
False
Attributes:
timestep : simtk.unit.Quantity

The timestep to use for HMC dynamics, which uses velocity Verlet following velocity randomization (time units).

n_steps : int

The number of dynamics steps to take before Metropolis acceptance/rejection.

context_cache : openmmtools.cache.ContextCache

The ContextCache to use for Context creation. If None, the global cache openmmtools.cache.global_context_cache is used.

Methods

apply(thermodynamic_state, sampler_state) Apply the MCMC move.
__init__(timestep=Quantity(value=1.0, unit=femtosecond), n_steps=1000, **kwargs)[source]

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

Methods

__init__([timestep, unit, n_steps]) x.__init__(…) initializes x; see help(type(x)) for signature
apply(thermodynamic_state, sampler_state) Apply the MCMC move.
apply(thermodynamic_state, sampler_state)[source]

Apply the MCMC move.

This modifies the given sampler_state.

Parameters:
thermodynamic_state : openmmtools.states.ThermodynamicState

The thermodynamic state to use when applying the MCMC move.

sampler_state : openmmtools.states.SamplerState

The sampler state to apply the move to. This is modified.