coniii.utils module
Utility functions for the coniii package.
This subpackage was split out of a single monolithic utils.py to
keep concerns separate while preserving the public API: every name in
__all__ is still importable as coniii.utils.<name> and via
from coniii.utils import *.
Internal modules (leading underscore — do not import directly):
coniii.utils._indexingCoordinate <-> index helpers, state generation (
bin_states,xbin_states,xpotts_states), base representation, matrix-vector shape helpers (vec2mat,mat2vec,replace_diag,zero_diag).coniii.utils._correlationspair_corr,k_corr,calc_de,calc_overlap,convert_corr,state_probs.coniii.utils._paramsconvert_params,ising_convert_params,split_concat_params, and thedefine_*_helper_functionsfactories.coniii.utils._graphadj,adj_sym,coarse_grain_with_func.
- coniii.utils.sub_to_ind(n, i, j)[source]
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.
- Parameters:
n (int) – Dimension of square array
i (int) – coordinates
j (int) – coordinates
- Return type:
int
- coniii.utils.ind_to_sub(n, ix)[source]
Convert index from flattened upper triangular matrix to pair subindex.
- Parameters:
n (int) – Dimension size of square array.
ix (int) – Index to convert.
- Returns:
subix – (i,j)
- Return type:
tuple
- coniii.utils.bin_states(n, sym=False)[source]
Generate all possible binary spin states.
- Parameters:
n (int) – Number of spins.
sym (bool) – If true, return states in {-1,1} basis.
- Returns:
v
- Return type:
ndarray
- coniii.utils.xbin_states(n, sym=False)[source]
Generator for iterating through all possible binary states.
- Parameters:
n (int) – Number of spins.
sym (bool) – If true, return states in {-1,1} basis.
- Return type:
generator
- coniii.utils.xpotts_states(n, k)[source]
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.
- Parameters:
n (int) – Number of spins.
k (int) – Number of distinct states. These are labeled by integers starting from 0 and must be <=36.
- Return type:
generator
- coniii.utils.base_repr(i, base)[source]
Return decimal number in given base as list.
- Parameters:
i (int)
base (int)
- Return type:
list
- coniii.utils.unique_rows(mat, return_inverse=False)[source]
Return unique rows indices of a numeric numpy array.
- Parameters:
mat (ndarray)
return_inverse (bool) – If True, return inverse that returns back indices of unique array that would return the original array
- Returns:
u (ndarray) – Unique elements of matrix.
idx (ndarray) – row indices of given mat that will give unique array
- coniii.utils.vec2mat(multipliers, separate_fields=False)[source]
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.
- Parameters:
multipliers (ndarray) – Vector of fields and couplings.
separate_fields (bool, False)
- Returns:
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.mat2vec(multipliers)[source]
Convert matrix form of Ising parameters to a vector.
This is specific to the Ising model.
- Parameters:
multipliers (ndarray) – Matrix of couplings with diagonal elements as fields.
- Returns:
Vector of fields and couplings, respectively.
- Return type:
ndarray
- coniii.utils.pair_corr(X, weights=None, concat=False, exclude_empty=False, subtract_mean=False, laplace_count=False)[source]
Calculate averages and pairwise correlations of spins.
- Parameters:
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
- Returns:
(si,sisj) or np.concatenate((si,sisj))
- Return type:
twople
- coniii.utils.k_corr(X, k, weights=None, exclude_empty=False)[source]
Calculate kth order correlations of spins.
- Parameters:
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.
- Returns:
Kth order correlations <s_{i_1} * s_{i_2} * … * s_{i_k}>.
- Return type:
ndarray
- coniii.utils.calc_de(s, i)[source]
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}.
- Parameters:
s (ndarray) – Two-dimensional vector of spins where each row is a state.
i (int)
- Returns:
dE – Derivative of hamiltonian with respect to ith parameter, i.e. the corresponding observable.
- Return type:
float
- coniii.utils.calc_overlap(sample, ignore_zeros=False)[source]
<si_a si_b> between all pairs of replicas a and b
Params:
sample ignore_zeros (bool=False)
Instead of normalizing by the number of spins, normalize by the minimum number of nonzero spins.
- coniii.utils.convert_corr(si, sisj, convert_to, concat=False)[source]
Convert single spin means and pairwise correlations between {0,1} and {-1,1} formulations.
- Parameters:
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.
- Returns:
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.state_probs(v, allstates=None, weights=None, normalized=True)[source]
Get probability of unique states. There is an option to allow for weighted counting.
- Parameters:
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.
- Returns:
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.convert_params(h, J, convert_to, concat=False)[source]
Convert Ising model fields and couplings from {0,1} basis to {-1,1} and vice versa.
- Parameters:
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.
- Returns:
ndarray – Mean bias h vector. Concatenated vector of h and J if concat is True.
ndarray, optional – Vector of J.
- coniii.utils.ising_convert_params(oparams, convert_to, concat=False)[source]
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.
- Parameters:
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)
- Returns:
params – New parameters in order of lowest to highest order interactions to mean biases. Can all be concatenated together if concat switch is True.
- Return type:
tuple of lists or list
- coniii.utils.split_concat_params(p, n)[source]
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.
- Parameters:
p (list-like)
- Returns:
Parameters increasing in order: (h, Jij, Kijk, … ).
- Return type:
list of list-like
- coniii.utils.define_ising_helper_functions()[source]
Functions for plugging into solvers for +/-1 Ising model with fields h_i and couplings J_ij.
- Returns:
function – calc_e
function – calc_observables
function – mch_approximation
- coniii.utils.define_ising_helper_functions_sym()[source]
Functions for plugging into solvers for +/-1 Ising model with couplings J_ij and no fields.
- Returns:
function – calc_e
function – calc_observables
function – mch_approximation
- coniii.utils.define_potts_helper_functions(k)[source]
Helper functions for calculating quantities in k-state Potts model.
- Parameters:
k (int) – Number of possible states.
- Returns:
function – calc_e
function – calc_observables
function – mch_approximation
- coniii.utils.define_pseudo_ising_helper_functions(N)[source]
Define helper functions for using Pseudo method on Ising model.
- Parameters:
N (int) – System size.
- Returns:
function – get_multipliers_r
function – calc_observables_r
- coniii.utils.define_pseudo_potts_helper_functions(n, k)[source]
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.
- Parameters:
n (int) – System size.
k (int) – Number of possible configurations in Potts model.
- Returns:
function – get_multipliers_r
function – calc_observables_r
- coniii.utils.adj(s, n_random_neighbors=0)[source]
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.
- Parameters:
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.
- Returns:
neighbors – Each row is a neighbor. s.size + n_random_neighbors are returned.
- Return type:
ndarray
- coniii.utils.adj_sym(s, n_random_neighbors=False)[source]
Symmetric version of adj() where spins are in {-1,1}.
- coniii.utils.replace_diag(mat, newdiag)[source]
Replace diagonal entries of square matrix.
- Parameters:
mat (ndarray)
newdiag (ndarray)
- Return type:
ndarray
- coniii.utils.zero_diag(mat)[source]
Replace diagonal entries of square matrix with zeros.
- Parameters:
mat (ndarray)
- Return type:
ndarray
- coniii.utils.coarse_grain_with_func(X, n_times, sim_func, coarse_func)[source]
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.
- Parameters:
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.
- Returns:
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.logsumexp(a, axis=None, b=None, keepdims=False, return_sign=False)[source]
Compute the log of the sum of exponentials of input elements.
- Parameters:
a (array_like) – Input array.
axis (None or int or tuple of ints, optional) –
Axis or axes over which the sum is taken. By default axis is None, and all elements are summed.
Added in version 0.11.0.
b (array-like, optional) –
Scaling factor for exp(a) must be of the same shape as a or broadcastable to a. These values may be negative in order to implement subtraction.
Added in version 0.12.0.
keepdims (bool, optional) –
If this is set to True, the axes which are reduced are left in the result as dimensions with size one. With this option, the result will broadcast correctly against the original array.
Added in version 0.15.0.
return_sign (bool, optional) –
If this is set to True, the result will be a pair containing sign information; if False, results that are negative will be returned as NaN. Default is False (no sign information).
Added in version 0.16.0.
- Returns:
res (ndarray) – The result,
np.log(np.sum(np.exp(a)))calculated in a numerically more stable way. If b is given thennp.log(np.sum(b*np.exp(a)))is returned. Ifreturn_signis True,rescontains the log of the absolute value of the argument.sgn (ndarray) – If
return_signis True, this will be an array of floating-point numbers matching res containing +1, 0, -1 (for real-valued inputs) or a complex phase (for complex inputs). This gives the sign of the argument of the logarithm inres. Ifreturn_signis False, only one result is returned.
See also
numpy.logaddexp,numpy.logaddexp2Notes
NumPy has a logaddexp function which is very similar to logsumexp, but only handles two arguments. logaddexp.reduce is similar to this function, but may be less stable.
The logarithm is a multivalued function: for each \(x\) there is an infinite number of \(z\) such that \(exp(z) = x\). The convention is to return the \(z\) whose imaginary part lies in \((-pi, pi]\).
Examples
>>> import numpy as np >>> from scipy.special import logsumexp >>> a = np.arange(10) >>> logsumexp(a) 9.4586297444267107 >>> np.log(np.sum(np.exp(a))) 9.4586297444267107
With weights
>>> a = np.arange(10) >>> b = np.arange(10, 0, -1) >>> logsumexp(a, b=b) 9.9170178533034665 >>> np.log(np.sum(b*np.exp(a))) 9.9170178533034647
Returning a sign flag
>>> logsumexp([1,2],b=[1,-1],return_sign=True) (1.5413248546129181, -1.0)
Notice that logsumexp does not directly support masked arrays. To use it on a masked array, convert the mask into zero weights:
>>> a = np.ma.array([np.log(2), 2, np.log(3)], ... mask=[False, True, False]) >>> b = (~a.mask).astype(int) >>> logsumexp(a.data, b=b), np.log(5) 1.6094379124341005, 1.6094379124341005