summaryrefslogtreecommitdiff
path: root/src/python/gudhi/representations/vector_methods.py
blob: d4449e7d9bf65de8be49d0d08cad6b28e70cb1b2 (plain)
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
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
# This file is part of the Gudhi Library - https://gudhi.inria.fr/ - which is released under MIT.
# See file LICENSE or go to https://gudhi.inria.fr/licensing/ for full license details.
# Author(s):       Mathieu Carrière, Martin Royer
#
# Copyright (C) 2018-2020 Inria
#
# Modification(s):
#   - 2020/06 Martin: ATOL integration

import numpy as np
from sklearn.base          import BaseEstimator, TransformerMixin
from sklearn.preprocessing import MinMaxScaler, MaxAbsScaler
from sklearn.neighbors     import DistanceMetric
from sklearn.metrics       import pairwise

from .preprocessing import DiagramScaler, BirthPersistenceTransform

#############################################
# Finite Vectorization methods ##############
#############################################

class PersistenceImage(BaseEstimator, TransformerMixin):
    """
    This is a class for computing persistence images from a list of persistence diagrams. A persistence image is a 2D function computed from a persistence diagram by convolving the diagram points with a weighted Gaussian kernel. The plane is then discretized into an image with pixels, which is flattened and returned as a vector. See http://jmlr.org/papers/v18/16-337.html for more details.
    """
    def __init__(self, bandwidth=1., weight=lambda x: 1, resolution=[20,20], im_range=[np.nan, np.nan, np.nan, np.nan]):
        """
        Constructor for the PersistenceImage class.

        Parameters:
            bandwidth (double): bandwidth of the Gaussian kernel (default 1.).
            weight (function): weight function for the persistence diagram points (default constant function, ie lambda x: 1). This function must be defined on 2D points, ie lists or numpy arrays of the form [p_x,p_y].
            resolution ([int,int]): size (in pixels) of the persistence image (default [20,20]).
            im_range ([double,double,double,double]): minimum and maximum of each axis of the persistence image, of the form [x_min, x_max, y_min, y_max] (default [numpy.nan, numpy.nan, numpy.nan, numpyp.nan]). If one of the values is numpy.nan, it can be computed from the persistence diagrams with the fit() method.
        """
        self.bandwidth, self.weight = bandwidth, weight
        self.resolution, self.im_range = resolution, im_range

    def fit(self, X, y=None):
        """
        Fit the PersistenceImage class on a list of persistence diagrams: if any of the values in **im_range** is numpy.nan, replace it with the corresponding value computed on the given list of persistence diagrams.

        Parameters:
            X (list of n x 2 numpy arrays): input persistence diagrams.
            y (n x 1 array): persistence diagram labels (unused).
        """
        if np.isnan(np.array(self.im_range)).any():
            new_X = BirthPersistenceTransform().fit_transform(X)
            pre = DiagramScaler(use=True, scalers=[([0], MinMaxScaler()), ([1], MinMaxScaler())]).fit(new_X,y)
            [mx,my],[Mx,My] = [pre.scalers[0][1].data_min_[0], pre.scalers[1][1].data_min_[0]], [pre.scalers[0][1].data_max_[0], pre.scalers[1][1].data_max_[0]]
            self.im_range = np.where(np.isnan(np.array(self.im_range)), np.array([mx, Mx, my, My]), np.array(self.im_range))
        return self

    def transform(self, X):
        """
        Compute the persistence image for each persistence diagram individually and store the results in a single numpy array.

        Parameters:
            X (list of n x 2 numpy arrays): input persistence diagrams.
    
        Returns:
            numpy array with shape (number of diagrams) x (number of pixels = **resolution[0]** x **resolution[1]**): output persistence images.
        """
        num_diag, Xfit = len(X), []
        new_X = BirthPersistenceTransform().fit_transform(X)

        for i in range(num_diag):

            diagram, num_pts_in_diag = new_X[i], X[i].shape[0]

            w = np.empty(num_pts_in_diag)
            for j in range(num_pts_in_diag):
                w[j] = self.weight(diagram[j,:])

            x_values, y_values = np.linspace(self.im_range[0], self.im_range[1], self.resolution[0]), np.linspace(self.im_range[2], self.im_range[3], self.resolution[1])
            Xs, Ys = np.tile((diagram[:,0][:,np.newaxis,np.newaxis]-x_values[np.newaxis,np.newaxis,:]),[1,self.resolution[1],1]), np.tile(diagram[:,1][:,np.newaxis,np.newaxis]-y_values[np.newaxis,:,np.newaxis],[1,1,self.resolution[0]])
            image = np.tensordot(w, np.exp((-np.square(Xs)-np.square(Ys))/(2*np.square(self.bandwidth)))/(np.square(self.bandwidth)*2*np.pi), 1)

            Xfit.append(image.flatten()[np.newaxis,:])

        Xfit = np.concatenate(Xfit,0)

        return Xfit

    def __call__(self, diag):
        """
        Apply PersistenceImage on a single persistence diagram and outputs the result.

        Parameters:
            diag (n x 2 numpy array): input persistence diagram.

        Returns:
            numpy array with shape (number of pixels = **resolution[0]** x **resolution[1]**):: output persistence image.
        """
        return self.fit_transform([diag])[0,:]

