blob: 9033331f8b62dc4d2a9fa9df3000a5aacf3951af (
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
|
Persistence graphical tools user manual
=======================================
Definition
----------
.. include:: persistence_graphical_tools_sum.rst
Show palette values
-------------------
This function is useful to show the color palette values of dimension:
.. testcode::
import gudhi
plt = gudhi.show_palette_values(alpha=1.0)
plt.show()
.. plot::
import gudhi
plt = gudhi.show_palette_values(alpha=1.0)
plt.show()
Show persistence as a barcode
-----------------------------
This function can display the persistence result as a barcode:
.. testcode::
import gudhi
periodic_cc = gudhi.PeriodicCubicalComplex(perseus_file=gudhi.__root_source_dir__ + \
'/data/bitmap/3d_torus.txt')
diag = periodic_cc.persistence()
plt = gudhi.plot_persistence_barcode(diag)
plt.show()
.. plot::
import gudhi
periodic_cc = gudhi.PeriodicCubicalComplex(perseus_file=gudhi.__root_source_dir__ + \
'/data/bitmap/3d_torus.txt')
diag = periodic_cc.persistence()
print("diag = ", diag)
plt = gudhi.plot_persistence_barcode(diag)
plt.show()
Show persistence as a diagram
-----------------------------
This function can display the persistence result as a diagram:
.. testcode::
import gudhi
rips_complex = gudhi.RipsComplex(off_file=gudhi.__root_source_dir__ + \
'/data/points/tore3D_1307.off', max_edge_length=0.2)
simplex_tree = rips_complex.create_simplex_tree(max_dimension=3)
diag = simplex_tree.persistence()
plt = gudhi.plot_persistence_diagram(diag, band_boot=0.13)
plt.show()
.. plot::
import gudhi
rips_complex = gudhi.RipsComplex(off_file=gudhi.__root_source_dir__ + \
'/data/points/tore3D_1307.off', max_edge_length=0.2)
simplex_tree = rips_complex.create_simplex_tree(max_dimension=3)
diag = simplex_tree.persistence()
plt = gudhi.plot_persistence_diagram(diag, band_boot=0.13)
plt.show()
|