RankedBlockState#

class graph_tool.inference.RankedBlockState(g, b=None, u=None, entropy_args={}, **kwargs)[source]#

Bases: MCMCState, MultiflipMCMCState, MultilevelMCMCState, GibbsMCMCState, DrawBlockState

Obtain the ordered partition of a network according to the ranked stochastic block model.

Parameters:
gGraph

Graph to be modelled.

bPropertyMap (optional, default: None)

Initial partition. If not supplied, a partition into a single group will be used.

uPropertyMap or iterable (optional, default: None)

Ordering of the group labels. It should contain a map from each group label to the unit interval \([0,1]\), inidicating how the groups should be ordered. If not supplied, the numeric values of the group lalbels will be used to initialize the ordering.

entropy_args: ``dict`` (optional, default: ``{}``)

Override default arguments for entropy() method and releated operations.

References

[peixoto-ordered-2022]

Tiago P. Peixoto, “Ordered community detection in directed networks”, Phys. Rev. E 106, 024305 (2022), DOI: 10.1103/PhysRevE.106.024305 [sci-hub, @tor], arXiv: 2203.16460

Methods

collect_vertex_marginals([p, b, update])

Collect the vertex marginal histogram, which counts the number of times a node was assigned to a given block.

copy(**kwargs)

Copies the state.

draw(**kwargs)

Convenience wrapper to graph_draw() that draws the state of the graph as colors on the vertices and edges.

entropy([adjacency, dl, partition_dl, ...])

Return the description length (negative joint log-likelihood).

get_B()

Returns the total number of blocks.

get_Be()

Returns the effective number of blocks, defined as \(e^{H}\), with \(H=-\sum_r\frac{n_r}{N}\ln \frac{n_r}{N}\), where \(n_r\) is the number of nodes in group r.

get_E()

Return the number of edges.

get_Es()

Return the number of dowstream, lateral, and upstream edges.

get_N()

Return the number of nodes.

get_block_order()

Returns an array indexed by the group label containing its rank order.

get_block_state([b, vweight, deg_corr])

Returns a BlockState corresponding to the block graph (i.e. the blocks of the current state become the nodes).

get_blocks()

Returns the property map which contains the block labels for each vertex.

get_edge_colors()

Return EdgePropertyMap containing the edge colors according to their rank direction: upstream (blue), downstream (red), lateral (grey).

get_edge_dir()

Return an edge PropertyMap containing the edge direction: -1 (downstream), 0 (lateral), +1 (upstream).

get_entropy_args()

Return the current default values for the parameters of the function entropy(), together with other operations that depend on them.

get_nonempty_B()

Alias to get_B().

get_state()

Alias to get_blocks().

get_vertex_order()

Returns a vertex PropertyMap with the rank order for every vertex.

get_vertex_position()

Returns a vertex PropertyMap with vertex positions in the unit interval \([0,1]\).

gibbs_sweep([beta, niter, entropy_args, ...])

Perform niter sweeps of a rejection-free Gibbs MCMC to sample network partitions.

mcmc_sweep([beta, c, d, niter, ...])

Perform niter sweeps of a Metropolis-Hastings acceptance-rejection MCMC to sample network partitions.

move_vertex(v, s)

Move vertex v to block s.

multiflip_mcmc_sweep([pmovelabel])