class Landscape(BaseEstimator, TransformerMixin):
    """
    This is a class for computing persistence landscapes from a list of persistence diagrams. A persistence landscape is a collection of 1D piecewise-linear functions computed from the rank function associated to the persistence diagram. These piecewise-linear functions are then sampled evenly on a given range and the corresponding vectors of samples are concatenated and returned. See http://jmlr.org/papers/v16/bubenik15a.html for more details.
    """
    def __init__(self, num_landscapes=5, resolution=100, sample_range=[np.nan, np.nan]):
        """
        Constructor for the Landscape class.

        Parameters:
            num_landscapes (int): number of piecewise-linear functions to output (default 5).
            resolution (int): number of sample for all piecewise-linear functions (default 100).
            sample_range ([double, double]): minimum and maximum of all piecewise-linear function domains, of the form [x_min, x_max] (default [numpy.nan, numpy.nan]). It is the interval on which samples will be drawn evenly. If one of the values is numpy.nan, it can be computed from the persistence diagrams with the fit() method.
        """
        self.num_landscapes, self.resolution, self.sample_range = num_landscapes, resolution, sample_range
        self.nan_in_range = np.isnan(np.array(self.sample_range))
        self.new_resolution = self.resolution + self.nan_in_range.sum()

    def fit(self, X, y=None):
        """
        Fit the Landscape class on a list of persistence diagrams: if any of the values in **sample_range** is numpy.nan, replace it with the corresponding value computed on the given list of persistence diagrams.

        Parameters:
            X (list of n x 2 numpy arrays): input persistence diagrams.
            y (n x 1 array): persistence diagram labels (unused).
        """
        if self.nan_in_range.any():
            pre = DiagramScaler(use=True, scalers=[([0], MinMaxScaler()), ([1], MinMaxScaler())]).fit(X,y)
            [mx,my],[Mx,My] = [pre.scalers[0][1].data_min_[0], pre.scalers[1][1].data_min_[0]], [pre.scalers[0][1].data_max_[0], pre.scalers[1][1].data_max_[0]]
            self.sample_range = np.where(self.nan_in_range, np.array([mx, My]), np.array(self.sample_range))
        return self

    def transform(self, X):
        """
        Compute the persistence landscape for each persistence diagram individually and concatenate the results.

        Parameters:
            X (list of n x 2 numpy arrays): input persistence diagrams.
    
        Returns:
            numpy array with shape (number of diagrams) x (number of samples = **num_landscapes** x **resolution**): output persistence landscapes.
        """
        num_diag, Xfit = len(X), []
        x_values = np.linspace(self.sample_range[0], self.sample_range[1], self.new_resolution)
        step_x = x_values[1] - x_values[0]

        for i in range(num_diag):

            diagram, num_pts_in_diag = X[i], X[i].shape[0]

            ls = np.zeros([self.num_landscapes, self.new_resolution])

            events = []
            for j in range(self.new_resolution):
                events.append([])

            for j in range(num_pts_in_diag):
                [px,py] = diagram[j,:2]
                min_idx = np.clip(np.ceil((px          - self.sample_range[0]) / step_x).astype(int), 0, self.new_resolution)
                mid_idx = np.clip(np.ceil((0.5*(py+px) - self.sample_range[0]) / step_x).astype(int), 0, self.new_resolution)
                max_idx = np.clip(np.ceil((py          - self.sample_range[0]) / step_x).astype(int), 0, self.new_resolution)

                if min_idx < self.new_resolution and max_idx > 0:

                    landscape_value = self.sample_range[0] + min_idx * step_x - px
                    for k in range(min_idx, mid_idx):
                        events[k].append(landscape_value)
                        landscape_value += step_x

                    landscape_value = py - self.sample_range[0] - mid_idx * step_x
                    for k in range(mid_idx, max_idx):
                        events[k].append(landscape_value)
                        landscape_value -= step_x

            for j in range(self.new_resolution):
                events[j].sort(reverse=True)
                for k in range( min(self.num_landscapes, len(events[j])) ):
                    ls[k,j] = events[j][k]

            if self.nan_in_range[0]:
                ls = ls[:,1:]
            if self.nan_in_range[1]:
                ls = ls[:,:-1]
            ls = np.sqrt(2)*np.reshape(ls,[1,-1])
            Xfit.append(ls)

        Xfit = np.concatenate(Xfit,0)

        return Xfit

    def __call__(self, diag):
        """
        Apply Landscape on a single persistence diagram and outputs the result.

        Parameters:
            diag (n x 2 numpy array): input persistence diagram.

        Returns:
            numpy array with shape (number of samples = **num_landscapes** x **resolution**): output persistence landscape.
        """
        return self.fit_transform([diag])[0,:]

