algorithms package

Submodules

algorithms.contagion module

algorithms.contagion.Gillespie_SIR(H, tau, gamma, transmission_function=<function threshold>, initial_infecteds=None, initial_recovereds=None, rho=None, tmin=0, tmax=inf, **args)[source]

A continuous-time SIR model for hypergraphs similar to the model in “The effect of heterogeneity on hypergraph contagion models” by Landry and Restrepo https://doi.org/10.1063/5.0020034 and implemented for networks in the EoN package by Joel C. Miller https://epidemicsonnetworks.readthedocs.io/en/latest/

Parameters:
  • H (HyperNetX Hypergraph object) –

  • tau (dictionary) – Edge sizes as keys (must account for all edge sizes present) and rates of infection for each size (float)

  • gamma (float) – The healing rate

  • transmission_function (lambda function, default: threshold) – A lambda function that has required arguments (node, status, edge) and optional arguments

  • initial_infecteds (list or numpy array, default: None) – Iterable of initially infected node uids

  • initial_recovereds (list or numpy array, default: None) – An iterable of initially recovered node uids

  • rho (float from 0 to 1, default: None) – The fraction of initially infected individuals. Both rho and initially infected cannot be specified.

  • tmin (float, default: 0) – Time at the start of the simulation

  • tmax (float, default: float('Inf')) – Time at which the simulation should be terminated if it hasn’t already.

  • return_full_data (bool, default: False) – This returns all the infection and recovery events at each time if True.

  • **args (Optional arguments to transmission function) – This allows user-defined transmission functions with extra parameters.

Returns:

t, S, I, R – time (t), number of susceptible (S), infected (I), and recovered (R) at each time.

Return type:

numpy arrays

Notes

Example:

>>> import hypernetx.algorithms.contagion as contagion
>>> import random
>>> import hypernetx as hnx
>>> n = 1000
>>> m = 10000
>>> hyperedgeList = [random.sample(range(n), k=random.choice([2,3])) for i in range(m)]
>>> H = hnx.Hypergraph(hyperedgeList)
>>> tau = {2:0.1, 3:0.1}
>>> gamma = 0.1
>>> tmax = 100
>>> t, S, I, R = contagion.Gillespie_SIR(H, tau, gamma, rho=0.1, tmin=0, tmax=tmax)
algorithms.contagion.Gillespie_SIS(H, tau, gamma, transmission_function=<function threshold>, initial_infecteds=None, rho=None, tmin=0, tmax=inf, return_full_data=False, sim_kwargs=None, **args)[source]

A continuous-time SIS model for hypergraphs similar to the model in “The effect of heterogeneity on hypergraph contagion models” by Landry and Restrepo https://doi.org/10.1063/5.0020034 and implemented for networks in the EoN package by Joel C. Miller https://epidemicsonnetworks.readthedocs.io/en/latest/

Parameters:
  • H (HyperNetX Hypergraph object) –

  • tau (dictionary) – Edge sizes as keys (must account for all edge sizes present) and rates of infection for each size (float)

  • gamma (float) – The healing rate

  • transmission_function (lambda function, default: threshold) – A lambda function that has required arguments (node, status, edge) and optional arguments

  • initial_infecteds (list or numpy array, default: None) – Iterable of initially infected node uids

  • rho (float from 0 to 1, default: None) – The fraction of initially infected individuals. Both rho and initially infected cannot be specified.

  • tmin (float, default: 0) – Time at the start of the simulation

  • tmax (float, default: 100) – Time at which the simulation should be terminated if it hasn’t already.

  • return_full_data (bool, default: False) – This returns all the infection and recovery events at each time if True.

  • **args (Optional arguments to transmission function) – This allows user-defined transmission functions with extra parameters.

Returns:

t, S, I – time (t), number of susceptible (S), and infected (I) at each time.

Return type:

numpy arrays

Notes

Example:

>>> import hypernetx.algorithms.contagion as contagion
>>> import random
>>> import hypernetx as hnx
>>> n = 1000
>>> m = 10000
>>> hyperedgeList = [random.sample(range(n), k=random.choice([2,3])) for i in range(m)]
>>> H = hnx.Hypergraph(hyperedgeList)
>>> tau = {2:0.1, 3:0.1}
>>> gamma = 0.1
>>> tmax = 100
>>> t, S, I = contagion.Gillespie_SIS(H, tau, gamma, rho=0.1, tmin=0, tmax=tmax)
algorithms.contagion.collective_contagion(node, status, edge)[source]

The collective contagion mechanism described in “The effect of heterogeneity on hypergraph contagion models” by Landry and Restrepo https://doi.org/10.1063/5.0020034

Parameters:
  • node (hashable) – the node uid to infect (If it doesn’t have status “S”, it will automatically return False)

  • status (dictionary) – The nodes are keys and the values are statuses (The infected state denoted with “I”)

  • edge (iterable) – Iterable of node ids (node must be in the edge or it will automatically return False)

Returns:

False if there is no potential to infect and True if there is.

Return type:

bool

Notes

Example:

>>> status = {0:"S", 1:"I", 2:"I", 3:"S", 4:"R"}
>>> collective_contagion(0, status, (0, 1, 2))
    True
>>> collective_contagion(1, status, (0, 1, 2))
    False
>>> collective_contagion(3, status, (0, 1, 2))
    False
algorithms.contagion.contagion_animation(fig, H, transition_events, node_state_color_dict, edge_state_color_dict, node_radius=1, fps=1)[source]

A function to animate discrete-time contagion models for hypergraphs. Currently only supports a circular layout.

Parameters:
  • fig (matplotlib Figure object) –

  • H (HyperNetX Hypergraph object) –

  • transition_events (dictionary) – The dictionary that is output from the discrete_SIS and discrete_SIR functions with return_full_data=True

  • node_state_color_dict (dictionary) – Dictionary which specifies the colors of each node state. All node states must be specified.

  • edge_state_color_dict (dictionary) – Dictionary with keys that are edge states and values which specify the colors of each edge state (can specify an alpha parameter). All edge-dependent transition states must be specified (most common is “I”) and there must be a a default “OFF” setting.

  • node_radius (float, default: 1) – The radius of the nodes to draw

  • fps (int > 0, default: 1) – Frames per second of the animation

Return type:

matplotlib Animation object

Notes

Example:

>>> import hypernetx.algorithms.contagion as contagion
>>> import random
>>> import hypernetx as hnx
>>> import matplotlib.pyplot as plt
>>> from IPython.display import HTML
>>> n = 1000
>>> m = 10000
>>> hyperedgeList = [random.sample(range(n), k=random.choice([2,3])) for i in range(m)]
>>> H = hnx.Hypergraph(hyperedgeList)
>>> tau = {2:0.1, 3:0.1}
>>> gamma = 0.1
>>> tmax = 100
>>> dt = 0.1
>>> transition_events = contagion.discrete_SIS(H, tau, gamma, rho=0.1, tmin=0, tmax=tmax, dt=dt, return_full_data=True)
>>> node_state_color_dict = {"S":"green", "I":"red", "R":"blue"}
>>> edge_state_color_dict = {"S":(0, 1, 0, 0.3), "I":(1, 0, 0, 0.3), "R":(0, 0, 1, 0.3), "OFF": (1, 1, 1, 0)}
>>> fps = 1
>>> fig = plt.figure()
>>> animation = contagion.contagion_animation(fig, H, transition_events, node_state_color_dict, edge_state_color_dict, node_radius=1, fps=fps)
>>> HTML(animation.to_jshtml())
algorithms.contagion.discrete_SIR(H, tau, gamma, transmission_function=<function threshold>, initial_infecteds=None, initial_recovereds=None, rho=None, tmin=0, tmax=inf, dt=1.0, return_full_data=False, **args)[source]

A discrete-time SIR model for hypergraphs similar to the construction described in “The effect of heterogeneity on hypergraph contagion models” by Landry and Restrepo https://doi.org/10.1063/5.0020034 and “Simplicial models of social contagion” by Iacopini et al. https://doi.org/10.1038/s41467-019-10431-6

Parameters:
  • H (HyperNetX Hypergraph object) –

  • tau (dictionary) – Edge sizes as keys (must account for all edge sizes present) and rates of infection for each size (float)

  • gamma (float) – The healing rate

  • transmission_function (lambda function, default: threshold) – A lambda function that has required arguments (node, status, edge) and optional arguments

  • initial_infecteds (list or numpy array, default: None) – Iterable of initially infected node uids

  • initial_recovereds (list or numpy array, default: None) – An iterable of initially recovered node uids

  • rho (float from 0 to 1, default: None) – The fraction of initially infected individuals. Both rho and initially infected cannot be specified.

  • tmin (float, default: 0) – Time at the start of the simulation

  • tmax (float, default: float('Inf')) – Time at which the simulation should be terminated if it hasn’t already.

  • dt (float > 0, default: 1.0) – Step forward in time that the simulation takes at each step.

  • return_full_data (bool, default: False) – This returns all the infection and recovery events at each time if True.

  • **args (Optional arguments to transmission function) – This allows user-defined transmission functions with extra parameters.

Returns:

  • if return_full_data

    dictionary

    Time as the keys and events that happen as the values.

  • else

    t, S, I, Rnumpy arrays

    time (t), number of susceptible (S), infected (I), and recovered (R) at each time.

Notes

Example:

>>> import hypernetx.algorithms.contagion as contagion
>>> import random
>>> import hypernetx as hnx
>>> n = 1000
>>> m = 10000
>>> hyperedgeList = [random.sample(range(n), k=random.choice([2,3])) for i in range(m)]
>>> H = hnx.Hypergraph(hyperedgeList)
>>> tau = {2:0.1, 3:0.1}
>>> gamma = 0.1
>>> tmax = 100
>>> dt = 0.1
>>> t, S, I, R = contagion.discrete_SIR(H, tau, gamma, rho=0.1, tmin=0, tmax=tmax, dt=dt)
algorithms.contagion.discrete_SIS(H, tau, gamma, transmission_function=<function threshold>, initial_infecteds=None, rho=None, tmin=0, tmax=100, dt=1.0, return_full_data=False, **args)[source]

A discrete-time SIS model for hypergraphs as implemented in “The effect of heterogeneity on hypergraph contagion models” by Landry and Restrepo https://doi.org/10.1063/5.0020034 and “Simplicial models of social contagion” by Iacopini et al. https://doi.org/10.1038/s41467-019-10431-6

