API

Exact optimal transport

PythonOT.emdFunction
emd(μ, ν, C; kwargs...)

Compute the optimal transport plan for the Monge-Kantorovich problem with source and target marginals μ and ν and cost matrix C of size (length(μ), length(ν)).

The optimal transport plan γ is of the same size as C and solves

\[\inf_{\gamma \in \Pi(\mu, \nu)} \langle \gamma, C \rangle.\]

This function is a wrapper of the function emd in the Python Optimal Transport package. Keyword arguments are listed in the documentation of the Python function.

Examples

julia> μ = [0.5, 0.2, 0.3];

julia> ν = [0.0, 1.0];

julia> C = [0.0 1.0; 2.0 0.0; 0.5 1.5];

julia> emd(μ, ν, C)
3×2 Matrix{Float64}:
 0.0  0.5
 0.0  0.2
 0.0  0.3

See also: emd2

source
PythonOT.emd2Function
emd2(μ, ν, C; kwargs...)

Compute the optimal transport cost for the Monge-Kantorovich problem with source and target marginals μ and ν and cost matrix C of size (length(μ), length(ν)).

The optimal transport cost is the scalar value

\[\inf_{\gamma \in \Pi(\mu, \nu)} \langle \gamma, C \rangle.\]

This function is a wrapper of the function emd2 in the Python Optimal Transport package. Keyword arguments are listed in the documentation of the Python function.

Examples

julia> μ = [0.5, 0.2, 0.3];

julia> ν = [0.0, 1.0];

julia> C = [0.0 1.0; 2.0 0.0; 0.5 1.5];

julia> emd2(μ, ν, C)
0.95

See also: emd

source
PythonOT.emd_1dFunction
emd_1d(xsource, xtarget; kwargs...)

Compute the optimal transport plan for the Monge-Kantorovich problem with univariate discrete measures with support xsource and xtarget as source and target marginals.

This function is a wrapper of the function emd_1d in the Python Optimal Transport package. Keyword arguments are listed in the documentation of the Python function.

Examples

julia> xsource = [0.2, 0.5];

julia> xtarget = [0.8, 0.3];

julia> emd_1d(xsource, xtarget)
2×2 Matrix{Float64}:
 0.0  0.5
 0.5  0.0

julia> histogram_source = [0.8, 0.2];

julia> histogram_target = [0.7, 0.3];

julia> emd_1d(xsource, xtarget; a=histogram_source, b=histogram_target)
2×2 Matrix{Float64}:
 0.5  0.3
 0.2  0.0

See also: emd, emd2_1d

source
PythonOT.emd2_1dFunction
emd2_1d(xsource, xtarget; kwargs...)

Compute the optimal transport cost for the Monge-Kantorovich problem with univariate discrete measures with support xsource and xtarget as source and target marginals.

This function is a wrapper of the function emd2_1d in the Python Optimal Transport package. Keyword arguments are listed in the documentation of the Python function.

Examples

julia> xsource = [0.2, 0.5];

julia> xtarget = [0.8, 0.3];

julia> round(emd2_1d(xsource, xtarget); sigdigits=6)
0.05

julia> histogram_source = [0.8, 0.2];

julia> histogram_target = [0.7, 0.3];

julia> round(emd2_1d(xsource, xtarget; a=histogram_source, b=histogram_target); sigdigits=6)
0.201

See also: emd2, emd2_1d

source

Regularized optimal transport

PythonOT.sinkhornFunction
sinkhorn(μ, ν, C, ε; kwargs...)

Compute the optimal transport plan for the entropic regularization optimal transport problem with source and target marginals μ and ν, cost matrix C of size (length(μ), length(ν)), and entropic regularization parameter ε.

The optimal transport plan γ is of the same size as C and solves

\[\inf_{\gamma \in \Pi(\mu, \nu)} \langle \gamma, C \rangle + \varepsilon \Omega(\gamma),\]

where $\Omega(\gamma) = \sum_{i,j} \gamma_{i,j} \log \gamma_{i,j}$ is the entropic regularization term.

This function is a wrapper of the function sinkhorn in the Python Optimal Transport package. Keyword arguments are listed in the documentation of the Python function.

Examples

julia> μ = [0.5, 0.2, 0.3];

julia> ν = [0.0, 1.0];