class Silhouette(BaseEstimator, TransformerMixin):
    """
    This is a class for computing persistence silhouettes from a list of persistence diagrams. A persistence silhouette is computed by taking a weighted average of the collection of 1D piecewise-linear functions given by the persistence landscapes, and then by evenly sampling this average on a given range. Finally, the corresponding vector of samples is returned. See https://arxiv.org/abs/1312.0308 for more details.
    """
    def __init__(self, weight=lambda x: 1, resolution=100, sample_range=[np.nan, np.nan]):
        """
        Constructor for the Silhouette class.

        Parameters:
            weight (function): weight function for the persistence diagram points (default constant function, ie lambda x: 1). This function must be defined on 2D points, ie on lists or numpy arrays of the form [p_x,p_y].
            resolution (int): number of samples for the weighted average (default 100).
            sample_range ([double, double]): minimum and maximum for the weighted average domain, of the form [x_min, x_max] (default [numpy.nan, numpy.nan]). It is the interval on which samples will be drawn evenly. If one of the values is numpy.nan, it can be computed from the persistence diagrams with the fit() method.
        """
        self.weight, self.resolution, self.sample_range = weight, resolution, sample_range

    def fit(self, X, y=None):
        """
        Fit the Silhouette class on a list of persistence diagrams: if any of the values in **sample_range** is numpy.nan, replace it with the corresponding value computed on the given list of persistence diagrams.

        Parameters:
            X (list of n x 2 numpy arrays): input persistence diagrams.
            y (n x 1 array): persistence diagram labels (unused).
        """
        if np.isnan(np.array(self.sample_range)).any():
            pre = DiagramScaler(use=True, scalers=[([0], MinMaxScaler()), ([1], MinMaxScaler())]).fit(X,y)
            [mx,my],[Mx,My] = [pre.scalers[0][1].data_min_[0], pre.scalers[1][1].data_min_[0]], [pre.scalers[0][1].data_max_[0], pre.scalers[1][1].data_max_[0]]
            self.sample_range = np.where(np.isnan(np.array(self.sample_range)), np.array([mx, My]), np.array(self.sample_range))
        return self

    def transform(self, X):
        """
        Compute the persistence silhouette for each persistence diagram individually and concatenate the results.

        Parameters:
            X (list of n x 2 numpy arrays): input persistence diagrams.
    
        Returns:
            numpy array with shape (number of diagrams) x (**resolution**): output persistence silhouettes.
        """
        num_diag, Xfit = len(X), []
        x_values = np.linspace(self.sample_range[0], self.sample_range[1], self.resolution)
        step_x = x_values[1] - x_values[0]

        for i in range(num_diag):

            diagram, num_pts_in_diag = X[i], X[i].shape[0]

            sh, weights = np.zeros(self.resolution), np.zeros(num_pts_in_diag)
            for j in range(num_pts_in_diag):
                weights[j] = self.weight(diagram[j,:])
            total_weight = np.sum(weights)

            for j in range(num_pts_in_diag):

                [px,py] = diagram[j,:2]
                weight  = weights[j] / total_weight
                min_idx = np.clip(np.ceil((px          - self.sample_range[0]) / step_x).astype(int), 0, self.resolution)
                mid_idx = np.clip(np.ceil((0.5*(py+px) - self.sample_range[0]) / step_x).astype(int), 0, self.resolution)
                max_idx = np.clip(np.ceil((py          - self.sample_range[0]) / step_x).astype(int), 0, self.resolution)

                if min_idx < self.resolution and max_idx > 0:

                    silhouette_value = self.sample_range[0] + min_idx * step_x - px
                    for k in range(min_idx, mid_idx):
                        sh[k] += weight * silhouette_value
                        silhouette_value += step_x

                    silhouette_value = py - self.sample_range[0] - mid_idx * step_x
                    for k in range(mid_idx, max_idx):
                        sh[k] += weight * silhouette_value
                        silhouette_value -= step_x

            Xfit.append(np.reshape(np.sqrt(2) * sh, [1,-1]))

        Xfit = np.concatenate(Xfit, 0)

        return Xfit 

    def __call__(self, diag):
        """
        Apply Silhouette on a single persistence diagram and outputs the result.

        Parameters:
            diag (n x 2 numpy array): input persistence diagram.

        Returns:
            numpy array with shape (**resolution**): output persistence silhouette.
        """
        return self.fit_transform([diag])[0,:]

