uncertain_blockmodel.py 43.1 KB
Newer Older
1 2 3 4 5
#! /usr/bin/env python
# -*- coding: utf-8 -*-
#
# graph_tool -- a general graph manipulation python module
#
Tiago Peixoto's avatar
Tiago Peixoto committed
6
# Copyright (C) 2006-2018 Tiago de Paula Peixoto <tiago@skewed.de>
7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25
#
# This program is free software: you can redistribute it and/or modify
# it under the terms of the GNU General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# This program is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
# GNU General Public License for more details.
#
# You should have received a copy of the GNU General Public License
# along with this program.  If not, see <http://www.gnu.org/licenses/>.

from __future__ import division, absolute_import, print_function
import sys
if sys.version_info < (3,):
    range = xrange

26 27
from .. import _degree, _prop, Graph, GraphView, libcore, _get_rng, \
    PropertyMap, edge_endpoint_property
28 29 30 31 32 33 34 35

from .. dl_import import dl_import
dl_import("from . import libgraph_tool_inference as libinference")

from . blockmodel import *
from . nested_blockmodel import *
from . blockmodel import _bm_test

36 37
import collections

38
def get_uentropy_args(kargs):
39
    ea = get_entropy_args(kargs, ignore=["latent_edges", "density"])
40 41
    uea = libinference.uentropy_args(ea)
    uea.latent_edges = kargs.get("latent_edges", True)
42
    uea.density = kargs.get("density", True)
43 44
    return uea

45
class UncertainBaseState(object):
46 47
    r"""Base state for uncertain network inference."""

48
    def __init__(self, g, nested=True, state_args={}, bstate=None,
49
                 self_loops=False, init_empty=False):
50 51 52 53

        self.g = g

        if bstate is None:
54 55 56
            if init_empty:
                self.u = Graph(directed=g.is_directed())
                self.u.add_vertex(g.num_vertices())
57 58 59 60 61
                self.eweight = self.u.new_ep("int", val=1)
            elif "g" in state_args:
                self.u = state_args.pop("g")
                self.eweight = state_args.pop("eweight",
                                              self.u.new_ep("int", val=1))
62 63
            else:
                self.u = g.copy()
64
                self.eweight = self.u.new_ep("int", val=1)
65 66 67 68 69 70
        else:
            self.u = bstate.g
            if nested:
                self.eweight = bstate.levels[0].eweight
            else:
                self.eweight = bstate.eweight
71
        self.u.set_fast_edge_removal()
72 73 74 75 76 77 78 79 80 81

        self.self_loops = self_loops
        N = self.u.num_vertices()
        if self.u.is_directed():
            if self_loops:
                M = N * N
            else:
                M = N * (N - 1)
        else:
            if self_loops:
82
                M = (N * (N + 1)) / 2
83
            else:
84
                M = (N * (N - 1)) / 2
85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107

        self.M = M

        if bstate is None:
            if nested:
                state_args["state_args"] = state_args.get("state_args", {})
                state_args["state_args"]["eweight"] = self.eweight
                state_args["sampling"] = True
                self.nbstate = NestedBlockState(self.u, **dict(state_args,
                                                               sampling=True))
                self.bstate = self.nbstate.levels[0]
            else:
                self.nbstate = None
                self.bstate = BlockState(self.u, eweight=self.eweight,
                                         **state_args)
        else:
            if nested:
                self.nbstate = bstate
                self.bstate = bstate.levels[0]
            else:
                self.nbstate = None
                self.bstate = bstate

108
        edges = self.g.get_edges()
109 110 111 112 113 114
        edges = numpy.concatenate((edges,
                                   numpy.ones(edges.shape,
                                              dtype=edges.dtype) * (N + 1)))
        self.slist = Vector_size_t(init=edges[:,0])
        self.tlist = Vector_size_t(init=edges[:,1])

115 116
        init_q_cache()

117
    def get_block_state(self):
118 119 120 121 122 123 124 125
        """Return the underlying block state, which can be either
        :class:`~graph_tool.inference.blockmodel.BlockState` or
        :class:`~graph_tool.inference.nested_blockmodel.NestedBlockState`.
        """
        if self.nbstate is None:
            return self.bstate
        else:
            return self.nbstate
126

127
    def entropy(self, latent_edges=True, density=True, **kwargs):
128
        """Return the entropy, i.e. negative log-likelihood."""
129
        S = self._state.entropy(latent_edges, density)
130
        if self.nbstate is None:
131
            S += self.bstate.entropy(**kwargs)
132
        else:
133
            S += self.nbstate.entropy(**kwargs)
134 135

        if kwargs.get("test", True) and _bm_test():
136
            args = kwargs.copy()
137 138 139 140 141 142 143 144 145 146 147
            assert not isnan(S) and not isinf(S), \
                "invalid entropy %g (%s) " % (S, str(args))
            args["test"] = False
            state_copy = self.copy()
            Salt = state_copy.entropy(latent_edges, density, **args)

            assert math.isclose(S, Salt, abs_tol=1e-8), \
                "entropy discrepancy after copying (%g %g %g)" % (S, Salt,
                                                                  S - Salt)
        return S

