coniii.samplers module
Monte Carlo samplers for Boltzmann-type (maxent) models.
This module provides the sampling machinery used by the iterative
solvers in coniii.solvers (notably MCH)
and exposed to users through the convenience function
sample_ising().
Public API (see __all__)
SamplerAbstract base class for all samplers.
MetropolisSingle-spin-flip Metropolis sampler for the Ising model. The workhorse; accelerated by the optional Boost C++ extension (
coniii.samplers_ext) when it is compiled, otherwise a numba-jitted pure-Python path is used.WolffIsingWolff cluster sampler.
ParallelTemperingReplica-exchange (parallel tempering) wrapper around Metropolis.
Potts3Three-state Potts variant of
Metropolis.sample_ising()One-call helper to draw Metropolis samples from an Ising model.
Work-in-progress samplers (SWIsing, HamiltonianMC,
Heisenberg3DSampler) live in coniii.experimental.samplers.
- class coniii.samplers.Sampler(n, theta, **kwargs)[source]
Bases:
objectBase class for MCMC sampling.
- class coniii.samplers.Metropolis(n, theta, calc_e=None, rng=None, boost=True, iprint=True)[source]
Bases:
Sampler- generate_sample_boost(sample_size, n_iters=1000, burn_in=None, systematic_iter=False)[source]
Generate Metropolis samples using C++ and boost.
- Parameters:
sample_size (int) – Number of samples.
n_iters (int, 1000) – Number of Metropolis iterations between samples.
burn_in (int, None) – If not set, will be the same value as n_iters.
systematic_iter (bool, False) – If True, iterate through each element of system by increment index by one.
- Returns:
Saved array of energies at each sampling step.
- Return type:
ndarray, optional
- generate_sample_py(sample_size, n_iters=1000, burn_in=None, systematic_iter=False, saveHistory=False, initial_sample=None)[source]
Generate Metropolis samples using a for loop.
- Parameters:
sample_size (int) – Number of samples.
n_iters (int, 1000) – Number of iterations to run the sampler floor.
burn_in (int, None) – If None, n_iters is used.
systematic_iter (bool, False) – If True, iterate through each element of system by increment index by one.
saveHistory (bool, False) – If True, also save the energy of each sample at each sampling step.
initial_sample (ndarray, None) – Start with this sample (i.e. to avoid warming up). Otherwise, self._samples is the initial sample.
- Returns:
Saved array of energies at each sampling step.
- Return type:
ndarray, optional
- generate_sample_parallel_boost(sample_size, n_iters=1000, burn_in=None, systematic_iter=False)[source]
Generate samples in parallel. Each replica in self._samples runs on its own thread and a sample is generated every n_iters.
In order to control the random number generator, we pass in seeds that are samples from the class instance’s rng.
- Parameters:
sample_size (int) – Number of samples.
n_iters (int, 1000) – Number of iterations between taking a random sample.
burn_in (int, None) – If None, n_iters is used.
systematic_iter (bool, False) – If True, iterate through spins systematically instead of choosing them randomly.
- generate_sample_parallel_py(sample_size, n_iters=1000, burn_in=None, initial_sample=None, systematic_iter=False)[source]
Generate samples in parallel. Each replica in self._samples runs on its own thread and a sample is generated every n_iters.
In order to control the random number generator, we pass in seeds that are samples from the class instance’s rng.
- Parameters:
sample_size (int) – Number of samples.
n_iters (int, 1000) – Number of iterations between taking a random sample.
burn_in (int, None) – If None, n_iters is used.
initial_sample (ndarray, None) – Starting set of replicas otherwise self._samples is used.
systematic_iter (bool, False) – If True, iterate through spins systematically instead of choosing them randomly.
- generate_cond_sample_boost(sample_size, fixed_subset, burn_in=1000, n_iters=1000, initial_sample=None, systematic_iter=False, parallel=True)[source]
Generate samples from conditional distribution (while a subset of the spins are held fixed).
- Parameters:
sample_size (int)
fixed_subset (list of duples) – Each duple is the index of the spin and the value to fix it at. These should be ordered by spin index.
burn_in (int, 1000) – Burn in.
n_iters (int, 1000)
initial_sample (ndarray, None) – Option to set initial random sample.
systematic_iter (bool, False) – Iterate through spins systematically instead of choosing them randomly.
parallel (bool, True) – If True, use parallelized routine.
- Returns:
ndarray – Samples from distribution.
ndarray – Energy of each sample.
- generate_cond_sample_py(sample_size, fixed_subset, burn_in=1000, n_iters=1000, initial_sample=None, systematic_iter=False, parallel=False)[source]
Generate samples from conditional distribution (while a subset of the spins are held fixed). Samples are generated in parallel.
NOTE: There is a bug with multiprocess where many calls to the parallel sampling routine in a row leads to increasingly slow evaluation of the code.
- Parameters:
sample_size (int)
fixed_subset (list of duples) – Each duple is the index of the spin and the value to fix it at. These should be ordered by spin index.
burn_in (int, 1000) – Burn in.
n_iters (int, 1000)
initial_sample (ndarray, None) – Option to set initial random sample.
systematic_iter (bool, False) – Iterate through spins systematically instead of choosing them randomly.
parallel (bool, False) – If True, use parallelized routine.
- Returns:
ndarray – Samples from distribution.
ndarray – Energy of each sample.
- sample_metropolis(sample0, E0, rng=None, flip_site=None, calc_e=None)[source]
Metropolis sampling given an arbitrary energy function.
- Parameters:
sample0 (ndarray) – Sample to start with. Passed by ref and changed.
E0 (ndarray) – Initial energy of state.
rng (np.random.RandomState) – Random number generator.
flip_site (int) – Site to flip.
calc_e (function) – If another function to calculate energy should be used
- Returns:
delta energy.
- Return type:
float
- class coniii.samplers.WolffIsing(J, h)[source]
Bases:
Sampler- generate_sample(samplesize, n_iters, initialSample=None, save_history=False)[source]
Generate samples by starting from random initial states.
- generate_sample_parallel(samplesize, n_iters, initialSample=None)[source]
Generate samples by starting from random or given initial states.
- one_step(state, initialsite=None)[source]
Run one iteration of the Wolff algorithm that involves finding a cluster and possibly flipping it.
- find_neighbors(state, site, alreadyMarked)[source]
Return neighbors of given site that need to be visited excluding sites that have already been visited. This is the implementation of the Wolff algorithm for finding neighbors such that detailed balance is satisfied. I have modified to include random fields such tha the probability of adding a neighbors depends both on its coupling with the current site and the neighbor’s magnetic field.
- Parameters:
state
site
alreadyMarked
- class coniii.samplers.ParallelTempering(n, theta, calc_e, n_replicas, Tbds=(1.0, 3.0), sample_size=1000, replica_burnin=None, rep_ex_burnin=None, rng=None)[source]
Bases:
Sampler- update_replica_parameters()[source]
Update parameters for each replica. Remember that the parameters include the factor of beta.
- setup_replicas()[source]
Initialise a set of replicas at different temperatures using the Metropolis algorithm and optimize the temperatures. Replicas are burned in and ready to sample.
- burn_in_replicas(pool=None, close_pool=True, n_iters=None)[source]
Run each replica separately.
- Parameters:
pool (multiprocess.Pool, None)
close_pool (bool, True) – If True, call pool.close() at end.
n_iters (int, None) – Default value is self.replicaBurnin.
- generate_sample(sample_size, save_exchange_trajectory=False)[source]
Burn in, run replica exchange simulation, then sample.
- Parameters:
sample_size (int) – Number of samples to take for each replica.
save_exchange_trajectory (bool, False) – If True, keep track of the location of each replica in beta space and return the history.
- Returns:
Trajectory of each replica through beta space. Each row is tells where each index is located in beta space.
- Return type:
ndarray, optional
- static iterate_beta(beta, acceptance_ratio)[source]
Apply algorithm from Hukushima but reversed to maintain one replica at T=1.
- Parameters:
beta (ndarray) – Inverse temperature.
acceptance_ratio (ndarray) – Estimate of acceptance ratio.
- Returns:
New beta.
- Return type:
ndarray
- optimize_beta(n_samples, n_iters, tol=0.01, max_iter=10)[source]
Find suitable temperature range for replicas. Sets self.beta.
- Parameters:
n_samples (int) – Number of samples to use to estimate acceptance ratio. Acceptance ratio is estimated as the average of these samples.
n_iters (int) – Number of sampling iterations for each replica.
tol (float, .1) – Average change in beta to reach before stopping.
max_iter (int, 10) – Number of times to iterate algorithm for beta. Each iteration involves sampling from replicas.
- class coniii.samplers.Potts3(n, theta, calc_e=None, n_cpus=None, rng=None, boost=True)[source]
Bases:
Metropolis- generate_sample_parallel_boost(sample_size, n_iters=1000, burn_in=None, systematic_iter=False)[source]
Generate samples in parallel. Each replica in self._samples runs on its own thread and a sample is generated every n_iters.
In order to control the random number generator, we pass in seeds that are samples from the class instance’s rng.
- Parameters:
sample_size (int) – Number of samples.
n_iters (int, 1000) – Number of iterations between taking a random sample.
burn_in (int, None) – If None, n_iters is used.
systematic_iter (bool, False) – If True, iterate through spins systematically instead of choosing them randomly.
- sample_metropolis(sample0, E0, rng=None, flip_site=None, calc_e=None)[source]
Metropolis sampling given an arbitrary sampling function.
- Parameters:
sample0 (ndarray) – Sample to start with. Passed by ref and changed.
E0 (ndarray) – Initial energy of state.
rng (np.random.RandomState) – Random number generator.
flip_site (int) – Site to flip.
calc_e (function) – If another function to calculate energy should be used
- Returns:
delta energy.
- Return type:
float
- coniii.samplers.sample_ising(multipliers, n_samples, seed=None, parallel=True, generate_sample_kw={})[source]
Easy way to Metropolis sample from Ising model.
- Parameters:
multipliers (ndarray or twople) – N individual fields followed by N(N-1)/2 pairwise couplings.
n_samples (int) – Number of samples to take.
seed (int, None) – Random number generator seed.
parallel (bool, True) – If True, use parallelized sampling routine.
generate_sample_kw (dict, {}) –
- Any extra arguments to send into Metropolis sampler. Default args are
n_iters=1000 systematic_iter=False saveHistory=False initial_sample=None
- Returns:
Matrix of dimensions (n_samples, n).
- Return type:
ndarray