Call MultiflipMCMCState.multiflip_mcmc_sweep() with default parameter pmovelabel=1.`

multilevel_mcmc_sweep([niter, beta, c, d, ...])

Perform niter sweeps of a multilevel agglomerative acceptance-rejection pseudo-MCMC (i.e. detailed balance is not preserved) to sample network partitions, that uses a bisection search on the number of groups, together with group merges and singe-node moves.

reset_entropy_args()

Reset the current default values for the parameters of the function entropy(), together with other operations that depend on them.

sample_graph(**kwargs)

Sample a new graph from the fitted model.

update_entropy_args(**kwargs)

Update the default values for the parameters of the function entropy() from the keyword arguments, in a stateful way, together with other operations that depend on them.

virtual_vertex_move(v, s, **kwargs)

Computes the entropy difference if vertex v is moved to block s.

collect_vertex_marginals(p=None, b=None, update=1)[source]#

Collect the vertex marginal histogram, which counts the number of times a node was assigned to a given block.

Parameters:
pVertexPropertyMap (optional, default: None)

Vertex property map with vector-type values, storing the previous block membership counts. If not provided, an empty histogram will be created.

bVertexPropertyMap (optional, default: None)

Vertex property map with group partition. If not provided, the state’s partition will be used.

updateint (optional, default: 1)

Each call increases the current count by the amount given by this parameter.

Returns:
pVertexPropertyMap

Vertex property map with vector-type values, storing the accumulated block membership counts.

copy(**kwargs)[source]#

Copies the state. The parameters override the state properties, and have the same meaning as in the constructor.

draw(**kwargs)[source]#

Convenience wrapper to graph_draw() that draws the state of the graph as colors on the vertices and edges.

entropy(adjacency=True, dl=True, partition_dl=True, degree_dl=True, degree_dl_kind='distributed', edges_dl=True, dense=False, multigraph=True, deg_entropy=True, recs=True, recs_dl=True, beta_dl=1.0, Bfield=True, exact=True, **kwargs)#

Return the description length (negative joint log-likelihood). See BlockState.entropy() for documentation.

Warning

The default arguments of this function are overriden by those obtained from get_entropy_args(). To update the defaults in a stateful way, update_entropy_args() should be called.

get_B()[source]#

Returns the total number of blocks.

get_Be()[source]#

Returns the effective number of blocks, defined as \(e^{H}\), with \(H=-\sum_r\frac{n_r}{N}\ln \frac{n_r}{N}\), where \(n_r\) is the number of nodes in group r.

get_E()[source]#

Return the number of edges.

get_Es()[source]#

Return the number of dowstream, lateral, and upstream edges.

get_N()[source]#

Return the number of nodes.

get_block_order()[source]#

Returns an array indexed by the group label containing its rank order.

get_block_state(b=None, vweight=None, deg_corr=False, **kwargs)[source]#

Returns a BlockState corresponding to the block graph (i.e. the blocks of the current state become the nodes).

get_blocks()[source]#

Returns the property map which contains the block labels for each vertex.

get_edge_colors()[source]#

Return EdgePropertyMap containing the edge colors according to their rank direction: upstream (blue), downstream (red), lateral (grey).

get_edge_dir()[source]#

Return an edge PropertyMap containing the edge direction: -1 (downstream), 0 (lateral), +1 (upstream).

get_entropy_args()#

Return the current default values for the parameters of the function entropy(), together with other operations that depend on them.

get_nonempty_B()[source]#

Alias to get_B().

get_state()[source]#

Alias to get_blocks().

get_vertex_order()[source]#

Returns a vertex PropertyMap with the rank order for every vertex.

get_vertex_position()[source]#

Returns a vertex PropertyMap with vertex positions in the unit interval \([0,1]\).

gibbs_sweep(beta=1.0, niter=1, entropy_args={}, allow_new_group=True, sequential=True, deterministic=False, vertices=None, verbose=False, **kwargs)#

Perform niter sweeps of a rejection-free Gibbs MCMC to sample network partitions.

Parameters:
betafloat (optional, default: 1.)

Inverse temperature.

niterint (optional, default: 1)

Number of sweeps to perform. During each sweep, a move attempt is made for each node.

entropy_argsdict (optional, default: {})

Entropy arguments, with the same meaning and defaults as in graph_tool.inference.BlockState.entropy().

allow_new_groupbool (optional, default: True)

Allow the number of groups to increase and decrease.

sequentialbool (optional, default: True)

If sequential == True each vertex move attempt is made sequentially, where vertices are visited in random order. Otherwise the moves are attempted by sampling vertices randomly, so that the same vertex can be moved more than once, before other vertices had the chance to move.

deterministicbool (optional, default: False)

If sequential == True and deterministic == True the vertices will be visited in deterministic order.

verticeslist of ints (optional, default: None)

If provided, this should be a list of vertices which will be moved. Otherwise, all vertices will.

verbosebool (optional, default: False)

If verbose == True, detailed information will be displayed.

Returns:
dSfloat

Entropy difference after the sweeps.

nattemptsint

Number of vertex moves attempted.

nmovesint

Number of vertices moved.

Notes

This algorithm has an \(O(E\times B)\) complexity, where \(B\) is the number of groups, and \(E\) is the number of edges.

mcmc_sweep(beta=1.0, c=0.5, d=0.01, niter=1, entropy_args={}, allow_vacate=True, sequential=True, deterministic=False, vertices=None, verbose=False, **kwargs)#

Perform niter sweeps of a Metropolis-Hastings acceptance-rejection MCMC to sample network partitions.

Parameters:
betafloat (optional, default: 1.)

Inverse temperature.

cfloat (optional, default: .5)

Sampling parameter c for move proposals: For \(c\to 0\) the blocks are sampled according to the local neighborhood of a given node and their block connections; for \(c\to\infty\) the blocks are sampled randomly. Note that only for \(c > 0\) the MCMC is guaranteed to be ergodic.

dfloat (optional, default: .01)

Probability of selecting a new (i.e. empty) group for a given move.

niterint (optional, default: 1)

Number of sweeps to perform. During each sweep, a move attempt is made for each node.

entropy_argsdict (optional, default: {})

Entropy arguments, with the same meaning and defaults as in graph_tool.inference.BlockState.entropy().

allow_vacatebool (optional, default: True)

Allow groups to be vacated.

sequentialbool (optional, default: True)

If sequential == True each vertex move attempt is made sequentially, where vertices are visited in random order. Otherwise the moves are attempted by sampling vertices randomly, so that the same vertex can be moved more than once, before other vertices had the chance to move.

deterministicbool (optional, default: False)

If sequential == True and deterministic == True the vertices will be visited in deterministic order.

verticeslist of ints (optional, default: None)

If provided, this should be a list of vertices which will be moved. Otherwise, all vertices will.

verbosebool (optional, default: False)

If verbose == True, detailed information will be displayed.

Returns:
dSfloat

Entropy difference after the sweeps.

nattemptsint

Number of vertex moves attempted.

nmovesint

Number of vertices moved.

Notes

This algorithm has an \(O(E)\) complexity, where \(E\) is the number of edges (independent of the number of groups).

References

[peixoto-efficient-2014]

Tiago P. Peixoto, “Efficient Monte Carlo and greedy heuristic for the inference of stochastic block models”, Phys. Rev. E 89, 012804 (2014), DOI: 10.1103/PhysRevE.89.012804 [sci-hub, @tor], arXiv: 1310.4378

move_vertex(v, s)[source]#

Move vertex v to block s.

multiflip_mcmc_sweep(pmovelabel=1, **kwargs)[source]#

Call MultiflipMCMCState.multiflip_mcmc_sweep() with default parameter pmovelabel=1.`