148 149 150 151 152 153 154 155 156
    def virtual_remove_edge(self, u, v, entropy_args={}):
        dentropy_args = dict(self.bstate._entropy_args, **entropy_args)
        entropy_args = get_uentropy_args(dentropy_args)
        return self._state.remove_edge_dS(int(u), int(v), entropy_args)

    def virtual_add_edge(self, u, v, entropy_args={}):
        dentropy_args = dict(self.bstate._entropy_args, **entropy_args)
        entropy_args = get_uentropy_args(dentropy_args)
        return self._state.add_edge_dS(int(u), int(v), entropy_args)
157 158

    def _algo_sweep(self, algo, r=.5, **kwargs):
159
        kwargs = kwargs.copy()
160 161
        beta = kwargs.get("beta", 1.)
        niter = kwargs.get("niter", 1)
162
        edges_only = kwargs.pop("edges_only", False)
163 164 165 166 167
        verbose = kwargs.get("verbose", False)
        slist = self.slist
        tlist = self.tlist
        dentropy_args = dict(self.bstate._entropy_args,
                             **kwargs.get("entropy_args", {}))
168
        entropy_args = get_uentropy_args(dentropy_args)
169 170
        kwargs.get("entropy_args", {}).pop("latent_edges", None)
        kwargs.get("entropy_args", {}).pop("density", None)
171
        state = self._state
172 173 174 175 176 177

        mcmc_state = DictState(dict(kwargs, **locals()))

        kwargs.pop("xlog", None)
        kwargs.pop("xstep", None)
        kwargs.pop("xdefault", None)
178 179 180 181

        if _bm_test():
            Si = self.entropy(**dentropy_args)

182 183 184 185 186 187 188 189 190
        if self.nbstate is None:
            self.bstate._clear_egroups()
        else:
            self.nbstate._clear_egroups()
        if numpy.random.random() < r:
            edges = True
            dS, nattempts, nmoves = self._mcmc_sweep(mcmc_state)
        else:
            edges = False
191
            if self.nbstate is None:
192
                dS, nattempts, nmoves = algo(self.bstate, **kwargs)
193
            else:
194 195 196 197 198 199 200
                dS, nattempts, nmoves = algo(self.nbstate, **kwargs)

        if _bm_test():
            Sf = self.entropy(**dentropy_args)
            assert math.isclose(dS, (Sf - Si), abs_tol=1e-8), \
                "inconsistent entropy delta %g (%g): %s %s" % (dS, Sf - Si, edges,
                                                               str(dentropy_args))
201 202 203 204 205 206 207 208

        return dS, nattempts, nmoves

    def mcmc_sweep(self, r=.5, **kwargs):
        r"""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
209
        :meth:`~graph_tool.inference.blockmodel.BlockState.mcmc_sweep`.
210 211 212 213 214 215 216 217 218 219
        """

        return self._algo_sweep(lambda s, **kw: s.mcmc_sweep(**kw),
                                r=r, **kwargs)

    def multiflip_mcmc_sweep(self, r=.5, **kwargs):
        r"""Perform sweeps of a multiflip 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
220
        will be passed to :meth:`~graph_tool.inference.blockmodel.BlockState.multiflip_mcmc_sweep`.
221 222 223 224 225 226 227 228
        """

        return self._algo_sweep(lambda s, **kw: s.multiflip_mcmc_sweep(**kw),
                                r=r, **kwargs)

    def get_edge_prob(self, u, v, entropy_args={}, epsilon=1e-8):
        r"""Return conditional posterior probability of edge :math:`(u,v)`."""
        entropy_args = dict(self.bstate._entropy_args, **entropy_args)
229
        ea = get_uentropy_args(entropy_args)
230
        return self._state.get_edge_prob(u, v, ea, epsilon)
231

232 233 234 235
    def get_edges_prob(self, elist, entropy_args={}, epsilon=1e-8):
        r"""Return conditional posterior probability of an edge list, with
        shape :math:`(E,2)`."""
        entropy_args = dict(self.bstate._entropy_args, **entropy_args)
236
        ea = get_uentropy_args(entropy_args)
237 238 239 240
        elist = numpy.asarray(elist, dtype="uint64")
        probs = numpy.zeros(elist.shape[0])
        self._state.get_edges_prob(elist, probs, ea, epsilon)
        return probs
241

242 243 244 245 246 247 248 249 250 251 252
    def get_graph(self):
        r"""Return the current inferred graph."""
        if self.self_loops:
            u = GraphView(self.u, efilt=self.eweight.fa > 0)
        else:
            es = edge_endpoint_property(self.u, self.u.vertex_index, "source")
            et = edge_endpoint_property(self.u, self.u.vertex_index, "target")
            u = GraphView(self.u, efilt=numpy.logical_and(self.eweight.fa > 0,
                                                          es.fa != et.fa))
        return u

253 254 255 256 257 258
    def collect_marginal(self, g=None):
        r"""Collect marginal inferred network during MCMC runs.

        Parameters
        ----------
        g : :class:`~graph_tool.Graph` (optional, default: ``None``)
259
            Previous marginal graph.
260 261 262 263 264 265 266 267 268

        Returns
        -------
        g : :class:`~graph_tool.Graph`
            New marginal graph, with internal edge :class:`~graph_tool.PropertyMap`
            ``"eprob"``, containing the marginal probabilities for each edge.

        Notes
        -----
269
        The posterior marginal probability of an edge :math:`(i,j)` is defined as
270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288

        .. math::

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

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

        """

        if g is None:
            g = Graph(directed=self.g.is_directed())
            g.add_vertex(self.g.num_vertices())
            g.gp.count = g.new_gp("int", 0)
            g.ep.count = g.new_ep("int")

        if "eprob" not in g.ep:
            g.ep.eprob = g.new_ep("double")

