TransportMaps.Likelihoods.LikelihoodBase

Module Contents

Classes

LogLikelihood

Abstract class for log-likelihood \(\log \pi({\bf y} \vert {\bf x})\)

AdditiveLogLikelihood

Log-likelihood \(\log \pi({\bf y} \vert {\bf x})=\log\pi({\bf y} - T({\bf x}))\)

AdditiveNormalLogLikelihood

Define the log-likelihood for the additive Gaussian model

AdditiveLinearNormalLogLikelihood

Define the log-likelihood for the additive linear Gaussian model

AdditiveConditionallyLinearNormalLogLikelihood

Define the log-likelihood for the additive linear Gaussian model

LogisticLogLikelihood

Class for log-likelihood of \({\bf y}\), where \(y_i \vert {\bf x}; {\bf f}_i \sim \mathrm{Bernoulli}(\mathrm{logit}^{-1}({\bf f}_i^\top {\bf x})\)

TwoParametersLogisticLogLikelihood

See Stan example

IndependentLogLikelihood

Handling likelihoods in the form \(\pi({\bf y}\vert{\bf x}) = \prod_{i=1}^{n}\pi_i({\bf y}_i\vert{\bf x}_i)\)

Attributes

AdditiveLinearGaussianLogLikelihood

AdditiveConditionallyLinearGaussianLogLikelihood

class TransportMaps.Likelihoods.LikelihoodBase.LogLikelihood(y, dim)[source]

Bases: TransportMaps.Maps.Functionals.Functional

Abstract class for log-likelihood \(\log \pi({\bf y} \vert {\bf x})\)

Note that \(\log\pi:\mathbb{R}^d \rightarrow \mathbb{R}\) is considered a function of \({\bf x}\), while the data \({\bf y}\) is fixed.

Parameters:
  • y (ndarray) – data

  • dim (int) – input dimension $d$

y[source]
_get_y()[source]
_set_y(y)[source]
abstract evaluate(x, *args, **kwargs)[source]

[Abstract] Evaluate \(\log\pi({\bf y} \vert {\bf x})\).

Parameters:

x (ndarray [\(m,d\)]) – evaluation points

Returns:

(ndarray [\(m,1\)]) – function evaluations

abstract grad_x(x, *args, **kwargs)[source]

[Abstract] Evaluate \(\nabla_{\bf x}\log\pi({\bf y} \vert {\bf x})\).

Parameters:

x (ndarray [\(m,d\)]) – evaluation points

Returns:

(ndarray [\(m,1,d\)]) – gradient evaluations

tuple_grad_x(x, cache=None, **kwargs)[source]

Evaluate \(\left(\log\pi({\bf y} \vert {\bf x}),\nabla_{\bf x}\log\pi({\bf y} \vert {\bf x})\right)\).

Parameters:

x (ndarray [\(m,d\)]) – evaluation points

Returns:

(tuple) –

\(\left(\log\pi({\bf y} \vert {\bf x}),\nabla_{\bf x}\log\pi({\bf y} \vert {\bf x})\right)\)

abstract hess_x(x, *args, **kwargs)[source]

[Abstract] Evaluate \(\nabla^2_{\bf x}\log\pi({\bf y} \vert {\bf x})\).

Parameters:

x (ndarray [\(m,d\)]) – evaluation points

Returns:

(ndarray [\(m,1,d,d\)]) – Hessian evaluations

abstract action_hess_x(x, dx, *args, **kwargs)[source]

[Abstract] Evaluate \(\langle \nabla^2_{\bf x}\log\pi({\bf y} \vert {\bf x}),\delta{\bf x}\rangle\).

Parameters:
  • x (ndarray [\(m,d\)]) – evaluation points

  • dx (ndarray [\(m,d\)]) – direction on which to evaluate the Hessian

Returns:

(ndarray [\(m,1,d\)]) – Hessian evaluations

class TransportMaps.Likelihoods.LikelihoodBase.AdditiveLogLikelihood(y, pi, T)[source]

Bases: LogLikelihood

Log-likelihood \(\log \pi({\bf y} \vert {\bf x})=\log\pi({\bf y} - T({\bf x}))\)

Parameters:
  • y (ndarray \([d_y]\)) – data

  • pi (Distribution) – distribution \(\nu_\pi\)

  • T (Map) – map \(T:\mathbb{R}^d\rightarrow\mathbb{R}^{d_y}\)

property isPiCond[source]
get_ncalls_tree(indent='')[source]
get_nevals_tree(indent='')[source]
get_teval_tree(indent='')[source]
update_ncalls_tree(obj)[source]
update_nevals_tree(obj)[source]
update_teval_tree(obj)[source]
reset_counters()[source]
evaluate(x, idxs_slice=slice(None, None, None), cache=None, **kwargs)[source]

Evaluate \(\log\pi({\bf y} \vert {\bf x})\).

Parameters:
  • x (ndarray [\(m,d\)]) – evaluation points

  • cache (dict) – cache

Returns:

(ndarray [\(m,1\)]) – function evaluations

grad_x(x, idxs_slice=slice(None, None, None), cache=None, **kwargs)[source]

Evaluate \(\nabla_{\bf x}\log\pi({\bf y} \vert {\bf x})\).

Parameters:
  • x (ndarray [\(m,d\)]) – evaluation points

  • cache (dict) – cache

Returns:

(ndarray [\(m,1,d\)]) – gradient evaluations

Todo

caching is not implemented

tuple_grad_x(x, idxs_slice=slice(None, None, None), cache=None, **kwargs)[source]

Evaluate \(\log\pi({\bf y} \vert {\bf x}), \nabla_{\bf x}\log\pi({\bf y} \vert {\bf x})\).

Parameters:
  • x (ndarray [\(m,d\)]) – evaluation points

  • cache (dict) – cache

Returns:

(tuple) – function and gradient evaluations

Todo

caching is not implemented

hess_x(x, idxs_slice=slice(None, None, None), cache=None, **kwargs)[source]

Evaluate \(\nabla^2_{\bf x}\log\pi({\bf y} \vert {\bf x})\).

Parameters:
  • x (ndarray [\(m,d\)]) – evaluation points

  • cache (dict) – cache

Returns:

(ndarray [\(m,1,d,d\)]) – Hessian evaluations

Todo

caching is not implemented

action_hess_x(x, dx, idxs_slice=slice(None, None, None), cache=None, **kwargs)[source]

Evaluate \(\langle\nabla^2_{\bf x}\log\pi({\bf y} \vert {\bf x}),\delta{\bf x}\rangle\).

Parameters:
  • x (ndarray [\(m,d\)]) – evaluation points

  • dx (ndarray [\(m,d\)]) – direction on which to evaluate the Hessian

  • cache (dict) – cache

Returns:

(ndarray [\(m,1,d\)]) – action of Hessian evaluations

Todo

caching is not implemented

class TransportMaps.Likelihoods.LikelihoodBase.AdditiveNormalLogLikelihood(y, T, mu, covariance=None, precision=None)[source]

Bases: AdditiveLogLikelihood

Define the log-likelihood for the additive Gaussian model

The model is

\[{\bf y} = {\bf T}({\bf x}) + \varepsilon \;, \quad \varepsilon \sim \mathcal{N}(\mu, \Sigma)\]

where \(T \in \mathbb{R}^{d_y \times d_x}\), \(\mu \in \mathbb{R}^{d_y}\) and \(\Sigma \in \mathbb{R}^{d_y \times d_y}\) is symmetric positve definite

__init__(y, T, mu, covariance=None, precision=None)[source]
Parameters:
  • y (ndarray [\(d_y\)]) – data

  • T (Map) – map \(T:\mathbb{R}^d\rightarrow\mathbb{R}^{d_y}\)

  • mu (ndarray [\(d_y\)] or Map) – noise mean or parametric map returning the mean

  • covariance (ndarray [\(d_y,d_y\)] or Map) – noise covariance or parametric map returning the covariance

  • precision (ndarray [\(d_y,d_y\)] or Map) – noise precision matrix or parametric map returning the precision matrix

class TransportMaps.Likelihoods.LikelihoodBase.AdditiveLinearNormalLogLikelihood(y, c, T, mu, covariance=None, precision=None)[source]

Bases: AdditiveNormalLogLikelihood

Define the log-likelihood for the additive linear Gaussian model

The model is

\[{\bf y} = {\bf c} + {\bf T}{\bf x} + \varepsilon \;, \quad \varepsilon \sim \mathcal{N}(\mu, \Sigma)\]

where \(T \in \mathbb{R}^{d_y \times d_x}\), \(\mu \in \mathbb{R}^{d_y}\) and \(\Sigma \in \mathbb{R}^{d_y \times d_y}\) is symmetric positve definite

Parameters:
  • y (ndarray [\(d_y\)]) – data

  • c (ndarray [\(d_y\)] or Map) – system constant or parametric map returning the constant

  • T (ndarray [\(d_y,d_x\)] or Map) – system matrix or parametric map returning the system matrix

  • mu (ndarray [\(d_y\)] or Map) – noise mean or parametric map returning the mean

  • covariance (ndarray [\(d_y,d_y\)] or Map) – noise covariance or parametric map returning the covariance

  • precision (ndarray [\(d_y,d_y\)] or Map) – noise precision matrix or parametric map returning the precision matrix

TransportMaps.Likelihoods.LikelihoodBase.AdditiveLinearGaussianLogLikelihood[source]
class TransportMaps.Likelihoods.LikelihoodBase.AdditiveConditionallyLinearNormalLogLikelihood(y, c, T, mu, covariance=None, precision=None, active_vars_system=[], active_vars_distribution=[], coeffs=None)[source]

Bases: AdditiveLogLikelihood

Define the log-likelihood for the additive linear Gaussian model

The model is

\[{\bf y} = {\bf c}(\theta) + {\bf T}(\theta){\bf x} + \varepsilon \;, \quad \varepsilon \sim \mathcal{N}(\mu(\theta), \Sigma(\theta))\]

where \(T \in \mathbb{R}^{d_y \times d_x}\), \(\mu \in \mathbb{R}^{d_y}\) and \(\Sigma \in \mathbb{R}^{d_y \times d_y}\) is symmetric positve definite

Parameters:
  • y (ndarray [\(d_y\)]) – data

  • c (ndarray [\(d_y\)] or Map) – system constant or parametric map returning the constant

  • T (ndarray [\(d_y,d_x\)] or Map) – system matrix or parametric map returning the system matrix

  • mu (ndarray [\(d_y\)] or Map) – noise mean or parametric map returning the mean

  • covariance (ndarray [\(d_y,d_y\)] or Map) – noise covariance or parametric map returning the covariance

  • precision (ndarray [\(d_y,d_y\)] or Map) – noise precision matrix or parametric map returning the precision matrix

  • active_vars_system (list of int) – active variables identifying the parameters for for \(c(\theta), T(\theta)\).

  • active_vars_distribution (list of int) – active variables identifying the parameters for for \(\mu(\theta), \Sigma(\theta)\).

  • coeffs (ndarray) – initialization coefficients

property n_coeffs[source]
property coeffs[source]
TransportMaps.Likelihoods.LikelihoodBase.AdditiveConditionallyLinearGaussianLogLikelihood[source]
class TransportMaps.Likelihoods.LikelihoodBase.LogisticLogLikelihood(y, F)[source]

Bases: LogLikelihood

Class for log-likelihood of \({\bf y}\), where \(y_i \vert {\bf x}; {\bf f}_i \sim \mathrm{Bernoulli}(\mathrm{logit}^{-1}({\bf f}_i^\top {\bf x})\)

This implements the vectorized log-likelihood:

\[\log\pi({\bf y}\vert{\bf x}; {\bf f}_{1:N}) = \sum_{i=1}^N y_i \frac{1}{1+\exp(-{\bf f}_i^\top {\bf x}} + (1 - y_i) \left[ 1 - \frac{1}{1+\exp(-{\bf f}_i^\top {\bf x}} \right]\]
Parameters:
  • y (ndarray of int [N]) – binary data

  • F (ndarray [N,d-1]) – factors \({\bf f}_{1:N}\)

property F[source]
logit(x, *args, **kwargs)[source]
inv_logit(x, *args, **kwargs)[source]
evaluate(x, *args, **kwargs)[source]
Parameters:

x (ndarray [m,d]) – evaluation points

grad_x(x, *args, **kwargs)[source]

[Abstract] Evaluate \(\nabla_{\bf x}\log\pi({\bf y} \vert {\bf x})\).

Parameters:

x (ndarray [\(m,d\)]) – evaluation points

Returns:

(ndarray [\(m,1,d\)]) – gradient evaluations

class TransportMaps.Likelihoods.LikelihoodBase.TwoParametersLogisticLogLikelihood(y)[source]

Bases: LogLikelihood

See Stan example

\(I\): number of words, \(J\): number of people

Parameters:

y (ndarray of int [I,J]) – binary data

inv_logit(alpha, beta, theta, *args, **kwargs)[source]
evaluate(x, *args, **kwargs)[source]

[Abstract] Evaluate \(\log\pi({\bf y} \vert {\bf x})\).

Parameters:

x (ndarray [\(m,d\)]) – evaluation points

Returns:

(ndarray [\(m,1\)]) – function evaluations

grad_x(x, *args, **kwargs)[source]

[Abstract] Evaluate \(\nabla_{\bf x}\log\pi({\bf y} \vert {\bf x})\).

Parameters:

x (ndarray [\(m,d\)]) – evaluation points

Returns:

(ndarray [\(m,1,d\)]) – gradient evaluations

class TransportMaps.Likelihoods.LikelihoodBase.IndependentLogLikelihood(factors)[source]

Bases: TransportMaps.Maps.Functionals.Functional

Handling likelihoods in the form \(\pi({\bf y}\vert{\bf x}) = \prod_{i=1}^{n}\pi_i({\bf y}_i\vert{\bf x}_i)\)

Parameters:

factors (list of tuple) – each tuple contains a log-likelihood (LogLikelihood) and the sub-set of variables of \({\bf x}\) on which it acts.

Example

Let \(\pi(y_0,y_2\vert x_0,x_1,x_2)= \pi_0(y_0\vert x_0) \pi_2(y_2,x_2)\).

>>> factors = [(ll0, [0]),
>>>            (ll2, [2])]
>>> ll = IndependentLogLikelihood(factors)
property y[source]
property n_factors[source]
get_ncalls_tree(indent='')[source]
get_nevals_tree(indent='')[source]
get_teval_tree(indent='')[source]
update_ncalls_tree(obj)[source]
update_nevals_tree(obj)[source]
update_teval_tree(obj)[source]
reset_counters()[source]
append(factor)[source]

Add a new factor to the likelihood

Parameters:

factors (tuple) – tuple containing a log-likelihood (LogLikelihood) and the sub-set of variables of \({\bf x}\) on which it acts.

Example

Let \(\pi(y_0,y_2\vert x_0,x_1,x_2)= \pi_0(y_0\vert x_0) \pi_2(y_2,x_2)\) and let’s add the factor \(\pi_1(y_1\vert x_1)\), obtaining:

\[\pi(y_0,y_1,y_2\vert x_0,x_1,x_2)= \pi_0(y_0\vert x_0) \pi_1(y_1\vert x_1) \pi_2(y_2,x_2)\]
>>> factor = (ll1, [1])
>>> ll.append(factor)
evaluate(x, idxs_slice=slice(None, None, None), cache=None, **kwargs)[source]

Evaluate \(\log\pi({\bf y} \vert {\bf x})\).

Parameters:
  • x (ndarray [\(m,d\)]) – evaluation points

  • cache (dict) – cache

Returns:

(ndarray [\(m,1\)]) – function evaluations

grad_x(x, idxs_slice=slice(None, None, None), cache=None, **kwargs)[source]

Evaluate \(\nabla_{\bf x}\log\pi({\bf y} \vert {\bf x})\).

Parameters:
  • x (ndarray [\(m,d\)]) – evaluation points

  • cache (dict) – cache

Returns:

(ndarray [\(m,1,d\)]) – gradient evaluations

tuple_grad_x(x, idxs_slice=slice(None, None, None), cache=None, **kwargs)[source]

Evaluate \(\log\pi({\bf y} \vert {\bf x}), \nabla_{\bf x}\log\pi({\bf y} \vert {\bf x})\).

Parameters:
  • x (ndarray [\(m,d\)]) – evaluation points

  • cache (dict) – cache

Returns:

(tuple) – function and gradient evaluations

hess_x(x, idxs_slice=slice(None, None, None), cache=None, **kwargs)[source]

Evaluate \(\nabla^2_{\bf x}\log\pi({\bf y} \vert {\bf x})\).

Parameters:
  • x (ndarray [\(m,d\)]) – evaluation points

  • cache (dict) – cache

Returns:

(ndarray [\(m,1,d,d\)]) – Hessian evaluations

action_hess_x(x, dx, idxs_slice=slice(None, None, None), cache=None, **kwargs)[source]

Evaluate \(\langle\nabla^2_{\bf x}\log\pi({\bf y} \vert {\bf x}),\delta{\bf x}\rangle\).

Parameters:
  • x (ndarray [\(m,d\)]) – evaluation points

  • dx (ndarray [\(m,d\)]) – direction on which to evaluate the Hessian

  • cache (dict) – cache

Returns:

(ndarray [\(m,1,d\)]) – Hessian evaluations