Alpha complex

Classes

class  Gudhi::alpha_complex::Alpha_complex< Kernel >
 Alpha complex data structure. More...
 

Detailed Description

Author
Vincent Rouvreau

Definition

Alpha_complex is a simplicial complex 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.

alpha_complex_representation.png
Alpha-complex representation

Alpha_complex is constructing a Delaunay Triangulation [13] from CGAL (the Computational Geometry Algorithms Library [19]) and is able to create a SimplicialComplexForAlpha.

The complex is a template class requiring an Epick_d dD Geometry Kernel [18] from CGAL as template parameter.

Remarks
When the simplicial complex is constructed with an infinite value of alpha, the complex is a Delaunay complex.

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 <gudhi/Alpha_complex.h>
// to construct a simplex_tree from alpha complex
#include <gudhi/Simplex_tree.h>
#include <CGAL/Epick_d.h>
#include <iostream>
#include <string>
#include <vector>
#include <limits> // for numeric limits
using Kernel = CGAL::Epick_d< CGAL::Dimension_tag<2> >;
using Point = Kernel::Point_d;
using Vector_of_points = std::vector<Point>;
void usage(int nbArgs, char * const progName) {
std::cerr << "Error: Number of arguments (" << nbArgs << ") is not correct\n";
std::cerr << "Usage: " << progName << " [alpha_square_max_value]\n";
std::cerr << " i.e.: " << progName << " 60.0\n";
exit(-1); // ----- >>
}
int main(int argc, char **argv) {
if ((argc != 1) && (argc != 2)) usage(argc, (argv[0] - 1));
// Delaunay complex if alpha_square_max_value is not given by the user.
double alpha_square_max_value {std::numeric_limits<double>::infinity()};
if (argc == 2)
alpha_square_max_value = atof(argv[1]);
// ----------------------------------------------------------------------------
// Init of a list of points
// ----------------------------------------------------------------------------
Vector_of_points points;
points.push_back(Point(1.0, 1.0));
points.push_back(Point(7.0, 0.0));
points.push_back(Point(4.0, 6.0));
points.push_back(Point(9.0, 6.0));
points.push_back(Point(0.0, 14.0));
points.push_back(Point(2.0, 19.0));
points.push_back(Point(9.0, 17.0));
// ----------------------------------------------------------------------------
// Init of an alpha complex from the list of points
// ----------------------------------------------------------------------------
Gudhi::alpha_complex::Alpha_complex<Kernel> alpha_complex_from_points(points);
if (alpha_complex_from_points.create_complex(simplex, alpha_square_max_value)) {
// ----------------------------------------------------------------------------
// Display information about the alpha complex
// ----------------------------------------------------------------------------
std::cout << "Alpha complex is of dimension " << simplex.dimension() <<
" - " << simplex.num_simplices() << " simplices - " <<
simplex.num_vertices() << " vertices." << std::endl;
std::cout << "Iterator on alpha complex simplices in the filtration order, with [filtration value]:" << std::endl;
for (auto f_simplex : simplex.filtration_simplex_range()) {
std::cout << " ( ";
for (auto vertex : simplex.simplex_vertex_range(f_simplex)) {
std::cout << vertex << " ";
}
std::cout << ") -> " << "[" << simplex.filtration(f_simplex) << "] ";
std::cout << std::endl;
}
}
return 0;
}

When launching:

$> ./Alpha_complex_example_from_points

the program output is:

Alpha complex is of dimension 2 - 25 simplices - 7 vertices.
Iterator on alpha complex simplices in the filtration order, with [filtration value]:
( 0 ) -> [0]
( 1 ) -> [0]
( 2 ) -> [0]
( 3 ) -> [0]
( 4 ) -> [0]
( 5 ) -> [0]
( 6 ) -> [0]
( 3 2 ) -> [6.25]
( 5 4 ) -> [7.25]
( 2 0 ) -> [8.5]
( 1 0 ) -> [9.25]
( 3 1 ) -> [10]
( 2 1 ) -> [11.25]
( 3 2 1 ) -> [12.5]
( 2 1 0 ) -> [12.9959]
( 6 5 ) -> [13.25]
( 4 2 ) -> [20]
( 6 4 ) -> [22.7367]
( 6 5 4 ) -> [22.7367]
( 6 3 ) -> [30.25]
( 6 2 ) -> [36.5]
( 6 3 2 ) -> [36.5]
( 6 4 2 ) -> [37.2449]
( 4 0 ) -> [59.7107]
( 4 2 0 ) -> [59.7107]

Create complex algorithm

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, :

alpha_complex_doc.png
Simplicial complex structure construction example

