graph_tool.inference.PPBlockState#

class graph_tool.inference.PPBlockState(g, b=None)[source]#

Bases: MCMCState, MultiflipMCMCState, MultilevelMCMCState, GibbsMCMCState, DrawBlockState

Obtain the partition of a network according to the Bayesian planted partition 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.

References

[lizhi-statistical-2020]

Lizhi Zhang, Tiago P. Peixoto, “Statistical inference of assortative community structures”, Phys. Rev. Research 2 043271 (2020), DOI: 10.1103/PhysRevResearch.2.043271 [sci-hub, @tor], arXiv: 2006.14493

Methods

copy([g, b])

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([uniform, degree_dl_kind])

Return the model entropy (negative 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_blocks()

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

get_nonempty_B()

Alias to get_B().

get_state()

Alias to get_blocks().

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([beta, c, psingle, ...])

Perform niter sweeps of a Metropolis-Hastings acceptance-rejection MCMC with multiple simultaneous moves (i.e. merges and splits) to sample network partitions.

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

Perform niter sweeps of a multilevel agglomerative acceptance-rejection MCMC to sample network partitions, that uses a bisection search on the number of groups, together with group merges and singe-node moves.

virtual_vertex_move(v, s, **kwargs)

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

copy(g=None, b=None)[source]#

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

draw(**kwargs)#

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

entropy(uniform=False, degree_dl_kind='distributed', **kwargs)[source]#

Return the model entropy (negative log-likelihood).

Parameters:
uniformbool (optional, default: False)

If True, the uniform planted partition model is used, otherwise a non-uniform version is used.

degree_dl_kindstr (optional, default: "distributed")

This specifies the prior used for the degree sequence. It must be one of: "uniform" or "distributed" (default).

Notes

The “entropy” of the state is the negative log-likelihood of the microcanonical SBM, that includes the generated graph \(\boldsymbol{A}\) and the model parameters \(e_{\text{in}}\), \(e_{\text{out}}\), \(\boldsymbol{k}\) and \(\boldsymbol{b}\),

\[\begin{split}\Sigma &= - \ln P(\boldsymbol{A},e_{\text{in}},e_{\text{out}},\boldsymbol{k},\boldsymbol{b}) \\ &= - \ln P(\boldsymbol{A}|e_{\text{in}},e_{\text{out}},\boldsymbol{k},\boldsymbol{b}) - \ln P(e_{\text{in}},e_{\text{out}},\boldsymbol{k},\boldsymbol{b}).\end{split}\]

This value is also called the description length of the data, and it corresponds to the amount of information required to describe it (in nats).

For the uniform version of the model, the likelihood is

\[P(\boldsymbol{A}|\boldsymbol{k},\boldsymbol{b}) = \frac{e_{\text{in}}!e_{\text{out}}!} {\left(\frac{B}{2}\right)^{e_{\text{in}}}{B\choose 2}^{e_{\text{out}}}(E+1)^{1-\delta_{B,1}}\prod_re_r!}\times \frac{\prod_ik_i!}{\prod_{i<j}A_{ij}!\prod_i A_{ii}!!}.\]

where \(e_{\text{in}}\) and \(e_{\text{out}}\) are the number of edges inside and outside communities, respectively, and \(e_r\) is the sum of degrees in group \(r\).

For the non-uniform model we have instead:

\[P(\boldsymbol{A}|\boldsymbol{k},\boldsymbol{b}) = \frac{e_{\text{out}}!\prod_re_{rr}!!} {{B\choose 2}^{e_{\text{out}}}(E+1)^{1-\delta_{B,1}}\prod_re_r!}\times{B + e_{\text{in}} - 1 \choose e_{\text{in}}}^{-1}\times \frac{\prod_ik_i!}{\prod_{i<j}A_{ij}!\prod_i A_{ii}!!}.\]

Here there are two options for the prior on the degrees:

  1. degree_dl_kind == "uniform"

    \[P(\boldsymbol{k}|\boldsymbol{e},\boldsymbol{b}) = \prod_r\left(\!\!{n_r\choose e_r}\!\!\right)^{-1}.\]

    This corresponds to a noninformative prior, where the degrees are sampled from a uniform distribution.

  2. degree_dl_kind == "distributed" (default)

    \[P(\boldsymbol{k}|\boldsymbol{e},\boldsymbol{b}) = \prod_r\frac{\prod_k\eta_k^r!}{n_r!} \prod_r q(e_r, n_r)^{-1}\]

    with \(\eta_k^r\) being the number of nodes with degree \(k\) in group \(r\), and \(q(n,m)\) being the number of partitions of integer \(n\) into at most \(m\) parts.

    This corresponds to a prior for the degree sequence conditioned on the degree frequencies, which are themselves sampled from a uniform hyperprior. This option should be preferred in most cases.

For the partition prior \(P(\boldsymbol{b})\) please refer to entropy().

References

[lizhi-statistical-2020]

Lizhi Zhang, Tiago P. Peixoto, “Statistical inference of assortative community structures”, Phys. Rev. Research 2 043271 (2020), DOI: 10.1103/PhysRevResearch.2.043271 [sci-hub, @tor], arXiv: 2006.14493

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_blocks()[source]#

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

get_nonempty_B()[source]#

Alias to get_B().

get_state()[source]#

Alias to get_blocks().

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(beta=1.0, c=0.5, psingle=None, psplit=1, pmerge=1, pmergesplit=1, pmovelabel=0, d=0.01, gibbs_sweeps=10, niter=1, entropy_args={}, accept_stats=None, verbose=False, **kwargs)#

Perform niter sweeps of a Metropolis-Hastings acceptance-rejection MCMC with multiple simultaneous moves (i.e. merges and splits) 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.

psinglefloat (optional, default: None)

Relative probability of proposing a single node move. If None, it will be selected as the number of nodes in the graph.

psplitfloat (optional, default: 1)

Relative probability of proposing a group split.

pmergesplitfloat (optional, default: 1)

Relative probability of proposing a marge-split move.

pmovelabelfloat (optional, default: 0)

Relative probability of proposing a group label move.

dfloat (optional, default: 1)

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

gibbs_sweepsint (optional, default: 10)

Number of sweeps of Gibbs sampling to be performed (i.e. each node is attempted once per sweep) to refine a split proposal.

niterint (optional, default: 1)

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

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)\) complexity, where \(E\) is the number of edges (independent of the number of groups).

References

[peixoto-merge-split-2020]

Tiago P. Peixoto, “Merge-split Markov chain Monte Carlo for community detection”, Phys. Rev. E 102, 012305 (2020), DOI: 10.1103/PhysRevE.102.012305 [sci-hub, @tor], arXiv: 2003.07070

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

Perform niter sweeps of a multilevel agglomerative acceptance-rejection MCMC 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: 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.

psinglefloat (optional, default: None)

Relative probability of proposing a single node move. If None, it will be selected as the number of nodes in the graph.

pmultilevelfloat (optional, default: 1)

Relative probability of proposing a multilevel move.

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_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.

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

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().