julia> C = [0.0 1.0; 2.0 0.0; 0.5 1.5];

julia> sinkhorn(μ, ν, C, 0.01)
3×2 Matrix{Float64}:
 0.0  0.5
 0.0  0.2
 0.0  0.3

It is possible to provide multiple target marginals as columns of a matrix. In this case the optimal transport costs are returned:

julia> ν = [0.0 0.5; 1.0 0.5];

julia> round.(sinkhorn(μ, ν, C, 0.01); sigdigits=6)
2-element Vector{Float64}:
 0.95
 0.45

See also: sinkhorn2

source
PythonOT.sinkhorn2Function
sinkhorn2(μ, ν, C, ε; kwargs...)

Compute the optimal transport cost for the entropic regularization optimal transport problem with source and target marginals μ and ν, cost matrix C of size (length(μ), length(ν)), and entropic regularization parameter ε.

The optimal transport cost is the scalar value

\[\inf_{\gamma \in \Pi(\mu, \nu)} \langle \gamma, C \rangle + \varepsilon \Omega(\gamma),\]

where $\Omega(\gamma) = \sum_{i,j} \gamma_{i,j} \log \gamma_{i,j}$ is the entropic regularization term.

This function is a wrapper of the function sinkhorn2 in the Python Optimal Transport package. Keyword arguments are listed in the documentation of the Python function.

Examples

julia> μ = [0.5, 0.2, 0.3];

julia> ν = [0.0, 1.0];

julia> C = [0.0 1.0; 2.0 0.0; 0.5 1.5];

julia> round(sinkhorn2(μ, ν, C, 0.01); sigdigits=6)
0.95

It is possible to provide multiple target marginals as columns of a matrix.

julia> ν = [0.0 0.5; 1.0 0.5];

julia> round.(sinkhorn2(μ, ν, C, 0.01); sigdigits=6)
2-element Vector{Float64}:
 0.95
 0.45

See also: sinkhorn

source
PythonOT.empirical_sinkhorn_divergenceFunction
empirical_sinkhorn_divergence(xsource, xtarget, ε; kwargs...)

Compute the Sinkhorn divergence from empirical data, where xsource and xtarget are arrays representing samples in the source domain and target domain, respectively, and ε is the regularization term.

This function is a wrapper of the function ot.bregman.empirical_sinkhorn_divergence in the Python Optimal Transport package. Keyword arguments are listed in the documentation of the Python function.

Examples

julia> xsource = [1.0];

julia> xtarget = [2.0, 3.0];

julia> ε = 0.01;

julia> empirical_sinkhorn_divergence(xsource, xtarget, ε) ≈
       sinkhorn2([1], [0.5, 0.5], [1.0 4.0], ε) -
       (
           sinkhorn2([1], [1], zeros(1, 1), ε) +
           sinkhorn2([0.5, 0.5], [0.5, 0.5], Float64[0 1; 1 0], ε)
       ) / 2
true

See also: sinkhorn2

source
PythonOT.barycenterFunction
barycenter(A, C, ε; kwargs...)

Compute the entropically regularized Wasserstein barycenter with histograms A, cost matrix C, and entropic regularization parameter ε.

The Wasserstein barycenter is a histogram and solves

\[\inf_{a} \sum_{i} W_{\varepsilon,C}(a, a_i),\]

where the histograms $a_i$ are columns of matrix A and $W_{\varepsilon,C}(a, a_i)$ is the optimal transport cost for the entropically regularized optimal transport problem with marginals $a$ and $a_i$, cost matrix $C$, and entropic regularization parameter $\varepsilon$. Optionally, weights of the histograms $a_i$ can be provided with the keyword argument weights.

This function is a wrapper of the function barycenter in the Python Optimal Transport package. Keyword arguments are listed in the documentation of the Python function.

Examples

julia> A = rand(10, 3);

julia> A ./= sum(A; dims=1);

julia> C = rand(10, 10);

julia> isapprox(sum(barycenter(A, C, 0.01; method="sinkhorn_stabilized")), 1; atol=1e-4)
true
source

The submodule Smooth contains a function for solving regularized optimal transport problems with L2- and entropic regularization using the dual formulation. You can load the submodule with

using PythonOT.Smooth
PythonOT.Smooth.smooth_ot_dualFunction
smooth_ot_dual(μ, ν, C, ε; reg_type="l2", kwargs...)

