openmmtools.multistate.ReplicaExchangeSampler

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

Replica-exchange simulation facility.

This MultiStateSampler class provides a general replica-exchange simulation facility, allowing any set of thermodynamic states to be specified, along with a set of initial positions to be assigned to the replicas in a round-robin fashion.

No distinction is made between one-dimensional and multidimensional replica layout. By default, the replica mixing scheme attempts to mix all replicas to minimize slow diffusion normally found in multidimensional replica exchange simulations (Modification of the ‘replica_mixing_scheme’ setting will allow the traditional ‘neighbor swaps only’ scheme to be used.)

Stored configurations, energies, swaps, and restart information are all written to a single output file using the platform portable, robust, and efficient NetCDF4 library.

Parameters
mcmc_movesMCMCMove or list of MCMCMove, optional

The MCMCMove used to propagate the states. If a list of MCMCMoves, they will be assigned to the correspondent thermodynamic state on creation. If None is provided, Langevin dynamics with 2fm timestep, 5.0/ps collision rate, and 500 steps per iteration will be used.

number_of_iterationsint or infinity, optional, default: 1

The number of iterations to perform. Both float('inf') and numpy.inf are accepted for infinity. If you set this to infinity, be sure to set also online_analysis_interval.

replica_mixing_scheme‘swap-all’, ‘swap-neighbors’ or None, Default: ‘swap-all’

The scheme used to swap thermodynamic states between replicas.

online_analysis_intervalNone or Int >= 1, optional, default None

Choose the interval at which to perform online analysis of the free energy.

After every interval, the simulation will be stopped and the free energy estimated.

If the error in the free energy estimate is at or below online_analysis_target_error, then the simulation will be considered completed.

online_analysis_target_errorfloat >= 0, optional, default 0.2

The target error for the online analysis measured in kT per phase.

Once the free energy is at or below this value, the phase will be considered complete.

If online_analysis_interval is None, this option does nothing.

online_analysis_minimum_iterationsint >= 0, optional, default 50

Set the minimum number of iterations which must pass before online analysis is carried out.

Since the initial samples likely not to yield a good estimate of free energy, save time and just skip them If online_analysis_interval is None, this does nothing

Examples

Parallel tempering simulation of alanine dipeptide in implicit solvent (replica exchange among temperatures). This is just an illustrative example; use ParallelTempering class for actual production parallel tempering simulations.

Create the system.

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

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.
>>> temperatures = [T_min + (T_max - T_min) * (math.exp(float(i) / float(n_replicas-1)) - 1.0) / (math.e - 1.0)
...                 for i in range(n_replicas)]
>>> thermodynamic_states = [states.ThermodynamicState(system=testsystem.system, temperature=T)
...                         for T in temperatures]

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

>>> move = mcmc.GHMCMove(timestep=2.0*unit.femtoseconds, n_steps=50)
>>> simulation = ReplicaExchangeSampler(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 = multistate.MultiStateReporter(storage_path, checkpoint_interval=1)
>>> simulation.create(thermodynamic_states=thermodynamic_states,
...                   sampler_states=states.SamplerState(testsystem.positions),
...                   storage=reporter)
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()  # This runs for a maximum of 2 iterations.
>>> simulation.iteration
2
>>> simulation.run(n_iterations=1)
>>> simulation.iteration
2

To resume a simulation from an existing storage file and extend it beyond the original number of iterations.

>>> del simulation
>>> simulation = ReplicaExchangeSampler.from_storage(reporter)
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.extend(n_iterations=1)
>>> simulation.iteration
3

You can extract several information from the NetCDF file using the Reporter class while the simulation is running. This reads the SamplerStates of every run iteration.

>>> reporter = multistate.MultiStateReporter(storage=storage_path, open_mode='r', checkpoint_interval=1)
>>> sampler_states = reporter.read_sampler_states(iteration=1)
>>> len(sampler_states)
3
>>> sampler_states[0].positions.shape  # Alanine dipeptide has 22 atoms.
(22, 3)

Clean up.

>>> os.remove(storage_path)
Parameters
  • number_of_iterations – Maximum number of integer iterations that will be run

  • replica_mixing_scheme – Scheme which describes how replicas are exchanged each iteration as string

  • online_analysis_interval – How frequently to carry out online analysis in number of iterations

  • online_analysis_target_error – Target free energy difference error float at which simulation will be stopped during online analysis, in dimensionless energy

  • online_analysis_minimum_iterations – Minimum number of iterations needed before online analysis is run as int

Attributes
n_replicas

The integer number of replicas (read-only).

iteration

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

mcmc_moves

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

sampler_states

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

metadata

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

is_completed

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

Methods

Status(iteration, target_error, is_completed)

Attributes

create(thermodynamic_states, sampler_states, ...)

Create new multistate sampler simulation.

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)[source]

Methods

__init__([replica_mixing_scheme])

create(thermodynamic_states, sampler_states, ...)

Create new multistate sampler simulation.

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.