coniii.utils module

coniii.utils.adj

Return one-flip neighbors and a set of random neighbors. This is written to be used with the solvers.MPF class. Use adj_sym() if symmetric spins in {-1,1} are needed.

NOTE: For random neighbors, there is no check to make sure neighbors don’t repeat but this shouldn’t be a problem as long as state space is large enough.

s : ndarray
State whose neighbors are found. One-dimensional vector of spins.
n_random_neighbors : int,0
If >0, return this many random neighbors. Neighbors are just random states, but they are called “neighbors” because of the terminology in MPF. They can provide coupling from s to states that are very different, increasing the equilibration rate.
neighbors : ndarray
Each row is a neighbor. s.size + n_random_neighbors are returned.
coniii.utils.adj_sym

Symmetric version of adj() where spins are in {-1,1}.

coniii.utils.base_repr

Return decimal number in given base as list.

i : int base : int

list

coniii.utils.bin_states(n, sym=False)

Generate all possible binary spin states.

n : int
Number of spins.
sym : bool
If true, return states in {-1,1} basis.

v : ndarray

coniii.utils.calc_de(s, i)

Calculate the derivative of the energy wrt parameters given the state and index of the parameter. In this case, the parameters are the concatenated vector of {h_i,J_ij}.

s : ndarray
Two-dimensional vector of spins where each row is a state.

i : int

dE : float
Derivative of hamiltonian with respect to ith parameter, i.e. the corresponding observable.
coniii.utils.calc_overlap(sample, ignore_zeros=False)

<si_a si_b> between all pairs of replicas a and b

sample ignore_zeros (bool=False)

Instead of normalizing by the number of spins, normalize by the minimum number of nonzero spins.
coniii.utils.coarse_grain_with_func(X, n_times, sim_func, coarse_func)

Iteratively coarse-grain X by combining pairs with the highest similarity. Both the function to measure similarity and to implement the coarse-graining must be supplied.

X : ndarray
Each col is a variable and each row is an observation (n_samples, n_system).
n_times : int
Number of times to coarse grain.
sim_func : function
Takes an array like X and returns a vector of ncol*(ncol-1)//2 pairwise similarities.
coarse_func : function
Takes a two col array and returns a single vector.
ndarray
Coarse-grained version of X.
list of lists of ints
Each list specifies which columns of X have been coarse-grained into each col of the coarse X.
coniii.utils.convert_corr(si, sisj, convert_to, concat=False, **kwargs)

Convert single spin means and pairwise correlations between {0,1} and {-1,1} formulations.

si : ndarray
Individual means.
sisj : ndarray
Pairwise correlations.
convert_to : str
‘11’ will convert {0,1} formulation to +/-1 and ‘01’ will convert +/-1 formulation to {0,1}
concat : bool, False
If True, return concatenation of means and pairwise correlations.
ndarray
Averages <si>. Converted to appropriate basis. Returns concatenated vector <si> and <sisj> if concat is True.
ndarray, optional
Pairwise correlations <si*sj>. Converted to appropriate basis.
coniii.utils.convert_params(h, J, convert_to, concat=False)

Convert Ising model fields and couplings from {0,1} basis to {-1,1} and vice versa.

h : ndarray
Fields.
J : ndarray
Couplings.
convert_to : str
Either ‘01’ or ‘11’.
concat : bool, False
If True, return a vector concatenating fields and couplings.
ndarray
Mean bias h vector. Concatenated vector of h and J if concat is True.
ndarray, optional
Vector of J.
coniii.utils.define_ising_helper_functions()

Functions for plugging into solvers for +/-1 Ising model with fields h_i and couplings J_ij.

function
calc_e
function
calc_observables
function
mch_approximation
coniii.utils.define_ising_helper_functions_sym()

Functions for plugging into solvers for +/-1 Ising model with couplings J_ij and no fields.

function
calc_e
function
calc_observables
function
mch_approximation
coniii.utils.define_potts_helper_functions(k)

Helper functions for calculating quantities in k-state Potts model.

k : int
Number of possible states.
function
calc_e
function
calc_observables
function
mch_approximation
coniii.utils.define_pseudo_ising_helper_functions(N)

Define helper functions for using Pseudo method on Ising model.

N : int
System size.
function
get_multipliers_r
function
calc_observables_r
coniii.utils.define_pseudo_potts_helper_functions(n, k)

Define helper functions for using Pseudo method on Potts model with simple form for couplings that are only nonzero when the spins are occupying the same state.

n : int
System size.
k : int
Number of possible configurations in Potts model.
function
get_multipliers_r
function
calc_observables_r
coniii.utils.define_ternary_helper_functions()
coniii.utils.define_triplet_helper_functions()
coniii.utils.ind_to_sub

Convert index from flattened upper triangular matrix to pair subindex.

n : int
Dimension size of square array.
ix : int
Index to convert.
subix : tuple
(i,j)
coniii.utils.ising_convert_params(oparams, convert_to, concat=False)

General conversion of parameters from 01 to 11 basis.

Take set of Ising model parameters up to nth order interactions in either {0,1} or {-1,1} basis and convert to other basis.

oparams : tuple of lists
Tuple of lists of interactions between spins starting with the lowest order interactions. Each list should consist of all interactions of that order such that the length of each list should be binomial(n,i) for all i starting with i>=1.