289
        u = self.get_graph()
290
        libinference.collect_marginal(g._Graph__graph,
291
                                      u._Graph__graph,
292 293 294 295 296 297
                                      _prop("e", g, g.ep.count))
        g.gp.count += 1
        g.ep.eprob.fa = g.ep.count.fa
        g.ep.eprob.fa /= g.gp.count
        return g

298
class UncertainBlockState(UncertainBaseState):
299 300 301 302 303 304
    r"""The inference state of an uncertain graph, using the stochastic block model
    as a prior.

    Parameters
    ----------
    g : :class:`~graph_tool.Graph`
305
        Measured graph.
306 307 308 309 310 311 312 313
    q : :class:`~graph_tool.PropertyMap`
        Edge probabilities in range :math:`[0,1]`.
    q_default : ``float`` (optional, default: ``0.``)
        Non-edge probability in range :math:`[0,1]`.
    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``)
314
        If ``True``, a :class:`~graph_tool.inference.nested_blockmodel.NestedBlockState`
315
        will be used, otherwise
316
        :class:`~graph_tool.inference.blockmodel.BlockState`.
317 318
    state_args : ``dict`` (optional, default: ``{}``)
        Arguments to be passed to
319 320 321
        :class:`~graph_tool.inference.nested_blockmodel.NestedBlockState` or
        :class:`~graph_tool.inference.blockmodel.BlockState`.
    bstate : :class:`~graph_tool.inference.nested_blockmodel.NestedBlockState` or :class:`~graph_tool.inference.blockmodel.BlockState`  (optional, default: ``None``)
322 323 324 325 326 327
        If passed, this will be used to initialize the block state
        directly.
    self_loops : bool (optional, default: ``False``)
        If ``True``, it is assumed that the uncertain graph can contain
        self-loops.

328 329 330 331 332
    References
    ----------
    .. [peixoto-reconstructing-2018] Tiago P. Peixoto, "Reconstructing networks
       with unknown and heterogeneous errors" :arxiv:`1806.07956`

333 334
    """

335
    def __init__(self, g, q, q_default=0., aE=numpy.nan, nested=True, state_args={},
336
                 bstate=None, self_loops=False, **kwargs):
337 338 339 340

        super(UncertainBlockState, self).__init__(g, nested=nested,
                                                  state_args=state_args,
                                                  bstate=bstate,
341 342
                                                  self_loops=self_loops,
                                                  **kwargs)
343 344 345
        self._q = q
        self._q_default = q_default

346
        self.p = (q.fa.sum() + (self.M - g.num_edges()) * q_default) / self.M
347 348

        self.q = self.g.new_ep("double", vals=log(q.fa) - log1p(-q.fa))
349
        self.q.fa -= log(self.p) - log1p(-self.p)
350 351 352 353 354 355 356 357 358 359
        if q_default > 0:
            self.q_default = log(q_default) - log1p(q_default)
            self.q_default -= log(self.p) - log1p(-self.p)
        else:
            self.q_default = -numpy.inf

        self.S_const = (log1p(-q.fa[q.fa<1]).sum() +
                        log1p(-q_default) * (self.M - self.g.num_edges())
                        - self.M * log1p(-self.p))

360 361
        self.aE = aE
        if numpy.isnan(aE):
362 363 364
            self.E_prior = False
        else:
            self.E_prior = True
365 366 367 368 369

        self._state = libinference.make_uncertain_state(self.bstate._state,
                                                        self)
    def __getstate__(self):
        return dict(g=self.g, q=self._q, q_default=self._q_default,
370
                    aE=self.aE, nested=self.nbstate is not None,
371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394
                    bstate=(self.nbstate.copy() if self.nbstate is not None else
                            self.bstate.copy()), self_loops=self.self_loops)

    def __setstate__(self, state):
        self.__init__(**state)

    def copy(self, **kwargs):
        """Return a copy of the state."""
        return UncertainBlockState(**dict(self.__getstate__(), **kwargs))

    def __copy__(self):
        return self.copy()

    def __repr__(self):
        return "<UncertainBlockState object with %s, at 0x%x>" % \
            (self.nbstate if self.nbstate is not None else self.bstate,
             id(self))

    def _mcmc_sweep(self, mcmc_state):
        return libinference.mcmc_uncertain_sweep(mcmc_state,
                                                 self._state,
                                                 _get_rng())

