summaryrefslogtreecommitdiff
path: root/src/Alpha_complex/doc/Intro_alpha_complex.h
blob: c068b268b1e6609d3bc5e6c8db3c074f8c664115 (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
/*    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):       Vincent Rouvreau
 *
 *    Copyright (C) 2015 Inria
 *
 *    Modification(s):
 *      - YYYY/MM Author: Description of the modification
 */

#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
 *
 * @{
 *
<div class="toc">
Table of Contents
<ul>
<li class="level1"><a href="#definition">Definition</a></li>
<li class="level1"><a href="#pointsexample">Example from points</a></li>
<li class="level1"><a href="#createcomplexalgorithm">Create complex algorithm</a></li>
<li class="level1"><a href="#weightedversion">Weighted specific version</a></li>
<li class="level1"><a href="#offexample">Example from OFF file</a></li>
<li class="level1"><a href="#weighted3dexample">3d specific version</a></li>
</ul>
</div>

 * \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 \f$ > \alpha^2 \f$ are removed from the Delaunay complex
 * when creating the simplicial complex if it is specified.
 *
 * \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-19b from <a target="_blank" href="http://www.cgal.org/">CGAL</a> (the Computational Geometry
 * Algorithms Library \cite cgal:eb-19b) 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-19b from CGAL as template parameter.
 *
 * \remark
 * - When the simplicial complex is constructed with an infinite value of \f$ \alpha^2 \f$, the complex is a Delaunay
 * complex with special filtration values. The Delaunay complex without filtration values is also available by passing
 * `default_filtration_value=true` to `Alpha_complex::create_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.
 * - Using the default `CGAL::Epeck_d` makes the construction safe. If you pass `exact=true` to
 * `Alpha_complex::create_complex`, the filtration values are the exact ones converted to the filtration value type of
 * the simplicial complex. This can be very slow. If you pass `exact=false` (the default), the filtration values are
 * only guaranteed to have a small multiplicative error compared to the exact value, see <code>
 * <a class="el" target="_blank" href="https://doc.cgal.org/latest/Number_types/classCGAL_1_1Lazy__exact__nt.html">
 * CGAL::Lazy_exact_nt<NT>::set_relative_precision_of_to_double</a></code> for details. A drawback, when computing
 * persistence, is that an empty exact interval [10^12,10^12] may become a non-empty approximate interval
 * [10^12,10^12+10^6]. Using `CGAL::Epick_d` makes the computations slightly faster, and the combinatorics are still
 * exact, but the computation of filtration values can exceptionally be arbitrarily bad. In all cases, we still
 * guarantee that the output is a valid filtration (faces have a filtration value no larger than their cofaces).
 * - For performances reasons, it is advised to use \ref eigen &ge; 3.3.5 and \ref cgal &ge; 5.2.0.
 *
 * \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 } \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 \f$ \alpha^2 \f$ 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 \f$ \alpha^2 \f$ value (cf.
 * `SimplicialComplexForAlpha::prune_above_filtration()`).
 * In the following example, the value is given by the user as argument of the program.
 *
 * \section weightedversion Weighted specific version
 * <b>Requires:</b> \ref eigen &ge; 3.1.0 and \ref cgal &ge; 5.1.0.
 * 
 * A weighted version for Alpha complex is available (cf. Alpha_complex). It is like a usual Alpha complex, but based
 * on a <a href="https://doc.cgal.org/latest/Triangulation/index.html#title20">CGAL regular triangulation</a> instead
 * of Delaunay.
 *
 * This example builds the CGAL weighted alpha shapes from a small molecule, and initializes the alpha complex with
 * it. This example is taken from <a href="https://doc.cgal.org/latest/Alpha_shapes_3/index.html#title13">CGAL 3d
 * weighted alpha shapes</a>.
 *
 * Then, it is asked to display information about the alpha complex.
 *
 * \include Alpha_complex/Weighted_alpha_complex_from_points.cpp
 *
 * When launching:
 *
 * \code $> ./Weighted_alpha_complex_example_from_points
 * \endcode
 *
 * the program output is:
 *
 * \include Alpha_complex/weightedalpha3dfrompoints_for_doc.txt
 *
 *
 * \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
 *
 *
 * \section weighted3dexample 3d specific version
 *
 * A specific module for Alpha complex is available in 3d (cf. Alpha_complex_3d) and allows to construct standard,
 * weighted, periodic or weighted and periodic versions of alpha complexes. Alpha values computation can be
 * Gudhi::alpha_complex::complexity::FAST, Gudhi::alpha_complex::complexity::SAFE (default value) or
 * Gudhi::alpha_complex::complexity::EXACT.
 *
 * This example builds the CGAL 3d weighted alpha shapes from a small molecule, and initializes the alpha complex with
 * it. This example is taken from <a href="https://doc.cgal.org/latest/Alpha_shapes_3/index.html#title13">CGAL 3d
 * weighted alpha shapes</a>.
 *
 * Then, it is asked to display information about the alpha complex.
 *
 * \include Alpha_complex/Weighted_alpha_complex_3d_from_points.cpp
 *
 * The results will be the same as in \ref weightedversion .
 *
 */
/** @} */  // end defgroup alpha_complex

}  // namespace alpha_complex

namespace alphacomplex = alpha_complex;

}  // namespace Gudhi

#endif  // DOC_ALPHA_COMPLEX_INTRO_ALPHA_COMPLEX_H_