LatentClosureBlockState#
- class graph_tool.inference.LatentClosureBlockState(g, L=1, b=None, aE=nan, nested=True, state_args={}, g_orig=None, ew=None, ex=None, entropy_args={}, **kwargs)[source]#
Bases:
LatentLayerBaseState
Inference state of the stochastic block model with latent triadic closure edges.
- Parameters:
- g
Graph
Observed graph.
- L
int
(optional, default:1
) Maximum number of triadic closure generations.
- b
VertexPropertyMap
(optional, default:None
) Inital partition (or hierarchical partition
nested=True
).- aE
float
(optional, default:NaN
) Expected total number of edges used in prior. If
NaN
, a flat prior will be used instead.- nested
boolean
(optional, default:True
) If
True
, aNestedBlockState
will be used, otherwiseBlockState
.- state_args
dict
(optional, default:{}
) Arguments to be passed to
NestedBlockState
orBlockState
.- g_orig
Graph
(optional, default:None
) Original graph, if
g
is used to initialize differently from a graph with no triadic closure edges.- ewlist of
EdgePropertyMap
(optional, default:None
) List of edge property maps of
g
, containing the initial weights (counts) at each triadic generation.- exlist of
EdgePropertyMap
(optional, default:None
) List of edge property maps of
g
, each containing a list of integers with the ego graph memberships of every edge, for every triadic generation.- entropy_args: ``dict`` (optional, default: ``{}``)
Override default arguments for
entropy()
method and releated operations.
- g
References
[peixoto-disentangling-2022]Tiago P. Peixoto, “Disentangling homophily, community structure and triadic closure in networks”, Phys. Rev. X 12, 011004 (2022), DOI: 10.1103/PhysRevX.12.011004 [sci-hub, @tor], arXiv: 2101.02510
Methods
collect_marginal
([gs])Collect marginal inferred network during MCMC runs.
Collect marginal latent multigraphs during MCMC runs.
copy
(**kwargs)Return a copy of the state.
entropy
([latent_edges, density, aE, sbm])Return the entropy, i.e. negative log-likelihood.
get_ec
([ew])Return edge property map with layer membership.
Return the current default values for the parameters of the function
entropy()
, together with other operations that depend on them.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()
withmultiflip=True
.Reset the current default values for the parameters of the function
entropy()
, together with other operations that depend on them.sample_graph
([sample_sbm, canonical_sbm, ...])Sample graph from inferred 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.- collect_marginal(gs=None)[source]#
Collect marginal inferred network during MCMC runs.
- Parameters:
- glist of
Graph
(optional, default:None
) Previous marginal graphs.
- glist of
- Returns:
- glist
Graph
New list of marginal graphs, each with internal
EdgePropertyMap
"eprob"
, containing the marginal probabilities for each edge, as well asVertexPropertyMap
"t"
,"m"
,"c"
, containing the average number of closures, open triads, and fraction of closed triads on each node.
- glist
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.
This function returns a list with the marginal graphs for every layer.
- collect_marginal_multigraph(gs=None)#
Collect marginal latent multigraphs during MCMC runs.
- Parameters:
- glist of
Graph
(optional, default:None
) Previous marginal multigraphs.
- glist of
- Returns:
- glist of
Graph
New marginal multigraphs, each with internal edge
EdgePropertyMap
"w"
and"wcount"
, containing the edge multiplicities and their respective counts.
- glist of
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.
This function returns a list with the marginal graphs for every layer.
- entropy(latent_edges=True, density=False, aE=1.0, sbm=True, **kwargs)#
Return the entropy, i.e. negative log-likelihood.
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_ec(ew=None)#
Return edge property map with layer membership.
- get_entropy_args()#
Return the current default values for the parameters of the function
entropy()
, together with other operations that depend on them.
- 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 tomcmc_sweep()
ormultiflip_mcmc_sweep()
, ifmultiflip=True
.
- multiflip_mcmc_sweep(**kwargs)#
Alias for
mcmc_sweep()
withmultiflip=True
.
- 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(sample_sbm=True, canonical_sbm=False, sample_params=True, canonical_closure=True)[source]#
Sample graph from inferred model.
- Parameters:
- sample_sbm
boolean
(optional, default:True
) If
True
, the substrate network will be sampled anew from the SBM parameters. Otherwise, it will be the same as the current posterior state.- canonical_sbm
boolean
(optional, default:False
) If
True
, the canonical SBM will be used, otherwise the microcanonical SBM will be used.- sample_params
bool
(optional, default:True
) If
True
, andcanonical_sbm == True
the count parameters (edges between groups and node degrees) will be sampled from their posterior distribution conditioned on the actual state. Otherwise, their maximum-likelihood values will be used.- canonical_closure
boolean
(optional, default:True
) If
True
, the canonical version of triadic clousre will be used (i.e. conditioned on a probability), otherwise the microcanonical version will be used (i.e. conditional on the count number).
- sample_sbm
- Returns:
- ulist
Graph
Sampled graph, with internal edge
EdgePropertyMap
"gen"
, containing the triadic generation of each edge.
- ulist