class MeasuredBlockState(UncertainBaseState):
395 396 397 398 399 400
    r"""The inference state of a measured graph, using the stochastic block model as
    a prior.

    Parameters
    ----------
    g : :class:`~graph_tool.Graph`
401
        Measured graph.
402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421
    n : :class:`~graph_tool.PropertyMap`
        Edge property map of type ``int``, containing the total number of
        measurements for each edge.
    x : :class:`~graph_tool.PropertyMap`
        Edge property map of type ``int``, containing the number of
        positive measurements for each edge.
    n_default : ``int`` (optional, default: ``1``)
        Total number of measurements for each non-edge.
    x_default : ``int`` (optional, default: ``1``)
        Total number of positive measurements for each non-edge.
    fn_params : ``dict`` (optional, default: ``dict(alpha=1, beta=1)``)
        Beta distribution hyperparameters for the probability of missing
        edges (false negatives).
    fp_params : ``dict`` (optional, default: ``dict(mu=1, nu=1)``)
        Beta distribution hyperparameters for the probability of spurious
        edges (false positives).
    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``)
422
        If ``True``, a :class:`~graph_tool.inference.nested_blockmodel.NestedBlockState`
423
        will be used, otherwise
424
        :class:`~graph_tool.inference.blockmodel.BlockState`.
425 426
    state_args : ``dict`` (optional, default: ``{}``)
        Arguments to be passed to
427 428 429
        :class:`~graph_tool.inference.nested_blockmodel.NestedBlockState` or
        :class:`~graph_tool.inference.blockmodel.BlockState`.
    bstate : :class:`~graph_tool.inference.nested_blockmodel.NestedBlockState` or :class:`~graph_tool.inference.blockmodel.BlockState`  (optional, default: ``None``)
430 431 432 433 434
        If passed, this will be used to initialize the block state
        directly.
    self_loops : bool (optional, default: ``False``)
        If ``True``, it is assumed that the uncertain graph can contain
        self-loops.
435 436 437 438 439

    References
    ----------
    .. [peixoto-reconstructing-2018] Tiago P. Peixoto, "Reconstructing networks
       with unknown and heterogeneous errors" :arxiv:`1806.07956`
440 441
    """

442 443
    def __init__(self, g, n, x, n_default=1, x_default=0,
                 fn_params=dict(alpha=1, beta=1), fp_params=dict(mu=1, nu=1),
444
                 aE=numpy.nan, nested=True, state_args={}, bstate=None,
445
                 self_loops=False, **kwargs):
446 447 448

        super(MeasuredBlockState, self).__init__(g, nested=nested,
                                                 state_args=state_args,
449
                                                 bstate=bstate, **kwargs)
450

451 452
        self.aE = aE
        if numpy.isnan(aE):
453 454 455
            self.E_prior = False
        else:
            self.E_prior = True
456 457 458 459 460

        self.n = n
        self.x = x
        self.n_default = n_default
        self.x_default = x_default
461 462 463 464
        self.alpha = fn_params.get("alpha", 1)
        self.beta = fn_params.get("beta", 1)
        self.mu = fp_params.get("mu", 1)
        self.nu = fp_params.get("nu", 1)
465 466 467 468 469 470 471

        self._state = libinference.make_measured_state(self.bstate._state,
                                                       self)

    def __getstate__(self):
        return dict(g=self.g, n=self.n, x=self.x, n_default=self.n_default,
                    x_default=self.x_default,
472
                    fn_params=dict(alpha=self.alpha, beta=self.beta),
473
                    fp_params=dict(mu=self.mu, nu=self.nu), aE=self.aE,
474
                    nested=self.nbstate is not None,
475 476
                    bstate=(self.nbstate if self.nbstate is not None
                            else self.bstate), self_loops=self.self_loops)
477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496

    def __setstate__(self, state):
        self.__init__(**state)

    def copy(self, **kwargs):
        """Return a copy of the state."""
        return MeasuredBlockState(**dict(self.__getstate__(), **kwargs))

    def __repr__(self):
        return "<MeasuredBlockState object with %s, at 0x%x>" % \
            (self.nbstate if self.nbstate is not None else self.bstate,
             id(self))

    def _mcmc_sweep(self, mcmc_state):
        return libinference.mcmc_measured_sweep(mcmc_state,
                                                self._state,
                                                _get_rng())

    def set_hparams(self, alpha, beta, mu, nu):
        """Set edge and non-edge hyperparameters."""
497
        self._state.set_hparams(alpha, beta, mu, nu)
498 499 500 501 502 503
        self.alpha = alpha
        self.beta = beta
        self.mu = mu
        self.nu = nu

    def get_p_posterior(self):
504
        """Get beta distribution parameters for the posterior probability of missing edges."""
505 506
        T = self._state.get_T()
        M = self._state.get_M()
507
        return M - T + self.alpha, T + self.beta
508 509

    def get_q_posterior(self):
510
        """Get beta distribution parameters for the posterior probability of spurious edges."""
511 512 513 514
        N = self._state.get_N()
        X = self._state.get_X()
        T = self._state.get_T()
        M = self._state.get_M()
515
        return X - T + self.mu, N - X - (M - T) + self.nu
516 517