class BettiCurve(BaseEstimator, TransformerMixin):
    """
    This is a class for computing Betti curves from a list of persistence diagrams. A Betti curve is a 1D piecewise-constant function obtained from the rank function. It is sampled evenly on a given range and the vector of samples is returned. See https://www.researchgate.net/publication/316604237_Time_Series_Classification_via_Topological_Data_Analysis for more details.
    """
    def __init__(self, resolution=100, sample_range=[np.nan, np.nan]):
        """
        Constructor for the BettiCurve class.

        Parameters:
            resolution (int): number of sample for the piecewise-constant function (default 100).
            sample_range ([double, double]): minimum and maximum of the piecewise-constant function domain, of the form [x_min, x_max] (default [numpy.nan, numpy.nan]). It is the interval on which samples will be drawn evenly. If one of the values is numpy.nan, it can be computed from the persistence diagrams with the fit() method.
        """
        self.resolution, self.sample_range = resolution, sample_range

    def fit(self, X, y=None):
        """
        Fit the BettiCurve class on a list of persistence diagrams: if any of the values in **sample_range** is numpy.nan, replace it with the corresponding value computed on the given list of persistence diagrams.

        Parameters:
            X (list of n x 2 numpy arrays): input persistence diagrams.
            y (n x 1 array): persistence diagram labels (unused).
        """
        if np.isnan(np.array(self.sample_range)).any():
            pre = DiagramScaler(use=True, scalers=[([0], MinMaxScaler()), ([1], MinMaxScaler())]).fit(X,y)
            [mx,my],[Mx,My] = [pre.scalers[0][1].data_min_[0], pre.scalers[1][1].data_min_[0]], [pre.scalers[0][1].data_max_[0], pre.scalers[1][1].data_max_[0]]
            self.sample_range = np.where(np.isnan(np.array(self.sample_range)), np.array([mx, My]), np.array(self.sample_range))
        return self

    def transform(self, X):
        """
        Compute the Betti curve for each persistence diagram individually and concatenate the results.

        Parameters:
            X (list of n x 2 numpy arrays): input persistence diagrams.
    
        Returns:
            numpy array with shape (number of diagrams) x (**resolution**): output Betti curves.
        """
        Xfit = []
        x_values = np.linspace(self.sample_range[0], self.sample_range[1], self.resolution)
        step_x = x_values[1] - x_values[0]

        for diagram in X:
            diagram_int = np.clip(np.ceil((diagram[:,:2] - self.sample_range[0]) / step_x), 0, self.resolution).astype(int)
            bc =  np.zeros(self.resolution)
            for interval in diagram_int:
                bc[interval[0]:interval[1]] += 1
            Xfit.append(np.reshape(bc,[1,-1]))

        Xfit = np.concatenate(Xfit, 0)

        return Xfit

    def __call__(self, diag):
        """
        Apply BettiCurve on a single persistence diagram and outputs the result.

        Parameters:
            diag (n x 2 numpy array): input persistence diagram.

        Returns:
            numpy array with shape (**resolution**): output Betti curve.
        """
        return self.fit_transform([diag])[0,:]

