graph_tool.inference.MeasuredBlockState#

class graph_tool.inference.MeasuredBlockState(g, n, x, n_default=1, x_default=0, lp=nan, lq=nan, fn_params={'alpha': 1, 'beta': 1}, fp_params={'mu': 1, 'nu': 1}, aE=nan, nested=True, state_args={}, bstate=None, self_loops=False, **kwargs)[source]#

Bases: UncertainBaseState

Inference state of a measured graph, using the stochastic block model as a prior.

Parameters:
gGraph

Measured graph.

nEdgePropertyMap

Edge property map of type int, containing the total number of measurements for each edge.

xEdgePropertyMap

Edge property map of type int, containing the number of positive measurements for each edge.

n_defaultint (optional, default: 1)

Total number of measurements for each non-edge.

x_defaultint (optional, default: 0)

Total number of positive measurements for each non-edge.

lpfloat (optional, default: NaN)

Log-probability of missing edges (false negatives). If given as NaN, it is assumed this is an unknown sampled from a Beta distribution, with hyperparameters given by fn_params`.  Otherwise the  values of ``fn_params are ignored.

lqfloat (optional, default: NaN)

Log-probability of spurious edges (false positives). If given as NaN, it is assumed this is an unknown sampled from a Beta distribution, with hyperparameters given by fp_params`. Otherwise the values of ``fp_params are ignored.

fn_paramsdict (optional, default: dict(alpha=1, beta=1))

Beta distribution hyperparameters for the probability of missing edges (false negatives).

fp_paramsdict (optional, default: dict(mu=1, nu=1))

Beta distribution hyperparameters for the probability of spurious edges (false positives).

aEfloat (optional, default: NaN)

Expected total number of edges used in prior. If NaN, a flat prior will be used instead.

nestedboolean (optional, default: True)

If True, a NestedBlockState will be used, otherwise BlockState.

state_argsdict (optional, default: {})

Arguments to be passed to NestedBlockState or BlockState.

bstateNestedBlockState or BlockState (optional, default: None)

If passed, this will be used to initialize the block state directly.

self_loopsbool (optional, default: False)

If True, it is assumed that the uncertain graph can contain self-loops.

References

[peixoto-reconstructing-2018]

Tiago P. Peixoto, “Reconstructing networks with unknown and heterogeneous errors”, Phys. Rev. X 8 041011 (2018). DOI: 10.1103/PhysRevX.8.041011 [sci-hub, @tor], arXiv: 1806.07956

Methods

collect_marginal([g])

Collect marginal inferred network during MCMC runs.

collect_marginal_multigraph([g])

Collect marginal latent multigraph during MCMC runs.

copy(**kwargs)

Return a copy of the state.

entropy([latent_edges, density])

Return the entropy, i.e. negative log-likelihood.

get_block_state()

Return the underlying block state, which can be either BlockState or NestedBlockState.

get_edge_prob(u, v[, entropy_args, epsilon])

Return conditional posterior log-probability of edge \((u,v)\).

get_edges_prob(elist[, entropy_args, epsilon])

Return conditional posterior log-probability of an edge list, with shape \((E,2)\).

get_graph()

Return the current inferred graph.

get_p_posterior()

Get beta distribution parameters for the posterior probability of missing edges.

get_q_posterior()

Get beta distribution parameters for the posterior probability of spurious edges.

mcmc_sweep([r, multiflip])

Perform sweeps of a Metropolis-Hastings acceptance-rejection sampling MCMC to sample network partitions and latent edges.

multiflip_mcmc_sweep(**kwargs)

Alias for mcmc_sweep() with multiflip=True.

set_hparams(alpha, beta, mu, nu)

Set edge and non-edge hyperparameters.

set_state(g, w)

virtual_add_edge(u, v[, entropy_args])

virtual_remove_edge(u, v[, entropy_args])

collect_marginal(g=None)#

Collect marginal inferred network during MCMC runs.

Parameters:
gGraph (optional, default: None)

Previous marginal graph.

Returns:
gGraph

New marginal graph, with internal edge EdgePropertyMap "eprob", containing the marginal probabilities for each edge.

Notes

The posterior marginal probability of an edge \((i,j)\) is defined as

\[\pi_{ij} = \sum_{\boldsymbol A}A_{ij}P(\boldsymbol A|\boldsymbol D)\]

where \(P(\boldsymbol A|\boldsymbol D)\) is the posterior probability given the data.

collect_marginal_multigraph(g=None)#

Collect marginal latent multigraph during MCMC runs.

Parameters:
gGraph (optional, default: None)

Previous marginal multigraph.

Returns:
gGraph

New marginal graph, with internal edge EdgePropertyMap "w" and "wcount", containing the edge multiplicities and their respective counts.

Notes

The mean posterior marginal multiplicity distribution of a multi-edge \((i,j)\) is defined as

\[\pi_{ij}(w) = \sum_{\boldsymbol G}\delta_{w,G_{ij}}P(\boldsymbol G|\boldsymbol D)\]

where \(P(\boldsymbol G|\boldsymbol D)\) is the posterior probability of a multigraph \(\boldsymbol G\) given the data.

copy(**kwargs)[source]#

Return a copy of the state.

entropy(latent_edges=True, density=True, **kwargs)#

Return the entropy, i.e. negative log-likelihood.

get_block_state()#

Return the underlying block state, which can be either BlockState or NestedBlockState.

get_edge_prob(u, v, entropy_args={}, epsilon=1e-08)#

Return conditional posterior log-probability of edge \((u,v)\).

get_edges_prob(elist, entropy_args={}, epsilon=1e-08)#

Return conditional posterior log-probability of an edge list, with shape \((E,2)\).

get_graph()#

Return the current inferred graph.

get_p_posterior()[source]#

Get beta distribution parameters for the posterior probability of missing edges.

get_q_posterior()[source]#

Get beta distribution parameters for the posterior probability of spurious edges.

mcmc_sweep(r=0.5, multiflip=True, **kwargs)#

Perform sweeps of a Metropolis-Hastings acceptance-rejection sampling MCMC to sample network partitions and latent edges. The parameter r controls the probability with which edge move will be attempted, instead of partition moves. The remaining keyword parameters will be passed to mcmc_sweep() or multiflip_mcmc_sweep(), if multiflip=True.

multiflip_mcmc_sweep(**kwargs)#

Alias for mcmc_sweep() with multiflip=True.

set_hparams(alpha, beta, mu, nu)[source]#

Set edge and non-edge hyperparameters.

set_state(g, w)#
virtual_add_edge(u, v, entropy_args={})#
virtual_remove_edge(u, v, entropy_args={})#