1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
|
# -*- coding: utf-8 -*-
"""
Solvers for the original linear program OT problem
"""
# Author: Remi Flamary <remi.flamary@unice.fr>
#
# License: MIT License
import multiprocessing
import numpy as np
from scipy.sparse import coo_matrix
from .import cvx
# import compiled emd
from .emd_wrap import emd_c, check_result, emd_1d_sorted
from ..utils import parmap
from .cvx import barycenter
from ..utils import dist
__all__=['emd', 'emd2', 'barycenter', 'free_support_barycenter', 'cvx',
'emd_1d', 'emd2_1d', 'wasserstein_1d']
def emd(a, b, M, numItermax=100000, log=False):
"""Solves the Earth Movers distance problem and returns the OT matrix
.. math::
\gamma = arg\min_\gamma <\gamma,M>_F
s.t. \gamma 1 = a
\gamma^T 1= b
\gamma\geq 0
where :
- M is the metric cost matrix
- a and b are the sample weights
Uses the algorithm proposed in [1]_
Parameters
----------
a : (ns,) ndarray, float64
Source histogram (uniform weigth if empty list)
b : (nt,) ndarray, float64
Target histogram (uniform weigth if empty list)
M : (ns,nt) ndarray, float64
loss matrix
numItermax : int, optional (default=100000)
The maximum number of iterations before stopping the optimization
algorithm if it has not converged.
log: boolean, optional (default=False)
If True, returns a dictionary containing the cost and dual
variables. Otherwise returns only the optimal transportation matrix.
Returns
-------
gamma: (ns x nt) ndarray
Optimal transportation matrix for the given parameters
log: dict
If input log is true, a dictionary containing the cost and dual
variables and exit status
Examples
--------
Simple example with obvious solution. The function emd accepts lists and
perform automatic conversion to numpy arrays
>>> import ot
>>> a=[.5,.5]
>>> b=[.5,.5]
>>> M=[[0.,1.],[1.,0.]]
>>> ot.emd(a,b,M)
array([[ 0.5, 0. ],
[ 0. , 0.5]])
References
----------
.. [1] Bonneel, N., Van De Panne, M., Paris, S., & Heidrich, W.
(2011, December). Displacement interpolation using Lagrangian mass
transport. In ACM Transactions on Graphics (TOG) (Vol. 30, No. 6, p.
158). ACM.
See Also
--------
ot.bregman.sinkhorn : Entropic regularized OT
ot.optim.cg : General regularized OT"""
a = np.asarray(a, dtype=np.float64)
b = np.asarray(b, dtype=np.float64)
M = np.asarray(M, dtype=np.float64)
# if empty array given then use uniform distributions
if len(a) == 0:
a = np.ones((M.shape[0],), dtype=np.float64) / M.shape[0]
if len(b) == 0:
b = np.ones((M.shape[1],), dtype=np.float64) / M.shape[1]
G, cost, u, v, result_code = emd_c(a, b, M, numItermax)
result_code_string = check_result(result_code)
if log:
log = {}
log['cost'] = cost
log['u'] = u
log['v'] = v
log['warning'] = result_code_string
log['result_code'] = result_code
return G, log
return G
def emd2(a, b, M, processes=multiprocessing.cpu_count(),
numItermax=100000, log=False, return_matrix=False):
"""Solves the Earth Movers distance problem and returns the loss
.. math::
\gamma = arg\min_\gamma <\gamma,M>_F
s.t. \gamma 1 = a
\gamma^T 1= b
\gamma\geq 0
where :
- M is the metric cost matrix
- a and b are the sample weights
Uses the algorithm proposed in [1]_
Parameters
----------
a : (ns,) ndarray, float64
Source histogram (uniform weigth if empty list)
b : (nt,) ndarray, float64
Target histogram (uniform weigth if empty list)
M : (ns,nt) ndarray, float64
loss matrix
numItermax : int, optional (default=100000)
The maximum number of iterations before stopping the optimization
algorithm if it has not converged.
log: boolean, optional (default=False)
If True, returns a dictionary containing the cost and dual
variables. Otherwise returns only the optimal transportation cost.
return_matrix: boolean, optional (default=False)
If True, returns the optimal transportation matrix in the log.
Returns
-------
gamma: (ns x nt) ndarray
Optimal transportation matrix for the given parameters
log: dict
If input log is true, a dictionary containing the cost and dual
variables and exit status
Examples
--------
Simple example with obvious solution. The function emd accepts lists and
perform automatic conversion to numpy arrays
>>> import ot
>>> a=[.5,.5]
>>> b=[.5,.5]
>>> M=[[0.,1.],[1.,0.]]
>>> ot.emd2(a,b,M)
0.0
References
----------
.. [1] Bonneel, N., Van De Panne, M., Paris, S., & Heidrich, W.
(2011, December). Displacement interpolation using Lagrangian mass
transport. In ACM Transactions on Graphics (TOG) (Vol. 30, No. 6, p.
158). ACM.
See Also
--------
ot.bregman.sinkhorn : Entropic regularized OT
ot.optim.cg : General regularized OT"""
a = np.asarray(a, dtype=np.float64)
b = np.asarray(b, dtype=np.float64)
M = np.asarray(M, dtype=np.float64)
# if empty array given then use uniform distributions
if len(a) == 0:
a = np.ones((M.shape[0],), dtype=np.float64) / M.shape[0]
if len(b) == 0:
b = np.ones((M.shape[1],), dtype=np.float64) / M.shape[1]
if log or return_matrix:
def f(b):
G, cost, u, v, resultCode = emd_c(a, b, M, numItermax)
result_code_string = check_result(resultCode)
log = {}
if return_matrix:
log['G'] = G
log['u'] = u
log['v'] = v
log['warning'] = result_code_string
log['result_code'] = resultCode
return [cost, log]
else:
def f(b):
G, cost, u, v, result_code = emd_c(a, b, M, numItermax)
check_result(result_code)
return cost
if len(b.shape) == 1:
return f(b)
nb = b.shape[1]
res = parmap(f, [b[:, i] for i in range(nb)], processes)
return res
def free_support_barycenter(measures_locations, measures_weights, X_init, b=None, weights=None, numItermax=100, stopThr=1e-7, verbose=False, log=None):
"""
Solves the free support (locations of the barycenters are optimized, not the weights) Wasserstein barycenter problem (i.e. the weighted Frechet mean for the 2-Wasserstein distance)
The function solves the Wasserstein barycenter problem when the barycenter measure is constrained to be supported on k atoms.
This problem is considered in [1] (Algorithm 2). There are two differences with the following codes:
- we do not optimize over the weights
- we do not do line search for the locations updates, we use i.e. theta = 1 in [1] (Algorithm 2). This can be seen as a discrete implementation of the fixed-point algorithm of [2] proposed in the continuous setting.
Parameters
----------
measures_locations : list of (k_i,d) np.ndarray
The discrete support of a measure supported on k_i locations of a d-dimensional space (k_i can be different for each element of the list)
measures_weights : list of (k_i,) np.ndarray
Numpy arrays where each numpy array has k_i non-negatives values summing to one representing the weights of each discrete input measure
X_init : (k,d) np.ndarray
Initialization of the support locations (on k atoms) of the barycenter
b : (k,) np.ndarray
Initialization of the weights of the barycenter (non-negatives, sum to 1)
weights : (k,) np.ndarray
Initialization of the coefficients of the barycenter (non-negatives, sum to 1)
numItermax : int, optional
Max number of iterations
stopThr : float, optional
Stop threshol on error (>0)
verbose : bool, optional
Print information along iterations
log : bool, optional
record log if True
Returns
-------
X : (k,d) np.ndarray
Support locations (on k atoms) of the barycenter
References
----------
.. [1] Cuturi, Marco, and Arnaud Doucet. "Fast computation of Wasserstein barycenters." International Conference on Machine Learning. 2014.
.. [2] Álvarez-Esteban, Pedro C., et al. "A fixed-point approach to barycenters in Wasserstein space." Journal of Mathematical Analysis and Applications 441.2 (2016): 744-762.
"""
iter_count = 0
N = len(measures_locations)
k = X_init.shape[0]
d = X_init.shape[1]
if b is None:
b = np.ones((k,))/k
if weights is None:
weights = np.ones((N,)) / N
X = X_init
log_dict = {}
displacement_square_norms = []
displacement_square_norm = stopThr + 1.
while ( displacement_square_norm > stopThr and iter_count < numItermax ):
T_sum = np.zeros((k, d))
for (measure_locations_i, measure_weights_i, weight_i) in zip(measures_locations, measures_weights, weights.tolist()):
M_i = dist(X, measure_locations_i)
T_i = emd(b, measure_weights_i, M_i)
T_sum = T_sum + weight_i * np.reshape(1. / b, (-1, 1)) * np.matmul(T_i, measure_locations_i)
displacement_square_norm = np.sum(np.square(T_sum-X))
if log:
displacement_square_norms.append(displacement_square_norm)
X = T_sum
if verbose:
print('iteration %d, displacement_square_norm=%f\n', iter_count, displacement_square_norm)
iter_count += 1
if log:
log_dict['displacement_square_norms'] = displacement_square_norms
return X, log_dict
else:
return X
def emd_1d(x_a, x_b, a=None, b=None, metric='sqeuclidean', p=1., dense=True,
log=False):
"""Solves the Earth Movers distance problem between 1d measures and returns
the OT matrix
.. math::
\gamma = arg\min_\gamma \sum_i \sum_j \gamma_{ij} d(x_a[i], x_b[j])
s.t. \gamma 1 = a,
\gamma^T 1= b,
\gamma\geq 0
where :
- d is the metric
- x_a and x_b are the samples
- a and b are the sample weights
When 'minkowski' is used as a metric, :math:`d(x, y) = |x - y|^p`.
Uses the algorithm detailed in [1]_
Parameters
----------
x_a : (ns,) or (ns, 1) ndarray, float64
Source dirac locations (on the real line)
x_b : (nt,) or (ns, 1) ndarray, float64
Target dirac locations (on the real line)
a : (ns,) ndarray, float64, optional
Source histogram (default is uniform weight)
b : (nt,) ndarray, float64, optional
Target histogram (default is uniform weight)
metric: str, optional (default='sqeuclidean')
Metric to be used. Only strings listed in :func:`ot.dist` are accepted.
Due to implementation details, this function runs faster when
`'sqeuclidean'`, `'cityblock'`, or `'euclidean'` metrics are used.
p: float, optional (default=1.0)
The p-norm to apply for if metric='minkowski'
dense: boolean, optional (default=True)
If True, returns math:`\gamma` as a dense ndarray of shape (ns, nt).
Otherwise returns a sparse representation using scipy's `coo_matrix`
format. Due to implementation details, this function runs faster when
`'sqeuclidean'`, `'minkowski'`, `'cityblock'`, or `'euclidean'` metrics
are used.
log: boolean, optional (default=False)
If True, returns a dictionary containing the cost.
Otherwise returns only the optimal transportation matrix.
Returns
-------
gamma: (ns, nt) ndarray
Optimal transportation matrix for the given parameters
log: dict
If input log is True, a dictionary containing the cost
Examples
--------
Simple example with obvious solution. The function emd_1d accepts lists and
performs automatic conversion to numpy arrays
>>> import ot
>>> a=[.5, .5]
>>> b=[.5, .5]
>>> x_a = [2., 0.]
>>> x_b = [0., 3.]
>>> ot.emd_1d(x_a, x_b, a, b)
array([[0. , 0.5],
[0.5, 0. ]])
>>> ot.emd_1d(x_a, x_b)
array([[0. , 0.5],
[0.5, 0. ]])
References
----------
.. [1] Peyré, G., & Cuturi, M. (2017). "Computational Optimal
Transport", 2018.
See Also
--------
ot.lp.emd : EMD for multidimensional distributions
ot.lp.emd2_1d : EMD for 1d distributions (returns cost instead of the
transportation matrix)
"""
a = np.asarray(a, dtype=np.float64)
b = np.asarray(b, dtype=np.float64)
x_a = np.asarray(x_a, dtype=np.float64)
x_b = np.asarray(x_b, dtype=np.float64)
assert (x_a.ndim == 1 or x_a.ndim == 2 and x_a.shape[1] == 1), \
"emd_1d should only be used with monodimensional data"
assert (x_b.ndim == 1 or x_b.ndim == 2 and x_b.shape[1] == 1), \
"emd_1d should only be used with monodimensional data"
# if empty array given then use uniform distributions
if a.ndim == 0 or len(a) == 0:
a = np.ones((x_a.shape[0],), dtype=np.float64) / x_a.shape[0]
if b.ndim == 0 or len(b) == 0:
b = np.ones((x_b.shape[0],), dtype=np.float64) / x_b.shape[0]
x_a_1d = x_a.reshape((-1, ))
x_b_1d = x_b.reshape((-1, ))
perm_a = np.argsort(x_a_1d)
perm_b = np.argsort(x_b_1d)
G_sorted, indices, cost = emd_1d_sorted(a, b,
x_a_1d[perm_a], x_b_1d[perm_b],
metric=metric, p=p)
G = coo_matrix((G_sorted, (perm_a[indices[:, 0]], perm_b[indices[:, 1]])),
shape=(a.shape[0], b.shape[0]))
if dense:
G = G.toarray()
if log:
log = {'cost': cost}
return G, log
return G
def emd2_1d(x_a, x_b, a=None, b=None, metric='sqeuclidean', p=1., dense=True,
log=False):
"""Solves the Earth Movers distance problem between 1d measures and returns
the loss
.. math::
\gamma = arg\min_\gamma \sum_i \sum_j \gamma_{ij} d(x_a[i], x_b[j])
s.t. \gamma 1 = a,
\gamma^T 1= b,
\gamma\geq 0
where :
- d is the metric
- x_a and x_b are the samples
- a and b are the sample weights
When 'minkowski' is used as a metric, :math:`d(x, y) = |x - y|^p`.
Uses the algorithm detailed in [1]_
Parameters
----------
x_a : (ns,) or (ns, 1) ndarray, float64
Source dirac locations (on the real line)
x_b : (nt,) or (ns, 1) ndarray, float64
Target dirac locations (on the real line)
a : (ns,) ndarray, float64, optional
Source histogram (default is uniform weight)
b : (nt,) ndarray, float64, optional
Target histogram (default is uniform weight)
metric: str, optional (default='sqeuclidean')
Metric to be used. Only strings listed in :func:`ot.dist` are accepted.
Due to implementation details, this function runs faster when
`'sqeuclidean'`, `'minkowski'`, `'cityblock'`, or `'euclidean'` metrics
are used.
p: float, optional (default=1.0)
The p-norm to apply for if metric='minkowski'
dense: boolean, optional (default=True)
If True, returns math:`\gamma` as a dense ndarray of shape (ns, nt).
Otherwise returns a sparse representation using scipy's `coo_matrix`
format. Only used if log is set to True. Due to implementation details,
this function runs faster when dense is set to False.
log: boolean, optional (default=False)
If True, returns a dictionary containing the transportation matrix.
Otherwise returns only the loss.
Returns
-------
loss: float
Cost associated to the optimal transportation
log: dict
If input log is True, a dictionary containing the Optimal transportation
matrix for the given parameters
Examples
--------
Simple example with obvious solution. The function emd2_1d accepts lists and
performs automatic conversion to numpy arrays
>>> import ot
>>> a=[.5, .5]
>>> b=[.5, .5]
>>> x_a = [2., 0.]
>>> x_b = [0., 3.]
>>> ot.emd2_1d(x_a, x_b, a, b)
0.5
>>> ot.emd2_1d(x_a, x_b)
0.5
References
----------
.. [1] Peyré, G., & Cuturi, M. (2017). "Computational Optimal
Transport", 2018.
See Also
--------
ot.lp.emd2 : EMD for multidimensional distributions
ot.lp.emd_1d : EMD for 1d distributions (returns the transportation matrix
instead of the cost)
"""
# If we do not return G (log==False), then we should not to cast it to dense
# (useless overhead)
G, log_emd = emd_1d(x_a=x_a, x_b=x_b, a=a, b=b, metric=metric, p=p,
dense=dense and log, log=True)
cost = log_emd['cost']
if log:
log_emd = {'G': G}
return cost, log_emd
return cost
def wasserstein_1d(x_a, x_b, a=None, b=None, p=1.):
"""Solves the p-Wasserstein distance problem between 1d measures and returns
the distance
.. math::
\gamma = arg\min_\gamma \left( \sum_i \sum_j \gamma_{ij}
|x_a[i] - x_b[j]|^p \\right)^{1/p}
s.t. \gamma 1 = a,
\gamma^T 1= b,
\gamma\geq 0
where :
- x_a and x_b are the samples
- a and b are the sample weights
Uses the algorithm detailed in [1]_
Parameters
----------
x_a : (ns,) or (ns, 1) ndarray, float64
Source dirac locations (on the real line)
x_b : (nt,) or (ns, 1) ndarray, float64
Target dirac locations (on the real line)
a : (ns,) ndarray, float64, optional
Source histogram (default is uniform weight)
b : (nt,) ndarray, float64, optional
Target histogram (default is uniform weight)
p: float, optional (default=1.0)
The order of the p-Wasserstein distance to be computed
Returns
-------
dist: float
p-Wasserstein distance
Examples
--------
Simple example with obvious solution. The function wasserstein_1d accepts
lists and performs automatic conversion to numpy arrays
>>> import ot
>>> a=[.5, .5]
>>> b=[.5, .5]
>>> x_a = [2., 0.]
>>> x_b = [0., 3.]
>>> ot.wasserstein_1d(x_a, x_b, a, b)
0.5
>>> ot.wasserstein_1d(x_a, x_b)
0.5
References
----------
.. [1] Peyré, G., & Cuturi, M. (2017). "Computational Optimal
Transport", 2018.
See Also
--------
ot.lp.emd_1d : EMD for 1d distributions
"""
cost_emd = emd2_1d(x_a=x_a, x_b=x_b, a=a, b=b, metric='minkowski', p=p,
dense=False, log=False)
return np.power(cost_emd, 1. / p)
|