class MixedMeasuredBlockState(UncertainBaseState):
518 519 520 521 522 523
    r"""The inference state of a measured graph with heterogeneous errors, using the
    stochastic block model as a prior.

    Parameters
    ----------
    g : :class:`~graph_tool.Graph`
524
        Measured graph.
525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544
    n : :class:`~graph_tool.PropertyMap`
        Edge property map of type ``int``, containing the total number of
        measurements for each edge.
    x : :class:`~graph_tool.PropertyMap`
        Edge property map of type ``int``, containing the number of
        positive measurements for each edge.
    n_default : ``int`` (optional, default: ``1``)
        Total number of measurements for each non-edge.
    x_default : ``int`` (optional, default: ``1``)
        Total number of positive measurements for each non-edge.
    fn_params : ``dict`` (optional, default: ``dict(alpha=1, beta=10)``)
        Beta distribution hyperparameters for the probability of missing
        edges (false negatives).
    fp_params : ``dict`` (optional, default: ``dict(mu=1, nu=10)``)
        Beta distribution hyperparameters for the probability of spurious
        edges (false positives).
    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``)
545
        If ``True``, a :class:`~graph_tool.inference.nested_blockmodel.NestedBlockState`
546
        will be used, otherwise
547
        :class:`~graph_tool.inference.blockmodel.BlockState`.
548 549
    state_args : ``dict`` (optional, default: ``{}``)
        Arguments to be passed to
550 551 552
        :class:`~graph_tool.inference.nested_blockmodel.NestedBlockState` or
        :class:`~graph_tool.inference.blockmodel.BlockState`.
    bstate : :class:`~graph_tool.inference.nested_blockmodel.NestedBlockState` or :class:`~graph_tool.inference.blockmodel.BlockState`  (optional, default: ``None``)
553 554 555 556 557 558
        If passed, this will be used to initialize the block state
        directly.
    self_loops : bool (optional, default: ``False``)
        If ``True``, it is assumed that the uncertain graph can contain
        self-loops.

559 560 561 562
    References
    ----------
    .. [peixoto-reconstructing-2018] Tiago P. Peixoto, "Reconstructing networks
       with unknown and heterogeneous errors" :arxiv:`1806.07956`
563 564
    """

565 566
    def __init__(self, g, n, x, n_default=1, x_default=0,
                 fn_params=dict(alpha=1, beta=10), fp_params=dict(mu=1, nu=10),
567
                 aE=numpy.nan, nested=True, state_args={}, bstate=None,
568
                 self_loops=False, **kwargs):
569 570 571

        super(MixedMeasuredBlockState, self).__init__(g, nested=nested,
                                                      state_args=state_args,
572
                                                      bstate=bstate, **kwargs)
573 574
        self.aE = aE
        if numpy.isnan(aE):
575 576 577 578 579 580 581 582 583
            self.E_prior = False
        else:
            self.E_prior = True

        self.n = n
        self.x = x
        self.n_default = n_default
        self.x_default = x_default
        self.alpha = fn_params.get("alpha", 1)
584
        self.beta = fn_params.get("beta", 10)
585
        self.mu = fp_params.get("mu", 1)
586
        self.nu = fp_params.get("nu", 10)
587 588 589 590 591 592 593 594 595 596 597

        self._state = None

        self.q = self.g.new_ep("double")
        self.sync_q()

        self._state = libinference.make_uncertain_state(self.bstate._state,
                                                        self)

    def sync_q(self):
        ra, rb = self.transform(self.n.fa, self.x.fa)
598
        self.q.fa = ra - rb
599
        dra, drb = self.transform(self.n_default, self.x_default)
600
        self.q_default = dra - drb
601

602
        self.S_const = (self.M - self.g.num_edges()) * drb + rb.sum()
603 604 605 606 607 608

        if self._state is not None:
            self._state.set_q_default(self.q_default)
            self._state.set_S_const(self.S_const)

    def transform(self, na, xa):
609 610 611 612
        ra = (scipy.special.betaln(na - xa + self.alpha, xa + self.beta) -
              scipy.special.betaln(self.alpha, self.beta))
        rb = (scipy.special.betaln(xa + self.mu, na - xa + self.nu) -
              scipy.special.betaln(self.mu, self.nu))
613 614 615 616 617 618 619 620 621 622 623 624 625
        return ra, rb

    def set_hparams(self, alpha, beta, mu, nu):
        """Set edge and non-edge hyperparameters."""
        self.alpha = alpha
        self.beta = beta
        self.mu = mu
        self.nu = nu
        self.sync_q()

    def __getstate__(self):
        return dict(g=self.g, n=self.n, x=self.x, n_default=self.n_default,
                    x_default=self.x_default,
626
                    fn_params=dict(alpha=self.alpha, beta=self.beta),
627
                    fp_params=dict(mu=self.mu, nu=self.nu), aE=self.aE,
628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649
                    nested=self.nbstate is not None,
                    bstate=(self.nbstate if self.nbstate is not None
                            else self.bstate), self_loops=self.self_loops)

    def __setstate__(self, state):
        self.__init__(**state)

    def copy(self, **kwargs):
        """Return a copy of the state."""
        return MixedMeasuredBlockState(**dict(self.__getstate__(), **kwargs))

    def __copy__(self):
        return self.copy()

    def __setstate__(self, state):
        self.__init__(**state)

    def __repr__(self):
        return "<MixedMeasuredBlockState object with %s, at 0x%x>" % \
            (self.nbstate if self.nbstate is not None else self.bstate,
             id(self))

650
    def _algo_sweep(self, algo, r=.5, h=.1, hstep=1, niter=1, **kwargs):
651
        if numpy.random.random() < h:
652 653 654 655
            dS = nt = nm = 0
            for i in range(niter):
                hs = [self.alpha, self.beta, self.mu, self.nu]
                j = numpy.random.randint(len(hs))
656

657 658
                f_dh = [max(0, hs[j] - hstep), hs[j] + hstep]
                pf = 1./(f_dh[1] - f_dh[0])
659

660 661
                old_hs = hs[j]
                hs[j] = f_dh[0] + numpy.random.random() * (f_dh[1] - f_dh[0])
662

663 664
                b_dh = [max(0, hs[j] - hstep), hs[j] + hstep]
                pb = 1./min(1, hs[j])
665

666 667
                latent_edges = kwargs.get("entropy_args", {}).get("latent_edges", True)
                density = False
668

669
                Sb = self._state.entropy(latent_edges, density)
670
                self.set_hparams(*hs)
671 672 673 674 675 676 677 678 679 680
                Sa = self._state.entropy(latent_edges, density)

                nt += 1
                if Sa < Sb or numpy.random.random() < exp(-(Sa-Sb) + log(pb) - log(pf)):
                    dS += Sa - Sb
                    nm +=1
                else:
                    hs[j] = old_hs
                    self.set_hparams(*hs)
            return (dS, nt, nm)
681 682 683 684 685 686 687
        else:
            return super(MixedMeasuredBlockState, self)._algo_sweep(algo, r, **kwargs)

    def _mcmc_sweep(self, mcmc_state):
        return libinference.mcmc_uncertain_sweep(mcmc_state,
                                                 self._state,
                                                 _get_rng())
688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097

class DynamicsBlockStateBase(UncertainBaseState):
    def __init__(self, g, s, t, x=None, aE=numpy.nan, nested=True,
                 state_args={}, bstate=None, self_loops=False,
                 **kwargs):
        super(DynamicsBlockStateBase, self).__init__(g, nested=nested,
                                                     state_args=state_args,
                                                     bstate=bstate,
                                                     self_loops=self_loops)
        self.s = [g.own_property(x) for x in s]
        self.t = [g.own_property(x) for x in t]
        if x is None:
            x = self.u.new_ep("double")
        else:
            x = self.u.copy_property(x, g=x.get_graph())
        self.x = x

        self.aE = aE
        if numpy.isnan(aE):
            self.E_prior = False
        else:
            self.E_prior = True

        for k in kwargs.keys():
            v = kwargs[k]
            if isinstance(v, PropertyMap):
                kwargs[k] = g.own_property(v)
            elif (isinstance(v, collections.Iterable) and len(v) > 0 and
                  isinstance(v[0], PropertyMap)):
                kwargs[k] = [g.own_property(x) for x in v]
        self.params = kwargs
        self.os = [ns._get_any() for ns in s]
        self.ot = [nt._get_any() for nt in t]
        self._state = self._make_state()

    def set_params(self, params):
        self.params = dict(self.params, **params)
        self._state.set_params(self.params)

    def __getstate__(self):
        return dict(g=self.g, s=self.s, t=self.t, x=self.x, aE=self.aE,
                    nested=self.nbstate is not None,
                    bstate=(self.nbstate.copy() if self.nbstate is not None else
                            self.bstate.copy()), self_loops=self.self_loops,
                    **self.params)

    def __setstate__(self, state):
        self.__init__(**state)

    def copy(self, **kwargs):
        """Return a copy of the state."""
        return type(self)(**dict(self.__getstate__(), **kwargs))

    def __copy__(self):
        return self.copy()

    def __repr__(self):
        return "<%s object with %s, at 0x%x>" % \
            (self.__class__.__name__,
             self.nbstate if self.nbstate is not None else self.bstate,
             id(self))

    def _move_proposal(self, name, beta, step, rg, transform, entropy_args):
        x = x_orig = self.params[name]

        if isinstance(x, collections.Iterable):
            idx = numpy.random.randint(len(x))
            x = x[idx]
        else:
            idx = None

        if transform is not None:
            x = transform[1](x)

        if rg is not None:
            mi = max(rg[0], x - step)
            ma = min(rg[1], x + step)
        else:
            mi = x - step
            ma = x + step
        nx = numpy.random.random() * (ma - mi) + mi

        a = 0
        if rg is not None:
            a -= -log(ma - mi)
            mi = max(rg[0], nx - step)
            ma = min(rg[1], nx + step)
            a += -log(ma - mi)

        if transform is not None:
            nx = transform[0](nx)

        latent_edges = entropy_args.get("latent_edges", True)
        density = False

        Sb = self._state.entropy(latent_edges, density)

        if idx is not None:
            y = self.params[name]
            y[idx] = nx
            nx = y

        self.set_params({name:nx})
        Sa = self._state.entropy(latent_edges, density)

        a += beta * (Sb - Sa)

        if a > 0 or numpy.random.random() < exp(a):
            self.set_params({name:nx})
            return Sa-Sb, 1, 1

        if idx is not None:
            y = self.params[name]
            y[idx] = x_orig
            x_orig = y
        self.set_params({name:x_orig})

        return 0, 0, 0

    def get_x(self):
        return x

    def get_edge_prob(self, u, v, x, entropy_args={}, epsilon=1e-8):
        r"""Return conditional posterior probability of edge :math:`(u,v)`."""
        entropy_args = dict(self.bstate._entropy_args, **entropy_args)
        ea = get_uentropy_args(entropy_args)
        return self._state.get_edge_prob(u, v, x, ea, epsilon)

    def collect_marginal(self, g=None):
        r"""Collect marginal inferred network during MCMC runs.

        Parameters
        ----------
        g : :class:`~graph_tool.Graph` (optional, default: ``None``)
            Previous marginal graph.

        Returns
        -------
        g : :class:`~graph_tool.Graph`
            New marginal graph, with internal edge :class:`~graph_tool.PropertyMap`
            ``"eprob"``, containing the marginal probabilities for each edge.

        Notes
        -----
        The posterior marginal probability of an edge :math:`(i,j)` is defined as

        .. math::

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

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

        """

        if g is None:
            g = Graph(directed=self.g.is_directed())
            g.add_vertex(self.g.num_vertices())
            g.gp.count = g.new_gp("int", 0)
            g.ep.count = g.new_ep("int")
            g.ep.eprob = g.new_ep("double")

        if "x" not in g.ep:
            g.ep.xsum = g.new_ep("double")
            g.ep.x2sum = g.new_ep("double")
            g.ep.x = g.new_ep("double")
            g.ep.xdev = g.new_ep("double")

        u = self.get_graph()
        x = self.get_x()
        libinference.collect_xmarginal(g._Graph__graph,
                                       u._Graph__graph,
                                       _prop("e", u, x),
                                       _prop("e", g, g.ep.count),
                                       _prop("e", g, g.ep.xsum),
                                       _prop("e", g, g.ep.x2sum))
        g.gp.count += 1
        g.ep.eprob.fa = g.ep.count.fa / g.gp.count
        g.ep.x.fa = g.ep.xsum.fa / g.gp.count
        g.ep.xdev.fa = sqrt(g.ep.x2sum.fa / g.gp.count - g.ep.x.fa ** 2)
        return g

