openmmtools.multistate.MultiStateSampler

class openmmtools.multistate.MultiStateSampler(mcmc_moves=None, number_of_iterations=1, online_analysis_interval=200, online_analysis_target_error=0.0, online_analysis_minimum_iterations=200, locality=None)[source]

Base class for samplers that sample multiple thermodynamic states using one or more replicas.

This base class provides a general simulation facility for multistate from multiple thermodynamic states, allowing any set of thermodynamic states to be specified. If instantiated on its own, the thermodynamic state indices associated with each state are specified and replica mixing does not change any thermodynamic states, meaning that each replica remains in its original thermodynamic state.

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 thermodynamic 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.

online_analysis_intervalNone or Int >= 1, optional, default: 200

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.

If set to None, then no online analysis is performed

online_analysis_target_errorfloat >= 0, optional, default 0.0

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.

Default is set to 0.0 since online analysis runs by default, but a finite number_of_iterations should also be set to ensure there is some stop condition. If target error is 0 and an infinite number of iterations is set, then the sampler will run until the user stop it manually.

online_analysis_minimum_iterationsint >= 0, optional, default 200

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

localityint > 0, optional, default None

If None, the energies at all states will be computed for every replica each iteration. If int > 0, energies will only be computed for states range(max(0, state-locality), min(n_states, state+locality)).

Examples

Sampling multiple states of an alanine dipeptide in implicit solvent system.

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

Create thermodynamic states

>>> n_replicas = 3
>>> 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)]
>>> 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 = MultiStateSampler(mcmc_moves=move, number_of_iterations=2)

Create simulation and store output in temporary file

>>> storage_path = tempfile.NamedTemporaryFile(delete=False).name + '.nc'
>>> reporter = MultiStateReporter(storage_path, checkpoint_interval=1)
>>> simulation.create(thermodynamic_states=thermodynamic_states,
...                   sampler_states=states.SamplerState(testsystem.positions), storage=reporter)

Optionally, specify unlimited context cache attributes using the fastest mixed precision platform

>>> from openmmtools.cache import ContextCache
>>> from openmmtools.utils import get_fastest_platform
>>> platform = get_fastest_platform(minimum_precision='mixed')
>>> simulation.energy_context_cache = ContextCache(capacity=None, time_to_live=None, platform=platform)
>>> simulation.sampler_context_cache = ContextCache(capacity=None, time_to_live=None, platform=platform)

Run the simulation >>> simulation.run()

Note that to avoid the repex cycling problem upon resuming a simulation, make sure to specify do the optional energy and sampler context caches.

>>> reporter = MultiStateReporter(reporter_file, checkpoint_interval=10)
>>> simulation = HybridRepexSampler.from_storage(reporter)
>>> simulation.energy_context_cache = cache.ContextCache(capacity=None, time_to_live=None, platform=platform)
>>> simulation.sampler_context_cache = cache.ContextCache(capacity=None, time_to_live=None, platform=platform)
>>> simulation.extend(n_iterations=1)
Attributes
n_replicas

The integer number of replicas (read-only).

n_states

The integer number of thermodynamic states (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)

energy_context_cacheopenmmtools.cache.ContextCache, default=openmmtools.cache.global_context_cache

Context cache to be used for energy computations. Defaults to using global context cache.

sampler_context_cacheopenmmtools.cache.ContextCache, default=openmmtools.cache.global_context_cache

Context cache to be used for propagation. Defaults to using global context cache.

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__(mcmc_moves=None, number_of_iterations=1, online_analysis_interval=200, online_analysis_target_error=0.0, online_analysis_minimum_iterations=200, locality=None)[source]

Methods

__init__([mcmc_moves, number_of_iterations, ...])

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).

online_analysis_interval

interval to carry out online analysis

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.