class Entropy(BaseEstimator, TransformerMixin):
    """
    This is a class for computing persistence entropy. Persistence entropy is a statistic for persistence diagrams inspired from Shannon entropy. This statistic can also be used to compute a feature vector, called the entropy summary function. See https://arxiv.org/pdf/1803.08304.pdf for more details. Note that a previous implementation was contributed by Manuel Soriano-Trigueros.
    """
    def __init__(self, mode="scalar", normalized=True, resolution=100, sample_range=[np.nan, np.nan]):
        """
        Constructor for the Entropy class.

        Parameters:
            mode (string): what entropy to compute: either "scalar" for computing the entropy statistics, or "vector" for computing the entropy summary functions (default "scalar").
            normalized (bool): whether to normalize the entropy summary function (default True). Used only if **mode** = "vector". 
            resolution (int): number of sample for the entropy summary function (default 100). Used only if **mode** = "vector".
            sample_range ([double, double]): minimum and maximum of the entropy summary function domain, of the form [x_min, x_max] (default [numpy.nan, numpy.nan]). It is the interval on which samples will be drawn evenly. If one of the values is numpy.nan, it can be computed from the persistence diagrams with the fit() method. Used only if **mode** = "vector".
        """
        self.mode, self.normalized, self.resolution, self.sample_range = mode, normalized, resolution, sample_range

    def fit(self, X, y=None):
        """
        Fit the Entropy class on a list of persistence diagrams.

        Parameters:
            X (list of n x 2 numpy arrays): input persistence diagrams.
            y (n x 1 array): persistence diagram labels (unused).
        """
        if np.isnan(np.array(self.sample_range)).any():
            pre = DiagramScaler(use=True, scalers=[([0], MinMaxScaler()), ([1], MinMaxScaler())]).fit(X,y)
            [mx,my],[Mx,My] = [pre.scalers[0][1].data_min_[0], pre.scalers[1][1].data_min_[0]], [pre.scalers[0][1].data_max_[0], pre.scalers[1][1].data_max_[0]]
            self.sample_range = np.where(np.isnan(np.array(self.sample_range)), np.array([mx, My]), np.array(self.sample_range))
        return self

    def transform(self, X):
        """
        Compute the entropy for each persistence diagram individually and concatenate the results.

        Parameters:
            X (list of n x 2 numpy arrays): input persistence diagrams.
    
        Returns:
            numpy array with shape (number of diagrams) x (1 if **mode** = "scalar" else **resolution**): output entropy.
        """
        num_diag, Xfit = len(X), []
        x_values = np.linspace(self.sample_range[0], self.sample_range[1], self.resolution)
        step_x = x_values[1] - x_values[0]
        new_X = BirthPersistenceTransform().fit_transform(X)        

        for i in range(num_diag):

            orig_diagram, diagram, num_pts_in_diag = X[i], new_X[i], X[i].shape[0]
            new_diagram = DiagramScaler(use=True, scalers=[([1], MaxAbsScaler())]).fit_transform([diagram])[0]

            if self.mode == "scalar":
                ent = - np.sum( np.multiply(new_diagram[:,1], np.log(new_diagram[:,1])) )
                Xfit.append(np.array([[ent]]))

            else:
                ent = np.zeros(self.resolution)
                for j in range(num_pts_in_diag):
                    [px,py] = orig_diagram[j,:2]
                    min_idx = np.clip(np.ceil((px - self.sample_range[0]) / step_x).astype(int), 0, self.resolution)
                    max_idx = np.clip(np.ceil((py - self.sample_range[0]) / step_x).astype(int), 0, self.resolution)
                    for k in range(min_idx, max_idx):
                        ent[k] += (-1) * new_diagram[j,1] * np.log(new_diagram[j,1])
                    if self.normalized:
                        ent = ent / np.linalg.norm(ent, ord=1)
                    Xfit.append(np.reshape(ent,[1,-1]))

        Xfit = np.concatenate(Xfit, 0)

        return Xfit

    def __call__(self, diag):
        """
        Apply Entropy on a single persistence diagram and outputs the result.

        Parameters:
            diag (n x 2 numpy array): input persistence diagram.

        Returns:
            numpy array with shape (1 if **mode** = "scalar" else **resolution**): output entropy.
        """
        return self.fit_transform([diag])[0,:]