multilevel_mcmc_sweep(niter=1, beta=inf, c=0.5, d=0.01, r=0.9, random_bisect=True, merge_sweeps=10, mh_sweeps=10, init_r=0.99, init_min_iter=5, init_beta=1.0, gibbs=False, B_min=1, B_max=18446744073709551615, b_min=None, b_max=None, M=None, cache_states=True, force_accept=False, parallel=False, entropy_args={}, verbose=False, **kwargs)#

Perform niter sweeps of a multilevel agglomerative acceptance-rejection pseudo-MCMC (i.e. detailed balance is not preserved) to sample network partitions, that uses a bisection search on the number of groups, together with group merges and singe-node moves.

Parameters:
niterint (optional, default: 1)

Number of sweeps to perform. During each sweep, a move attempt is made for each node, on average.

betafloat (optional, default: numpy.inf)

Inverse temperature.

cfloat (optional, default: .5)

Sampling parameter c for move proposals: For \(c\to 0\) the blocks are sampled according to the local neighborhood of a given node and their block connections; for \(c\to\infty\) the blocks are sampled randomly. Note that only for \(c > 0\) the MCMC is guaranteed to be ergodic.

dfloat (optional, default: .01)

Probability of selecting a new (i.e. empty) group for a given single-node move.

rfloat (optional, default: 0.9)