Compute the optimal transport plan for a regularized optimal transport problem with source and target marginals μ and ν, cost matrix C of size (length(μ), length(ν)), and regularization parameter ε.

The optimal transport map γ is of the same size as C and solves

\[\inf_{\gamma \in \Pi(\mu, \nu)} \langle \gamma, C \rangle + \varepsilon \Omega(\gamma),\]

where $\Omega(\gamma)$ is the L2-regularization term $\Omega(\gamma) = \|\gamma\|_F^2/2$ if reg_type="l2" (the default) or the entropic regularization term $\Omega(\gamma) = \sum_{i,j} \gamma_{i,j} \log \gamma_{i,j}$ if reg_type="kl".

The function solves the dual formulation[BSR2018]

\[\max_{\alpha, \beta} \mu^{\mathsf{T}} \alpha + \nu^{\mathsf{T}} \beta − \sum_{j} \delta_{\Omega}(\alpha + \beta_j - C_j),\]

where $C_j$ is the $j$th column of the cost matrix and $\delta_{\Omega}$ is the conjugate of the regularization term $\Omega$.

This function is a wrapper of the function smooth_ot_dual in the Python Optimal Transport package. Keyword arguments are listed in the documentation of the Python function.

Examples

julia> μ = [0.5, 0.2, 0.3];

julia> ν = [0.0, 1.0];

julia> C = [0.0 1.0; 2.0 0.0; 0.5 1.5];

julia> smooth_ot_dual(μ, ν, C, 0.01)
3×2 Matrix{Float64}:
 0.0  0.5
 0.0  0.2
 0.0  0.300001
source

Unbalanced optimal transport

PythonOT.sinkhorn_unbalancedFunction
sinkhorn_unbalanced(μ, ν, C, ε, λ; kwargs...)

Compute the optimal transport plan for the unbalanced entropic regularization optimal transport problem with source and target marginals μ and ν, cost matrix C of size (length(μ), length(ν)), entropic regularization parameter ε, and marginal relaxation term λ.

The optimal transport plan γ is of the same size as C and solves

\[\inf_{\gamma} \langle \gamma, C \rangle + \varepsilon \Omega(\gamma) + \lambda \mathrm{KL}(\gamma 1, \mu) + \lambda \mathrm{KL}(\gamma^{\mathsf{T}} 1, \nu),\]

where $\Omega(\gamma) = \sum_{i,j} \gamma_{i,j} \log \gamma_{i,j}$ is the entropic regularization term and $\mathrm{KL}$ is the Kullback-Leibler divergence.

This function is a wrapper of the function sinkhorn_unbalanced in the Python Optimal Transport package. Keyword arguments are listed in the documentation of the Python function.

Examples

julia> μ = [0.5, 0.2, 0.3];

julia> ν = [0.0, 1.0];

julia> C = [0.0 1.0; 2.0 0.0; 0.5 1.5];

julia> round.(sinkhorn_unbalanced(μ, ν, C, 0.01, 1_000); sigdigits=4)
3×2 Matrix{Float64}:
 0.0  0.5
 0.0  0.2002
 0.0  0.2998

It is possible to provide multiple target marginals as columns of a matrix. In this case the optimal transport costs are returned:

julia> ν = [0.0 0.5; 1.0 0.5];

julia> round.(sinkhorn_unbalanced(μ, ν, C, 0.01, 1_000); sigdigits=4)
2-element Vector{Float64}:
 0.9497
 0.4494

See also: sinkhorn_unbalanced2

source
PythonOT.sinkhorn_unbalanced2Function
sinkhorn_unbalanced2(μ, ν, C, ε, λ; kwargs...)

Compute the optimal transport cost for the unbalanced entropic regularization optimal transport problem with source and target marginals μ and ν, cost matrix C of size (length(μ), length(ν)), entropic regularization parameter ε, and marginal relaxation term λ.

The optimal transport cost is the scalar value

\[\inf_{\gamma} \langle \gamma, C \rangle + \varepsilon \Omega(\gamma) + \lambda \mathrm{KL}(\gamma 1, \mu) + \lambda \mathrm{KL}(\gamma^{\mathsf{T}} 1, \nu),\]