class TopologicalVector(BaseEstimator, TransformerMixin):
    """
    This is a class for computing topological vectors from a list of persistence diagrams. The topological vector associated to a persistence diagram is the sorted vector of a slight modification of the pairwise distances between the persistence diagram points. See https://diglib.eg.org/handle/10.1111/cgf12692 for more details.
    """
    def __init__(self, threshold=10):
        """
        Constructor for the TopologicalVector class.

        Parameters:
            threshold (int): number of distances to keep (default 10). This is the dimension of the topological vector. If -1, this threshold is computed from the list of persistence diagrams by considering the one with the largest number of points and using the dimension of its corresponding topological vector as threshold. 
        """
        self.threshold = threshold

    def fit(self, X, y=None):
        """
        Fit the TopologicalVector class on a list of persistence diagrams (this function actually does nothing but is useful when TopologicalVector is included in a scikit-learn Pipeline).

        Parameters:
            X (list of n x 2 or n x 1 numpy arrays): input persistence diagrams.
            y (n x 1 array): persistence diagram labels (unused).
        """
        return self

    def transform(self, X):
        """
        Compute the topological vector for each persistence diagram individually and concatenate the results.

        Parameters:
            X (list of n x 2 numpy arrays): input persistence diagrams.
    
        Returns:
            numpy array with shape (number of diagrams) x (**threshold**): output topological vectors.
        """
        if self.threshold == -1:
            thresh = np.array([X[i].shape[0] for i in range(len(X))]).max()
        else:
            thresh = self.threshold

        num_diag = len(X)
        Xfit = np.zeros([num_diag, thresh])

        for i in range(num_diag):

            diagram, num_pts_in_diag = X[i], X[i].shape[0]
            pers = 0.5 * (diagram[:,1]-diagram[:,0])
            min_pers = np.minimum(pers,np.transpose(pers))
            distances = DistanceMetric.get_metric("chebyshev").pairwise(diagram)
            vect = np.flip(np.sort(np.triu(np.minimum(distances, min_pers)), axis=None), 0)
            dim = min(len(vect), thresh)
            Xfit[i, :dim] = vect[:dim]

        return Xfit

    def __call__(self, diag):
        """
        Apply TopologicalVector on a single persistence diagram and outputs the result.

        Parameters:
            diag (n x 2 numpy array): input persistence diagram.

        Returns:
            numpy array with shape (**threshold**): output topological vector.
        """
        return self.fit_transform([diag])[0,:]

