openmmtools.mcmc.GHMCMove

class openmmtools.mcmc.GHMCMove(timestep=Quantity(value=1.0, unit=femtosecond), collision_rate=Quantity(value=20.0, unit=/picosecond), n_steps=1000, **kwargs)[source]

Generalized hybrid Monte Carlo (GHMC) Markov chain Monte Carlo move.

This move uses generalized Hybrid Monte Carlo (GHMC), a form of Metropolized Langevin dynamics, to propagate the system.

Parameters:
timestepopenmm.unit.Quantity, optional

The timestep to use for Langevin integration (time units, default is 1*openmm.unit.femtoseconds).

collision_rateopenmm.unit.Quantity, optional

The collision rate with fictitious bath particles (1/time units, default is 20/openmm.unit.picoseconds).

n_stepsint, optional

The number of integration timesteps to take each time the move is applied (default is 1000).

References

Lelievre T, Stoltz G, Rousset M. Free energy computations: A mathematical perspective. World Scientific, 2010.

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 openmm 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 = GHMCMove()

or create a GHMC move with specified parameters.

>>> move = GHMCMove(timestep=0.5*unit.femtoseconds,
...                 collision_rate=20.0/unit.picoseconds, 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,context_cache=context_cache)
>>> 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,context_cache=context_cache)
>>> np.allclose(sampler_state.positions, test.positions)
False
Attributes:
timestepopenmm.unit.Quantity

The timestep to use for Langevin integration (time units).

collision_rateopenmm.unit.Quantity

The collision rate with fictitious bath particles (1/time units).

n_stepsint

The number of integration timesteps to take each time the move is applied.

n_acceptedint

The number of accepted steps.

n_proposedint

The number of attempted steps.

fraction_accepted

Ratio between accepted over attempted moves (read-only).

Methods

apply(thermodynamic_state, sampler_state[, ...])

Apply the GHMC MCMC move.

reset_statistics()

Reset the internal statistics of number of accepted and attempted moves.

__init__(timestep=Quantity(value=1.0, unit=femtosecond), collision_rate=Quantity(value=20.0, unit=/picosecond), n_steps=1000, **kwargs)[source]

Methods

__init__([timestep, collision_rate, n_steps])

apply(thermodynamic_state, sampler_state[, ...])

Apply the GHMC MCMC move.

reset_statistics()

Reset the internal statistics of number of accepted and attempted moves.

Attributes

context_cache

fraction_accepted

Ratio between accepted over attempted moves (read-only).

statistics

The acceptance statistics as a dictionary.