Group shrink ratio. The number of groups is reduced by this fraction at each merge sweep.

random_bisectbool (optional, default: True)

If True, bisections are done at randomly chosen intervals. Otherwise a Fibonacci sequence is used.

merge_sweepsint (optional, default: 10)

Number of sweeps spent to find good merge proposals.

mh_sweepsint (optional, default: 10)

Number of single-node Metropolis-Hastings sweeps between merge splits.

init_rdouble (optional, default: 0.99)

Stopping criterion for the intialization phase, after each node is put in their own group, to set the initial upper bound of the bisection search. A number of single-node Metropolis-Hastings sweeps is done until the number of groups is shrunk by a factor that is larger than this parameter.

init_min_iterint (optional, default: 5)

Minimum number of iterations at the intialization phase.

init_betafloat (optional, default: 1.)

Inverse temperature to be used for the very first sweep of the initialization phase.

gibbsbool (optional, default: False)

If True, the single node moves use (slower) Gibbs sampling, rather than Metropolis-Hastings.

B_minint (optional, default: 1)

Minimum number of groups to be considered in the search.

b_minVertexPropertyMap (optional, default: None)

If provided, this will be used for the partition corresponding to B_min.

B_maxint (optional, default: 1)

Maximum number of groups to be considered in the search.

b_maxVertexPropertyMap (optional, default: None)

If provided, this will be used for the partition corresponding to B_max.

Mint (optional, default: None)

Maximum number of groups to select for the multilevel move. If None is provided, then all groups are always elected.

cache_statesbool (optional, default: True)

If True, intermediary states will be cached during the bisection search.

force_acceptbool (optional, default: False)

If True, new state will be accepted even if it does not improve the objective function.

parallelbool (optional, default: False)

If True, the algorithm will run in parallel (if enabled during compilation).

entropy_argsdict (optional, default: {})

Entropy arguments, with the same meaning and defaults as in graph_tool.inference.BlockState.entropy().

verbosebool (optional, default: False)

If verbose == True, detailed information will be displayed.

Returns:
dSfloat

Entropy difference after the sweeps.

nattemptsint

Number of vertex moves attempted.

nmovesint

Number of vertices moved.

Notes

This algorithm has an \(O(E\ln^2 N)\) complexity, where \(E\) is the number of edges and \(N\) is the number of nodes (independently of the number of groups).

References

[peixoto-efficient-2014]

Tiago P. Peixoto, “Efficient Monte Carlo and greedy heuristic for the inference of stochastic block models”, Phys. Rev. E 89, 012804 (2014), DOI: 10.1103/PhysRevE.89.012804 [sci-hub, @tor], arXiv: 1310.4378

reset_entropy_args()#

Reset the current default values for the parameters of the function entropy(), together with other operations that depend on them.

sample_graph(**kwargs)[source]#

Sample a new graph from the fitted model. See BlockState.sample_graph() for documentation.

update_entropy_args(**kwargs)#

Update the default values for the parameters of the function entropy() from the keyword arguments, in a stateful way, together with other operations that depend on them.

Values updated in this manner are preserved by the copying or pickling of the state.

virtual_vertex_move(v, s, **kwargs)[source]#

Computes the entropy difference if vertex v is moved to block s. The remaining parameters are the same as in entropy().