openmmtools.mcmc.MCMCSampler

class openmmtools.mcmc.MCMCSampler(thermodynamic_state, sampler_state, move)[source]

Basic Markov chain Monte Carlo sampler.

Parameters:
thermodynamic_stateopenmmtools.states.ThermodynamicState

Initial thermodynamic state.

sampler_stateopenmmtools.states.SamplerState

Initial sampler state.

move_setcontainer of MarkovChainMonteCarloMove objects

Moves to attempt during MCMC run. If list or tuple, will run all moves each iteration in specified sequence (e.g. [move1, move2, move3]). If dict, will use specified unnormalized weights (e.g. { move1 : 0.3, move2 : 0.5, move3, 0.9 })

Examples

>>> import numpy as np
>>> from openmm import unit
>>> from openmmtools import testsystems
>>> from openmmtools.states import ThermodynamicState, SamplerState

Create and run an alanine dipeptide simulation with a weighted move.

>>> test = testsystems.AlanineDipeptideVacuum()
>>> thermodynamic_state = ThermodynamicState(system=test.system,
...                                          temperature=298*unit.kelvin)
>>> sampler_state = SamplerState(positions=test.positions)
>>> # Create a move set specifying probabilities fo each type of move.
>>> move = WeightedMove([(HMCMove(n_steps=10), 0.5), (LangevinDynamicsMove(n_steps=10), 0.5)])
>>> # Create an MCMC sampler instance and run 10 iterations of the simulation.
>>> sampler = MCMCSampler(thermodynamic_state, sampler_state, move=move)
>>> sampler.run(n_iterations=2)
>>> np.allclose(sampler.sampler_state.positions, test.positions)
False

NPT ensemble simulation of a Lennard Jones fluid with a sequence of moves.

>>> test = testsystems.LennardJonesFluid(nparticles=200)
>>> thermodynamic_state = ThermodynamicState(system=test.system, temperature=298*unit.kelvin,
...                                          pressure=1*unit.atmospheres)
>>> sampler_state = SamplerState(positions=test.positions)
>>> # Create a move set that includes a Monte Carlo barostat move.
>>> move = SequenceMove([GHMCMove(n_steps=50), MonteCarloBarostatMove(n_attempts=5)])
>>> # Create an MCMC sampler instance and run 5 iterations of the simulation.
>>> sampler = MCMCSampler(thermodynamic_state, sampler_state, move=move)
>>> sampler.run(n_iterations=2)
>>> np.allclose(sampler.sampler_state.positions, test.positions)
False
Attributes:
thermodynamic_stateopenmmtools.states.ThermodynamicState

Current thermodynamic state.

sampler_stateopenmmtools.states.SamplerState

Current sampler state.

move_setcontainer of MarkovChainMonteCarloMove objects

Moves to attempt during MCMC run. If list or tuple, will run all moves each iteration in specified sequence (e.g. [move1, move2, move3]). If dict, will use specified unnormalized weights (e.g. { move1 : 0.3, move2 : 0.5, move3, 0.9 })

Methods

minimize([tolerance, max_iterations, ...])

Minimize the current configuration.

run([n_iterations, context_cache])

Run the sampler for a specified number of iterations.

__init__(thermodynamic_state, sampler_state, move)[source]

Methods

__init__(thermodynamic_state, sampler_state, ...)

minimize([tolerance, max_iterations, ...])

Minimize the current configuration.

run([n_iterations, context_cache])

Run the sampler for a specified number of iterations.