summaryrefslogtreecommitdiff
path: root/cython/doc/rips_complex_user.rst
blob: 96ba9944f9556cfb092af3fd29e748ba614d847f (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
Rips complex user manual
=========================
Definition
----------

=======================================================  =====================================  =====================================
:Authors: Clément Maria, Pawel Dlotko, Vincent Rouvreau  :Introduced in: GUDHI 2.0.0            :Copyright: GPL v3
=======================================================  =====================================  =====================================

+-------------------------------------------+----------------------------------------------------------------------+
| :doc:`rips_complex_user`                  | :doc:`rips_complex_ref`                                              |
+-------------------------------------------+----------------------------------------------------------------------+

`Rips complex <https://en.wikipedia.org/wiki/Vietoris%E2%80%93Rips_complex>`_ is a one skeleton graph that allows to
construct a simplicial complex from it. The input can be a point cloud with a given distance function, or a distance
matrix.

The filtration value of each edge is computed from a user-given distance function, or directly from the distance
matrix.

All edges that have a filtration value strictly greater than a given threshold value are not inserted into the complex.

When creating a simplicial complex from this one skeleton graph, Rips inserts the one skeleton graph into the data
structure, and then expands the simplicial complex when required.

Vertex name correspond to the index of the point in the given range (aka. the point cloud).

.. figure::
    ../../doc/Rips_complex/rips_complex_representation.png
    :align: center

    Rips-complex one skeleton graph representation

On this example, as edges (4,5), (4,6) and (5,6) are in the complex, simplex (4,5,6) is added with the filtration value
set with :math:`max(filtration(4,5), filtration(4,6), filtration(5,6))`. And so on for simplex (0,1,2,3).

If the Rips_complex interfaces are not detailed enough for your need, please refer to rips_persistence_step_by_step.cpp
example, where the graph construction over the Simplex_tree is more detailed.

Point cloud
-----------

Example from a point cloud
^^^^^^^^^^^^^^^^^^^^^^^^^^

This example builds the one skeleton graph from the given points, and max_edge_length value.
Then it creates a :doc:`Simplex_tree <simplex_tree_ref>` with it.

Finally, it is asked to display information about the simplicial complex.

.. testcode::

    import gudhi
    rips_complex = gudhi.RipsComplex(points=[[1, 1], [7, 0], [4, 6], [9, 6], [0, 14], [2, 19], [9, 17]],
        max_edge_length=12.0)

    simplex_tree = rips_complex.create_simplex_tree(max_dimension=1)
    result_str = 'Rips complex is of dimension ' + repr(simplex_tree.dimension()) + ' - ' + \
        repr(simplex_tree.num_simplices()) + ' simplices - ' + \
        repr(simplex_tree.num_vertices()) + ' vertices.'
    print(result_str)
    fmt = '%s -> %.2f'
    for filtered_value in simplex_tree.get_filtration():
        print(fmt % tuple(filtered_value))

When launching (Rips maximal distance between 2 points is 12.0, is expanded
until dimension 1 - one skeleton graph in other words), the output is:

.. testoutput::

    Rips complex is of dimension 1 - 18 simplices - 7 vertices.
    [0] -> 0.00
    [1] -> 0.00
    [2] -> 0.00
    [3] -> 0.00
    [4] -> 0.00
    [5] -> 0.00
    [6] -> 0.00
    [2, 3] -> 5.00
    [4, 5] -> 5.39
    [0, 2] -> 5.83
    [0, 1] -> 6.08
    [1, 3] -> 6.32
    [1, 2] -> 6.71
    [5, 6] -> 7.28
    [2, 4] -> 8.94
    [0, 3] -> 9.43
    [4, 6] -> 9.49
    [3, 6] -> 11.00

Example from OFF file
^^^^^^^^^^^^^^^^^^^^^

This example builds the :doc:`Rips_complex <rips_complex_ref>` from the given
points in an OFF file, and max_edge_length value.
Then it creates a :doc:`Simplex_tree <simplex_tree_ref>` with it.

Finally, it is asked to display information about the Rips complex.


.. testcode::

    import gudhi
    rips_complex = gudhi.RipsComplex(off_file=gudhi.__root_source_dir__ + \
        '/data/points/alphacomplexdoc.off', max_edge_length=12.0)
    simplex_tree = rips_complex.create_simplex_tree(max_dimension=1)
    result_str = 'Rips complex is of dimension ' + repr(simplex_tree.dimension()) + ' - ' + \
        repr(simplex_tree.num_simplices()) + ' simplices - ' + \
        repr(simplex_tree.num_vertices()) + ' vertices.'
    print(result_str)
    fmt = '%s -> %.2f'
    for filtered_value in simplex_tree.get_filtration():
        print(fmt % tuple(filtered_value))

the program output is:

.. testoutput::

    Rips complex is of dimension 1 - 18 simplices - 7 vertices.
    [0] -> 0.00
    [1] -> 0.00
    [2] -> 0.00
    [3] -> 0.00
    [4] -> 0.00
    [5] -> 0.00
    [6] -> 0.00
    [2, 3] -> 5.00
    [4, 5] -> 5.39
    [0, 2] -> 5.83
    [0, 1] -> 6.08
    [1, 3] -> 6.32
    [1, 2] -> 6.71
    [5, 6] -> 7.28
    [2, 4] -> 8.94
    [0, 3] -> 9.43
    [4, 6] -> 9.49
    [3, 6] -> 11.00

Distance matrix
---------------

Example from a distance matrix
^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^

This example builds the one skeleton graph from the given distance matrix, and max_edge_length value.
Then it creates a :doc:`Simplex_tree <simplex_tree_ref>` with it.

Finally, it is asked to display information about the simplicial complex.

.. testcode::

    import gudhi
    rips_complex = gudhi.RipsComplex(distance_matrix=[[],
                                                      [6.0827625303],
                                                      [5.8309518948, 6.7082039325],
                                                      [9.4339811321, 6.3245553203, 5],
                                                      [13.0384048104, 15.6524758425, 8.94427191, 12.0415945788],
                                                      [18.0277563773, 19.6468827044, 13.152946438, 14.7648230602, 5.3851648071],
                                                      [17.88854382, 17.1172427686, 12.0830459736, 11, 9.4868329805, 7.2801098893]],
                                     max_edge_length=12.0)

    simplex_tree = rips_complex.create_simplex_tree(max_dimension=1)
    result_str = 'Rips complex is of dimension ' + repr(simplex_tree.dimension()) + ' - ' + \
        repr(simplex_tree.num_simplices()) + ' simplices - ' + \
        repr(simplex_tree.num_vertices()) + ' vertices.'
    print(result_str)
    fmt = '%s -> %.2f'
    for filtered_value in simplex_tree.get_filtration():
        print(fmt % tuple(filtered_value))

When launching (Rips maximal distance between 2 points is 12.0, is expanded
until dimension 1 - one skeleton graph in other words), the output is:

.. testoutput::

    Rips complex is of dimension 1 - 18 simplices - 7 vertices.
    [0] -> 0.00
    [1] -> 0.00
    [2] -> 0.00
    [3] -> 0.00
    [4] -> 0.00
    [5] -> 0.00
    [6] -> 0.00
    [2, 3] -> 5.00
    [4, 5] -> 5.39
    [0, 2] -> 5.83
    [0, 1] -> 6.08
    [1, 3] -> 6.32
    [1, 2] -> 6.71
    [5, 6] -> 7.28
    [2, 4] -> 8.94
    [0, 3] -> 9.43
    [4, 6] -> 9.49
    [3, 6] -> 11.00

Example from csv file
^^^^^^^^^^^^^^^^^^^^^

This example builds the :doc:`Rips_complex <rips_complex_ref>` from the given
points in an OFF file, and max_edge_length value.
Then it creates a :doc:`Simplex_tree <simplex_tree_ref>` with it.

Finally, it is asked to display information about the Rips complex.


.. testcode::

    import gudhi
    rips_complex = gudhi.RipsComplex(csv_file=gudhi.__root_source_dir__ + \
        '/data/distance_matrix/full_square_distance_matrix.csv', max_edge_length=12.0)
    simplex_tree = rips_complex.create_simplex_tree(max_dimension=1)
    result_str = 'Rips complex is of dimension ' + repr(simplex_tree.dimension()) + ' - ' + \
        repr(simplex_tree.num_simplices()) + ' simplices - ' + \
        repr(simplex_tree.num_vertices()) + ' vertices.'
    print(result_str)
    fmt = '%s -> %.2f'
    for filtered_value in simplex_tree.get_filtration():
        print(fmt % tuple(filtered_value))

the program output is:

.. testoutput::

    Rips complex is of dimension 1 - 18 simplices - 7 vertices.
    [0] -> 0.00
    [1] -> 0.00
    [2] -> 0.00
    [3] -> 0.00
    [4] -> 0.00
    [5] -> 0.00
    [6] -> 0.00
    [2, 3] -> 5.00
    [4, 5] -> 5.39
    [0, 2] -> 5.83
    [0, 1] -> 6.08
    [1, 3] -> 6.32
    [1, 2] -> 6.71
    [5, 6] -> 7.28
    [2, 4] -> 8.94
    [0, 3] -> 9.43
    [4, 6] -> 9.49
    [3, 6] -> 11.00