NestedBlockState#
- class graph_tool.inference.NestedBlockState(g, bs=None, base_type=<class 'graph_tool.inference.blockmodel.BlockState'>, state_args={}, hstate_args={}, hentropy_args={}, **kwargs)[source]#
Bases:
objectThe nested stochastic block model state of a given graph.
- Parameters:
- g
Graph Graph to be modeled.
- bs
listofVertexPropertyMapornumpy.ndarray(optional, default:None) Hierarchical node partition. If not provided it will correspond to a single-group hierarchy of length \(\lceil\log_2(N)\rceil\).
- base_type
type(optional, default:BlockState) State type for lowermost level (e.g.
BlockState,OverlapBlockStateorLayeredBlockState)- hstate_args
dict(optional, default: {}) Keyword arguments to be passed to the constructor of the higher-level states.
- hentropy_args
dict(optional, default: {}) Keyword arguments to be passed to the
entropy()method of the higher-level states.- state_args
dict(optional, default:{}) Keyword arguments to be passed to base type constructor.
- **kwargskeyword arguments
Keyword arguments to be passed to base type constructor. The
state_argsparameter overrides this.
- g
Methods
add_vertex(v, r)Add vertex
vto blockr.collect_partition_histogram([h, update])Collect a histogram of partitions.
copy(**kwargs)Copies the block state.
draw(**kwargs)Convenience wrapper to
draw_hierarchy()that draws the hierarchical state.entropy(**kwargs)Obtain the description length (i.e. negative joint log-likelihood) for the hierarchical partition.
get_bs()Get hierarchy levels as a list of
numpy.ndarrayobjects with the group memberships at each level.Return the nested levels as individual graphs.
get_clabel(l)Get clabel for level
l.get_edges_prob(missing[, spurious, entropy_args])Compute the joint log-probability of the missing and spurious edges given by
missingandspurious(a list of(source, target)tuples, orEdge()instances), together with the observed edges.Get hierarchy levels as a list of
BlockStateinstances.Alias to
get_bs().gibbs_sweep(**kwargs)Perform
nitersweeps of a rejection-free Gibbs sampling MCMC to sample network partitions.level_entropy(l[, bstate])Compute the entropy of level
l.mcmc_sweep(**kwargs)Perform
nitersweeps of a Metropolis-Hastings acceptance-rejection MCMC to sample hierarchical network partitions.move_vertex(v, s)Move vertex
vto blocks.multicanonical_sweep(m_state, **kwargs)Perform
nitersweeps of a non-Markovian multicanonical sampling using the Wang-Landau algorithm.multiflip_mcmc_sweep(**kwargs)Perform
nitersweeps of a Metropolis-Hastings acceptance-rejection MCMC with multiple moves to sample hierarchical network partitions.multilevel_mcmc_sweep(**kwargs)Perform
nitersweeps of a Metropolis-Hastings acceptance-rejection MCMC with multilevel moves to sample hierarchical network partitions.Print a hierarchy summary.
Project the partition at level
lonto the lowest level, and return the corresponding state.project_partition(j, l)Project partition of level
jonto levell, and return it.Project base clabel to level
l.Remove vertex
vfrom its current group.set_state(bs)Sets the internal nested partition of the state.
- add_vertex(v, r)[source]#
Add vertex
vto blockr.This optionally accepts a list of vertices and blocks to add.
Warning
This can leave the state in an inconsistent state if a vertex is added twice to the same group.
- collect_partition_histogram(h=None, update=1)[source]#
Collect a histogram of partitions.
This should be called multiple times, e.g. after repeated runs of the
graph_tool.inference.NestedBlockState.mcmc_sweep()function.- Parameters:
- h
PartitionHist(optional, default:None) Partition histogram. If not provided, an empty histogram will be created.
- updatefloat (optional, default:
1) Each call increases the current count by the amount given by this parameter.
- h
- Returns:
- h
PartitionHist(optional, default:None) Updated Partition histogram.
- h
- copy(**kwargs)[source]#
Copies the block state. The parameters override the state properties, and have the same meaning as in the constructor.
- draw(**kwargs)[source]#
Convenience wrapper to
draw_hierarchy()that draws the hierarchical state.
- entropy(**kwargs)[source]#
Obtain the description length (i.e. negative joint log-likelihood) for the hierarchical partition.
The keyword arguments are passed to the
entropy()method of the underlying state objects (e.g.graph_tool.inference.BlockState.entropy,graph_tool.inference.OverlapBlockState.entropy, orgraph_tool.inference.LayeredBlockState.entropy).
- get_bs()[source]#
Get hierarchy levels as a list of
numpy.ndarrayobjects with the group memberships at each level.
- get_bstack()[source]#
Return the nested levels as individual graphs.
This returns a list of
Graphinstances representing the inferred hierarchy at each level. Each graph has two internal vertex and edge property maps named “count” which correspond to the vertex and edge counts at the lower level, respectively. Additionally, an internal vertex property map named “b” specifies the block partition.
- get_edges_prob(missing, spurious=[], entropy_args={})[source]#
Compute the joint log-probability of the missing and spurious edges given by
missingandspurious(a list of(source, target)tuples, orEdge()instances), together with the observed edges.More precisely, the log-likelihood returned is
\[\ln \frac{P(\boldsymbol G + \delta \boldsymbol G | \boldsymbol b)}{P(\boldsymbol G| \boldsymbol b)}\]where \(\boldsymbol G + \delta \boldsymbol G\) is the modified graph (with missing edges added and spurious edges deleted).
The values in
entropy_argsare passed tograph_tool.inference.BlockState.entropy()to calculate the log-probability.
- get_levels()[source]#
Get hierarchy levels as a list of
BlockStateinstances.
- gibbs_sweep(**kwargs)[source]#
Perform
nitersweeps of a rejection-free Gibbs sampling MCMC to sample network partitions.The arguments accepted are the same as in
graph_tool.inference.BlockState.gibbs_sweep().Warning
This function performs
nitersweeps at each hierarchical level once. This means that in order for the chain to equilibrate, we need to call this function several times, i.e. it is not enough to call it once with a large value ofniter.
- mcmc_sweep(**kwargs)[source]#
Perform
nitersweeps of a Metropolis-Hastings acceptance-rejection MCMC to sample hierarchical network partitions.The arguments accepted are the same as in
graph_tool.inference.BlockState.mcmc_sweep().If the parameter
cis a scalar, the values used at each level arec * 2 ** lforlin the range[0, L-1]. Optionally, a list of values may be passed instead, which specifies the value ofc[l]to be used at each level.Warning
This function performs
nitersweeps at each hierarchical level once. This means that in order for the chain to equilibrate, we need to call this function several times, i.e. it is not enough to call it once with a large value ofniter.
- multicanonical_sweep(m_state, **kwargs)[source]#
Perform
nitersweeps of a non-Markovian multicanonical sampling using the Wang-Landau algorithm.The arguments accepted are the same as in
graph_tool.inference.BlockState.multicanonical_sweep().
- multiflip_mcmc_sweep(**kwargs)[source]#
Perform
nitersweeps of a Metropolis-Hastings acceptance-rejection MCMC with multiple moves to sample hierarchical network partitions.The arguments accepted are the same as in
graph_tool.inference.BlockState.multiflip_mcmc_sweep().If the parameter
cis a scalar, the values used at each level arec * 2 ** lforlin the range[0, L-1]. Optionally, a list of values may be passed instead, which specifies the value ofc[l]to be used at each level.Warning
This function performs
nitersweeps at each hierarchical level once. This means that in order for the chain to equilibrate, we need to call this function several times, i.e. it is not enough to call it once with a large value ofniter.
- multilevel_mcmc_sweep(**kwargs)[source]#
Perform
nitersweeps of a Metropolis-Hastings acceptance-rejection MCMC with multilevel moves to sample hierarchical network partitions.The arguments accepted are the same as in
graph_tool.inference.BlockState.multilevel_mcmc_sweep().If the parameter
cis a scalar, the values used at each level arec * 2 ** lforlin the range[0, L-1]. Optionally, a list of values may be passed instead, which specifies the value ofc[l]to be used at each level.Warning
This function performs
nitersweeps at each hierarchical level once. This means that in order for the chain to equilibrate, we need to call this function several times, i.e. it is not enough to call it once with a large value ofniter.
- project_level(l)[source]#
Project the partition at level
lonto the lowest level, and return the corresponding state.