perses.samplers.ExpandedEnsembleSampler

class perses.samplers.ExpandedEnsembleSampler(sampler, topology, state_key, proposal_engine, geometry_engine, log_weights=None, options=None, platform=None, envname=None, storage=None, ncmc_write_interval=1)[source]

Method of expanded ensembles sampling engine.

The acceptance criteria is given in the reference document. Roughly, the proposal scheme is:

  • Draw a proposed chemical state k’, and calculate reverse proposal probability
  • Conditioned on k’ and the current positions x, generate new positions with the GeometryEngine
  • With new positions, jump to a hybrid system at lambda=0
  • Anneal from lambda=0 to lambda=1, accumulating work
  • Jump from the hybrid system at lambda=1 to the k’ system, and compute reverse GeometryEngine proposal
  • Add weight of chemical states k and k’ to acceptance probabilities

References

[1] Lyubartsev AP, Martsinovski AA, Shevkunov SV, and Vorontsov-Velyaminov PN. New approach to Monte Carlo calculation of the free energy: Method of expanded ensembles. JCP 96:1776, 1992 http://dx.doi.org/10.1063/1.462133

Examples

>>> # Create a test system
>>> test = testsystems.AlanineDipeptideVacuum()
>>> # Create a SystemGenerator and rebuild the System.
>>> from perses.rjmc.topology_proposal import SystemGenerator
>>> system_generator = SystemGenerator(['amber99sbildn.xml'], forcefield_kwargs={ 'nonbondedMethod' : app.NoCutoff, 'implicitSolvent' : None, 'constraints' : None })
>>> test.system = system_generator.build_system(test.topology)
>>> # Create a sampler state.
>>> sampler_state = SamplerState(system=test.system, positions=test.positions)
>>> # Create a thermodynamic state.
>>> thermodynamic_state = ThermodynamicState(system=test.system, temperature=298.0*unit.kelvin)
>>> # Create an MCMC sampler
>>> mcmc_sampler = MCMCSampler(thermodynamic_state, sampler_state)
>>> # Turn off verbosity
>>> mcmc_sampler.verbose = False
>>> # Create an Expanded Ensemble sampler
>>> from perses.rjmc.topology_proposal import PointMutationEngine
>>> from perses.rjmc.geometry import FFAllAngleGeometryEngine
>>> geometry_engine = FFAllAngleGeometryEngine(metadata={})
>>> allowed_mutations = [[('2','ALA')],[('2','VAL'),('2','LEU')]]
>>> proposal_engine = PointMutationEngine(test.topology, system_generator, max_point_mutants=1, chain_id='1', proposal_metadata=None, allowed_mutations=allowed_mutations)
>>> exen_sampler = ExpandedEnsembleSampler(mcmc_sampler, test.topology, 'ACE-ALA-NME', proposal_engine, geometry_engine)
>>> # Run the sampler
>>> exen_sampler.run()
Attributes:
state_keys

Methods

get_log_weight(state_key) Get the log weight of the specified state.
run([niterations]) Run the sampler for the specified number of iterations
update() Update the sampler with one step of sampling.
update_positions([n_iterations]) Sample new positions.
update_state() Sample the thermodynamic state.
update_statistics() Update sampler statistics.
__init__(sampler, topology, state_key, proposal_engine, geometry_engine, log_weights=None, options=None, platform=None, envname=None, storage=None, ncmc_write_interval=1)[source]

Create an expanded ensemble sampler.

p(x,k) propto exp[-u_k(x) + g_k]

where g_k is the log weight.

Parameters:
sampler : MCMCSampler

MCMCSampler initialized with current SamplerState

topology : simtk.openmm.app.Topology

Current topology

state : hashable object

Current chemical state

proposal_engine : ProposalEngine

ProposalEngine to use for proposing new chemical states

geometry_engine : GeometryEngine

GeometryEngine to use for dimension matching

log_weights : dict of object

Log weights to use for expanded ensemble biases.

options : dict, optional, default=dict()

Options for initializing switching scheme, such as ‘timestep’, ‘nsteps’, ‘functions’ for NCMC

platform : simtk.openmm.Platform, optional, default=None

Platform to use for NCMC switching. If None, default (fastest) platform is used.

storage : NetCDFStorageView, optional, default=None

If specified, use this storage layer.

ncmc_write_interval : int, default 1

How frequently to write out NCMC protocol steps.

Methods

__init__(sampler, topology, state_key, …) Create an expanded ensemble sampler.
get_log_weight(state_key) Get the log weight of the specified state.
run([niterations]) Run the sampler for the specified number of iterations
update() Update the sampler with one step of sampling.
update_positions([n_iterations]) Sample new positions.
update_state() Sample the thermodynamic state.
update_statistics() Update sampler statistics.

Attributes

state_keys
get_log_weight(state_key)[source]

Get the log weight of the specified state.

Parameters:
state_key : hashable object

The state key (e.g. chemical state key) to look up.

Returns:
log_weight : float

The log weight of the provided state key.

run(niterations=1)[source]

Run the sampler for the specified number of iterations

Parameters:
niterations : int, optional, default=1

Number of iterations to run the sampler for.

update()[source]

Update the sampler with one step of sampling.

update_positions(n_iterations=1)[source]

Sample new positions.

update_state()[source]

Sample the thermodynamic state.

update_statistics()[source]

Update sampler statistics.