where $\Omega(\gamma) = \sum_{i,j} \gamma_{i,j} \log \gamma_{i,j}$ is the entropic regularization term and $\mathrm{KL}$ is the Kullback-Leibler divergence.

This function is a wrapper of the function sinkhorn_unbalanced2 in the Python Optimal Transport package. Keyword arguments are listed in the documentation of the Python function.

Examples

julia> μ = [0.5, 0.2, 0.3];

julia> ν = [0.0, 1.0];

julia> C = [0.0 1.0; 2.0 0.0; 0.5 1.5];

julia> round.(sinkhorn_unbalanced2(μ, ν, C, 0.01, 1_000); sigdigits=4)
0.9497

It is possible to provide multiple target marginals as columns of a matrix:

julia> ν = [0.0 0.5; 1.0 0.5];

julia> round.(sinkhorn_unbalanced2(μ, ν, C, 0.01, 1_000); sigdigits=4)
2-element Vector{Float64}:
 0.9497
 0.4494

See also: sinkhorn_unbalanced

source
PythonOT.barycenter_unbalancedFunction
barycenter_unbalanced(A, C, ε, λ; kwargs...)

Compute the entropically regularized unbalanced Wasserstein barycenter with histograms A, cost matrix C, entropic regularization parameter ε and marginal relaxation parameter λ.

The Wasserstein barycenter is a histogram and solves

\[\inf_{a} \sum_{i} W_{\varepsilon,C,\lambda}(a, a_i),\]

where the histograms $a_i$ are columns of matrix A and $W_{\varepsilon,C,\lambda}(a, a_i)}$ is the optimal transport cost for the entropically regularized optimal transport problem with marginals $a$ and $a_i$, cost matrix $C$, entropic regularization parameter $\varepsilon$ and marginal relaxation parameter $\lambda$. Optionally, weights of the histograms $a_i$ can be provided with the keyword argument weights.

This function is a wrapper of the function barycenter_unbalanced in the Python Optimal Transport package. Keyword arguments are listed in the documentation of the Python function.

Examples

julia> A = rand(10, 3);

julia> A ./= sum(A; dims=1);

julia> C = rand(10, 10);

julia> isapprox(sum(barycenter_unbalanced(A, C, 0.01, 1; method="sinkhorn_stabilized")), 1; atol=1e-4)
false

julia> isapprox(sum(barycenter_unbalanced(
           A, C, 0.01, 10_000; method="sinkhorn_stabilized", numItermax=5_000
       )), 1; atol=1e-4)
true

See also: barycenter

source
PythonOT.mm_unbalancedFunction
mm_unbalanced(a, b, M, reg_m; reg=0, c=a*b', kwargs...)

Solve the unbalanced optimal transport problem and return the OT plan. The function solves the following optimization problem:

\[W = \min_{\gamma \geq 0} \langle \gamma, M \rangle_F + \mathrm{reg_{m1}} \cdot \operatorname{div}(\gamma \mathbf{1}, a) + \mathrm{reg_{m2}} \cdot \operatorname{div}(\gamma^\mathsf{T} \mathbf{1}, b) + \mathrm{reg} \cdot \operatorname{div}(\gamma, c)\]

where

  • M is the metric cost matrix,
  • a and b are source and target unbalanced distributions,
  • c is a reference distribution for the regularization,
  • reg_m is the marginal relaxation term (if it is a scalar or an indexable object of length 1, then the same term is applied to both marginal relaxations), and
  • reg is a regularization term.

This function is a wrapper of the function mm_unbalanced in the Python Optimal Transport package. Keyword arguments are listed in the documentation of the Python function.

Examples

julia> a=[.5, .5];

julia> b=[.5, .5];

julia> M=[1. 36.; 9. 4.];

julia> round.(mm_unbalanced(a, b, M, 5, div="kl"), digits=2)
2×2 Matrix{Float64}:
 0.45  0.0
 0.0   0.34

julia> round.(mm_unbalanced(a, b, M, 5, div="l2"), digits=2)
2×2 Matrix{Float64}:
 0.4  0.0
 0.0  0.1
source
  • BSR2018Blondel, M., Seguy, V., & Rolet, A. (2018). Smooth and Sparse Optimal Transport. In Proceedings of the Twenty-First International Conference on Artificial Intelligence and Statistics (AISTATS).