perses.samplers.SAMSSampler

class perses.samplers.SAMSSampler(sampler, logZ=None, log_target_probabilities=None, update_method='two-stage', storage=None, second_stage_start=1000)[source]

Self-adjusted mixture sampling engine.

References

[1] Tan, Z. (2015) Optimally adjusted mixture sampling and locally weighted histogram analysis, Journal of Computational and Graphical Statistics, to appear. (Supplement) http://www.stat.rutgers.edu/home/ztan/Publication/SAMS_redo4.pdf

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
>>> from perses.rjmc.geometry import FFAllAngleGeometryEngine
>>> geometry_engine = FFAllAngleGeometryEngine(metadata={})
>>> # Create an Expanded Ensemble sampler
>>> from perses.rjmc.topology_proposal import PointMutationEngine
>>> 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)
>>> # Create a SAMS sampler
>>> sams_sampler = SAMSSampler(exen_sampler)
>>> # Run the sampler
>>> sams_sampler.run() # doctest: +ELLIPSIS
...
Attributes:
state_keys

Methods

run([niterations]) Run the sampler for the specified number of iterations
update() Update the sampler with one step of sampling.
update_logZ_estimates() Update the logZ estimates according to self.update_method.
update_sampler() Update the underlying expanded ensembles sampler.
__init__(sampler, logZ=None, log_target_probabilities=None, update_method='two-stage', storage=None, second_stage_start=1000)[source]

Create a SAMS Sampler.

Parameters:
sampler : ExpandedEnsembleSampler

The expanded ensemble sampler used to sample both configurations and discrete thermodynamic states.

logZ : dict of key

If specified, the log partition functions for each state will be initialized to the specified dictionary.

log_target_probabilities : dict of key

If specified, unnormalized target probabilities; default is all 0.

update_method : str, optional, default=’default’

SAMS update algorithm

storage : NetCDFStorageView, optional, default=None
second_state_start : int, optional, default None

At what iteration number to switch to the optimal gain decay

Methods

__init__(sampler[, logZ, …]) Create a SAMS Sampler.
run([niterations]) Run the sampler for the specified number of iterations
update() Update the sampler with one step of sampling.
update_logZ_estimates() Update the logZ estimates according to self.update_method.
update_sampler() Update the underlying expanded ensembles sampler.

Attributes

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

Update the logZ estimates according to self.update_method.

update_sampler()[source]

Update the underlying expanded ensembles sampler.