graph_tool.inference.EpidemicsBlockState#

class graph_tool.inference.EpidemicsBlockState(g, s, beta, r, r_v=None, global_beta=None, active=None, t=[], exposed=False, aE=nan, nested=True, state_args={}, bstate=None, self_loops=False, **kwargs)[source]#

Bases: DynamicsBlockStateBase

Inference state for network reconstruction based on epidemic dynamics, using the stochastic block model as a prior.

Parameters:
gGraph

Initial graph state.

slist of VertexPropertyMap

Collection of time-series with node states over time. Each entry in this list must be a VertexPropertyMap with type vector<int> containing the states of each node in each time step. A value of 1 means infected and 0 susceptible. Other values are allowed (e.g. for recovered), but their actual value is unimportant for reconstruction.

If the parameter t below is given, each property map value for a given node should contain only the states for the same points in time given by that parameter.

betafloat or EdgePropertyMap

Initial value of the global or local transmission probability for each edge.

rfloat

Spontaneous infection probability.

r_vVertexPropertyMap (optional, default: None)

If given, this will set the initial spontaneous infection probability for each node, and trigger the use of a model where this quantity is in principle different for each node.

global_betafloat (optional, default: None)

If provided, and beta is None this will trigger the use of a model where all transmission probabilities on edges are the same, and given (initially) by this value.

tlist of VertexPropertyMap (optional, default: [])

If nonempty, this allows for a compressed representation of the time-series parameter s, corresponding only to points in time where the state of each node changes. Each entry in this list must be a VertexPropertyMap with type vector<int> containing the points in time where the state of each node change. The corresponding state of the nodes at these times are given by parameter s.

activelist of VertexPropertyMap (optional, default: None)

If given, this specifies the points in time where each node is “active”, and prepared to change its state according to the state of its neighbors. Each entry in this list must be a VertexPropertyMap with type vector<int> containing the states of each node in each time step. A value of 1 means active and 0 inactive.

exposedboolean (optional, default: False)

If True, the data is supposed to come from a SEI, SEIR, etc. model, where a susceptible node (valued 0) first transits to an exposed state (valued -1) upon transmission, before transiting to the infective state (valued 1).

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 inferred graph can contain self-loops.

References

[peixoto-network-2019]

Tiago P. Peixoto, “Network reconstruction and community detection from dynamics”, Phys. Rev. Lett. 123 128301 (2019), DOI: 10.1103/PhysRevLett.123.128301 [sci-hub, @tor], arXiv: 1903.10833

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, x[, 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_x()

Return latent edge transmission probabilities.

mcmc_sweep([r, p, pstep, h, hstep, xstep, ...])

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_params(params)

Sets the model parameters via the dictionary params.

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, x, 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_x()[source]#

Return latent edge transmission probabilities.

mcmc_sweep(r=0.5, p=0.1, pstep=0.1, h=0.1, hstep=1, xstep=0.1, multiflip=True, **kwargs)[source]#

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 parameter h controls the relative probability with which moves for the parameters r_v will be attempted, and hstep is the size of the step. The parameter p controls the relative probability with which moves for the parameters global_beta and r will be attempted, and pstep is the size of the step. The paramter xstep determines the size of the attempted steps for the edge transmission probabilities.

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

Sets the model parameters via the dictionary params.

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