Parameters:
  • H (HyperNetX Hypergraph object) –

  • tau (dictionary) – Edge sizes as keys (must account for all edge sizes present) and rates of infection for each size (float)

  • gamma (float) – The healing rate

  • transmission_function (lambda function, default: threshold) – A lambda function that has required arguments (node, status, edge) and optional arguments

  • initial_infecteds (list or numpy array, default: None) – Iterable of initially infected node uids

  • rho (float from 0 to 1, default: None) – The fraction of initially infected individuals. Both rho and initially infected cannot be specified.

  • tmin (float, default: 0) – Time at the start of the simulation

  • tmax (float, default: 100) – Time at which the simulation should be terminated if it hasn’t already.

  • dt (float > 0, default: 1.0) – Step forward in time that the simulation takes at each step.

  • return_full_data (bool, default: False) – This returns all the infection and recovery events at each time if True.

  • **args (Optional arguments to transmission function) – This allows user-defined transmission functions with extra parameters.

Returns:

  • if return_full_data

    dictionary

    Time as the keys and events that happen as the values.

  • else

    t, S, Inumpy arrays

    time (t), number of susceptible (S), and infected (I) at each time.

Notes

Example:

>>> import hypernetx.algorithms.contagion as contagion
>>> import random
>>> import hypernetx as hnx
>>> n = 1000
>>> m = 10000
>>> hyperedgeList = [random.sample(range(n), k=random.choice([2,3])) for i in range(m)]
>>> H = hnx.Hypergraph(hyperedgeList)
>>> tau = {2:0.1, 3:0.1}
>>> gamma = 0.1
>>> tmax = 100
>>> dt = 0.1
>>> t, S, I = contagion.discrete_SIS(H, tau, gamma, rho=0.1, tmin=0, tmax=tmax, dt=dt)
algorithms.contagion.individual_contagion(node, status, edge)[source]

The individual contagion mechanism described in “The effect of heterogeneity on hypergraph contagion models” by Landry and Restrepo https://doi.org/10.1063/5.0020034

Parameters:
  • node (hashable) – The node uid to infect (If it doesn’t have status “S”, it will automatically return False)

  • status (dictionary) – The nodes are keys and the values are statuses (The infected state denoted with “I”)

  • edge (iterable) – Iterable of node ids (node must be in the edge or it will automatically return False)

Returns:

False if there is no potential to infect and True if there is.

Return type:

bool

Notes

Example:

>>> status = {0:"S", 1:"I", 2:"I", 3:"S", 4:"R"}
>>> individual_contagion(0, status, (0, 1, 3))
    True
>>> individual_contagion(1, status, (0, 1, 2))
    False
>>> collective_contagion(3, status, (0, 3, 4))
    False
algorithms.contagion.majority_vote(node, status, edge)[source]

The majority vote contagion mechanism. If a majority of neighbors are contagious, it is possible for an individual to change their opinion. If opinions are divided equally, choose randomly.

