summaryrefslogtreecommitdiff
path: root/doc/Alpha_complex/Intro_alpha_complex.h
blob: 7a375c9f4fdd0046814acd5ccc02812155dcad45 (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
/*    This file is part of the Gudhi Library. The Gudhi library
 *    (Geometric Understanding in Higher Dimensions) is a generic C++
 *    library for computational topology.
 *
 *    Author(s):       Vincent Rouvreau
 *
 *    Copyright (C) 2015 Inria
 *
 *    This program is free software: you can redistribute it and/or modify
 *    it under the terms of the GNU General Public License as published by
 *    the Free Software Foundation, either version 3 of the License, or
 *    (at your option) any later version.
 *
 *    This program is distributed in the hope that it will be useful,
 *    but WITHOUT ANY WARRANTY; without even the implied warranty of
 *    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 *    GNU General Public License for more details.
 *
 *    You should have received a copy of the GNU General Public License
 *    along with this program.  If not, see <http://www.gnu.org/licenses/>.
 */

#ifndef DOC_ALPHA_COMPLEX_INTRO_ALPHA_COMPLEX_H_
#define DOC_ALPHA_COMPLEX_INTRO_ALPHA_COMPLEX_H_

// needs namespace for Doxygen to link on classes
namespace Gudhi {
// needs namespace for Doxygen to link on classes
namespace alpha_complex {

/**  \defgroup alpha_complex Alpha complex
 * 
 * \author    Vincent Rouvreau
 *
 * @{
 * 
 * \section definition Definition
 * 
 * Alpha_complex is a <a target="_blank" href="https://en.wikipedia.org/wiki/Simplicial_complex">simplicial complex</a>
 * constructed from the finite cells of a Delaunay Triangulation.
 * 
 * The filtration value of each simplex is computed as the square of the circumradius of the simplex if the
 * circumsphere is empty (the simplex is then said to be Gabriel), and as the minimum of the filtration
 * values of the codimension 1 cofaces that make it not Gabriel otherwise.
 * 
 * All simplices that have a filtration value strictly greater than a given alpha squared value are not inserted into
 * the complex.
 * 
 * \image html "alpha_complex_representation.png" "Alpha-complex representation"
 * 
 * Alpha_complex is constructing a <a target="_blank"
 * href="http://doc.cgal.org/latest/Triangulation/index.html#Chapter_Triangulations">Delaunay Triangulation</a>
 * \cite cgal:hdj-t-15b from <a target="_blank" href="http://www.cgal.org/">CGAL</a> (the Computational Geometry
 * Algorithms Library \cite cgal:eb-15b) and is able to create a `SimplicialComplexForAlpha`.
 * 
 * The complex is a template class requiring an Epick_d <a target="_blank"
 * href="http://doc.cgal.org/latest/Kernel_d/index.html#Chapter_dD_Geometry_Kernel">dD Geometry Kernel</a>
 * \cite cgal:s-gkd-15b from CGAL as template parameter.
 * 
 * \remark
 * - When the simplicial complex is constructed with an infinite value of alpha, the complex is a Delaunay
 * complex.
 * - For people only interested in the topology of the \ref alpha_complex (for instance persistence),
 * \ref alpha_complex is equivalent to the \ref cech_complex and much smaller if you do not bound the radii.
 * \ref cech_complex can still make sense in higher dimension precisely because you can bound the radii.
 *
 * \section pointsexample Example from points
 * 
 * This example builds the Delaunay triangulation from the given points in a 2D static kernel, and creates a
 * `Simplex_tree` with it.
 * 
 * Then, it is asked to display information about the simplicial complex.
 * 
 * \include Alpha_complex/Alpha_complex_from_points.cpp
 * 
 * When launching:
 * 
 * \code $> ./Alpha_complex_example_from_points
 * \endcode
 *
 * the program output is:
 * 
 * \include Alpha_complex/alphaoffreader_for_doc_60.txt
 * 
 * \section createcomplexalgorithm Create complex algorithm
 * 
 * \subsection datastructure Data structure
 * 
 * In order to create the simplicial complex, first, it is built from the cells of the Delaunay Triangulation.
 * The filtration values are set to NaN, which stands for unknown value.
 * 
 * In example, :
 * \image html "alpha_complex_doc.png" "Simplicial complex structure construction example"
 *
 * \subsection filtrationcomputation Filtration value computation algorithm
 * <br>
 * \f$
 * \textbf{for } \text{i : dimension } \rightarrow 0 \textbf{ do}\\
 * \quad \textbf{for all } \sigma \text{ of dimension i}\\
 * \quad\quad \textbf{if } \text{filtration(} \sigma ) \text{ is NaN} \textbf{ then}\\
 * \quad\quad\quad \text{filtration(} \sigma ) = \alpha^2( \sigma )\\
 * \quad\quad \textbf{end if}\\
 * \quad\quad \textbf{for all } \tau \text{ face of } \sigma \textbf{ do}\quad\quad
 * \textit{// propagate alpha filtration value}\\
 * \quad\quad\quad \textbf{if } \text{filtration(} \tau ) \text{ is not NaN} \textbf{ then}\\
 * \quad\quad\quad\quad \text{filtration(} \tau \text{) = min( filtration(} \tau \text{), filtration(} \sigma
 * \text{) )}\\
 * \quad\quad\quad \textbf{else}\\
 * \quad\quad\quad\quad \textbf{if } \textbf{if } \tau \text{ is not Gabriel for } \sigma \textbf{ then}\\
 * \quad\quad\quad\quad\quad \text{filtration(} \tau \text{) = filtration(} \sigma \text{)}\\
 * \quad\quad\quad\quad \textbf{end if}\\
 * \quad\quad\quad \textbf{end if}\\
 * \quad\quad \textbf{end for}\\
 * \quad \textbf{end for}\\
 * \textbf{end for}\\
 * \text{make_filtration_non_decreasing()}\\
 * \text{prune_above_filtration()}\\
 * \f$
 *
 * \subsubsection dimension2 Dimension 2
 * 
 * From the example above, it means the algorithm looks into each triangle ([0,1,2], [0,2,4], [1,2,3], ...),
 * computes the filtration value of the triangle, and then propagates the filtration value as described
 * here :
 * \image html "alpha_complex_doc_420.png" "Filtration value propagation example"
 * 
 * \subsubsection dimension1 Dimension 1
 * 
 * Then, the algorithm looks into each edge ([0,1], [0,2], [1,2], ...),
 * computes the filtration value of the edge (in this case, propagation will have no effect).
 * 
 * \subsubsection dimension0 Dimension 0
 * 
 * Finally, the algorithm looks into each vertex ([0], [1], [2], [3], [4], [5] and [6]) and
 * sets the filtration value (0 in case of a vertex - propagation will have no effect).
 * 
 * \subsubsection nondecreasing Non decreasing filtration values
 * 
 * As the squared radii computed by CGAL are an approximation, it might happen that these alpha squared values do not
 * quite define a proper filtration (i.e. non-decreasing with respect to inclusion).
 * We fix that up by calling `SimplicialComplexForAlpha::make_filtration_non_decreasing()`.
 * 
 * \subsubsection pruneabove Prune above given filtration value
 * 
 * The simplex tree is pruned from the given maximum alpha squared value (cf.
 * `SimplicialComplexForAlpha::prune_above_filtration()`).
 * In the following example, the value is given by the user as argument of the program.
 * 
 * 
 * \section offexample Example from OFF file
 * 
 * This example builds the Delaunay triangulation in a dynamic kernel, and initializes the alpha complex with it.
 * 
 * 
 * Then, it is asked to display information about the alpha complex.
 * 
 * \include Alpha_complex/Alpha_complex_from_off.cpp
 * 
 * When launching:
 * 
 * \code $> ./Alpha_complex_example_from_off ../../data/points/alphacomplexdoc.off 32.0
 * \endcode
 *
 * the program output is:
 * 
 * \include Alpha_complex/alphaoffreader_for_doc_32.txt
 * 
 */
/** @} */  // end defgroup alpha_complex

}  // namespace alpha_complex

namespace alphacomplex = alpha_complex;

}  // namespace Gudhi

#endif  // DOC_ALPHA_COMPLEX_INTRO_ALPHA_COMPLEX_H_