convert_to : str concat : bool,False

params : tuple of lists or list
New parameters in order of lowest to highest order interactions to mean biases. Can all be concatenated together if concat switch is True.
coniii.utils.k_corr(X, k, weights=None, exclude_empty=False)

Calculate kth order correlations of spins.

X : ndarray
Dimensions (n_samples, n_dim).
k : int
Order of correlation function <s_{i_1} * s_{i_2} * … * s_{i_k}>.
weights : np.ndarray, None :
Calculate single and pairwise means given fractional weights for each state in the data such that a state only appears with some weight, typically less than one.
exclude_empty : bool, False
When using with {-1,1} basis, you can leave entries with 0 and those will not be counted for any pair. If True, the weights option doesn’t do anything.
ndarray
Kth order correlations <s_{i_1} * s_{i_2} * … * s_{i_k}>.
coniii.utils.mat2vec(multipliers)

Convert matrix form of Ising parameters to a vector.

This is specific to the Ising model.

multipliers : ndarray
Matrix of couplings with diagonal elements as fields.
ndarray
Vector of fields and couplings, respectively.
coniii.utils.multinomial(*args)
coniii.utils.pair_corr(X, weights=None, concat=False, exclude_empty=False, subtract_mean=False, laplace_count=False)

Calculate averages and pairwise correlations of spins.

X : ndarray
Dimensions (n_samples,n_dim).
weights : float or np.ndarray or twople, None
If an array is passed, it must be the length of the data and each data point will be given the corresponding weight. Otherwise, the two element tuple should contain the normalization for each mean and each pairwise correlation, in that order. In other words, the first array should be length {s_i} and the second length {si*s_j}.
concat : bool, False
Return means concatenated with the pairwise correlations into one array.
exclude_empty : bool, False
When using with {-1,1} basis, you can leave entries with 0 and those will not be counted for any pair. If True, the weights option doesn’t do anything.
subtract_mean : bool, False
If True, return pairwise correlations with product of individual means subtracted.

laplace_count :

twople
(si,sisj) or np.concatenate((si,sisj))
coniii.utils.replace_diag(mat, newdiag)

Replace diagonal entries of square matrix.

mat : ndarray newdiag : ndarray

ndarray

coniii.utils.split_concat_params(p, n)

Split parameters for Ising model that have all been concatenated together into a single list into separate lists. Assumes that the parameters are increasing in order of interaction and that all parameters are present.

p : list-like

list of list-like
Parameters increasing in order: (h, Jij, Kijk, … ).
coniii.utils.state_probs(v, allstates=None, weights=None, normalized=True)

Get probability of unique states. There is an option to allow for weighted counting.

states : ndarray
Sample of states on which to extract probabilities of unique configurations with dimensions (n_samples,n_dimension).
allstates : ndarray, None
Unique configurations to look for with dimensions (n_samples, n_dimension).
weights : vector, None
For weighted counting of each state given in allstate kwarg.
normalized : bool, True
If True, return probability distribution instead of frequency count.
ndarray
Vector of the probabilities of each state.
ndarray
All unique states found in the data. Each state is a row. Only returned if allstates kwarg is not provided.
coniii.utils.sub_to_ind

Convert pair of coordinates of a symmetric square array into consecutive index of flattened upper triangle. This is slimmed down so it won’t throw errors like if i>n or j>n or if they’re negative. Only checking for if the returned index is negative which could be problematic with wrapped indices.

n : int
Dimension of square array
i,j : int
coordinates

int

coniii.utils.unique_rows(mat, return_inverse=False)

Return unique rows indices of a numeric numpy array.

mat : ndarray return_inverse : bool

If True, return inverse that returns back indices of unique array that would return the original array
u : ndarray
Unique elements of matrix.
idx : ndarray
row indices of given mat that will give unique array
coniii.utils.unravel_index(ijk, n)

Unravel multi-dimensional index to flattened index but specifically for multi-dimensional analog of an upper triangular array (lower triangle indices are not counted).

ijk : tuple
Raveled index to unravel.
n : int
System size.
ix : int
Unraveled index.
coniii.utils.vec2mat(multipliers, separate_fields=False)

Convert vector of parameters containing fields and couplings to a matrix where the diagonal elements are the fields and the remaining elements are the couplings. Fields can be returned separately with the separate_fields keyword argument.

This is specific to the Ising model.

multipliers : ndarray
Vector of fields and couplings.

separate_fields : bool, False

ndarray
n x n matrix. Diagonal elements are fields unless separate_fields keyword argument is True, in which case the diagonal elements are 0.
ndarray (optional)
Fields if separate_fields keyword argument is True.
coniii.utils.xbin_states(n, sym=False)

Generator for iterating through all possible binary states.

n : int
Number of spins.
sym : bool
If true, return states in {-1,1} basis.

generator

coniii.utils.xpotts_states

Generator for iterating through all states for Potts model with k distinct states. This is a faster version of calling xbin_states(n, False) except with strings returned as elements instead of integers.

n : int
Number of spins.
k : int
Number of distinct states. These are labeled by integers starting from 0 and must be <=36.

generator

coniii.utils.zero_diag(mat)

Replace diagonal entries of square matrix with zeros.

mat : ndarray

ndarray