class ComplexPolynomial(BaseEstimator, TransformerMixin):
    """
    This is a class for computing complex polynomials from a list of persistence diagrams. The persistence diagram points are seen as the roots of some complex polynomial, whose coefficients are returned in a complex vector. See https://link.springer.com/chapter/10.1007%2F978-3-319-23231-7_27 for more details.
    """
    def __init__(self, polynomial_type="R", threshold=10):
        """
        Constructor for the ComplexPolynomial class.

        Parameters:
           polynomial_type (char): either "R", "S" or "T" (default "R"). Type of complex polynomial that is going to be computed (explained in https://link.springer.com/chapter/10.1007%2F978-3-319-23231-7_27).
           threshold (int): number of coefficients (default 10). This is the dimension of the complex vector of coefficients, i.e. the number of coefficients corresponding to the largest degree terms of the polynomial. If -1, this threshold is computed from the list of persistence diagrams by considering the one with the largest number of points and using the dimension of its corresponding complex vector of coefficients as threshold. 
        """
        self.threshold, self.polynomial_type = threshold, polynomial_type

    def fit(self, X, y=None):
        """
        Fit the ComplexPolynomial class on a list of persistence diagrams (this function actually does nothing but is useful when ComplexPolynomial is included in a scikit-learn Pipeline).

        Parameters:
            X (list of n x 2 or n x 1 numpy arrays): input persistence diagrams.
            y (n x 1 array): persistence diagram labels (unused).
        """
        return self

    def transform(self, X):
        """
        Compute the complex vector of coefficients for each persistence diagram individually and concatenate the results.

        Parameters:
            X (list of n x 2 numpy arrays): input persistence diagrams.
    
        Returns:
            numpy array with shape (number of diagrams) x (**threshold**): output complex vectors of coefficients.
        """
        if self.threshold == -1:
            thresh = np.array([X[i].shape[0] for i in range(len(X))]).max()
        else:
            thresh = self.threshold

        Xfit = np.zeros([len(X), thresh]) + 1j * np.zeros([len(X), thresh])
        for d in range(len(X)):
            D, N = X[d], X[d].shape[0]
            if self.polynomial_type == "R":
                roots = D[:,0] + 1j * D[:,1]
            elif self.polynomial_type == "S":
                alpha = np.linalg.norm(D, axis=1)
                alpha = np.where(alpha==0, np.ones(N), alpha)
                roots = np.multiply( np.multiply(  (D[:,0]+1j*D[:,1]), (D[:,1]-D[:,0])  ), 1./(np.sqrt(2)*alpha) )
            elif self.polynomial_type == "T":
                alpha = np.linalg.norm(D, axis=1)
                roots = np.multiply(  (D[:,1]-D[:,0])/2, np.cos(alpha) - np.sin(alpha) + 1j * (np.cos(alpha) + np.sin(alpha))  )
            coeff = [0] * (N+1)
            coeff[N] = 1
            for i in range(1, N+1): 
                for j in range(N-i-1, N): 
                    coeff[j] += ((-1) * roots[i-1] * coeff[j+1])   
            coeff = np.array(coeff[::-1])[1:]
            Xfit[d, :min(thresh, coeff.shape[0])] = coeff[:min(thresh, coeff.shape[0])]
        return Xfit

    def __call__(self, diag):
        """
        Apply ComplexPolynomial on a single persistence diagram and outputs the result.

        Parameters:
            diag (n x 2 numpy array): input persistence diagram.

        Returns:
            numpy array with shape (**threshold**): output complex vector of coefficients.
        """
        return self.fit_transform([diag])[0,:]

def _lapl_contrast(measure, centers, inertias):
    """contrast function for vectorising `measure` in ATOL"""
    return np.exp(-pairwise.pairwise_distances(measure, Y=centers) / inertias)

def _gaus_contrast(measure, centers, inertias):
    """contrast function for vectorising `measure` in ATOL"""
    return np.exp(-pairwise.pairwise_distances(measure, Y=centers, squared=True) / inertias**2)

def _indicator_contrast(diags, centers, inertias):
    """contrast function for vectorising `measure` in ATOL"""
    robe_curve = np.clip(2-pairwise.pairwise_distances(diags, Y=centers)/inertias, 0, 1)
    return robe_curve

def _cloud_weighting(measure):
    """automatic uniform weighting with mass 1 for `measure` in ATOL"""
    return np.ones(shape=measure.shape[0])

def _iidproba_weighting(measure):
    """automatic uniform weighting with mass 1/N for `measure` in ATOL"""
    return np.ones(shape=measure.shape[0]) / measure.shape[0]