class EpidemicsBlockState(DynamicsBlockStateBase):
    def __init__(self, g, s, beta, r, r_v=None, global_beta=None, active=None,
                 t=[], exposed=False, aE=numpy.nan, nested=True, state_args={},
                 bstate=None, self_loops=False, **kwargs):
        if beta is None:
            beta = global_beta
        x = kwargs.pop("x", None)
        if x is None:
            if bstate is not None:
                if nested:
                    u = bstate.levels[0].g
                else:
                    u = bstate.g
            else:
                u = g
            if isinstance(beta, PropertyMap):
                x = u.own_property(beta.copy())
            else:
                x = u.new_ep("double", val=beta)
            x.fa = log1p(-x.fa)
        super(EpidemicsBlockState, self).__init__(g, s=s, t=t,
                                                  global_beta=global_beta,
                                                  active=active,
                                                  x=x, r=r, r_v=r_v,
                                                  exposed=exposed, aE=aE,
                                                  nested=nested,
                                                  state_args=state_args,
                                                  bstate=bstate,
                                                  self_loops=self_loops,
                                                  **kwargs)
    def _make_state(self):
        return libinference.make_epidemics_state(self.bstate._state, self)

    def __setstate__(self, state):
        beta = state.pop("beta", None)
        if "x" not in state:
            state["global_beta"] = beta
        self.__init__(**state, beta=beta)

    def set_params(self, params):
        self.params = dict(self.params, **params)
        self._state.set_params(self.params)
        beta = self.params["global_beta"]
        if beta is not None:
            self.x.fa = log1p(-beta)

    def get_x(self):
        x = self.x.copy()
        x.fa = 1-exp(x.fa)
        return x

    def _algo_sweep(self, algo, r=.5, p=.1, pstep=.1, h=.1, hstep=1,
                    xstep=.1, niter=1, **kwargs):
        if numpy.random.random() < p:
            if self.params["r_v"] is None:
                h = 0
            if numpy.random.random() < h:
                mcmc_state = DictState(dict(beta=kwargs.get("beta", 1),
                                            hstep=hstep,
                                            verbose=kwargs.get("verbose", False),
                                            niter=niter,
                                            state=self._state))
                return libinference.mcmc_epidemics_sweep_r(mcmc_state,
                                                           self._state,
                                                           _get_rng())
            else:
                dS = nt = nm = 0
                for i in range(niter):
                    if self.params["global_beta"] is not None:
                        ret = self._move_proposal("global_beta",
                                                  kwargs.get("beta", 1),
                                                  pstep, (0, 1),
                                                  (lambda beta: log1p(-beta),
                                                   lambda x: 1-exp(x)),
                                                  kwargs.get("entropy_args", {}))
                        dS += ret[0]
                        nt += ret[1]
                        nm += ret[2]
                    ret = self._move_proposal("r",
                                              kwargs.get("beta", 1),
                                              pstep, (0, 1),
                                              None,
                                              kwargs.get("entropy_args", {}))
                    dS += ret[0]
                    nt += ret[1]
                    nm += ret[2]
                return (dS, nt, nm)
        else:
            beta = self.params["global_beta"]
            if beta is not None:
                xdefault = log1p(-beta)
            else:
                xdefault = 0
            return super(EpidemicsBlockState, self)._algo_sweep(algo, r,
                                                                xlog=True,
                                                                xstep=xstep,
                                                                xdefault=xdefault,
                                                                **kwargs)

    def _mcmc_sweep(self, mcmc_state):
        return libinference.mcmc_epidemics_sweep(mcmc_state, self._state,
                                                 _get_rng())