Filtration value computation algorithm

\begin{algorithm} \caption{Filtration value computation algorithm}\label{alpha} \begin{algorithmic} \For{i : dimension $\rightarrow$ 0} \ForAll{$\sigma$ of dimension i} \If {filtration($\sigma$) is NaN} \State filtration($\sigma$) = $\alpha^2(\sigma)$ \EndIf \ForAll{$\tau$ face of $\sigma$} \Comment{propagate alpha filtration value} \If {filtration($\tau$) is not NaN} \State filtration($\tau$) = min (filtration($\tau$), filtration($\sigma$)) \Else \If {$\tau$ is not Gabriel for $\sigma$} \State filtration($\tau$) = filtration($\sigma$) \EndIf \EndIf \EndFor \EndFor \EndFor \State make\_filtration\_non\_decreasing() \State prune\_above\_filtration() \end{algorithmic} \end{algorithm}

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 :

alpha_complex_doc_420.png
Filtration value propagation example

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).

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).

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().

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.

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 <gudhi/Alpha_complex.h>
// to construct a simplex_tree from alpha complex
#include <gudhi/Simplex_tree.h>
#include <CGAL/Epick_d.h>
#include <iostream>
#include <string>
void usage(int nbArgs, char * const progName) {
std::cerr << "Error: Number of arguments (" << nbArgs << ") is not correct\n";
std::cerr << "Usage: " << progName << " filename.off alpha_square_max_value [ouput_file.txt]\n";
std::cerr << " i.e.: " << progName << " ../../data/points/alphacomplexdoc.off 60.0\n";
exit(-1); // ----- >>
}
int main(int argc, char **argv) {
if ((argc != 3) && (argc != 4)) usage(argc, (argv[0] - 1));
std::string off_file_name {argv[1]};
double alpha_square_max_value {atof(argv[2])};
// ----------------------------------------------------------------------------
// Init of an alpha complex from an OFF file
// ----------------------------------------------------------------------------
using Kernel = CGAL::Epick_d< CGAL::Dynamic_dimension_tag >;
Gudhi::alpha_complex::Alpha_complex<Kernel> alpha_complex_from_file(off_file_name);
std::streambuf* streambufffer;
std::ofstream ouput_file_stream;
if (argc == 4) {
ouput_file_stream.open(std::string(argv[3]));
streambufffer = ouput_file_stream.rdbuf();
} else {
streambufffer = std::cout.rdbuf();
}
if (alpha_complex_from_file.create_complex(simplex, alpha_square_max_value)) {
std::ostream output_stream(streambufffer);
// ----------------------------------------------------------------------------
// Display information about the alpha complex
// ----------------------------------------------------------------------------
output_stream << "Alpha complex is of dimension " << simplex.dimension() <<
" - " << simplex.num_simplices() << " simplices - " <<
simplex.num_vertices() << " vertices." << std::endl;
output_stream << "Iterator on alpha complex simplices in the filtration order, with [filtration value]:" <<
std::endl;
for (auto f_simplex : simplex.filtration_simplex_range()) {
output_stream << " ( ";
for (auto vertex : simplex.simplex_vertex_range(f_simplex)) {
output_stream << vertex << " ";
}
output_stream << ") -> " << "[" << simplex.filtration(f_simplex) << "] ";
output_stream << std::endl;
}
}
ouput_file_stream.close();
return 0;
}

When launching:

$> ./Alpha_complex_example_from_off ../../data/points/alphacomplexdoc.off 32.0

the program output is:

Alpha complex is of dimension 2 - 20 simplices - 7 vertices.
Iterator on alpha complex simplices in the filtration order, with [filtration value]:
( 0 ) -> [0]
( 1 ) -> [0]
( 2 ) -> [0]
( 3 ) -> [0]
( 4 ) -> [0]
( 5 ) -> [0]
( 6 ) -> [0]
( 3 2 ) -> [6.25]
( 5 4 ) -> [7.25]
( 2 0 ) -> [8.5]
( 1 0 ) -> [9.25]
( 3 1 ) -> [10]
( 2 1 ) -> [11.25]
( 3 2 1 ) -> [12.5]
( 2 1 0 ) -> [12.9959]
( 6 5 ) -> [13.25]
( 4 2 ) -> [20]
( 6 4 ) -> [22.7367]
( 6 5 4 ) -> [22.7367]
( 6 3 ) -> [30.25]
GUDHI  Version 2.0.0  - C++ library for Topological Data Analysis (TDA) and Higher Dimensional Geometry Understanding. Generated on Wed Apr 19 2017 22:26:16 for GUDHI by doxygen 1.8.11