class Atol(BaseEstimator, TransformerMixin):
    """
    This class allows to vectorise measures (e.g. point clouds, persistence diagrams, etc) after a quantisation step.

    ATOL paper: :cite:`royer2019atol`

    Example
    --------
    >>> from sklearn.cluster import KMeans
    >>> from gudhi.representations.vector_methods import Atol
    >>> import numpy as np
    >>> a = np.array([[1, 2, 4], [1, 4, 0], [1, 0, 4]])
    >>> b = np.array([[4, 2, 0], [4, 4, 0], [4, 0, 2]])
    >>> c = np.array([[3, 2, -1], [1, 2, -1]])
    >>> atol_vectoriser = Atol(quantiser=KMeans(n_clusters=2, random_state=202006))
    >>> atol_vectoriser.fit(X=[a, b, c]).centers
    >>> # array([[ 2.        ,  0.66666667,  3.33333333],
    >>> #        [ 2.6       ,  2.8       , -0.4       ]])
    >>> atol_vectoriser(a)
    >>> # array([1.18168665, 0.42375966])
    >>> atol_vectoriser(c)
    >>> # array([0.02062512, 1.25157463])
    >>> atol_vectoriser.transform(X=[a, b, c])
    >>> # array([[1.18168665, 0.42375966],
    >>> #        [0.29861028, 1.06330156],
    >>> #        [0.02062512, 1.25157463]])
    """
    def __init__(self, quantiser, weighting_method="cloud", contrast="gaussian"):
        """
        Constructor for the Atol measure vectorisation class.

        Parameters:
            quantiser (Object): Object with `fit` (sklearn API consistent) and `cluster_centers` and `n_clusters`
                attributes, e.g. sklearn.cluster.KMeans. It will be fitted when the Atol object function `fit` is called.
            weighting_method (string): constant generic function for weighting the measure points
                choose from {"cloud", "iidproba"}
                (default: constant function, i.e. the measure is seen as a point cloud by default).
                This will have no impact if weights are provided along with measures all the way: `fit` and `transform`.
            contrast (string): constant function for evaluating proximity of a measure with respect to centers
                choose from {"gaussian", "laplacian", "indicator"}
                (default: gaussian contrast function, see page 3 in the ATOL paper).
        """
        self.quantiser = quantiser
        self.contrast = {
            "gaussian": _gaus_contrast,
            "laplacian": _lapl_contrast,
            "indicator": _indicator_contrast,
        }.get(contrast, _gaus_contrast)
        self.weighting_method = {
            "cloud"   : _cloud_weighting,
            "iidproba": _iidproba_weighting,
        }.get(weighting_method, _cloud_weighting)

    def fit(self, X, y=None, sample_weight=None):
        """
        Calibration step: fit centers to the sample measures and derive inertias between centers.

        Parameters:
            X (list N x d numpy arrays): input measures in R^d from which to learn center locations and inertias
                (measures can have different N).
            y: Ignored, present for API consistency by convention.
            sample_weight (list of numpy arrays): weights for each measure point in X, optional.
                If None, the object's weighting_method will be used.

        Returns:
            self
        """
        if not hasattr(self.quantiser, 'fit'):
            raise TypeError("quantiser %s has no `fit` attribute." % (self.quantiser))
        if sample_weight is None:
            sample_weight = np.concatenate([self.weighting_method(measure) for measure in X])

        measures_concat = np.concatenate(X)
        self.quantiser.fit(X=measures_concat, sample_weight=sample_weight)
        self.centers = self.quantiser.cluster_centers_
        if self.quantiser.n_clusters == 1:
            dist_centers = pairwise.pairwise_distances(measures_concat)
            np.fill_diagonal(dist_centers, 0)
            self.inertias = np.array([np.max(dist_centers)/2])
        else:
            dist_centers = pairwise.pairwise_distances(self.centers)
            dist_centers[dist_centers == 0] = np.inf
            self.inertias = np.min(dist_centers, axis=0)/2
        return self

    def __call__(self, measure, sample_weight=None):
        """
        Apply measure vectorisation on a single measure.

        Parameters:
            measure (n x d numpy array): input measure in R^d.

        Returns:
            numpy array in R^self.quantiser.n_clusters.
        """
        if sample_weight is None:
            sample_weight = self.weighting_method(measure)
        return np.sum(sample_weight * self.contrast(measure, self.centers, self.inertias.T).T, axis=1)

    def transform(self, X, sample_weight=None):
        """
        Apply measure vectorisation on a list of measures.

        Parameters:
            X (list N x d numpy arrays): input measures in R^d from which to learn center locations and inertias
                (measures can have different N).
            sample_weight (list of numpy arrays): weights for each measure point in X, optional.
                If None, the object's weighting_method will be used.

        Returns:
            numpy array with shape (number of measures) x (self.quantiser.n_clusters).
        """
        if sample_weight is None:
            sample_weight = [self.weighting_method(measure) for measure in X]
        return np.stack([self(measure, sample_weight=weight) for measure, weight in zip(X, sample_weight)])