class IsingBaseBlockState(DynamicsBlockStateBase):
    def __init__(self, g, s, beta, x=None, t=None, h=None, aE=numpy.nan,
                 nested=True, state_args={}, bstate=None, self_loops=False,
                 has_zero=False, **kwargs):
        if isinstance(s, PropertyMap):
            s = [s]
        if isinstance(t, PropertyMap):
            t = [t]
        if t is None:
            t = []
        if h is None:
            h = g.new_vp("double")
        if bstate is not None:
            if nested:
                u = bstate.levels[0].g
            else:
                u = bstate.g
        else:
            u = g
        if x == 1:
            x = u.new_ep("double", 1)
        elif x is None:
            x = u.new_ep("double", numpy.random.random(g.num_edges()) - 1/2)
        super(IsingBaseBlockState, self).__init__(g, s=s, t=t, beta=beta, x=x,
                                                  aE=aE, nested=nested,
                                                  state_args=state_args,
                                                  bstate=bstate,
                                                  self_loops=self_loops, h=h,
                                                  has_zero=has_zero, **kwargs)
    def __setstate__(self, state):
        beta = state.pop("beta", None)
        if "x" not in state:
            state["x"] = 1
        self.__init__(**state, beta=beta)

    def get_x(self):
        x = self.x.copy()
        x.fa *= self.params["beta"]
        return x

    def _algo_sweep(self, algo, r=.5, p=.1, pstep=1, h=.5, hstep=1, niter=1,
                    xstep=1, **kwargs):
        if numpy.random.random() < p:
            if numpy.random.random() < h:
                mcmc_state = DictState(dict(beta=kwargs.get("beta", 1),
                                            hstep=hstep,
                                            verbose=kwargs.get("verbose", False),
                                            niter=niter,
                                            n=numpy.random.randint(len(self.s)),
                                            state=self._state))
                return self._mcmc_sweep_h(mcmc_state)
            else:
                dS = nt = nm = 0
                for i in range(niter):
                    ret = self._move_proposal("beta",
                                              kwargs.get("beta", 1),
                                              pstep, None, None,
                                              kwargs.get("entropy_args", {}))
                    dS += ret[0]
                    nt += ret[1]
                    nm += ret[2]
                return (dS, nt, nm)
        else:
            return super(IsingBaseBlockState, self)._algo_sweep(algo, r,
                                                                xlog=False,
                                                                xstep=xstep,
                                                                xdefault=1,
                                                                **kwargs)

class PseudoIsingBlockState(IsingBaseBlockState):
    def __init__(self, *args, **kwargs):
        super(PseudoIsingBlockState, self).__init__(*args, **kwargs)

    def _make_state(self):
        return libinference.make_pseudo_ising_state(self.bstate._state, self)

    def _mcmc_sweep(self, mcmc_state):
        return libinference.mcmc_pseudo_ising_sweep(mcmc_state, self._state,
                                                    _get_rng())
    def _mcmc_sweep_h(self, mcmc_state):
        return libinference.mcmc_pseudo_ising_sweep_h(mcmc_state, self._state,
                                                      _get_rng())

class IsingGlauberBlockState(IsingBaseBlockState):
    def __init__(self, *args, **kwargs):
        super(IsingGlauberBlockState, self).__init__(*args, **kwargs)

    def _make_state(self):
        return libinference.make_ising_glauber_state(self.bstate._state, self)

    def _mcmc_sweep(self, mcmc_state):
        return libinference.mcmc_ising_glauber_sweep(mcmc_state, self._state,
                                                     _get_rng())

    def _mcmc_sweep_h(self, mcmc_state):
        return libinference.mcmc_ising_glauber_sweep_h(mcmc_state, self._state,
                                                       _get_rng())

class PseudoCIsingBlockState(IsingBaseBlockState):
    def __init__(self, *args, **kwargs):
        super(PseudoCIsingBlockState, self).__init__(*args, **kwargs)

    def _make_state(self):
        return libinference.make_pseudo_cising_state(self.bstate._state, self)

    def _mcmc_sweep(self, mcmc_state):
        return libinference.mcmc_pseudo_cising_sweep(mcmc_state, self._state,
                                                     _get_rng())
    def _mcmc_sweep_h(self, mcmc_state):
        return libinference.mcmc_pseudo_cising_sweep_h(mcmc_state, self._state,
                                                       _get_rng())

class CIsingGlauberBlockState(IsingBaseBlockState):
    def __init__(self, *args, **kwargs):
        super(CIsingGlauberBlockState, self).__init__(*args, **kwargs)

    def _make_state(self):
        return libinference.make_cising_glauber_state(self.bstate._state, self)

    def _mcmc_sweep(self, mcmc_state):
        return libinference.mcmc_cising_glauber_sweep(mcmc_state, self._state,
                                                      _get_rng())
    def _mcmc_sweep_h(self, mcmc_state):
        return libinference.mcmc_cising_glauber_sweep_h(mcmc_state, self._state,
                                                        _get_rng())