summaryrefslogtreecommitdiff
path: root/src/Alpha_complex
diff options
context:
space:
mode:
Diffstat (limited to 'src/Alpha_complex')
-rw-r--r--src/Alpha_complex/include/gudhi/Alpha_complex.h39
-rw-r--r--src/Alpha_complex/test/Alpha_complex_dim3_unit_test.cpp117
-rw-r--r--src/Alpha_complex/test/Alpha_complex_unit_test.cpp97
-rw-r--r--src/Alpha_complex/test/CMakeLists.txt6
-rw-r--r--src/Alpha_complex/utilities/alpha_complex_3d_persistence.cpp2
-rw-r--r--src/Alpha_complex/utilities/alpha_complex_persistence.cpp2
6 files changed, 154 insertions, 109 deletions
diff --git a/src/Alpha_complex/include/gudhi/Alpha_complex.h b/src/Alpha_complex/include/gudhi/Alpha_complex.h
index aec8c1b1..a7372f19 100644
--- a/src/Alpha_complex/include/gudhi/Alpha_complex.h
+++ b/src/Alpha_complex/include/gudhi/Alpha_complex.h
@@ -17,8 +17,7 @@
// to construct Alpha_complex from a OFF file of points
#include <gudhi/Points_off_io.h>
-#include <stdlib.h>
-#include <math.h> // isnan, fmax
+#include <cmath> // isnan, fmax
#include <memory> // for std::unique_ptr
#include <cstddef> // for std::size_t
@@ -45,6 +44,7 @@
#include <utility> // std::pair
#include <stdexcept>
#include <numeric> // for std::iota
+#include <algorithm> // for std::sort
// Make compilation fail - required for external projects - https://github.com/GUDHI/gudhi-devel/issues/10
#if CGAL_VERSION_NR < 1041101000
@@ -101,13 +101,17 @@ template<typename D> struct Is_Epeck_D<CGAL::Epeck_d<D>> { static const bool val
*/
template<class Kernel = CGAL::Epeck_d<CGAL::Dynamic_dimension_tag>, bool Weighted = false>
class Alpha_complex {
+ private:
+ // Vertex_handle internal type (required by triangulation_ and vertices_).
+ using Internal_vertex_handle = std::ptrdiff_t;
+
public:
/** \brief Geometric traits class that provides the geometric types and predicates needed by the triangulations.*/
using Geom_traits = std::conditional_t<Weighted, CGAL::Regular_triangulation_traits_adapter<Kernel>, Kernel>;
// Add an int in TDS to save point index in the structure
using TDS = CGAL::Triangulation_data_structure<typename Geom_traits::Dimension,
- CGAL::Triangulation_vertex<Geom_traits, std::ptrdiff_t>,
+ CGAL::Triangulation_vertex<Geom_traits, Internal_vertex_handle>,
CGAL::Triangulation_full_cell<Geom_traits> >;
/** \brief A (Weighted or not) Delaunay triangulation of a set of points in \f$ \mathbb{R}^D\f$.*/
@@ -132,9 +136,6 @@ class Alpha_complex {
// Vertex_iterator type from CGAL.
using CGAL_vertex_iterator = typename Triangulation::Vertex_iterator;
- // size_type type from CGAL.
- using size_type = typename Triangulation::size_type;
-
// Structure to switch from simplex tree vertex handle to CGAL vertex iterator.
using Vector_vertex_iterator = std::vector< CGAL_vertex_iterator >;
@@ -146,6 +147,10 @@ class Alpha_complex {
std::unique_ptr<Triangulation> triangulation_;
/** \brief Kernel for triangulation_ functions access.*/
A_kernel_d kernel_;
+ /** \brief Vertices to be inserted first by the create_complex method to avoid quadratic complexity.
+ * It isn't just [0, n) if some points have multiplicity (only one copy appears in the complex).
+ */
+ std::vector<Internal_vertex_handle> vertices_;
/** \brief Cache for geometric constructions: circumcenter and squared radius of a simplex.*/
std::vector<Sphere> cache_, old_cache_;
@@ -257,11 +262,11 @@ class Alpha_complex {
std::vector<Point_d> point_cloud(first, last);
// Creates a vector {0, 1, ..., N-1}
- std::vector<std::ptrdiff_t> indices(boost::counting_iterator<std::ptrdiff_t>(0),
- boost::counting_iterator<std::ptrdiff_t>(point_cloud.size()));
+ std::vector<Internal_vertex_handle> indices(boost::counting_iterator<Internal_vertex_handle>(0),
+ boost::counting_iterator<Internal_vertex_handle>(point_cloud.size()));
using Point_property_map = boost::iterator_property_map<typename std::vector<Point_d>::iterator,
- CGAL::Identity_property_map<std::ptrdiff_t>>;
+ CGAL::Identity_property_map<Internal_vertex_handle>>;
using Search_traits_d = CGAL::Spatial_sort_traits_adapter_d<Geom_traits, Point_property_map>;
CGAL::spatial_sort(indices.begin(), indices.end(), Search_traits_d(std::begin(point_cloud)));
@@ -279,6 +284,9 @@ class Alpha_complex {
// structure to retrieve CGAL points from vertex handle - one vertex handle per point.
// Needs to be constructed before as vertex handles arrives in no particular order.
vertex_handle_to_iterator_.resize(point_cloud.size());
+ // List of sorted unique vertices in the triangulation. We take advantage of the existing loop to construct it
+ // Vertices list avoids quadratic complexity with the Simplex_tree. We should not fill it up with Toplex_map e.g.
+ vertices_.reserve(triangulation_->number_of_vertices());
// Loop on triangulation vertices list
for (CGAL_vertex_iterator vit = triangulation_->vertices_begin(); vit != triangulation_->vertices_end(); ++vit) {
if (!triangulation_->is_infinite(*vit)) {
@@ -286,8 +294,10 @@ class Alpha_complex {
std::clog << "Vertex insertion - " << vit->data() << " -> " << vit->point() << std::endl;
#endif // DEBUG_TRACES
vertex_handle_to_iterator_[vit->data()] = vit;
+ vertices_.push_back(vit->data());
}
}
+ std::sort(vertices_.begin(), vertices_.end());
// --------------------------------------------------------------------------------------------
}
}
@@ -384,12 +394,21 @@ class Alpha_complex {
// --------------------------------------------------------------------------------------------
// Simplex_tree construction from loop on triangulation finite full cells list
if (num_vertices() > 0) {
+ std::vector<Vertex_handle> one_vertex(1);
+ for (auto vertex : vertices_) {
+#ifdef DEBUG_TRACES
+ std::clog << "SimplicialComplex insertion " << vertex << std::endl;
+#endif // DEBUG_TRACES
+ one_vertex[0] = vertex;
+ complex.insert_simplex_and_subfaces(one_vertex, std::numeric_limits<double>::quiet_NaN());
+ }
+
for (auto cit = triangulation_->finite_full_cells_begin();
cit != triangulation_->finite_full_cells_end();
++cit) {
Vector_vertex vertexVector;
#ifdef DEBUG_TRACES
- std::clog << "Simplex_tree insertion ";
+ std::clog << "SimplicialComplex insertion ";
#endif // DEBUG_TRACES
for (auto vit = cit->vertices_begin(); vit != cit->vertices_end(); ++vit) {
if (*vit != nullptr) {
diff --git a/src/Alpha_complex/test/Alpha_complex_dim3_unit_test.cpp b/src/Alpha_complex/test/Alpha_complex_dim3_unit_test.cpp
new file mode 100644
index 00000000..e7c261f1
--- /dev/null
+++ b/src/Alpha_complex/test/Alpha_complex_dim3_unit_test.cpp
@@ -0,0 +1,117 @@
+/* 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
+ */
+
+#define BOOST_TEST_DYN_LINK
+#define BOOST_TEST_MODULE "alpha_complex_dim3"
+#include <boost/test/unit_test.hpp>
+#include <boost/mpl/list.hpp>
+
+#include <CGAL/Epick_d.h>
+#include <CGAL/Epeck_d.h>
+
+#include <stdexcept> // std::out_of_range
+#include <string>
+#include <vector>
+
+#include <gudhi/Alpha_complex.h>
+#include <gudhi/Simplex_tree.h>
+
+// Use dynamic_dimension_tag for the user to be able to set dimension
+typedef CGAL::Epeck_d< CGAL::Dynamic_dimension_tag > Exact_kernel_d;
+// Use static dimension_tag for the user not to be able to set dimension
+typedef CGAL::Epeck_d< CGAL::Dimension_tag<3> > Exact_kernel_s;
+// Use dynamic_dimension_tag for the user to be able to set dimension
+typedef CGAL::Epick_d< CGAL::Dynamic_dimension_tag > Inexact_kernel_d;
+// Use static dimension_tag for the user not to be able to set dimension
+typedef CGAL::Epick_d< CGAL::Dimension_tag<3> > Inexact_kernel_s;
+// The triangulation uses the default instantiation of the TriangulationDataStructure template parameter
+
+typedef boost::mpl::list<Exact_kernel_d, Exact_kernel_s, Inexact_kernel_d, Inexact_kernel_s> list_of_kernel_variants;
+
+BOOST_AUTO_TEST_CASE_TEMPLATE(Alpha_complex_from_OFF_file, TestedKernel, list_of_kernel_variants) {
+ // ----------------------------------------------------------------------------
+ //
+ // Init of an alpha-complex from a OFF file
+ //
+ // ----------------------------------------------------------------------------
+ std::string off_file_name("alphacomplexdoc.off");
+ double max_alpha_square_value = 60.0;
+ std::clog << "========== OFF FILE NAME = " << off_file_name << " - alpha²=" <<
+ max_alpha_square_value << "==========" << std::endl;
+
+ Gudhi::alpha_complex::Alpha_complex<TestedKernel> alpha_complex_from_file(off_file_name);
+
+ Gudhi::Simplex_tree<> simplex_tree_60;
+ BOOST_CHECK(alpha_complex_from_file.create_complex(simplex_tree_60, max_alpha_square_value));
+
+ std::clog << "alpha_complex_from_file.num_vertices()=" << alpha_complex_from_file.num_vertices() << std::endl;
+ BOOST_CHECK(alpha_complex_from_file.num_vertices() == 7);
+
+ std::clog << "simplex_tree_60.dimension()=" << simplex_tree_60.dimension() << std::endl;
+ BOOST_CHECK(simplex_tree_60.dimension() == 2);
+
+ std::clog << "simplex_tree_60.num_vertices()=" << simplex_tree_60.num_vertices() << std::endl;
+ BOOST_CHECK(simplex_tree_60.num_vertices() == 7);
+
+ std::clog << "simplex_tree_60.num_simplices()=" << simplex_tree_60.num_simplices() << std::endl;
+ BOOST_CHECK(simplex_tree_60.num_simplices() == 25);
+
+ max_alpha_square_value = 59.0;
+ std::clog << "========== OFF FILE NAME = " << off_file_name << " - alpha²=" <<
+ max_alpha_square_value << "==========" << std::endl;
+
+ Gudhi::Simplex_tree<> simplex_tree_59;
+ BOOST_CHECK(alpha_complex_from_file.create_complex(simplex_tree_59, max_alpha_square_value));
+
+ std::clog << "alpha_complex_from_file.num_vertices()=" << alpha_complex_from_file.num_vertices() << std::endl;
+ BOOST_CHECK(alpha_complex_from_file.num_vertices() == 7);
+
+ std::clog << "simplex_tree_59.dimension()=" << simplex_tree_59.dimension() << std::endl;
+ BOOST_CHECK(simplex_tree_59.dimension() == 2);
+
+ std::clog << "simplex_tree_59.num_vertices()=" << simplex_tree_59.num_vertices() << std::endl;
+ BOOST_CHECK(simplex_tree_59.num_vertices() == 7);
+
+ std::clog << "simplex_tree_59.num_simplices()=" << simplex_tree_59.num_simplices() << std::endl;
+ BOOST_CHECK(simplex_tree_59.num_simplices() == 23);
+}
+
+
+BOOST_AUTO_TEST_CASE_TEMPLATE(Alpha_complex_from_empty_points, TestedKernel, list_of_kernel_variants) {
+ std::clog << "========== Alpha_complex_from_empty_points ==========" << std::endl;
+
+ // ----------------------------------------------------------------------------
+ // Init of an empty list of points
+ // ----------------------------------------------------------------------------
+ std::vector<typename TestedKernel::Point_d> points;
+
+ // ----------------------------------------------------------------------------
+ // Init of an alpha complex from the list of points
+ // ----------------------------------------------------------------------------
+ Gudhi::alpha_complex::Alpha_complex<TestedKernel> alpha_complex_from_points(points);
+
+ std::clog << "alpha_complex_from_points.num_vertices()=" << alpha_complex_from_points.num_vertices() << std::endl;
+ BOOST_CHECK(alpha_complex_from_points.num_vertices() == points.size());
+
+ // Test to the limit
+ BOOST_CHECK_THROW (alpha_complex_from_points.get_point(0), std::out_of_range);
+
+ Gudhi::Simplex_tree<> simplex_tree;
+ BOOST_CHECK(!alpha_complex_from_points.create_complex(simplex_tree));
+
+ std::clog << "simplex_tree.num_simplices()=" << simplex_tree.num_simplices() << std::endl;
+ BOOST_CHECK(simplex_tree.num_simplices() == 0);
+
+ std::clog << "simplex_tree.dimension()=" << simplex_tree.dimension() << std::endl;
+ BOOST_CHECK(simplex_tree.dimension() == -1);
+
+ std::clog << "simplex_tree.num_vertices()=" << simplex_tree.num_vertices() << std::endl;
+ BOOST_CHECK(simplex_tree.num_vertices() == points.size());
+}
diff --git a/src/Alpha_complex/test/Alpha_complex_unit_test.cpp b/src/Alpha_complex/test/Alpha_complex_unit_test.cpp
index f74ad217..b474917f 100644
--- a/src/Alpha_complex/test/Alpha_complex_unit_test.cpp
+++ b/src/Alpha_complex/test/Alpha_complex_unit_test.cpp
@@ -13,81 +13,17 @@
#include <boost/test/unit_test.hpp>
#include <boost/mpl/list.hpp>
-#include <CGAL/Delaunay_triangulation.h>
#include <CGAL/Epick_d.h>
#include <CGAL/Epeck_d.h>
-#include <cmath> // float comparison
-#include <limits>
+#include <stdexcept> // std::out_of_range
#include <string>
#include <vector>
#include <gudhi/Alpha_complex.h>
-// to construct a simplex_tree from Delaunay_triangulation
-#include <gudhi/graph_simplicial_complex.h>
#include <gudhi/Simplex_tree.h>
#include <gudhi/Unitary_tests_utils.h>
-// Use dynamic_dimension_tag for the user to be able to set dimension
-typedef CGAL::Epeck_d< CGAL::Dynamic_dimension_tag > Exact_kernel_d;
-// Use static dimension_tag for the user not to be able to set dimension
-typedef CGAL::Epeck_d< CGAL::Dimension_tag<3> > Exact_kernel_s;
-// Use dynamic_dimension_tag for the user to be able to set dimension
-typedef CGAL::Epick_d< CGAL::Dynamic_dimension_tag > Inexact_kernel_d;
-// Use static dimension_tag for the user not to be able to set dimension
-typedef CGAL::Epick_d< CGAL::Dimension_tag<3> > Inexact_kernel_s;
-// The triangulation uses the default instantiation of the TriangulationDataStructure template parameter
-
-typedef boost::mpl::list<Exact_kernel_d, Exact_kernel_s, Inexact_kernel_d, Inexact_kernel_s> list_of_kernel_variants;
-
-BOOST_AUTO_TEST_CASE_TEMPLATE(Alpha_complex_from_OFF_file, TestedKernel, list_of_kernel_variants) {
- // ----------------------------------------------------------------------------
- //
- // Init of an alpha-complex from a OFF file
- //
- // ----------------------------------------------------------------------------
- std::string off_file_name("alphacomplexdoc.off");
- double max_alpha_square_value = 60.0;
- std::clog << "========== OFF FILE NAME = " << off_file_name << " - alpha²=" <<
- max_alpha_square_value << "==========" << std::endl;
-
- Gudhi::alpha_complex::Alpha_complex<TestedKernel> alpha_complex_from_file(off_file_name);
-
- Gudhi::Simplex_tree<> simplex_tree_60;
- BOOST_CHECK(alpha_complex_from_file.create_complex(simplex_tree_60, max_alpha_square_value));
-
- std::clog << "alpha_complex_from_file.num_vertices()=" << alpha_complex_from_file.num_vertices() << std::endl;
- BOOST_CHECK(alpha_complex_from_file.num_vertices() == 7);
-
- std::clog << "simplex_tree_60.dimension()=" << simplex_tree_60.dimension() << std::endl;
- BOOST_CHECK(simplex_tree_60.dimension() == 2);
-
- std::clog << "simplex_tree_60.num_vertices()=" << simplex_tree_60.num_vertices() << std::endl;
- BOOST_CHECK(simplex_tree_60.num_vertices() == 7);
-
- std::clog << "simplex_tree_60.num_simplices()=" << simplex_tree_60.num_simplices() << std::endl;
- BOOST_CHECK(simplex_tree_60.num_simplices() == 25);
-
- max_alpha_square_value = 59.0;
- std::clog << "========== OFF FILE NAME = " << off_file_name << " - alpha²=" <<
- max_alpha_square_value << "==========" << std::endl;
-
- Gudhi::Simplex_tree<> simplex_tree_59;
- BOOST_CHECK(alpha_complex_from_file.create_complex(simplex_tree_59, max_alpha_square_value));
-
- std::clog << "alpha_complex_from_file.num_vertices()=" << alpha_complex_from_file.num_vertices() << std::endl;
- BOOST_CHECK(alpha_complex_from_file.num_vertices() == 7);
-
- std::clog << "simplex_tree_59.dimension()=" << simplex_tree_59.dimension() << std::endl;
- BOOST_CHECK(simplex_tree_59.dimension() == 2);
-
- std::clog << "simplex_tree_59.num_vertices()=" << simplex_tree_59.num_vertices() << std::endl;
- BOOST_CHECK(simplex_tree_59.num_vertices() == 7);
-
- std::clog << "simplex_tree_59.num_simplices()=" << simplex_tree_59.num_simplices() << std::endl;
- BOOST_CHECK(simplex_tree_59.num_simplices() == 23);
-}
-
// Use static dimension_tag for the user not to be able to set dimension
typedef CGAL::Epeck_d< CGAL::Dimension_tag<4> > Kernel_4;
typedef Kernel_4::Point_d Point_4;
@@ -236,37 +172,6 @@ BOOST_AUTO_TEST_CASE(Alpha_complex_from_points) {
}
-BOOST_AUTO_TEST_CASE_TEMPLATE(Alpha_complex_from_empty_points, TestedKernel, list_of_kernel_variants) {
- std::clog << "========== Alpha_complex_from_empty_points ==========" << std::endl;
-
- // ----------------------------------------------------------------------------
- // Init of an empty list of points
- // ----------------------------------------------------------------------------
- std::vector<typename TestedKernel::Point_d> points;
-
- // ----------------------------------------------------------------------------
- // Init of an alpha complex from the list of points
- // ----------------------------------------------------------------------------
- Gudhi::alpha_complex::Alpha_complex<TestedKernel> alpha_complex_from_points(points);
-
- std::clog << "alpha_complex_from_points.num_vertices()=" << alpha_complex_from_points.num_vertices() << std::endl;
- BOOST_CHECK(alpha_complex_from_points.num_vertices() == points.size());
-
- // Test to the limit
- BOOST_CHECK_THROW (alpha_complex_from_points.get_point(0), std::out_of_range);
-
- Gudhi::Simplex_tree<> simplex_tree;
- BOOST_CHECK(!alpha_complex_from_points.create_complex(simplex_tree));
-
- std::clog << "simplex_tree.num_simplices()=" << simplex_tree.num_simplices() << std::endl;
- BOOST_CHECK(simplex_tree.num_simplices() == 0);
-
- std::clog << "simplex_tree.dimension()=" << simplex_tree.dimension() << std::endl;
- BOOST_CHECK(simplex_tree.dimension() == -1);
-
- std::clog << "simplex_tree.num_vertices()=" << simplex_tree.num_vertices() << std::endl;
- BOOST_CHECK(simplex_tree.num_vertices() == points.size());
-}
using Inexact_kernel_2 = CGAL::Epick_d< CGAL::Dimension_tag<2> >;
using Exact_kernel_2 = CGAL::Epeck_d< CGAL::Dimension_tag<2> >;
diff --git a/src/Alpha_complex/test/CMakeLists.txt b/src/Alpha_complex/test/CMakeLists.txt
index 0595ca92..dd2c235f 100644
--- a/src/Alpha_complex/test/CMakeLists.txt
+++ b/src/Alpha_complex/test/CMakeLists.txt
@@ -8,14 +8,18 @@ if (NOT CGAL_WITH_EIGEN3_VERSION VERSION_LESS 4.11.0)
add_executable ( Alpha_complex_test_unit Alpha_complex_unit_test.cpp )
target_link_libraries(Alpha_complex_test_unit ${CGAL_LIBRARY})
+ add_executable ( Alpha_complex_dim3_test_unit Alpha_complex_dim3_unit_test.cpp )
+ target_link_libraries(Alpha_complex_dim3_test_unit ${CGAL_LIBRARY})
add_executable ( Delaunay_complex_test_unit Delaunay_complex_unit_test.cpp )
target_link_libraries(Delaunay_complex_test_unit ${CGAL_LIBRARY})
if (TBB_FOUND)
target_link_libraries(Alpha_complex_test_unit ${TBB_LIBRARIES})
+ target_link_libraries(Alpha_complex_dim3_test_unit ${TBB_LIBRARIES})
target_link_libraries(Delaunay_complex_test_unit ${TBB_LIBRARIES})
endif()
gudhi_add_boost_test(Alpha_complex_test_unit)
+ gudhi_add_boost_test(Alpha_complex_dim3_test_unit)
gudhi_add_boost_test(Delaunay_complex_test_unit)
endif (NOT CGAL_WITH_EIGEN3_VERSION VERSION_LESS 4.11.0)
@@ -73,4 +77,4 @@ if (NOT CGAL_WITH_EIGEN3_VERSION VERSION_LESS 5.1.0)
endif()
gudhi_add_boost_test(Zero_weighted_alpha_complex_test_unit)
-endif (NOT CGAL_WITH_EIGEN3_VERSION VERSION_LESS 5.1.0) \ No newline at end of file
+endif (NOT CGAL_WITH_EIGEN3_VERSION VERSION_LESS 5.1.0)
diff --git a/src/Alpha_complex/utilities/alpha_complex_3d_persistence.cpp b/src/Alpha_complex/utilities/alpha_complex_3d_persistence.cpp
index 91899040..e65d8c6f 100644
--- a/src/Alpha_complex/utilities/alpha_complex_3d_persistence.cpp
+++ b/src/Alpha_complex/utilities/alpha_complex_3d_persistence.cpp
@@ -263,7 +263,7 @@ void program_options(int argc, char *argv[], std::string &off_file_points, bool
"cuboid-file,c", po::value<std::string>(&cuboid_file),
"Name of file describing the periodic domain. Format is:\n min_hx min_hy min_hz\n max_hx max_hy max_hz")(
"output-file,o", po::value<std::string>(&output_file_diag)->default_value(std::string()),
- "Name of file in which the persistence diagram is written. Default print in std::clog")(
+ "Name of file in which the persistence diagram is written. Default print in standard output")(
"max-alpha-square-value,r",
po::value<Filtration_value>(&alpha_square_max_value)
->default_value(std::numeric_limits<Filtration_value>::infinity()),
diff --git a/src/Alpha_complex/utilities/alpha_complex_persistence.cpp b/src/Alpha_complex/utilities/alpha_complex_persistence.cpp
index e86b34e2..29edbd8e 100644
--- a/src/Alpha_complex/utilities/alpha_complex_persistence.cpp
+++ b/src/Alpha_complex/utilities/alpha_complex_persistence.cpp
@@ -163,7 +163,7 @@ void program_options(int argc, char *argv[], std::string &off_file_points, bool
"weight-file,w", po::value<std::string>(&weight_file)->default_value(std::string()),
"Name of file containing a point weights. Format is one weight per line:\n W1\n ...\n Wn ")(
"output-file,o", po::value<std::string>(&output_file_diag)->default_value(std::string()),
- "Name of file in which the persistence diagram is written. Default print in std::clog")(
+ "Name of file in which the persistence diagram is written. Default print in standard output")(
"max-alpha-square-value,r", po::value<Filtration_value>(&alpha_square_max_value)
->default_value(std::numeric_limits<Filtration_value>::infinity()),
"Maximal alpha square value for the Alpha complex construction.")(