summaryrefslogtreecommitdiff
path: root/src/Alpha_complex/include
diff options
context:
space:
mode:
authorvrouvrea <vrouvrea@636b058d-ea47-450e-bf9e-a15bfbe3eedb>2018-11-27 09:10:33 +0000
committervrouvrea <vrouvrea@636b058d-ea47-450e-bf9e-a15bfbe3eedb>2018-11-27 09:10:33 +0000
commit954d62e67a1126d206ed5c14a3c3b5f3c81235f6 (patch)
tree51461185f148f2f5fab5318ead9cffaac4c84f7e /src/Alpha_complex/include
parent2ead33fd01bbe97e3b88070d169647c68c3db359 (diff)
parentbfd989dae36a22450c1da3dc21cea57bb7d2e96e (diff)
Merge alpha_complex_3d_module_vincent after second review round
git-svn-id: svn+ssh://scm.gforge.inria.fr/svnroot/gudhi/trunk@4014 636b058d-ea47-450e-bf9e-a15bfbe3eedb Former-commit-id: d2f4f38b34ef898ed3b1f267278aaa2ae4cab58f
Diffstat (limited to 'src/Alpha_complex/include')
-rw-r--r--src/Alpha_complex/include/gudhi/Alpha_complex_3d.h227
1 files changed, 107 insertions, 120 deletions
diff --git a/src/Alpha_complex/include/gudhi/Alpha_complex_3d.h b/src/Alpha_complex/include/gudhi/Alpha_complex_3d.h
index 00a47d5c..19445637 100644
--- a/src/Alpha_complex/include/gudhi/Alpha_complex_3d.h
+++ b/src/Alpha_complex/include/gudhi/Alpha_complex_3d.h
@@ -26,10 +26,6 @@
#include <boost/version.hpp>
#include <boost/variant.hpp>
-#if BOOST_VERSION >= 105400
-#include <boost/container/static_vector.hpp>
-#endif
-
#include <gudhi/Debug_utils.h>
#include <gudhi/Alpha_complex_options.h>
@@ -50,6 +46,8 @@
#include <CGAL/iterator.h>
#include <CGAL/version.h>
+#include <boost/container/static_vector.hpp>
+
#include <iostream>
#include <vector>
#include <unordered_map>
@@ -60,7 +58,7 @@
#if CGAL_VERSION_NR < 1041101000
// Make compilation fail - required for external projects - https://gitlab.inria.fr/GUDHI/gudhi-devel/issues/10
-static_assert(false, "Alpha_complex_3d is only available for CGAL >= 4.11");
+# error Alpha_complex_3d is only available for CGAL >= 4.11
#endif
namespace Gudhi {
@@ -74,62 +72,22 @@ thread_local
// Value_from_iterator returns the filtration value from an iterator on alpha shapes values
//
-// FAST SAFE EXACT
-// not weighted and *iterator Specific case due to CGAL CGAL::to_double(iterator->exact())
-// not periodic issue # 3153
-//
-// otherwise *iterator CGAL::to_double(*iterator) CGAL::to_double(iterator->exact())
+// FAST SAFE EXACT
+// CGAL::to_double(*iterator) CGAL::to_double(*iterator) CGAL::to_double(iterator->exact())
-template <complexity Complexity, bool Weighted_or_periodic>
+template <complexity Complexity>
struct Value_from_iterator {
template <typename Iterator>
static double perform(Iterator it) {
- // Default behaviour is to return the value pointed by the given iterator
- return *it;
- }
-};
-
-template <>
-struct Value_from_iterator<complexity::SAFE, true> {
- template <typename Iterator>
- static double perform(Iterator it) {
- // In SAFE mode, we are with Epick with EXACT value set to CGAL::Tag_true.
+ // Default behaviour
return CGAL::to_double(*it);
}
};
template <>
-struct Value_from_iterator<complexity::SAFE, false> {
+struct Value_from_iterator<complexity::EXACT> {
template <typename Iterator>
static double perform(Iterator it) {
- // In SAFE mode, we are with Epeck with EXACT value set to CGAL::Tag_true.
- // Specific case due to CGAL issue https://github.com/CGAL/cgal/issues/3153
- auto approx = it->approx();
- double r;
- if (CGAL::fit_in_double(approx, r)) return r;
-
- // If it's precise enough, then OK.
- if (CGAL::has_smaller_relative_precision(approx, RELATIVE_PRECISION_OF_TO_DOUBLE)) return CGAL::to_double(approx);
-
- it->exact();
- return CGAL::to_double(it->approx());
- }
-};
-
-template <>
-struct Value_from_iterator<complexity::EXACT, true> {
- template <typename Iterator>
- static double perform(Iterator it) {
- // In EXACT mode, we are with Epeck or Epick with EXACT value set to CGAL::Tag_true.
- return CGAL::to_double(it->exact());
- }
-};
-
-template <>
-struct Value_from_iterator<complexity::EXACT, false> {
- template <typename Iterator>
- static double perform(Iterator it) {
- // In EXACT mode, we are with Epeck or Epick with EXACT value set to CGAL::Tag_true.
return CGAL::to_double(it->exact());
}
};
@@ -145,8 +103,8 @@ struct Value_from_iterator<complexity::EXACT, false> {
* Shapes</a> from a range of points (can be read from an OFF file, cf. Points_off_reader).
* Duplicate points are inserted once in the Alpha_complex. This is the reason why the vertices may be not contiguous.
*
- * \tparam Complexity shall be `Gudhi::alpha_complex::complexity`. Default value is
- * `Gudhi::alpha_complex::complexity::FAST`.
+ * \tparam Complexity shall be `Gudhi::alpha_complex::complexity` type. Default value is
+ * `Gudhi::alpha_complex::complexity::SAFE`.
*
* \tparam Weighted Boolean used to set/unset the weighted version of Alpha_complex_3d. Default value is false.
*
@@ -169,22 +127,22 @@ struct Value_from_iterator<complexity::EXACT, false> {
* 3d Delaunay complex.
*
*/
-template <complexity Complexity = complexity::FAST, bool Weighted = false, bool Periodic = false>
+template <complexity Complexity = complexity::SAFE, bool Weighted = false, bool Periodic = false>
class Alpha_complex_3d {
// Epick = Exact_predicates_inexact_constructions_kernel
// Epeck = Exact_predicates_exact_constructions_kernel
- // ExactAlphaComparisonTag = exact version of CGAL Alpha_shape_3 and of its objects (Alpha_shape_vertex_base_3 and
+ // Exact_alpha_comparison_tag = exact version of CGAL Alpha_shape_3 and of its objects (Alpha_shape_vertex_base_3 and
// Alpha_shape_cell_base_3). Not available if weighted or periodic.
- // Can be CGAL::Tag_false or CGAL::Tag_true
+ // Can be CGAL::Tag_false or CGAL::Tag_true. Default is False.
// cf. https://doc.cgal.org/latest/Alpha_shapes_3/classCGAL_1_1Alpha__shape__3.html
//
+ // We could use Epick + CGAL::Tag_true for not weighted nor periodic, but during benchmark, we found a bug
+ // https://github.com/CGAL/cgal/issues/3460
+ // This is the reason we only use Epick + CGAL::Tag_false, or Epeck
//
- // FAST SAFE EXACT
- // not weighted and Epick + CGAL::Tag_false Epick + CGAL::Tag_true Epick + CGAL::Tag_true
- // not periodic
- //
- // otherwise Epick + CGAL::Tag_false Epeck Epeck
- using Predicates = typename std::conditional<((!Weighted && !Periodic) || (Complexity == complexity::FAST)),
+ // FAST SAFE EXACT
+ // Epick + CGAL::Tag_false Epeck Epeck
+ using Predicates = typename std::conditional<(Complexity == complexity::FAST),
CGAL::Exact_predicates_inexact_constructions_kernel,
CGAL::Exact_predicates_exact_constructions_kernel>::type;
@@ -192,14 +150,11 @@ class Alpha_complex_3d {
template <typename Predicates, bool Weighted_version, bool Periodic_version>
struct Kernel_3 {};
- template <typename Predicates>
- struct Kernel_3<Predicates, false, false> {
- using Kernel = Predicates;
- };
- template <typename Predicates>
- struct Kernel_3<Predicates, true, false> {
+ template <typename Predicates, bool Is_periodic>
+ struct Kernel_3<Predicates, Is_periodic, false> {
using Kernel = Predicates;
};
+
template <typename Predicates>
struct Kernel_3<Predicates, false, true> {
using Kernel = CGAL::Periodic_3_Delaunay_triangulation_traits_3<Predicates>;
@@ -211,15 +166,13 @@ class Alpha_complex_3d {
using Kernel = typename Kernel_3<Predicates, Weighted, Periodic>::Kernel;
- using Exact_tag = typename std::conditional<(Complexity == complexity::FAST), CGAL::Tag_false, CGAL::Tag_true>::type;
-
using TdsVb = typename std::conditional<Periodic, CGAL::Periodic_3_triangulation_ds_vertex_base_3<>,
CGAL::Triangulation_ds_vertex_base_3<>>::type;
using Tvb = typename std::conditional<Weighted, CGAL::Regular_triangulation_vertex_base_3<Kernel, TdsVb>,
CGAL::Triangulation_vertex_base_3<Kernel, TdsVb>>::type;
- using Vb = CGAL::Alpha_shape_vertex_base_3<Kernel, Tvb, Exact_tag>;
+ using Vb = CGAL::Alpha_shape_vertex_base_3<Kernel, Tvb>;
using TdsCb = typename std::conditional<Periodic, CGAL::Periodic_3_triangulation_ds_cell_base_3<>,
CGAL::Triangulation_ds_cell_base_3<>>::type;
@@ -227,64 +180,101 @@ class Alpha_complex_3d {
using Tcb = typename std::conditional<Weighted, CGAL::Regular_triangulation_cell_base_3<Kernel, TdsCb>,
CGAL::Triangulation_cell_base_3<Kernel, TdsCb>>::type;
- using Cb = CGAL::Alpha_shape_cell_base_3<Kernel, Tcb, Exact_tag>;
+ using Cb = CGAL::Alpha_shape_cell_base_3<Kernel, Tcb>;
using Tds = CGAL::Triangulation_data_structure_3<Vb, Cb>;
// The other way to do a conditional type. Here there 4 possibilities, cannot use std::conditional
template <typename Kernel, typename Tds, bool Weighted_version, bool Periodic_version>
- struct Triangulation {};
+ struct Triangulation_3 {};
template <typename Kernel, typename Tds>
- struct Triangulation<Kernel, Tds, false, false> {
- using Triangulation_3 = CGAL::Delaunay_triangulation_3<Kernel, Tds>;
+ struct Triangulation_3<Kernel, Tds, false, false> {
+ using Dt = CGAL::Delaunay_triangulation_3<Kernel, Tds>;
+ using Weighted_point_3 = void;
};
template <typename Kernel, typename Tds>
- struct Triangulation<Kernel, Tds, true, false> {
- using Triangulation_3 = CGAL::Regular_triangulation_3<Kernel, Tds>;
+ struct Triangulation_3<Kernel, Tds, true, false> {
+ using Dt = CGAL::Regular_triangulation_3<Kernel, Tds>;
+ using Weighted_point_3 = typename Dt::Weighted_point;
};
template <typename Kernel, typename Tds>
- struct Triangulation<Kernel, Tds, false, true> {
- using Triangulation_3 = CGAL::Periodic_3_Delaunay_triangulation_3<Kernel, Tds>;
+ struct Triangulation_3<Kernel, Tds, false, true> {
+ using Dt = CGAL::Periodic_3_Delaunay_triangulation_3<Kernel, Tds>;
+ using Weighted_point_3 = void;
};
template <typename Kernel, typename Tds>
- struct Triangulation<Kernel, Tds, true, true> {
- using Triangulation_3 = CGAL::Periodic_3_regular_triangulation_3<Kernel, Tds>;
+ struct Triangulation_3<Kernel, Tds, true, true> {
+ using Dt = CGAL::Periodic_3_regular_triangulation_3<Kernel, Tds>;
+ using Weighted_point_3 = typename Dt::Weighted_point;
};
- public:
- using Triangulation_3 = typename Triangulation<Kernel, Tds, Weighted, Periodic>::Triangulation_3;
-
- using Alpha_shape_3 = CGAL::Alpha_shape_3<Triangulation_3, Exact_tag>;
+ /** \brief Is either Delaunay_triangulation_3 (Weighted = false and Periodic = false),
+ * Regular_triangulation_3 (Weighted = true and Periodic = false),
+ * Periodic_3_Delaunay_triangulation_3 (Weighted = false and Periodic = true)
+ * or Periodic_3_regular_triangulation_3 (Weighted = true and Periodic = true).
+ *
+ * This type is required by `Gudhi::alpha_complex::Alpha_complex_3d::Alpha_shape_3`.
+ * */
+ using Dt = typename Triangulation_3<Kernel, Tds, Weighted, Periodic>::Dt;
+public:
+ /** \brief The <a href="https://doc.cgal.org/latest/Alpha_shapes_3/classCGAL_1_1Alpha__shape__3.html">CGAL 3D Alpha
+ * Shapes</a> type.
+ *
+ * The `Gudhi::alpha_complex::Alpha_complex_3d` is a wrapper on top of this class to ease the standard, weighted
+ * and/or periodic build of the Alpha complex 3d.*/
+ using Alpha_shape_3 = CGAL::Alpha_shape_3<Dt>;
+
+ /** \brief The alpha values type.
+ * Must be compatible with double. */
+ using FT = typename Alpha_shape_3::FT;
+
+ /** \brief Gives public access to the Point_3 type. Here is a Point_3 constructor example:
+\code{.cpp}
+using Alpha_complex_3d = Gudhi::alpha_complex::Alpha_complex_3d<Gudhi::alpha_complex::complexity::SAFE, false, false>;
+
+// x0 = 1., y0 = -1.1, z0 = -1..
+Alpha_complex_3d::Point_3 p0(1., -1.1, -1.);
+\endcode
+ * */
using Point_3 = typename Kernel::Point_3;
- private:
- using Alpha_value_type = typename Alpha_shape_3::FT;
+ /** \brief Gives public access to the Weighted_point_3 type. A Weighted point can be constructed as follows:
+\code{.cpp}
+using Weighted_alpha_complex_3d =
+ Gudhi::alpha_complex::Alpha_complex_3d<Gudhi::alpha_complex::complexity::SAFE, true, false>;
+
+// x0 = 1., y0 = -1.1, z0 = -1., weight = 4.
+Weighted_alpha_complex_3d::Weighted_point_3 wp0(Weighted_alpha_complex_3d::Point_3(1., -1.1, -1.), 4.);
+\endcode
+ *
+ * Note: This type is defined to void if Alpha complex is not weighted.
+ *
+ * */
+ using Weighted_point_3 = typename Triangulation_3<Kernel, Tds, Weighted, Periodic>::Weighted_point_3;
+
+private:
using Dispatch =
- CGAL::Dispatch_output_iterator<CGAL::cpp11::tuple<CGAL::Object, Alpha_value_type>,
+ CGAL::Dispatch_output_iterator<CGAL::cpp11::tuple<CGAL::Object, FT>,
CGAL::cpp11::tuple<std::back_insert_iterator<std::vector<CGAL::Object>>,
- std::back_insert_iterator<std::vector<Alpha_value_type>>>>;
+ std::back_insert_iterator<std::vector<FT>>>>;
using Cell_handle = typename Alpha_shape_3::Cell_handle;
using Facet = typename Alpha_shape_3::Facet;
using Edge = typename Alpha_shape_3::Edge;
using Alpha_vertex_handle = typename Alpha_shape_3::Vertex_handle;
-#if BOOST_VERSION >= 105400
using Vertex_list = boost::container::static_vector<Alpha_vertex_handle, 4>;
-#else
- using Vertex_list = std::vector<Alpha_vertex_handle>;
-#endif
public:
/** \brief Alpha_complex constructor from a list of points.
*
* @param[in] points Range of points to triangulate. Points must be in `Alpha_complex_3d::Point_3` or
- * `Alpha_complex_3d::Triangulation_3::Weighted_point`.
+ * `Alpha_complex_3d::Weighted_point_3`.
*
* @pre Available if Alpha_complex_3d is not Periodic.
*
* The type InputPointRange must be a range for which std::begin and std::end return input iterators on a
- * `Alpha_complex_3d::Point_3` or a `Alpha_complex_3d::Triangulation_3::Weighted_point`.
+ * `Alpha_complex_3d::Point_3` or a `Alpha_complex_3d::Weighted_point_3`.
*/
template <typename InputPointRange>
Alpha_complex_3d(const InputPointRange& points) {
@@ -298,15 +288,15 @@ class Alpha_complex_3d {
*
* @exception std::invalid_argument In debug mode, if points and weights do not have the same size.
*
- * @param[in] points Range of points to triangulate. Points must be in `Alpha_complex_3d::Point_3`
- * @param[in] weights Range of weights on points. Weights shall be in `Alpha_complex_3d::Alpha_shape_3::FT`
+ * @param[in] points Range of points to triangulate. Points must be in `Alpha_complex_3d::Point_3`.
+ * @param[in] weights Range of weights on points. Weights shall be in double.
*
* @pre Available if Alpha_complex_3d is Weighted and not Periodic.
*
* The type InputPointRange must be a range for which std::begin and
* std::end return input iterators on a `Alpha_complex_3d::Point_3`.
* The type WeightRange must be a range for which std::begin and
- * std::end return an input iterator on a `Alpha_complex_3d::Alpha_shape_3::FT`.
+ * std::end return an input iterator on a double.
*/
template <typename InputPointRange, typename WeightRange>
Alpha_complex_3d(const InputPointRange& points, WeightRange weights) {
@@ -315,7 +305,6 @@ class Alpha_complex_3d {
GUDHI_CHECK((weights.size() == points.size()),
std::invalid_argument("Points number in range different from weights range number"));
- using Weighted_point_3 = typename Triangulation_3::Weighted_point;
std::vector<Weighted_point_3> weighted_points_3;
std::size_t index = 0;
@@ -334,7 +323,7 @@ class Alpha_complex_3d {
* @exception std::invalid_argument In debug mode, if the size of the cuboid in every directions is not the same.
*
* @param[in] points Range of points to triangulate. Points must be in `Alpha_complex_3d::Point_3` or
- * `Alpha_complex_3d::Triangulation_3::Weighted_point`.
+ * `Alpha_complex_3d::Weighted_point_3`.
* @param[in] x_min Iso-oriented cuboid x_min.
* @param[in] y_min Iso-oriented cuboid y_min.
* @param[in] z_min Iso-oriented cuboid z_min.
@@ -345,14 +334,14 @@ class Alpha_complex_3d {
* @pre Available if Alpha_complex_3d is Periodic.
*
* The type InputPointRange must be a range for which std::begin and std::end return input iterators on a
- * `Alpha_complex_3d::Point_3` or a `Alpha_complex_3d::Triangulation_3::Weighted_point`.
+ * `Alpha_complex_3d::Point_3` or a `Alpha_complex_3d::Weighted_point_3`.
*
* @note In weighted version, please check weights are greater than zero, and lower than 1/64*cuboid length
* squared.
*/
template <typename InputPointRange>
- Alpha_complex_3d(const InputPointRange& points, Alpha_value_type x_min, Alpha_value_type y_min,
- Alpha_value_type z_min, Alpha_value_type x_max, Alpha_value_type y_max, Alpha_value_type z_max) {
+ Alpha_complex_3d(const InputPointRange& points, FT x_min, FT y_min,
+ FT z_min, FT x_max, FT y_max, FT z_max) {
static_assert(Periodic, "This constructor is not available for non-periodic versions of Alpha_complex_3d");
// Checking if the cuboid is the same in x,y and z direction. If not, CGAL will not process it.
GUDHI_CHECK(
@@ -360,7 +349,7 @@ class Alpha_complex_3d {
std::invalid_argument("The size of the cuboid in every directions is not the same."));
// Define the periodic cube
- Triangulation_3 pdt(typename Kernel::Iso_cuboid_3(x_min, y_min, z_min, x_max, y_max, z_max));
+ Dt pdt(typename Kernel::Iso_cuboid_3(x_min, y_min, z_min, x_max, y_max, z_max));
// Heuristic for inserting large point sets (if pts is reasonably large)
pdt.insert(std::begin(points), std::end(points), true);
// As pdt won't be modified anymore switch to 1-sheeted cover if possible
@@ -381,8 +370,8 @@ class Alpha_complex_3d {
* @exception std::invalid_argument In debug mode, if a weight is negative, zero, or greater than 1/64*cuboid length
* squared.
*
- * @param[in] points Range of points to triangulate. Points must be in `Alpha_complex_3d::Point_3`
- * @param[in] weights Range of weights on points. Weights shall be in `Alpha_complex_3d::Alpha_shape_3::FT`
+ * @param[in] points Range of points to triangulate. Points must be in `Alpha_complex_3d::Point_3`.
+ * @param[in] weights Range of weights on points. Weights shall be in double.
* @param[in] x_min Iso-oriented cuboid x_min.
* @param[in] y_min Iso-oriented cuboid y_min.
* @param[in] z_min Iso-oriented cuboid z_min.
@@ -395,12 +384,12 @@ class Alpha_complex_3d {
* The type InputPointRange must be a range for which std::begin and
* std::end return input iterators on a `Alpha_complex_3d::Point_3`.
* The type WeightRange must be a range for which std::begin and
- * std::end return an input iterator on a `Alpha_complex_3d::Alpha_shape_3::FT`.
- * The type of x_min, y_min, z_min, x_max, y_max and z_max is `Alpha_complex_3d::Alpha_shape_3::FT`.
+ * std::end return an input iterator on a double.
+ * The type of x_min, y_min, z_min, x_max, y_max and z_max must be a double.
*/
template <typename InputPointRange, typename WeightRange>
- Alpha_complex_3d(const InputPointRange& points, WeightRange weights, Alpha_value_type x_min, Alpha_value_type y_min,
- Alpha_value_type z_min, Alpha_value_type x_max, Alpha_value_type y_max, Alpha_value_type z_max) {
+ Alpha_complex_3d(const InputPointRange& points, WeightRange weights, FT x_min, FT y_min,
+ FT z_min, FT x_max, FT y_max, FT z_max) {
static_assert(Weighted, "This constructor is not available for non-weighted versions of Alpha_complex_3d");
static_assert(Periodic, "This constructor is not available for non-periodic versions of Alpha_complex_3d");
GUDHI_CHECK((weights.size() == points.size()),
@@ -410,7 +399,6 @@ class Alpha_complex_3d {
(x_max - x_min == y_max - y_min) && (x_max - x_min == z_max - z_min) && (z_max - z_min == y_max - y_min),
std::invalid_argument("The size of the cuboid in every directions is not the same."));
- using Weighted_point_3 = typename Triangulation_3::Weighted_point;
std::vector<Weighted_point_3> weighted_points_3;
std::size_t index = 0;
@@ -418,7 +406,7 @@ class Alpha_complex_3d {
#ifdef GUDHI_DEBUG
// Defined in GUDHI_DEBUG to avoid unused variable warning for GUDHI_CHECK
- Alpha_value_type maximal_possible_weight = 0.015625 * (x_max - x_min) * (x_max - x_min);
+ FT maximal_possible_weight = 0.015625 * (x_max - x_min) * (x_max - x_min);
#endif
while ((index < weights.size()) && (index < points.size())) {
@@ -431,7 +419,7 @@ class Alpha_complex_3d {
}
// Define the periodic cube
- Triangulation_3 pdt(typename Kernel::Iso_cuboid_3(x_min, y_min, z_min, x_max, y_max, z_max));
+ Dt pdt(typename Kernel::Iso_cuboid_3(x_min, y_min, z_min, x_max, y_max, z_max));
// Heuristic for inserting large point sets (if pts is reasonably large)
pdt.insert(std::begin(weighted_points_3), std::end(weighted_points_3), true);
// As pdt won't be modified anymore switch to 1-sheeted cover if possible
@@ -451,13 +439,12 @@ class Alpha_complex_3d {
* \tparam SimplicialComplexForAlpha3d must meet `SimplicialComplexForAlpha3d` concept.
*
* @param[in] complex SimplicialComplexForAlpha3d to be created.
- * @param[in] max_alpha_square maximum for alpha square value. Default value is +\f$\infty\f$.
+ * @param[in] max_alpha_square maximum for alpha square value. Default value is +\f$\infty\f$, and there is very
+ * little point using anything else since it does not save time.
*
* @return true if creation succeeds, false otherwise.
*
- * @pre The simplicial complex must be empty (no vertices)
- *
- * Initialization can be launched once.
+ * @pre The simplicial complex must be empty (no vertices).
*
*/
template <typename SimplicialComplexForAlpha3d,
@@ -481,9 +468,9 @@ class Alpha_complex_3d {
std::size_t count_cells = 0;
#endif // DEBUG_TRACES
std::vector<CGAL::Object> objects;
- std::vector<Alpha_value_type> alpha_values;
+ std::vector<FT> alpha_values;
- Dispatch dispatcher = CGAL::dispatch_output<CGAL::Object, Alpha_value_type>(std::back_inserter(objects),
+ Dispatch dispatcher = CGAL::dispatch_output<CGAL::Object, FT>(std::back_inserter(objects),
std::back_inserter(alpha_values));
alpha_shape_3_ptr_->filtration_with_alpha_values(dispatcher);
@@ -492,7 +479,7 @@ class Alpha_complex_3d {
#endif // DEBUG_TRACES
Alpha_shape_simplex_tree_map map_cgal_simplex_tree;
- using Alpha_value_iterator = typename std::vector<Alpha_value_type>::const_iterator;
+ using Alpha_value_iterator = typename std::vector<FT>::const_iterator;
Alpha_value_iterator alpha_value_iterator = alpha_values.begin();
for (auto object_iterator : objects) {
Vertex_list vertex_list;
@@ -559,7 +546,7 @@ class Alpha_complex_3d {
}
}
// Construction of the simplex_tree
- Filtration_value filtr = Value_from_iterator<Complexity, (Weighted || Periodic)>::perform(alpha_value_iterator);
+ Filtration_value filtr = Value_from_iterator<Complexity>::perform(alpha_value_iterator);
#ifdef DEBUG_TRACES
std::cout << "filtration = " << filtr << std::endl;