openmmtools.multistate.ParallelTemperingSampler

class openmmtools.multistate.ParallelTemperingSampler(replica_mixing_scheme='swap-all', **kwargs)[source]

Parallel tempering simulation facility.

This class provides a facility for parallel tempering simulations. It is a subclass of ReplicaExchange, but provides efficiency improvements for parallel tempering simulations, so should be preferred for this type of simulation. In particular, this makes use of the fact that the reduced potentials are linear in inverse temperature.

Examples

Create the system.

>>> from openmm import unit
>>> from openmmtools import testsystems, states, mcmc
>>> import tempfile
>>> testsystem = testsystems.AlanineDipeptideImplicit()
>>> import os

Create thermodynamic states for parallel tempering with exponentially-spaced schedule.

>>> n_replicas = 3  # Number of temperature replicas.
>>> T_min = 298.0 * unit.kelvin  # Minimum temperature.
>>> T_max = 600.0 * unit.kelvin  # Maximum temperature.
>>> reference_state = states.ThermodynamicState(system=testsystem.system, temperature=T_min)

Initialize simulation object with options. Run with a GHMC integrator.

>>> move = mcmc.GHMCMove(timestep=2.0*unit.femtoseconds, n_steps=50)
>>> simulation = ParallelTemperingSampler(mcmc_moves=move, number_of_iterations=2)

Create simulation with its storage file (in a temporary directory) and run.

>>> storage_path = tempfile.NamedTemporaryFile(delete=False).name + '.nc'
>>> reporter = MultiStateReporter(storage_path, checkpoint_interval=10)
>>> simulation.create(reference_state,
...                   states.SamplerState(testsystem.positions),
...                   reporter, min_temperature=T_min,
...                   max_temperature=T_max, n_temperatures=n_replicas)
Please cite the following:

        Friedrichs MS, Eastman P, Vaidyanathan V, Houston M, LeGrand S, Beberg AL, Ensign DL, Bruns CM, and Pande VS. Accelerating molecular dynamic simulations on graphics processing unit. J. Comput. Chem. 30:864, 2009. DOI: 10.1002/jcc.21209
        Eastman P and Pande VS. OpenMM: A hardware-independent framework for molecular simulations. Comput. Sci. Eng. 12:34, 2010. DOI: 10.1109/MCSE.2010.27
        Eastman P and Pande VS. Efficient nonbonded interactions for molecular dynamics on a graphics processing unit. J. Comput. Chem. 31:1268, 2010. DOI: 10.1002/jcc.21413
        Eastman P and Pande VS. Constant constraint matrix approximation: A robust, parallelizable constraint method for molecular simulations. J. Chem. Theor. Comput. 6:434, 2010. DOI: 10.1021/ct900463w
        Chodera JD and Shirts MR. Replica exchange and expanded ensemble simulations as Gibbs multistate: Simple improvements for enhanced mixing. J. Chem. Phys., 135:194110, 2011. DOI:10.1063/1.3660669

>>> simulation.run(n_iterations=1)

Clean up.

>>> os.remove(storage_path)
Attributes
is_completed

Check if we have reached any of the stop target criteria (read-only)

is_periodic

Return True if system is periodic, False if not, and None if not initialized

iteration

The integer current iteration of the simulation (read-only).

mcmc_moves

A copy of the MCMCMoves list used to propagate the simulation.

metadata

A copy of the metadata dictionary passed on creation (read-only).

n_replicas

The integer number of replicas (read-only).

n_states

The integer number of thermodynamic states (read-only).

options

dict of all class options (keyword arguments for __init__ for class and superclasses)

sampler_states

A copy of the sampler states list at the current iteration.

Methods

Status(iteration, target_error, is_completed)

Attributes

create(thermodynamic_state, sampler_states, ...)

Initialize a parallel tempering simulation object.

default_options()

dict of all default class options (keyword arguments for __init__ for class and superclasses)

equilibrate(n_iterations[, mcmc_moves])

Equilibrate all replicas.

extend(n_iterations)

Extend the simulation by the given number of iterations.

from_storage(storage)

Constructor from an existing storage file.

minimize([tolerance, max_iterations])

Minimize all replicas.

read_status(storage)

Read the status of the calculation from the storage file.

run([n_iterations])

Run the replica-exchange simulation.

__init__(replica_mixing_scheme='swap-all', **kwargs)

Methods

__init__([replica_mixing_scheme])

create(thermodynamic_state, sampler_states, ...)

Initialize a parallel tempering simulation object.

default_options()

dict of all default class options (keyword arguments for __init__ for class and superclasses)

equilibrate(n_iterations[, mcmc_moves])

Equilibrate all replicas.

extend(n_iterations)

Extend the simulation by the given number of iterations.

from_storage(storage)

Constructor from an existing storage file.

minimize([tolerance, max_iterations])

Minimize all replicas.

read_status(storage)

Read the status of the calculation from the storage file.

run([n_iterations])

Run the replica-exchange simulation.

Attributes

is_completed

Check if we have reached any of the stop target criteria (read-only)

is_periodic

Return True if system is periodic, False if not, and None if not initialized

iteration

The integer current iteration of the simulation (read-only).

mcmc_moves

A copy of the MCMCMoves list used to propagate the simulation.

metadata

A copy of the metadata dictionary passed on creation (read-only).

n_replicas

The integer number of replicas (read-only).

n_states

The integer number of thermodynamic states (read-only).

options

dict of all class options (keyword arguments for __init__ for class and superclasses)

sampler_states

A copy of the sampler states list at the current iteration.