Parameters:
  • node (hashable) – The node uid to infect (If it doesn’t have status “S”, it will automatically return False)

  • status (dictionary) – The nodes are keys and the values are statuses (The infected state denoted with “I”)

  • edge (iterable) – Iterable of node ids (node must be in the edge or it will automatically return False

Returns:

False if there is no potential to infect and True if there is.

Return type:

bool

Notes

Example:

>>> status = {0:"S", 1:"I", 2:"I", 3:"S", 4:"R"}
>>> majority_vote(0, status, (0, 1, 2))
    True
>>> majority_vote(0, status, (0, 1, 2, 3))
    True
>>> majority_vote(1, status, (0, 1, 2))
    False
>>> majority_vote(3, status, (0, 1, 2))
    False
algorithms.contagion.threshold(node, status, edge, tau=0.1)[source]

The threshold contagion mechanism

Parameters:
  • node (hashable) – The node uid to infect (If it doesn’t have status “S”, it will automatically return False)

  • status (dictionary) – The nodes are keys and the values are statuses (The infected state denoted with “I”)

  • edge (iterable) – Iterable of node ids (node must be in the edge or it will automatically return False)

  • tau (float between 0 and 1, default: 0.1) – The fraction of nodes in an edge that must be infected for the edge to be able to transmit to the node

Returns:

False if there is no potential to infect and True if there is.

Return type:

bool

Notes

Example:

>>> status = {0:"S", 1:"I", 2:"I", 3:"S", 4:"R"}
>>> threshold(0, status, (0, 2, 3, 4), tau=0.2)
    True
>>> threshold(0, status, (0, 2, 3, 4), tau=0.5)
    False
>>> threshold(3, status, (1, 2, 3), tau=1)
    False

algorithms.generative_models module

algorithms.generative_models.chung_lu_hypergraph(k1, k2)[source]

A function to generate an extension of Chung-Lu hypergraph as implemented by Mirah Shi and described for bipartite networks by Aksoy et al. in https://doi.org/10.1093/comnet/cnx001

Parameters:
  • k1 (dictionary) – This a dictionary where the keys are node ids and the values are node degrees.

  • k2 (dictionary) – This a dictionary where the keys are edge ids and the values are edge degrees also known as edge sizes.

Return type:

HyperNetX Hypergraph object

Notes

The sums of k1 and k2 should be roughly the same. If they are not the same, this function returns a warning but still runs. The output currently is a static Hypergraph object. Dynamic hypergraphs are not currently supported.

Example:

>>> import hypernetx.algorithms.generative_models as gm
>>> import random
>>> n = 100
>>> k1 = {i : random.randint(1, 100) for i in range(n)}
>>> k2 = {i : sorted(k1.values())[i] for i in range(n)}
>>> H = gm.chung_lu_hypergraph(k1, k2)
algorithms.generative_models.dcsbm_hypergraph(k1, k2, g1, g2, omega)[source]

A function to generate an extension of DCSBM hypergraph as implemented by Mirah Shi and described for bipartite networks by Larremore et al. in https://doi.org/10.1103/PhysRevE.90.012805

Parameters:
  • k1 (dictionary) – This a dictionary where the keys are node ids and the values are node degrees.

  • k2 (dictionary) – This a dictionary where the keys are edge ids and the values are edge degrees also known as edge sizes.

  • g1 (dictionary) – This a dictionary where the keys are node ids and the values are the group ids to which the node belongs. The keys must match the keys of k1.

  • g2 (dictionary) – This a dictionary where the keys are edge ids and the values are the group ids to which the edge belongs. The keys must match the keys of k2.

  • omega (2D numpy array) – This is a matrix with entries which specify the number of edges between a given node community and edge community. The number of rows must match the number of node communities and the number of columns must match the number of edge communities.

Return type:

HyperNetX Hypergraph object

Notes

The sums of k1 and k2 should be the same. If they are not the same, this function returns a warning but still runs. The sum of k1 (and k2) and omega should be the same. If they are not the same, this function returns a warning but still runs and the number of entries in the incidence matrix is determined by the omega matrix.

The output currently is a static Hypergraph object. Dynamic hypergraphs are not currently supported.

Example:

>>> n = 100
>>> k1 = {i : random.randint(1, 100) for i in range(n)}
>>> k2 = {i : sorted(k1.values())[i] for i in range(n)}
>>> g1 = {i : random.choice([0, 1]) for i in range(n)}
>>> g2 = {i : random.choice([0, 1]) for i in range(n)}
>>> omega = np.array([[100, 10], [10, 100]])
>>> H = gm.dcsbm_hypergraph(k1, k2, g1, g2, omega)
algorithms.generative_models.erdos_renyi_hypergraph(n, m, p, node_labels=None, edge_labels=None)[source]

A function to generate an Erdos-Renyi hypergraph as implemented by Mirah Shi and described for bipartite networks by Aksoy et al. in https://doi.org/10.1093/comnet/cnx001

Parameters:
  • n (int) – Number of nodes

  • m (int) – Number of edges

  • p (float) – The probability that a bipartite edge is created

  • node_labels (list, default=None) – Vertex labels

  • edge_labels (list, default=None) – Hyperedge labels

Return type:

HyperNetX Hypergraph object

Example:

>>> import hypernetx.algorithms.generative_models as gm
>>> n = 1000
>>> m = n
>>> p = 0.01
>>> H = gm.erdos_renyi_hypergraph(n, m, p)

algorithms.homology_mod2 module

Homology and Smith Normal Form

The purpose of computing the Homology groups for data generated hypergraphs is to identify data sources that correspond to interesting features in the topology of the hypergraph.

The elements of one of these Homology groups are generated by $k$ dimensional cycles of relationships in the original data that are not bound together by higher order relationships. Ideally, we want the briefest description of these cycles; we want a minimal set of relationships exhibiting interesting cyclic behavior. This minimal set will be a bases for the Homology group.

The cyclic relationships in the data are discovered using a boundary map represented as a matrix. To discover the bases we compute the Smith Normal Form of the boundary map.

Homology Mod2

This module computes the homology groups for data represented as an abstract simplicial complex with chain groups ${C_k}$ and $Z_2$ additions. The boundary matrices are represented as rectangular matrices over $Z_2$. These matrices are diagonalized and represented in Smith Normal Form. The kernel and image bases are computed and the Betti numbers and homology bases are returned.

Methods for obtaining SNF for Z/2Z are based on Ferrario’s work: http://www.dlfer.xyz/post/2016-10-27-smith-normal-form/

algorithms.homology_mod2.add_to_column(M, i, j)[source]

Replaces column i (of M) with logical xor between column i and j

Parameters:
  • M (np.array) – matrix

  • i (int) – index of column being altered

  • j (int) – index of column being added to altered

Returns:

N

Return type:

np.array

algorithms.homology_mod2.add_to_row(M, i, j)[source]

Replaces row i with logical xor between row i and j

Parameters:
  • M (np.array) –

  • i (int) – index of row being altered

  • j (int) – index of row being added to altered

Returns:

N

Return type:

np.array

algorithms.homology_mod2.betti(bd, k=None)[source]

Generate the kth-betti numbers for a chain complex with boundary matrices given by bd

Parameters:
  • bd (dict of k-boundary matrices keyed on dimension of domain) –

  • k (int, list or tuple, optional, default=None) – list must be min value and max value of k values inclusive if None, then all betti numbers for dimensions of existing cells will be computed.

Returns:

betti – Description

Return type:

dict

algorithms.homology_mod2.betti_numbers(h, k=None)[source]

Return the kth betti numbers for the simplicial homology of the ASC associated to h

Parameters:
  • h (hnx.Hypergraph) – Hypergraph to compute the betti numbers from

  • k (int or list, optional, default=None) – list must be min value and max value of k values inclusive if None, then all betti numbers for dimensions of existing cells will be computed.

Returns:

betti – A dictionary of betti numbers keyed by dimension

Return type:

dict

algorithms.homology_mod2.bkMatrix(km1basis, kbasis)[source]

Compute the boundary map from $C_{k-1}$-basis to $C_k$ basis with respect to $Z_2$

Parameters:
  • km1basis (indexable iterable) – Ordered list of $k-1$ dimensional cell

  • kbasis (indexable iterable) – Ordered list of $k$ dimensional cells

Returns:

bk – boundary matrix in $Z_2$ stored as boolean

Return type:

np.array

algorithms.homology_mod2.boundary_group(image_basis)[source]

Returns a csr_matrix with rows corresponding to the elements of the group generated by image basis over $mathbb{Z}_2$

Parameters:

image_basis (numpy.ndarray or scipy.sparse.csr_matrix) – 2d-array of basis elements

Return type:

scipy.sparse.csr_matrix

algorithms.homology_mod2.chain_complex(h, k=None)[source]

Compute the k-chains and k-boundary maps required to compute homology for all values in k

Parameters:
  • h (hnx.Hypergraph) –

  • k (int or list of length 2, optional, default=None) – k must be an integer greater than 0 or a list of length 2 indicating min and max dimensions to be computed. eg. if k = [1,2] then 0,1,2,3-chains and boundary maps for k=1,2,3 will be returned, if None than k = [1,max dimension of edge in h]

Returns:

C, bd – C is a dictionary of lists bd is a dictionary of numpy arrays

Return type:

dict

algorithms.homology_mod2.homology_basis(bd, k=None, boundary=False, **kwargs)[source]

Compute a basis for the kth-simplicial homology group, $H_k$, defined by a chain complex $C$ with boundary maps given by bd $= {k:partial_k }$

Parameters:
  • bd (dict) – dict of boundary matrices on k-chains to k-1 chains keyed on k if krange is a tuple then all boundary matrices k in [krange[0],..,krange[1]] inclusive must be in the dictionary

  • k (int or list of ints, optional, default=None) – k must be a positive integer or a list of 2 integers indicating min and max dimensions to be computed, if none given all homology groups will be computed from available boundary matrices in bd

  • boundary (bool) – option to return a basis for the boundary group from each dimension. Needed to compute the shortest generators in the homology group.

Returns:

  • basis (dict) – dict of generators as 0-1 tuples keyed by dim basis for dimension k will be returned only if bd[k] and bd[k+1] have been provided.

  • im (dict) – dict of boundary group generators keyed by dim

algorithms.homology_mod2.hypergraph_homology_basis(h, k=None, shortest=False, interpreted=True)[source]

Computes the kth-homology groups mod 2 for the ASC associated with the hypergraph h for k in krange inclusive

Parameters:
  • h (hnx.Hypergraph) –

  • k (int or list of length 2, optional, default = None) – k must be an integer greater than 0 or a list of length 2 indicating min and max dimensions to be computed

  • shortest (bool, optional, default=False) – option to look for shortest representative for each coset in the homology group, only good for relatively small examples

  • interpreted (bool, optional, default = True) – if True will return an explicit basis in terms of the k-chains

Returns:

  • basis (list) – list of generators as k-chains as boolean vectors

  • interpreted_basis – lists of kchains in basis

algorithms.homology_mod2.interpret(Ck, arr, labels=None)[source]

Returns the data as represented in Ck associated with the arr

Parameters:
  • Ck (list) – a list of k-cells being referenced by arr

  • arr (np.array) – array of 0-1 vectors

  • labels (dict, optional) – dictionary of labels to associate to the nodes in the cells

Returns:

list of k-cells referenced by data in Ck

Return type:

list

algorithms.homology_mod2.kchainbasis(h, k)[source]

Compute the set of k dimensional cells in the abstract simplicial complex associated with the hypergraph.

Parameters:
  • h (hnx.Hypergraph) –

  • k (int) – dimension of cell

Returns:

an ordered list of kchains represented as tuples of length k+1

Return type:

list

See also

hnx.hypergraph.toplexes

Notes

  • Method works best if h is simple [Berge], i.e. no edge contains another and there are no duplicate edges (toplexes).

  • Hypergraph node uids must be sortable.

algorithms.homology_mod2.logical_dot(ar1, ar2)[source]

Returns the boolean equivalent of the dot product mod 2 on two 1-d arrays of the same length.

Parameters:
  • ar1 (numpy.ndarray) – 1-d array

  • ar2 (numpy.ndarray) – 1-d array

Returns:

boolean value associated with dot product mod 2

Return type:

bool

Raises:

HyperNetXError – If arrays are not of the same length an error will be raised.

algorithms.homology_mod2.logical_matadd(mat1, mat2)[source]

Returns the boolean equivalent of matrix addition mod 2 on two binary arrays stored as type boolean

Parameters:
  • mat1 (np.ndarray) – 2-d array of boolean values

  • mat2 (np.ndarray) – 2-d array of boolean values

Returns:

mat – boolean matrix equivalent to the mod 2 matrix addition of the matrices as matrices over Z/2Z

Return type:

np.ndarray

Raises:

HyperNetXError – If dimensions are not equal an error will be raised.

algorithms.homology_mod2.logical_matmul(mat1, mat2)[source]

Returns the boolean equivalent of matrix multiplication mod 2 on two binary arrays stored as type boolean

Parameters:
  • mat1 (np.ndarray) – 2-d array of boolean values

  • mat2 (np.ndarray) – 2-d array of boolean values

Returns:

mat – boolean matrix equivalent to the mod 2 matrix multiplication of the matrices as matrices over Z/2Z

Return type:

np.ndarray

Raises:

HyperNetXError – If inner dimensions are not equal an error will be raised.

algorithms.homology_mod2.matmulreduce(arr, reverse=False)[source]

Recursively applies a ‘logical multiplication’ to a list of boolean arrays.

For arr = [arr[0],arr[1],arr[2]…arr[n]] returns product arr[0]arr[1]…arr[n] If reverse = True, returns product arr[n]arr[n-1]…arr[0]

Parameters:
  • arr (list of np.array) – list of nxm matrices represented as np.array

  • reverse (bool, optional) – order to multiply the matrices

Returns:

P – Product of matrices in the list

Return type:

np.array

algorithms.homology_mod2.reduced_row_echelon_form_mod2(M)[source]

Computes the invertible transformation matrices needed to compute the reduced row echelon form of M modulo 2

Parameters:

M (np.array) – a rectangular matrix with elements in $Z_2$

Returns:

L, S, Linv – LM = S where S is the reduced echelon form of M and M = LinvS

Return type:

np.arrays

algorithms.homology_mod2.smith_normal_form_mod2(M)[source]

Computes the invertible transformation matrices needed to compute the Smith Normal Form of M modulo 2

Parameters:
  • M (np.array) – a rectangular matrix with data type bool

  • track (bool) – if track=True will print out the transformation as Z/2Z matrix as it discovers L[i] and R[j]

Returns:

L, R, S, Linv – LMR = S is the Smith Normal Form of the matrix M.

Return type:

np.arrays

Note

Given a mxn matrix $M$ with entries in $Z_2$ we start with the equation: $L M R = S$, where $L = I_m$, and $R=I_n$ are identity matrices and $S = M$. We repeatedly apply actions to the left and right side of the equation to transform S into a diagonal matrix. For each action applied to the left side we apply its inverse action to the right side of I_m to generate $L^{-1}$. Finally we verify: $L M R = S$ and $LLinv = I_m$.

algorithms.homology_mod2.swap_columns(i, j, *args)[source]

Swaps ith and jth column of each matrix in args Returns a list of new matrices

Parameters:
  • i (int) –

  • j (int) –

  • args (np.arrays) –

Returns:

list of copies of args with ith and jth row swapped

Return type:

list

algorithms.homology_mod2.swap_rows(i, j, *args)[source]

Swaps ith and jth row of each matrix in args Returns a list of new matrices

Parameters:
  • i (int) –

  • j (int) –

  • args (np.arrays) –

Returns:

list of copies of args with ith and jth row swapped

Return type:

list

algorithms.hypergraph_modularity module

Hypergraph_Modularity

Modularity and clustering for hypergraphs using HyperNetX. Adapted from F. Théberge’s GitHub repository: Hypergraph Clustering See Tutorial 13 in the tutorials folder for library usage.

References

algorithms.hypergraph_modularity.conductance(H, A)[source]

Computes conductance [4] of hypergraph HG with respect to partition A.

Parameters:
  • H (Hypergraph) – The hypergraph

  • A (set) – Partition of the vertices in H

Returns:

The conductance function for partition A on H

Return type:

float

algorithms.hypergraph_modularity.dict2part(D)[source]

Given a dictionary mapping the part for each vertex, return a partition as a list of sets; inverse function to part2dict

Parameters:

D (dict) – Dictionary keyed by vertices with values equal to integer index of the partition the vertex belongs to

Returns:

List of sets; one set for each part in the partition

Return type:

list

algorithms.hypergraph_modularity.kumar(HG, delta=0.01, verbose=False)[source]

Compute a partition of the vertices in hypergraph HG as per Kumar’s algorithm [1]

Parameters:
  • HG (Hypergraph) –

  • delta (float, optional) – convergence stopping criterion

Returns:

A partition of the vertices in HG

Return type:

list of sets

algorithms.hypergraph_modularity.last_step(HG, A, wdc=<function linear>, delta=0.01, verbose=False)[source]

Given some initial partition L, compute a new partition of the vertices in HG as per Last-Step algorithm [2]

Note

This is a very simple algorithm that tries moving nodes between communities to improve hypergraph modularity. It requires an initial non-trivial partition which can be obtained for example via graph clustering on the 2-section of HG, or via Kumar’s algorithm.

Parameters:
  • HG (Hypergraph) –

  • L (list of sets) – some initial partition of the vertices in HG

  • wdc (func, optional) – Hyperparameter for hypergraph modularity [2]

  • delta (float, optional) – convergence stopping criterion

  • verbose (boolean, optional) – If set, also returns progress after each pass through the vertices

Returns:

A new partition for the vertices in HG

Return type:

list of sets

algorithms.hypergraph_modularity.linear(d, c)[source]

Hyperparameter for hypergraph modularity [2] for d-edge with c vertices in the majority class. This is the default choice for modularity() and last_step() functions.

Parameters:
  • d (int) – Number of vertices in an edge

  • c (int) – Number of vertices in the majority class

Returns:

c/d if c>d/2 else 0

Return type:

float

algorithms.hypergraph_modularity.majority(d, c)[source]

Hyperparameter for hypergraph modularity [2] for d-edge with c vertices in the majority class. This corresponds to the majority rule [3]

Parameters:
  • d (int) – Number of vertices in an edge

  • c (int) – Number of vertices in the majority class

Returns:

1 if c>d/2 else 0

Return type:

bool

algorithms.hypergraph_modularity.modularity(HG, A, wdc=<function linear>)[source]

Computes modularity of hypergraph HG with respect to partition A.

Parameters:
  • HG (Hypergraph) – The hypergraph with some precomputed attributes via: precompute_attributes(HG)

  • A (list of sets) – Partition of the vertices in HG

  • wdc (func, optional) – Hyperparameter for hypergraph modularity [2]

Note

For ‘wdc’, any function of the format w(d,c) that returns 0 when c <= d/2 and value in [0,1] otherwise can be used. Default is ‘linear’; other supplied choices are ‘majority’ and ‘strict’.

Returns:

The modularity function for partition A on HG

Return type:

float

algorithms.hypergraph_modularity.part2dict(A)[source]

Given a partition (list of sets), returns a dictionary mapping the part for each vertex; inverse function to dict2part

Parameters:

A (list of sets) – a partition of the vertices

Returns:

a dictionary with {vertex: partition index}

Return type:

dict

algorithms.hypergraph_modularity.strict(d, c)[source]

Hyperparameter for hypergraph modularity [2] for d-edge with c vertices in the majority class. This corresponds to the strict rule [3]

Parameters:
  • d (int) – Number of vertices in an edge

  • c (int) – Number of vertices in the majority class

Returns:

1 if c==d else 0

Return type:

bool

algorithms.hypergraph_modularity.two_section(HG)[source]

Creates a random walk based [1] 2-section igraph Graph with transition weights defined by the weights of the hyperedges.

Parameters:

HG (Hypergraph) –

Returns:

The 2-section graph built from HG

Return type:

igraph.Graph

algorithms.laplacians_clustering module

Hypergraph Probability Transition Matrices, Laplacians, and Clustering

We contruct hypergraph random walks utilizing optional “edge-dependent vertex weights”, which are weights associated with each vertex-hyperedge pair (i.e. cell weights on the incidence matrix). The probability transition matrix of this random walk is used to construct a normalized Laplacian matrix for the hypergraph. That normalized Laplacian then serves as the input for a spectral clustering algorithm. This spectral clustering algorithm, as well as the normalized Laplacian and other details of this methodology are described in

K. Hayashi, S. Aksoy, C. Park, H. Park, “Hypergraph random walks, Laplacians, and clustering”, Proceedings of the 29th ACM International Conference on Information & Knowledge Management. 2020. https://doi.org/10.1145/3340531.3412034

Please direct any inquiries concerning the clustering module to Sinan Aksoy, sinan.aksoy@pnnl.gov

algorithms.laplacians_clustering.get_pi(P)[source]

Returns the eigenvector corresponding to the largest eigenvalue (in magnitude), normalized so its entries sum to 1. Intended for the probability transition matrix of a random walk on a (connected) hypergraph, in which case the output can be interpreted as the stationary distribution.

Parameters:

P (csr matrix) – Probability transition matrix

Returns:

pi – Stationary distribution of random walk defined by P

Return type:

numpy.ndarray

algorithms.laplacians_clustering.norm_lap(H, weights=False, index=True)[source]

Normalized Laplacian matrix of the hypergraph. Symmetrizes the probability transition matrix of a hypergraph random walk using the stationary distribution, using the digraph Laplacian defined in:

Chung, Fan. “Laplacians and the Cheeger inequality for directed graphs.” Annals of Combinatorics 9.1 (2005): 1-19.

and studied in the context of hypergraphs in:

Hayashi, K., Aksoy, S. G., Park, C. H., & Park, H. Hypergraph random walks, laplacians, and clustering. In Proceedings of CIKM 2020, (2020): 495-504.

Parameters:
  • H (hnx.Hypergraph) – The hypergraph must be connected, meaning there is a path linking any two vertices

  • weight (bool, optional, default : False) – Uses cell_weights, if False, uniform weights are utilized.

  • index (bool, optional) – Whether to return matrix-index to vertex-label mapping

Returns:

  • P (scipy.sparse.csr.csr_matrix) – Probability transition matrix of the random walk on the hypergraph

  • id (list) – contains list of index of node ids for rows

algorithms.laplacians_clustering.prob_trans(H, weights=False, index=True, check_connected=True)[source]

The probability transition matrix of a random walk on the vertices of a hypergraph. At each step in the walk, the next vertex is chosen by:

  1. Selecting a hyperedge e containing the vertex with probability proportional to w(e)

  2. Selecting a vertex v within e with probability proportional to a gamma(v,e)

If weights are not specified, then all weights are uniform and the walk is equivalent to a simple random walk. If weights are specified, the hyperedge weights w(e) are determined from the weights gamma(v,e).

Parameters:
  • H (hnx.Hypergraph) – The hypergraph must be connected, meaning there is a path linking any two vertices

  • weights (bool, optional, default : False) – Use the cell_weights associated with the hypergraph If False, uniform weights are utilized.

  • index (bool, optional) – Whether to return matrix index to vertex label mapping

Returns:

  • P (scipy.sparse.csr.csr_matrix) – Probability transition matrix of the random walk on the hypergraph

  • index (list) – contains list of index of node ids for rows

algorithms.laplacians_clustering.spec_clus(H, k, existing_lap=None, weights=False)[source]

Hypergraph spectral clustering of the vertex set into k disjoint clusters using the normalized hypergraph Laplacian. Equivalent to the “RDC-Spec” Algorithm 1 in:

Hayashi, K., Aksoy, S. G., Park, C. H., & Park, H. Hypergraph random walks, laplacians, and clustering. In Proceedings of CIKM 2020, (2020): 495-504.

Parameters:
  • H (hnx.Hypergraph) – The hypergraph must be connected, meaning there is a path linking any two vertices

  • k (int) – Number of clusters

  • existing_lap (csr matrix, optional) – Whether to use an existing Laplacian; otherwise, normalized hypergraph Laplacian will be utilized

  • weights (bool, optional) – Use the cell_weights of the hypergraph. If False uniform weights are used.

Returns:

clusters – Vertex cluster dictionary, keyed by integers 0,…,k-1, with lists of vertices as values.

Return type:

dict

algorithms.s_centrality_measures module

S-Centrality Measures

We generalize graph metrics to s-metrics for a hypergraph by using its s-connected components. This is accomplished by computing the s edge-adjacency matrix and constructing the corresponding graph of the matrix. We then use existing graph metrics on this representation of the hypergraph. In essence we construct an s-line graph corresponding to the hypergraph on which to apply our methods.

S-Metrics for hypergraphs are discussed in depth in: Aksoy, S.G., Joslyn, C., Ortiz Marrero, C. et al. Hypernetwork science via high-order hypergraph walks. EPJ Data Sci. 9, 16 (2020). https://doi.org/10.1140/epjds/s13688-020-00231-0

algorithms.s_centrality_measures.s_betweenness_centrality(H, s=1, edges=True, normalized=True, return_singletons=True)[source]

A centrality measure for an s-edge(node) subgraph of H based on shortest paths. Equals the betweenness centrality of vertices in the edge(node) s-linegraph.

In a graph (2-uniform hypergraph) the betweenness centrality of a vertex $v$ is the ratio of the number of non-trivial shortest paths between any pair of vertices in the graph that pass through $v$ divided by the total number of non-trivial shortest paths in the graph.

The centrality of edge to all shortest s-edge paths $V$ = the set of vertices in the linegraph. $sigma(s,t)$ = the number of shortest paths between vertices $s$ and $t$. $sigma(s,t|v)$ = the number of those paths that pass through vertex $v$.

\[c_B(v) = \sum_{s \neq t \in V} \frac{\sigma(s, t|v)}{\sigma(s,t)}\]
Parameters:
  • H (hnx.Hypergraph) –

  • s (int) – s connectedness requirement

  • edges (bool, optional) – determines if edge or node linegraph

  • normalized – bool, default=False, If true the betweenness values are normalized by 2/((n-1)(n-2)), where n is the number of edges in H

  • return_singletons (bool, optional) – if False will ignore singleton components of linegraph

Returns:

A dictionary of s-betweenness centrality value of the edges.

Return type:

dict

algorithms.s_centrality_measures.s_closeness_centrality(H, s=1, edges=True, return_singletons=True, source=None)[source]

In a connected component the reciprocal of the sum of the distance between an edge(node) and all other edges(nodes) in the component times the number of edges(nodes) in the component minus 1.

$V$ = the set of vertices in the linegraph. $n = |V|$ $d$ = shortest path distance

\[C(u) = \frac{n - 1}{\sum_{v \neq u \in V} d(v, u)}\]
Parameters:
  • H (hnx.Hypergraph) –

  • s (int, optional) –

  • edges (bool, optional) – Indicates if method should compute edge linegraph (default) or node linegraph.

  • return_singletons (bool, optional) – Indicates if method should return values for singleton components.

  • source (str, optional) – Identifier of node or edge of interest for computing centrality

Returns:

returns the s-closeness centrality value of the edges(nodes). If source=None a dictionary of values for each s-edge in H is returned. If source then a single value is returned.

Return type:

dict or float

algorithms.s_centrality_measures.s_eccentricity(H, s=1, edges=True, source=None, return_singletons=True)[source]

The length of the longest shortest path from a vertex $u$ to every other vertex in the s-linegraph. $V$ = set of vertices in the s-linegraph $d$ = shortest path distance

\[\text{s-ecc}(u) = \text{max}\{d(u,v): v \in V\}\]
Parameters:
  • H (hnx.Hypergraph) –

  • s (int, optional) –

  • edges (bool, optional) – Indicates if method should compute edge linegraph (default) or node linegraph.

  • return_singletons (bool, optional) – Indicates if method should return values for singleton components.

  • source (str, optional) – Identifier of node or edge of interest for computing centrality

Returns:

returns the s-eccentricity value of the edges(nodes). If source=None a dictionary of values for each s-edge in H is returned. If source then a single value is returned. If the s-linegraph is disconnected, np.inf is returned.

Return type:

dict or float

algorithms.s_centrality_measures.s_harmonic_centrality(H, s=1, edges=True, source=None, normalized=False, return_singletons=True)[source]

A centrality measure for an s-edge subgraph of H. A value equal to 1 means the s-edge intersects every other s-edge in H. All values range between 0 and 1. Edges of size less than s return 0. If H contains only one s-edge a 0 is returned.

The denormalized reciprocal of the harmonic mean of all distances from $u$ to all other vertices. $V$ = the set of vertices in the linegraph. $d$ = shortest path distance

\[C(u) = \sum_{v \neq u \in V} \frac{1}{d(v, u)}\]

Normalized this becomes: $$C(u) = sum_{v neq u in V} frac{1}{d(v, u)}cdotfrac{2}{(n-1)(n-2)}$$ where $n$ is the number vertices.

Parameters:
  • H (hnx.Hypergraph) –

  • s (int, optional) –

  • edges (bool, optional) – Indicates if method should compute edge linegraph (default) or node linegraph.

  • return_singletons (bool, optional) – Indicates if method should return values for singleton components.

  • source (str, optional) – Identifier of node or edge of interest for computing centrality

Returns:

returns the s-harmonic closeness centrality value of the edges, a number between 0 and 1 inclusive. If source=None a dictionary of values for each s-edge in H is returned. If source then a single value is returned.

Return type:

dict or float

algorithms.s_centrality_measures.s_harmonic_closeness_centrality(H, s=1, edge=None)[source]

Module contents

algorithms.Gillespie_SIR(H, tau, gamma, transmission_function=<function threshold>, initial_infecteds=None, initial_recovereds=None, rho=None, tmin=0, tmax=inf, **args)[source]

A continuous-time SIR model for hypergraphs similar to the model in “The effect of heterogeneity on hypergraph contagion models” by Landry and Restrepo https://doi.org/10.1063/5.0020034 and implemented for networks in the EoN package by Joel C. Miller https://epidemicsonnetworks.readthedocs.io/en/latest/

Parameters:
  • H (HyperNetX Hypergraph object) –

  • tau (dictionary) – Edge sizes as keys (must account for all edge sizes present) and rates of infection for each size (float)

  • gamma (float) – The healing rate

  • transmission_function (lambda function, default: threshold) – A lambda function that has required arguments (node, status, edge) and optional arguments

  • initial_infecteds (list or numpy array, default: None) – Iterable of initially infected node uids

  • initial_recovereds (list or numpy array, default: None) – An iterable of initially recovered node uids

  • rho (float from 0 to 1, default: None) – The fraction of initially infected individuals. Both rho and initially infected cannot be specified.

  • tmin (float, default: 0) – Time at the start of the simulation

  • tmax (float, default: float('Inf')) – Time at which the simulation should be terminated if it hasn’t already.

  • return_full_data (bool, default: False) – This returns all the infection and recovery events at each time if True.

  • **args (Optional arguments to transmission function) – This allows user-defined transmission functions with extra parameters.

Returns:

t, S, I, R – time (t), number of susceptible (S), infected (I), and recovered (R) at each time.

Return type:

numpy arrays

Notes

Example:

>>> import hypernetx.algorithms.contagion as contagion
>>> import random
>>> import hypernetx as hnx
>>> n = 1000
>>> m = 10000
>>> hyperedgeList = [random.sample(range(n), k=random.choice([2,3])) for i in range(m)]
>>> H = hnx.Hypergraph(hyperedgeList)
>>> tau = {2:0.1, 3:0.1}
>>> gamma = 0.1
>>> tmax = 100
>>> t, S, I, R = contagion.Gillespie_SIR(H, tau, gamma, rho=0.1, tmin=0, tmax=tmax)
algorithms.Gillespie_SIS(H, tau, gamma, transmission_function=<function threshold>, initial_infecteds=None, rho=None, tmin=0, tmax=inf, return_full_data=False, sim_kwargs=None, **args)[source]

A continuous-time SIS model for hypergraphs similar to the model in “The effect of heterogeneity on hypergraph contagion models” by Landry and Restrepo https://doi.org/10.1063/5.0020034 and implemented for networks in the EoN package by Joel C. Miller https://epidemicsonnetworks.readthedocs.io/en/latest/

Parameters:
  • H (HyperNetX Hypergraph object) –

  • tau (dictionary) – Edge sizes as keys (must account for all edge sizes present) and rates of infection for each size (float)

  • gamma (float) – The healing rate

  • transmission_function (lambda function, default: threshold) – A lambda function that has required arguments (node, status, edge) and optional arguments

  • initial_infecteds (list or numpy array, default: None) – Iterable of initially infected node uids

  • rho (float from 0 to 1, default: None) – The fraction of initially infected individuals. Both rho and initially infected cannot be specified.

  • tmin (float, default: 0) – Time at the start of the simulation

  • tmax (float, default: 100) – Time at which the simulation should be terminated if it hasn’t already.

  • return_full_data (bool, default: False) – This returns all the infection and recovery events at each time if True.

  • **args (Optional arguments to transmission function) – This allows user-defined transmission functions with extra parameters.

Returns:

t, S, I – time (t), number of susceptible (S), and infected (I) at each time.

Return type:

numpy arrays

Notes

Example:

>>> import hypernetx.algorithms.contagion as contagion
>>> import random
>>> import hypernetx as hnx
>>> n = 1000
>>> m = 10000
>>> hyperedgeList = [random.sample(range(n), k=random.choice([2,3])) for i in range(m)]
>>> H = hnx.Hypergraph(hyperedgeList)
>>> tau = {2:0.1, 3:0.1}
>>> gamma = 0.1
>>> tmax = 100
>>> t, S, I = contagion.Gillespie_SIS(H, tau, gamma, rho=0.1, tmin=0, tmax=tmax)
algorithms.add_to_column(M, i, j)[source]

Replaces column i (of M) with logical xor between column i and j

Parameters:
  • M (np.array) – matrix

  • i (int) – index of column being altered

  • j (int) – index of column being added to altered

Returns:

N

Return type:

np.array

algorithms.add_to_row(M, i, j)[source]

Replaces row i with logical xor between row i and j

Parameters:
  • M (np.array) –

  • i (int) – index of row being altered

  • j (int) – index of row being added to altered

Returns:

N

Return type:

np.array

algorithms.betti(bd, k=None)[source]

Generate the kth-betti numbers for a chain complex with boundary matrices given by bd

Parameters:
  • bd (dict of k-boundary matrices keyed on dimension of domain) –

  • k (int, list or tuple, optional, default=None) – list must be min value and max value of k values inclusive if None, then all betti numbers for dimensions of existing cells will be computed.

Returns:

betti – Description

Return type:

dict

algorithms.betti_numbers(h, k=None)[source]

Return the kth betti numbers for the simplicial homology of the ASC associated to h

Parameters:
  • h (hnx.Hypergraph) – Hypergraph to compute the betti numbers from

  • k (int or list, optional, default=None) – list must be min value and max value of k values inclusive if None, then all betti numbers for dimensions of existing cells will be computed.

Returns:

betti – A dictionary of betti numbers keyed by dimension

Return type:

dict

algorithms.bkMatrix(km1basis, kbasis)[source]

Compute the boundary map from $C_{k-1}$-basis to $C_k$ basis with respect to $Z_2$

Parameters:
  • km1basis (indexable iterable) – Ordered list of $k-1$ dimensional cell

  • kbasis (indexable iterable) – Ordered list of $k$ dimensional cells

Returns:

bk – boundary matrix in $Z_2$ stored as boolean

Return type:

np.array

algorithms.boundary_group(image_basis)[source]

Returns a csr_matrix with rows corresponding to the elements of the group generated by image basis over $mathbb{Z}_2$

Parameters:

image_basis (numpy.ndarray or scipy.sparse.csr_matrix) – 2d-array of basis elements

Return type:

scipy.sparse.csr_matrix

algorithms.chain_complex(h, k=None)[source]

Compute the k-chains and k-boundary maps required to compute homology for all values in k

Parameters:
  • h (hnx.Hypergraph) –

  • k (int or list of length 2, optional, default=None) – k must be an integer greater than 0 or a list of length 2 indicating min and max dimensions to be computed. eg. if k = [1,2] then 0,1,2,3-chains and boundary maps for k=1,2,3 will be returned, if None than k = [1,max dimension of edge in h]

Returns:

C, bd – C is a dictionary of lists bd is a dictionary of numpy arrays

Return type:

dict

algorithms.chung_lu_hypergraph(k1, k2)[source]

A function to generate an extension of Chung-Lu hypergraph as implemented by Mirah Shi and described for bipartite networks by Aksoy et al. in https://doi.org/10.1093/comnet/cnx001

Parameters:
  • k1 (dictionary) – This a dictionary where the keys are node ids and the values are node degrees.

  • k2 (dictionary) – This a dictionary where the keys are edge ids and the values are edge degrees also known as edge sizes.

Return type:

HyperNetX Hypergraph object

Notes

The sums of k1 and k2 should be roughly the same. If they are not the same, this function returns a warning but still runs. The output currently is a static Hypergraph object. Dynamic hypergraphs are not currently supported.

Example:

>>> import hypernetx.algorithms.generative_models as gm
>>> import random
>>> n = 100
>>> k1 = {i : random.randint(1, 100) for i in range(n)}
>>> k2 = {i : sorted(k1.values())[i] for i in range(n)}
>>> H = gm.chung_lu_hypergraph(k1, k2)
algorithms.collective_contagion(node, status, edge)[source]

The collective contagion mechanism described in “The effect of heterogeneity on hypergraph contagion models” by Landry and Restrepo https://doi.org/10.1063/5.0020034

Parameters:
  • node (hashable) – the node uid to infect (If it doesn’t have status “S”, it will automatically return False)

  • status (dictionary) – The nodes are keys and the values are statuses (The infected state denoted with “I”)

  • edge (iterable) – Iterable of node ids (node must be in the edge or it will automatically return False)

Returns:

False if there is no potential to infect and True if there is.

Return type:

bool

Notes

Example:

>>> status = {0:"S", 1:"I", 2:"I", 3:"S", 4:"R"}
>>> collective_contagion(0, status, (0, 1, 2))
    True
>>> collective_contagion(1, status, (0, 1, 2))
    False
>>> collective_contagion(3, status, (0, 1, 2))
    False
algorithms.contagion_animation(fig, H, transition_events, node_state_color_dict, edge_state_color_dict, node_radius=1, fps=1)[source]

A function to animate discrete-time contagion models for hypergraphs. Currently only supports a circular layout.

Parameters:
  • fig (matplotlib Figure object) –

  • H (HyperNetX Hypergraph object) –

  • transition_events (dictionary) – The dictionary that is output from the discrete_SIS and discrete_SIR functions with return_full_data=True

  • node_state_color_dict (dictionary) – Dictionary which specifies the colors of each node state. All node states must be specified.

  • edge_state_color_dict (dictionary) – Dictionary with keys that are edge states and values which specify the colors of each edge state (can specify an alpha parameter). All edge-dependent transition states must be specified (most common is “I”) and there must be a a default “OFF” setting.

  • node_radius (float, default: 1) – The radius of the nodes to draw

  • fps (int > 0, default: 1) – Frames per second of the animation

Return type:

matplotlib Animation object

Notes

Example:

>>> import hypernetx.algorithms.contagion as contagion
>>> import random
>>> import hypernetx as hnx
>>> import matplotlib.pyplot as plt
>>> from IPython.display import HTML
>>> n = 1000
>>> m = 10000
>>> hyperedgeList = [random.sample(range(n), k=random.choice([2,3])) for i in range(m)]
>>> H = hnx.Hypergraph(hyperedgeList)
>>> tau = {2:0.1, 3:0.1}
>>> gamma = 0.1
>>> tmax = 100
>>> dt = 0.1
>>> transition_events = contagion.discrete_SIS(H, tau, gamma, rho=0.1, tmin=0, tmax=tmax, dt=dt, return_full_data=True)
>>> node_state_color_dict = {"S":"green", "I":"red", "R":"blue"}
>>> edge_state_color_dict = {"S":(0, 1, 0, 0.3), "I":(1, 0, 0, 0.3), "R":(0, 0, 1, 0.3), "OFF": (1, 1, 1, 0)}
>>> fps = 1
>>> fig = plt.figure()
>>> animation = contagion.contagion_animation(fig, H, transition_events, node_state_color_dict, edge_state_color_dict, node_radius=1, fps=fps)
>>> HTML(animation.to_jshtml())
algorithms.dcsbm_hypergraph(k1, k2, g1, g2, omega)[source]

A function to generate an extension of DCSBM hypergraph as implemented by Mirah Shi and described for bipartite networks by Larremore et al. in https://doi.org/10.1103/PhysRevE.90.012805

Parameters:
  • k1 (dictionary) – This a dictionary where the keys are node ids and the values are node degrees.

  • k2 (dictionary) – This a dictionary where the keys are edge ids and the values are edge degrees also known as edge sizes.

  • g1 (dictionary) – This a dictionary where the keys are node ids and the values are the group ids to which the node belongs. The keys must match the keys of k1.

  • g2 (dictionary) – This a dictionary where the keys are edge ids and the values are the group ids to which the edge belongs. The keys must match the keys of k2.

  • omega (2D numpy array) – This is a matrix with entries which specify the number of edges between a given node community and edge community. The number of rows must match the number of node communities and the number of columns must match the number of edge communities.

Return type:

HyperNetX Hypergraph object

Notes

The sums of k1 and k2 should be the same. If they are not the same, this function returns a warning but still runs. The sum of k1 (and k2) and omega should be the same. If they are not the same, this function returns a warning but still runs and the number of entries in the incidence matrix is determined by the omega matrix.

The output currently is a static Hypergraph object. Dynamic hypergraphs are not currently supported.

Example:

>>> n = 100
>>> k1 = {i : random.randint(1, 100) for i in range(n)}
>>> k2 = {i : sorted(k1.values())[i] for i in range(n)}
>>> g1 = {i : random.choice([0, 1]) for i in range(n)}
>>> g2 = {i : random.choice([0, 1]) for i in range(n)}
>>> omega = np.array([[100, 10], [10, 100]])
>>> H = gm.dcsbm_hypergraph(k1, k2, g1, g2, omega)
algorithms.dict2part(D)[source]

Given a dictionary mapping the part for each vertex, return a partition as a list of sets; inverse function to part2dict

Parameters:

D (dict) – Dictionary keyed by vertices with values equal to integer index of the partition the vertex belongs to

Returns:

List of sets; one set for each part in the partition

Return type:

list

algorithms.discrete_SIR(H, tau, gamma, transmission_function=<function threshold>, initial_infecteds=None, initial_recovereds=None, rho=None, tmin=0, tmax=inf, dt=1.0, return_full_data=False, **args)[source]

A discrete-time SIR model for hypergraphs similar to the construction described in “The effect of heterogeneity on hypergraph contagion models” by Landry and Restrepo https://doi.org/10.1063/5.0020034 and “Simplicial models of social contagion” by Iacopini et al. https://doi.org/10.1038/s41467-019-10431-6

Parameters:
  • H (HyperNetX Hypergraph object) –

  • tau (dictionary) – Edge sizes as keys (must account for all edge sizes present) and rates of infection for each size (float)

  • gamma (float) – The healing rate

  • transmission_function (lambda function, default: threshold) – A lambda function that has required arguments (node, status, edge) and optional arguments

  • initial_infecteds (list or numpy array, default: None) – Iterable of initially infected node uids

  • initial_recovereds (list or numpy array, default: None) – An iterable of initially recovered node uids

  • rho (float from 0 to 1, default: None) – The fraction of initially infected individuals. Both rho and initially infected cannot be specified.

  • tmin (float, default: 0) – Time at the start of the simulation

  • tmax (float, default: float('Inf')) – Time at which the simulation should be terminated if it hasn’t already.

  • dt (float > 0, default: 1.0) – Step forward in time that the simulation takes at each step.

  • return_full_data (bool, default: False) – This returns all the infection and recovery events at each time if True.

  • **args (Optional arguments to transmission function) – This allows user-defined transmission functions with extra parameters.

Returns:

  • if return_full_data

    dictionary

    Time as the keys and events that happen as the values.

  • else

    t, S, I, Rnumpy arrays

    time (t), number of susceptible (S), infected (I), and recovered (R) at each time.

Notes

Example:

>>> import hypernetx.algorithms.contagion as contagion
>>> import random
>>> import hypernetx as hnx
>>> n = 1000
>>> m = 10000
>>> hyperedgeList = [random.sample(range(n), k=random.choice([2,3])) for i in range(m)]
>>> H = hnx.Hypergraph(hyperedgeList)
>>> tau = {2:0.1, 3:0.1}
>>> gamma = 0.1
>>> tmax = 100
>>> dt = 0.1
>>> t, S, I, R = contagion.discrete_SIR(H, tau, gamma, rho=0.1, tmin=0, tmax=tmax, dt=dt)
algorithms.discrete_SIS(H, tau, gamma, transmission_function=<function threshold>, initial_infecteds=None, rho=None, tmin=0, tmax=100, dt=1.0, return_full_data=False, **args)[source]

A discrete-time SIS model for hypergraphs as implemented in “The effect of heterogeneity on hypergraph contagion models” by Landry and Restrepo https://doi.org/10.1063/5.0020034 and “Simplicial models of social contagion” by Iacopini et al. https://doi.org/10.1038/s41467-019-10431-6

Parameters:
  • H (HyperNetX Hypergraph object) –

  • tau (dictionary) – Edge sizes as keys (must account for all edge sizes present) and rates of infection for each size (float)

  • gamma (float) – The healing rate

  • transmission_function (lambda function, default: threshold) – A lambda function that has required arguments (node, status, edge) and optional arguments

  • initial_infecteds (list or numpy array, default: None) – Iterable of initially infected node uids

  • rho (float from 0 to 1, default: None) – The fraction of initially infected individuals. Both rho and initially infected cannot be specified.

  • tmin (float, default: 0) – Time at the start of the simulation

  • tmax (float, default: 100) – Time at which the simulation should be terminated if it hasn’t already.

  • dt (float > 0, default: 1.0) – Step forward in time that the simulation takes at each step.

  • return_full_data (bool, default: False) – This returns all the infection and recovery events at each time if True.

  • **args (Optional arguments to transmission function) – This allows user-defined transmission functions with extra parameters.

Returns:

  • if return_full_data

    dictionary

    Time as the keys and events that happen as the values.

  • else

    t, S, Inumpy arrays

    time (t), number of susceptible (S), and infected (I) at each time.

Notes

Example:

>>> import hypernetx.algorithms.contagion as contagion
>>> import random
>>> import hypernetx as hnx
>>> n = 1000
>>> m = 10000
>>> hyperedgeList = [random.sample(range(n), k=random.choice([2,3])) for i in range(m)]
>>> H = hnx.Hypergraph(hyperedgeList)
>>> tau = {2:0.1, 3:0.1}
>>> gamma = 0.1
>>> tmax = 100
>>> dt = 0.1
>>> t, S, I = contagion.discrete_SIS(H, tau, gamma, rho=0.1, tmin=0, tmax=tmax, dt=dt)
algorithms.erdos_renyi_hypergraph(n, m, p, node_labels=None, edge_labels=None)[source]

A function to generate an Erdos-Renyi hypergraph as implemented by Mirah Shi and described for bipartite networks by Aksoy et al. in https://doi.org/10.1093/comnet/cnx001

Parameters:
  • n (int) – Number of nodes

  • m (int) – Number of edges

  • p (float) – The probability that a bipartite edge is created

  • node_labels (list, default=None) – Vertex labels

  • edge_labels (list, default=None) – Hyperedge labels

Return type:

HyperNetX Hypergraph object

Example:

>>> import hypernetx.algorithms.generative_models as gm
>>> n = 1000
>>> m = n
>>> p = 0.01
>>> H = gm.erdos_renyi_hypergraph(n, m, p)
algorithms.get_pi(P)[source]

Returns the eigenvector corresponding to the largest eigenvalue (in magnitude), normalized so its entries sum to 1. Intended for the probability transition matrix of a random walk on a (connected) hypergraph, in which case the output can be interpreted as the stationary distribution.

Parameters:

P (csr matrix) – Probability transition matrix

Returns:

pi – Stationary distribution of random walk defined by P

Return type:

numpy.ndarray

algorithms.homology_basis(bd, k=None, boundary=False, **kwargs)[source]

Compute a basis for the kth-simplicial homology group, $H_k$, defined by a chain complex $C$ with boundary maps given by bd $= {k:partial_k }$

Parameters:
  • bd (dict) – dict of boundary matrices on k-chains to k-1 chains keyed on k if krange is a tuple then all boundary matrices k in [krange[0],..,krange[1]] inclusive must be in the dictionary

  • k (int or list of ints, optional, default=None) – k must be a positive integer or a list of 2 integers indicating min and max dimensions to be computed, if none given all homology groups will be computed from available boundary matrices in bd

  • boundary (bool) – option to return a basis for the boundary group from each dimension. Needed to compute the shortest generators in the homology group.

Returns:

  • basis (dict) – dict of generators as 0-1 tuples keyed by dim basis for dimension k will be returned only if bd[k] and bd[k+1] have been provided.

  • im (dict) – dict of boundary group generators keyed by dim

algorithms.hypergraph_homology_basis(h, k=None, shortest=False, interpreted=True)[source]

Computes the kth-homology groups mod 2 for the ASC associated with the hypergraph h for k in krange inclusive

Parameters:
  • h (hnx.Hypergraph) –

  • k (int or list of length 2, optional, default = None) – k must be an integer greater than 0 or a list of length 2 indicating min and max dimensions to be computed

  • shortest (bool, optional, default=False) – option to look for shortest representative for each coset in the homology group, only good for relatively small examples

  • interpreted (bool, optional, default = True) – if True will return an explicit basis in terms of the k-chains

Returns:

  • basis (list) – list of generators as k-chains as boolean vectors

  • interpreted_basis – lists of kchains in basis

algorithms.individual_contagion(node, status, edge)[source]

The individual contagion mechanism described in “The effect of heterogeneity on hypergraph contagion models” by Landry and Restrepo https://doi.org/10.1063/5.0020034

Parameters:
  • node (hashable) – The node uid to infect (If it doesn’t have status “S”, it will automatically return False)

  • status (dictionary) – The nodes are keys and the values are statuses (The infected state denoted with “I”)

  • edge (iterable) – Iterable of node ids (node must be in the edge or it will automatically return False)

Returns:

False if there is no potential to infect and True if there is.

Return type:

bool

Notes

Example:

>>> status = {0:"S", 1:"I", 2:"I", 3:"S", 4:"R"}
>>> individual_contagion(0, status, (0, 1, 3))
    True
>>> individual_contagion(1, status, (0, 1, 2))
    False
>>> collective_contagion(3, status, (0, 3, 4))
    False
algorithms.interpret(Ck, arr, labels=None)[source]

Returns the data as represented in Ck associated with the arr

Parameters:
  • Ck (list) – a list of k-cells being referenced by arr

  • arr (np.array) – array of 0-1 vectors

  • labels (dict, optional) – dictionary of labels to associate to the nodes in the cells

Returns:

list of k-cells referenced by data in Ck

Return type:

list

algorithms.kchainbasis(h, k)[source]

Compute the set of k dimensional cells in the abstract simplicial complex associated with the hypergraph.

Parameters:
  • h (hnx.Hypergraph) –

  • k (int) – dimension of cell

Returns:

an ordered list of kchains represented as tuples of length k+1

Return type:

list

See also

hnx.hypergraph.toplexes

Notes

  • Method works best if h is simple [Berge], i.e. no edge contains another and there are no duplicate edges (toplexes).

  • Hypergraph node uids must be sortable.

algorithms.kumar(HG, delta=0.01, verbose=False)[source]

Compute a partition of the vertices in hypergraph HG as per Kumar’s algorithm [1]

Parameters:
  • HG (Hypergraph) –

  • delta (float, optional) – convergence stopping criterion

Returns:

A partition of the vertices in HG

Return type:

list of sets

algorithms.last_step(HG, A, wdc=<function linear>, delta=0.01, verbose=False)[source]

Given some initial partition L, compute a new partition of the vertices in HG as per Last-Step algorithm [2]

Note

This is a very simple algorithm that tries moving nodes between communities to improve hypergraph modularity. It requires an initial non-trivial partition which can be obtained for example via graph clustering on the 2-section of HG, or via Kumar’s algorithm.

Parameters:
  • HG (Hypergraph) –

  • L (list of sets) – some initial partition of the vertices in HG

  • wdc (func, optional) – Hyperparameter for hypergraph modularity [2]

  • delta (float, optional) – convergence stopping criterion

  • verbose (boolean, optional) – If set, also returns progress after each pass through the vertices

Returns:

A new partition for the vertices in HG

Return type:

list of sets

algorithms.linear(d, c)[source]

Hyperparameter for hypergraph modularity [2] for d-edge with c vertices in the majority class. This is the default choice for modularity() and last_step() functions.

Parameters:
  • d (int) – Number of vertices in an edge

  • c (int) – Number of vertices in the majority class

Returns:

c/d if c>d/2 else 0

Return type:

float

algorithms.logical_dot(ar1, ar2)[source]

Returns the boolean equivalent of the dot product mod 2 on two 1-d arrays of the same length.

Parameters:
  • ar1 (numpy.ndarray) – 1-d array

  • ar2 (numpy.ndarray) – 1-d array

Returns:

boolean value associated with dot product mod 2

Return type:

bool

Raises:

HyperNetXError – If arrays are not of the same length an error will be raised.

algorithms.logical_matadd(mat1, mat2)[source]

Returns the boolean equivalent of matrix addition mod 2 on two binary arrays stored as type boolean

Parameters:
  • mat1 (np.ndarray) – 2-d array of boolean values

  • mat2 (np.ndarray) – 2-d array of boolean values

Returns:

mat – boolean matrix equivalent to the mod 2 matrix addition of the matrices as matrices over Z/2Z

Return type:

np.ndarray

Raises:

HyperNetXError – If dimensions are not equal an error will be raised.

algorithms.logical_matmul(mat1, mat2)[source]

Returns the boolean equivalent of matrix multiplication mod 2 on two binary arrays stored as type boolean

Parameters:
  • mat1 (np.ndarray) – 2-d array of boolean values

  • mat2 (np.ndarray) – 2-d array of boolean values

Returns:

mat – boolean matrix equivalent to the mod 2 matrix multiplication of the matrices as matrices over Z/2Z

Return type:

np.ndarray

Raises:

HyperNetXError – If inner dimensions are not equal an error will be raised.

algorithms.majority(d, c)[source]

Hyperparameter for hypergraph modularity [2] for d-edge with c vertices in the majority class. This corresponds to the majority rule [3]

Parameters:
  • d (int) – Number of vertices in an edge

  • c (int) – Number of vertices in the majority class

Returns:

1 if c>d/2 else 0

Return type:

bool

algorithms.majority_vote(node, status, edge)[source]

The majority vote contagion mechanism. If a majority of neighbors are contagious, it is possible for an individual to change their opinion. If opinions are divided equally, choose randomly.

Parameters:
  • node (hashable) – The node uid to infect (If it doesn’t have status “S”, it will automatically return False)

  • status (dictionary) – The nodes are keys and the values are statuses (The infected state denoted with “I”)

  • edge (iterable) – Iterable of node ids (node must be in the edge or it will automatically return False

Returns:

False if there is no potential to infect and True if there is.

Return type:

bool

Notes

Example:

>>> status = {0:"S", 1:"I", 2:"I", 3:"S", 4:"R"}
>>> majority_vote(0, status, (0, 1, 2))
    True
>>> majority_vote(0, status, (0, 1, 2, 3))
    True
>>> majority_vote(1, status, (0, 1, 2))
    False
>>> majority_vote(3, status, (0, 1, 2))
    False
algorithms.matmulreduce(arr, reverse=False)[source]

Recursively applies a ‘logical multiplication’ to a list of boolean arrays.

For arr = [arr[0],arr[1],arr[2]…arr[n]] returns product arr[0]arr[1]…arr[n] If reverse = True, returns product arr[n]arr[n-1]…arr[0]

Parameters:
  • arr (list of np.array) – list of nxm matrices represented as np.array

  • reverse (bool, optional) – order to multiply the matrices

Returns:

P – Product of matrices in the list

Return type:

np.array

algorithms.modularity(HG, A, wdc=<function linear>)[source]

Computes modularity of hypergraph HG with respect to partition A.

Parameters:
  • HG (Hypergraph) – The hypergraph with some precomputed attributes via: precompute_attributes(HG)

  • A (list of sets) – Partition of the vertices in HG

  • wdc (func, optional) – Hyperparameter for hypergraph modularity [2]

Note

For ‘wdc’, any function of the format w(d,c) that returns 0 when c <= d/2 and value in [0,1] otherwise can be used. Default is ‘linear’; other supplied choices are ‘majority’ and ‘strict’.

Returns:

The modularity function for partition A on HG

Return type:

float

algorithms.norm_lap(H, weights=False, index=True)[source]

Normalized Laplacian matrix of the hypergraph. Symmetrizes the probability transition matrix of a hypergraph random walk using the stationary distribution, using the digraph Laplacian defined in:

Chung, Fan. “Laplacians and the Cheeger inequality for directed graphs.” Annals of Combinatorics 9.1 (2005): 1-19.

and studied in the context of hypergraphs in:

Hayashi, K., Aksoy, S. G., Park, C. H., & Park, H. Hypergraph random walks, laplacians, and clustering. In Proceedings of CIKM 2020, (2020): 495-504.

Parameters:
  • H (hnx.Hypergraph) – The hypergraph must be connected, meaning there is a path linking any two vertices

  • weight (bool, optional, default : False) – Uses cell_weights, if False, uniform weights are utilized.

  • index (bool, optional) – Whether to return matrix-index to vertex-label mapping

Returns:

  • P (scipy.sparse.csr.csr_matrix) – Probability transition matrix of the random walk on the hypergraph

  • id (list) – contains list of index of node ids for rows

algorithms.part2dict(A)[source]

Given a partition (list of sets), returns a dictionary mapping the part for each vertex; inverse function to dict2part

Parameters:

A (list of sets) – a partition of the vertices

Returns:

a dictionary with {vertex: partition index}

Return type:

dict

algorithms.prob_trans(H, weights=False, index=True, check_connected=True)[source]

The probability transition matrix of a random walk on the vertices of a hypergraph. At each step in the walk, the next vertex is chosen by:

  1. Selecting a hyperedge e containing the vertex with probability proportional to w(e)

  2. Selecting a vertex v within e with probability proportional to a gamma(v,e)

If weights are not specified, then all weights are uniform and the walk is equivalent to a simple random walk. If weights are specified, the hyperedge weights w(e) are determined from the weights gamma(v,e).

Parameters:
  • H (hnx.Hypergraph) – The hypergraph must be connected, meaning there is a path linking any two vertices

  • weights (bool, optional, default : False) – Use the cell_weights associated with the hypergraph If False, uniform weights are utilized.

  • index (bool, optional) – Whether to return matrix index to vertex label mapping

Returns:

  • P (scipy.sparse.csr.csr_matrix) – Probability transition matrix of the random walk on the hypergraph

  • index (list) – contains list of index of node ids for rows

algorithms.reduced_row_echelon_form_mod2(M)[source]

Computes the invertible transformation matrices needed to compute the reduced row echelon form of M modulo 2

Parameters:

M (np.array) – a rectangular matrix with elements in $Z_2$

Returns:

L, S, Linv – LM = S where S is the reduced echelon form of M and M = LinvS

Return type:

np.arrays

algorithms.s_betweenness_centrality(H, s=1, edges=True, normalized=True, return_singletons=True)[source]

A centrality measure for an s-edge(node) subgraph of H based on shortest paths. Equals the betweenness centrality of vertices in the edge(node) s-linegraph.

In a graph (2-uniform hypergraph) the betweenness centrality of a vertex $v$ is the ratio of the number of non-trivial shortest paths between any pair of vertices in the graph that pass through $v$ divided by the total number of non-trivial shortest paths in the graph.

The centrality of edge to all shortest s-edge paths $V$ = the set of vertices in the linegraph. $sigma(s,t)$ = the number of shortest paths between vertices $s$ and $t$. $sigma(s,t|v)$ = the number of those paths that pass through vertex $v$.

\[c_B(v) = \sum_{s \neq t \in V} \frac{\sigma(s, t|v)}{\sigma(s,t)}\]
Parameters:
  • H (hnx.Hypergraph) –

  • s (int) – s connectedness requirement

  • edges (bool, optional) – determines if edge or node linegraph

  • normalized – bool, default=False, If true the betweenness values are normalized by 2/((n-1)(n-2)), where n is the number of edges in H

  • return_singletons (bool, optional) – if False will ignore singleton components of linegraph

Returns:

A dictionary of s-betweenness centrality value of the edges.

Return type:

dict

algorithms.s_closeness_centrality(H, s=1, edges=True, return_singletons=True, source=None)[source]

In a connected component the reciprocal of the sum of the distance between an edge(node) and all other edges(nodes) in the component times the number of edges(nodes) in the component minus 1.

$V$ = the set of vertices in the linegraph. $n = |V|$ $d$ = shortest path distance

\[C(u) = \frac{n - 1}{\sum_{v \neq u \in V} d(v, u)}\]
Parameters:
  • H (hnx.Hypergraph) –

  • s (int, optional) –

  • edges (bool, optional) – Indicates if method should compute edge linegraph (default) or node linegraph.

  • return_singletons (bool, optional) – Indicates if method should return values for singleton components.

  • source (str, optional) – Identifier of node or edge of interest for computing centrality

Returns:

returns the s-closeness centrality value of the edges(nodes). If source=None a dictionary of values for each s-edge in H is returned. If source then a single value is returned.

Return type:

dict or float

algorithms.s_eccentricity(H, s=1, edges=True, source=None, return_singletons=True)[source]

The length of the longest shortest path from a vertex $u$ to every other vertex in the s-linegraph. $V$ = set of vertices in the s-linegraph $d$ = shortest path distance

\[\text{s-ecc}(u) = \text{max}\{d(u,v): v \in V\}\]
Parameters:
  • H (hnx.Hypergraph) –

  • s (int, optional) –

  • edges (bool, optional) – Indicates if method should compute edge linegraph (default) or node linegraph.

  • return_singletons (bool, optional) – Indicates if method should return values for singleton components.

  • source (str, optional) – Identifier of node or edge of interest for computing centrality

Returns:

returns the s-eccentricity value of the edges(nodes). If source=None a dictionary of values for each s-edge in H is returned. If source then a single value is returned. If the s-linegraph is disconnected, np.inf is returned.

Return type:

dict or float

algorithms.s_harmonic_centrality(H, s=1, edges=True, source=None, normalized=False, return_singletons=True)[source]

A centrality measure for an s-edge subgraph of H. A value equal to 1 means the s-edge intersects every other s-edge in H. All values range between 0 and 1. Edges of size less than s return 0. If H contains only one s-edge a 0 is returned.

The denormalized reciprocal of the harmonic mean of all distances from $u$ to all other vertices. $V$ = the set of vertices in the linegraph. $d$ = shortest path distance

\[C(u) = \sum_{v \neq u \in V} \frac{1}{d(v, u)}\]

Normalized this becomes: $$C(u) = sum_{v neq u in V} frac{1}{d(v, u)}cdotfrac{2}{(n-1)(n-2)}$$ where $n$ is the number vertices.

Parameters:
  • H (hnx.Hypergraph) –

  • s (int, optional) –

  • edges (bool, optional) – Indicates if method should compute edge linegraph (default) or node linegraph.

  • return_singletons (bool, optional) – Indicates if method should return values for singleton components.

  • source (str, optional) – Identifier of node or edge of interest for computing centrality

Returns:

returns the s-harmonic closeness centrality value of the edges, a number between 0 and 1 inclusive. If source=None a dictionary of values for each s-edge in H is returned. If source then a single value is returned.

Return type:

dict or float

algorithms.s_harmonic_closeness_centrality(H, s=1, edge=None)[source]
algorithms.smith_normal_form_mod2(M)[source]

Computes the invertible transformation matrices needed to compute the Smith Normal Form of M modulo 2

Parameters:
  • M (np.array) – a rectangular matrix with data type bool

  • track (bool) – if track=True will print out the transformation as Z/2Z matrix as it discovers L[i] and R[j]

Returns:

L, R, S, Linv – LMR = S is the Smith Normal Form of the matrix M.

Return type:

np.arrays

Note

Given a mxn matrix $M$ with entries in $Z_2$ we start with the equation: $L M R = S$, where $L = I_m$, and $R=I_n$ are identity matrices and $S = M$. We repeatedly apply actions to the left and right side of the equation to transform S into a diagonal matrix. For each action applied to the left side we apply its inverse action to the right side of I_m to generate $L^{-1}$. Finally we verify: $L M R = S$ and $LLinv = I_m$.

algorithms.spec_clus(H, k, existing_lap=None, weights=False)[source]

Hypergraph spectral clustering of the vertex set into k disjoint clusters using the normalized hypergraph Laplacian. Equivalent to the “RDC-Spec” Algorithm 1 in:

Hayashi, K., Aksoy, S. G., Park, C. H., & Park, H. Hypergraph random walks, laplacians, and clustering. In Proceedings of CIKM 2020, (2020): 495-504.

Parameters:
  • H (hnx.Hypergraph) – The hypergraph must be connected, meaning there is a path linking any two vertices

  • k (int) – Number of clusters

  • existing_lap (csr matrix, optional) – Whether to use an existing Laplacian; otherwise, normalized hypergraph Laplacian will be utilized

  • weights (bool, optional) – Use the cell_weights of the hypergraph. If False uniform weights are used.

Returns:

clusters – Vertex cluster dictionary, keyed by integers 0,…,k-1, with lists of vertices as values.

Return type:

dict

algorithms.strict(d, c)[source]

Hyperparameter for hypergraph modularity [2] for d-edge with c vertices in the majority class. This corresponds to the strict rule [3]

Parameters:
  • d (int) – Number of vertices in an edge

  • c (int) – Number of vertices in the majority class

Returns:

1 if c==d else 0

Return type:

bool

algorithms.swap_columns(i, j, *args)[source]

Swaps ith and jth column of each matrix in args Returns a list of new matrices

Parameters:
  • i (int) –

  • j (int) –

  • args (np.arrays) –

Returns:

list of copies of args with ith and jth row swapped

Return type:

list

algorithms.swap_rows(i, j, *args)[source]

Swaps ith and jth row of each matrix in args Returns a list of new matrices

Parameters:
  • i (int) –

  • j (int) –

  • args (np.arrays) –

Returns:

list of copies of args with ith and jth row swapped

Return type:

list

algorithms.threshold(node, status, edge, tau=0.1)[source]

The threshold contagion mechanism

Parameters:
  • node (hashable) – The node uid to infect (If it doesn’t have status “S”, it will automatically return False)

  • status (dictionary) – The nodes are keys and the values are statuses (The infected state denoted with “I”)

  • edge (iterable) – Iterable of node ids (node must be in the edge or it will automatically return False)

  • tau (float between 0 and 1, default: 0.1) – The fraction of nodes in an edge that must be infected for the edge to be able to transmit to the node

Returns:

False if there is no potential to infect and True if there is.

Return type:

bool

Notes

Example:

>>> status = {0:"S", 1:"I", 2:"I", 3:"S", 4:"R"}
>>> threshold(0, status, (0, 2, 3, 4), tau=0.2)
    True
>>> threshold(0, status, (0, 2, 3, 4), tau=0.5)
    False
>>> threshold(3, status, (1, 2, 3), tau=1)
    False
algorithms.two_section(HG)[source]

Creates a random walk based [1] 2-section igraph Graph with transition weights defined by the weights of the hyperedges.

Parameters:

HG (Hypergraph) –

Returns:

The 2-section graph built from HG

Return type:

igraph.Graph