diff --git a/Alpha_shapes_2/doc/Alpha_shapes_2/Alpha_shapes_2.txt b/Alpha_shapes_2/doc/Alpha_shapes_2/Alpha_shapes_2.txt
index ba0a7c8e3db..abf2ca27b5b 100644
--- a/Alpha_shapes_2/doc/Alpha_shapes_2/Alpha_shapes_2.txt
+++ b/Alpha_shapes_2/doc/Alpha_shapes_2/Alpha_shapes_2.txt
@@ -108,7 +108,7 @@ such that the \f$ \alpha\f$-shape satisfies the following two properties,
or at least the second one if there is no such \f$ \alpha\f$ that both
are satisfied:
-- The number of components equals a number of your choice, and
+
- the number of components equals a number of your choice, and
- all data points are either on the boundary or in the interior of the regularized version of
the \f$ \alpha\f$-shape (no singular edges).
diff --git a/Alpha_shapes_2/doc/Alpha_shapes_2/CGAL/Alpha_shape_2.h b/Alpha_shapes_2/doc/Alpha_shapes_2/CGAL/Alpha_shape_2.h
index 8aabfbd533e..cd7d82058ed 100644
--- a/Alpha_shapes_2/doc/Alpha_shapes_2/CGAL/Alpha_shape_2.h
+++ b/Alpha_shapes_2/doc/Alpha_shapes_2/CGAL/Alpha_shape_2.h
@@ -371,8 +371,11 @@ size_type number_of_solid_components(const FT& alpha = get_alpha()) const;
Returns an iterator pointing to the first element with \f$ \alpha\f$-value
such that the alpha shape satisfies the following two properties:
-- `nb_components` equals the number of solid components and
-- all data points are either on the boundary or in the interior of the regularized version of the alpha shape.
+- All data points are either on the boundary or in the interior
+of the regularized version of the alpha shape.
+
+- The number of solid component of the alpha shape is equal to or
+smaller than `nb_components`.
If no such value is found, the iterator points to the first element with
\f$ \alpha\f$-value such that the alpha shape satisfies the second property.
diff --git a/Alpha_shapes_3/doc/Alpha_shapes_3/Alpha_shapes_3.txt b/Alpha_shapes_3/doc/Alpha_shapes_3/Alpha_shapes_3.txt
index 40edd85b84e..d440d65382d 100644
--- a/Alpha_shapes_3/doc/Alpha_shapes_3/Alpha_shapes_3.txt
+++ b/Alpha_shapes_3/doc/Alpha_shapes_3/Alpha_shapes_3.txt
@@ -145,9 +145,9 @@ order of apparition in the alpha complex when alpha increases.
Finally, it provides a function to determine
the smallest value \f$ \alpha\f$
-such that the alpha shape satisfies the following two properties
+such that the alpha shape satisfies the following two properties:
-- all data points are either on the boundary or in the interior
+
- All data points are either on the boundary or in the interior
of the regularized version of the alpha shape (no singular faces).
- The number of components is equal or less than a given number.
@@ -198,10 +198,13 @@ We currently do not specify concepts for the underlying triangulation
type. Models that work for a family of alpha-shapes are the instantiations
of the classes `Delaunay_triangulation_3` and
`Periodic_3_Delaunay_triangulation_3` (see
-example \ref AlphaShape_3DExampleForPeriodicAlphaShapes). A model that works for a fixed alpha-shape are the instantiations
+example \ref AlphaShape_3DExampleForPeriodicAlphaShapes).
+A model that works for a fixed alpha-shape are the instantiations
of the class `Delaunay_triangulation_3`.
-A model that works for a weighted alpha-shape is
-the class `Regular_triangulation_3`. The triangulation needs a geometric traits class
+Models that work for a weighted alpha-shape are the instantiations
+of the classes `Regular_triangulation_3` and
+`Periodic_3_regular_triangulation_3`.
+The triangulation needs a geometric traits class
and a triangulation data structure as template parameters.
\subsection AlphaShape3D_ConceptAndModelsAlphaShapes Alpha Shapes
@@ -234,19 +237,20 @@ and `Fixed_alpha_shape_cell_base_3`, respectively.
\subsection AlphaShape3D_ConceptAndModelsTDS Triangulation data structure
-Additionally requirements are put when using `Regular_triangulation_3` or
-`Periodic_3_Delaunay_triangulation_3` as underlying triangulations:
+Additionally requirements are put when using weighted or
+periodic triangulations as underlying triangulation:
-- When using `Regular_triangulation_3` as underlying triangulation, the vertex
+
- When using a weighted triangulation (`Regular_triangulation_3` or
+`Periodic_3_regular_triangulation_3`), the vertex
and cell classes must be models to both `AlphaShapeVertex_3` and
`RegularTriangulationVertexBase_3`, as well as
`AlphaShapeCell_3` and `RegularTriangulationCellBase_3` respectively
-(see example \ref AlphaShape_3DExampleforWeightedAlphaShapes).
-
- When using `Periodic_3_Delaunay_triangulation_3` as underlying
-triangulation the vertex and cell classes need to be models to both
-`AlphaShapeVertex_3` and `Periodic_3TriangulationDSVertexBase_3`, as well as
-`AlphaShapeCell_3` and `Periodic_3TriangulationDSCellBase_3`
-(see example \ref AlphaShape_3DExampleForPeriodicAlphaShapes).
+(see example: \ref AlphaShape_3DExampleforWeightedAlphaShapes).
+
- When using a periodic triangulation (`Periodic_3_Delaunay_triangulation_3`
+or `Periodic_3_regular_triangulation_3`), the vertex and cell classes must
+be models to both `AlphaShapeVertex_3` and `Periodic_3TriangulationDSVertexBase_3`,
+as well as `AlphaShapeCell_3` and `Periodic_3TriangulationDSCellBase_3`
+(see example: \ref AlphaShape_3DExampleForPeriodicAlphaShapes).
\section Alpha_shapes_3AlphaShape3OrFixedAlphaShape3 Alpha_shape_3 vs. Fixed_alpha_shape_3
@@ -320,14 +324,15 @@ them with a traits with inexact constructions, the tag
\subsection AlphaShape_3DExampleForPeriodicAlphaShapes Example for Periodic Alpha Shapes
-The following example shows how to use the periodic Delaunay
+The following example shows how to use a periodic Delaunay
triangulation (Chapter \ref Chapter_3D_Periodic_Triangulations "3D Periodic Triangulations") as underlying
-triangulation for the alpha shape computation.
+triangulation for the alpha shape computation. Usage of a weighted Delaunay periodic
+triangulation is presented in the example: \ref Alpha_shapes_3/ex_weighted_periodic_alpha_shapes_3.cpp "ex_weighted_periodic_alpha_shapes_3.cpp".
In order to define the original domain and to benefit from the
-built-in heuristic optimizations of the periodic Delaunay
-triangulation computation, it is recommended to first construct the
-triangulation and then construct the alpha shape from it. The alpha
+built-in heuristic optimizations of the periodic triangulation computation,
+it is recommended to first construct the triangulation and
+then construct the alpha shape from it. The alpha
shape constructor that takes a point range can be used as well but in
this case the original domain cannot be specified and the default unit
cube will be chosen and no optimizations will be used.
diff --git a/Alpha_shapes_3/doc/Alpha_shapes_3/CGAL/Alpha_shape_3.h b/Alpha_shapes_3/doc/Alpha_shapes_3/CGAL/Alpha_shape_3.h
index 79a508644c0..50734a63de7 100644
--- a/Alpha_shapes_3/doc/Alpha_shapes_3/CGAL/Alpha_shape_3.h
+++ b/Alpha_shapes_3/doc/Alpha_shapes_3/CGAL/Alpha_shape_3.h
@@ -17,8 +17,9 @@ The modifying functions `insert` and `remove` will overwrite
the one inherited from the underlying triangulation class `Dt`.
At the moment, only the static version is implemented.
-\tparam Dt must be either `Delaunay_triangulation_3`, `Regular_triangulation_3`
-or `Periodic_3_triangulation_3`. Note that `Dt::Geom_traits`, `Dt::Vertex`, and `Dt::Face`
+\tparam Dt must be either `Delaunay_triangulation_3`, `Regular_triangulation_3`,
+`Periodic_3_Delaunay_triangulation_3` or `Periodic_3_regular_triangulation_3`.
+Note that `Dt::Geom_traits`, `Dt::Vertex`, and `Dt::Face`
must be model the concepts `AlphaShapeTraits_3`,
`AlphaShapeVertex_3` and `AlphaShapeCell_3`, respectively.
@@ -438,13 +439,13 @@ size_type number_of_solid_components(const FT& alpha = get_alpha()) const;
/*!
Returns an iterator pointing to smallest \f$ \alpha\f$ value
-such that the alpha shape satisfies the following two properties:
+such that the alpha shape satisfies the following two properties:
-all data points are either on the boundary or in the interior
-of the regularized version of the alpha shape.
+- All data points are either on the boundary or in the interior
+of the regularized version of the alpha shape.
-The number of solid component of the alpha shape is equal to or
-smaller than `nb_components`.
+- The number of solid component of the alpha shape is equal to or
+smaller than `nb_components`.
*/
Alpha_iterator find_optimal_alpha(size_type nb_components) const;
diff --git a/Alpha_shapes_3/doc/Alpha_shapes_3/CGAL/Fixed_alpha_shape_3.h b/Alpha_shapes_3/doc/Alpha_shapes_3/CGAL/Fixed_alpha_shape_3.h
index 8d60d020c27..71b7394b4ef 100644
--- a/Alpha_shapes_3/doc/Alpha_shapes_3/CGAL/Fixed_alpha_shape_3.h
+++ b/Alpha_shapes_3/doc/Alpha_shapes_3/CGAL/Fixed_alpha_shape_3.h
@@ -12,8 +12,9 @@ represents connectivity and order among its faces. Each
\f$ k\f$-dimensional face of the `Dt` is associated with
a classification that specifies its status in the alpha complex, alpha being fixed.
-\tparam Dt must be either `Delaunay_triangulation_3`, `Regular_triangulation_3`
-or `Periodic_3_triangulation_3`. Note that `Dt::Geom_traits`, `Dt::Vertex`, and `Dt::Face`
+\tparam Dt must be either `Delaunay_triangulation_3`, `Regular_triangulation_3`,
+`Periodic_3_Delaunay_triangulation_3` or `Periodic_3_regular_triangulation_3`.
+Note that `Dt::Geom_traits`, `Dt::Vertex`, and `Dt::Face`
must be model the concepts `AlphaShapeTraits_3`,
`AlphaShapeVertex_3` and `AlphaShapeFace_3`, respectively.
diff --git a/Alpha_shapes_3/doc/Alpha_shapes_3/examples.txt b/Alpha_shapes_3/doc/Alpha_shapes_3/examples.txt
index ab726f84c6d..1ab05e5c8a9 100644
--- a/Alpha_shapes_3/doc/Alpha_shapes_3/examples.txt
+++ b/Alpha_shapes_3/doc/Alpha_shapes_3/examples.txt
@@ -5,5 +5,6 @@
\example Alpha_shapes_3/ex_fixed_weighted_alpha_shapes_3.cpp
\example Alpha_shapes_3/ex_periodic_alpha_shapes_3.cpp
\example Alpha_shapes_3/ex_weighted_alpha_shapes_3.cpp
+\example Alpha_shapes_3/ex_weighted_periodic_alpha_shapes_3.cpp
\example Alpha_shapes_2/ex_alpha_projection_traits.cpp
*/
diff --git a/Alpha_shapes_3/examples/Alpha_shapes_3/ex_weighted_periodic_alpha_shapes_3.cpp b/Alpha_shapes_3/examples/Alpha_shapes_3/ex_weighted_periodic_alpha_shapes_3.cpp
new file mode 100644
index 00000000000..a971aede00c
--- /dev/null
+++ b/Alpha_shapes_3/examples/Alpha_shapes_3/ex_weighted_periodic_alpha_shapes_3.cpp
@@ -0,0 +1,77 @@
+#include
+
+#include
+#include
+#include
+
+#include
+
+#include
+#include
+
+// Traits
+typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
+typedef CGAL::Periodic_3_regular_triangulation_traits_3 PK;
+
+// Vertex type
+typedef CGAL::Periodic_3_triangulation_ds_vertex_base_3<> DsVb;
+typedef CGAL::Regular_triangulation_vertex_base_3 Vb;
+typedef CGAL::Alpha_shape_vertex_base_3 AsVb;
+// Cell type
+typedef CGAL::Periodic_3_triangulation_ds_cell_base_3<> DsCb;
+typedef CGAL::Regular_triangulation_cell_base_3 Cb;
+typedef CGAL::Alpha_shape_cell_base_3 AsCb;
+
+typedef CGAL::Triangulation_data_structure_3 Tds;
+typedef CGAL::Periodic_3_regular_triangulation_3 P3RT3;
+typedef CGAL::Alpha_shape_3 Alpha_shape_3;
+
+typedef P3RT3::Bare_point Bare_point;
+typedef P3RT3::Weighted_point Weighted_point;
+
+int main()
+{
+ CGAL::Random random(8);
+ std::list pts;
+
+ // read input
+ std::ifstream is("./data/bunny_1000");
+ int n;
+ is >> n;
+ std::cout << "Reading " << n << " points " << std::endl;
+ Bare_point bp;
+ for( ; n>0 ; n--) {
+ is >> bp;
+ Weighted_point p(bp, 0.0001 * random.get_double(0., 0.015625)); // arbitrary weights
+ pts.push_back(p);
+ }
+
+ // Define the periodic cube
+ P3RT3 prt(PK::Iso_cuboid_3(-0.1,0.,-0.1, 0.1,0.2,0.1));
+
+ // Heuristic for inserting large point sets (if pts is reasonably large)
+ prt.insert(pts.begin(), pts.end(), true);
+
+ // As prt won't be modified anymore switch to 1-sheeted cover if possible
+ if(prt.is_triangulation_in_1_sheet())
+ {
+ std::cout << "Switching to 1-sheeted covering" << std::endl;
+ prt.convert_to_1_sheeted_covering();
+ }
+
+ std::cout << "Periodic Regular computed." << std::endl;
+
+ // compute alpha shape
+ Alpha_shape_3 as(prt);
+ std::cout << "Alpha shape computed in REGULARIZED mode by default." << std::endl;
+
+ // find optimal alpha values
+ Alpha_shape_3::NT alpha_solid = as.find_alpha_solid();
+ Alpha_shape_3::Alpha_iterator opt = as.find_optimal_alpha(1);
+ std::cout << "Smallest alpha value to get a solid through data points is " << alpha_solid << std::endl;
+ std::cout << "Optimal alpha value to get one connected component is " << *opt << std::endl;
+ as.set_alpha(*opt);
+ assert(as.number_of_solid_components() == 1);
+
+ return 0;
+}
diff --git a/Alpha_shapes_3/include/CGAL/Alpha_shape_3.h b/Alpha_shapes_3/include/CGAL/Alpha_shape_3.h
index 24ecf925f53..61b31cb82fc 100644
--- a/Alpha_shapes_3/include/CGAL/Alpha_shape_3.h
+++ b/Alpha_shapes_3/include/CGAL/Alpha_shape_3.h
@@ -148,12 +148,6 @@ public:
using Dt::finite_vertices_end;
using Dt::finite_cells_begin;
using Dt::finite_cells_end;
- using Dt::VERTEX;
- using Dt::EDGE;
- using Dt::FACET;
- using Dt::CELL;
- using Dt::OUTSIDE_CONVEX_HULL;
- using Dt::OUTSIDE_AFFINE_HULL;
using Dt::vertex_triple_index;
using Dt::is_infinite;
using Dt::is_Gabriel;
@@ -161,15 +155,22 @@ public:
using Dt::incident_vertices;
using Dt::incident_facets;
using Dt::locate;
+ using Dt::point;
+
+ using Dt::VERTEX;
+ using Dt::EDGE;
+ using Dt::FACET;
+ using Dt::CELL;
+ using Dt::OUTSIDE_CONVEX_HULL;
+ using Dt::OUTSIDE_AFFINE_HULL;
+
+ enum Classification_type {EXTERIOR,
+ SINGULAR,
+ REGULAR,
+ INTERIOR};
- enum Classification_type {EXTERIOR,
- SINGULAR,
- REGULAR,
- INTERIOR};
-
enum Mode {GENERAL, REGULARIZED};
-
typedef CGAL::Alpha_status< NT > Alpha_status;
typedef Compact_container Alpha_status_container;
typedef typename Alpha_status_container::const_iterator
@@ -744,44 +745,38 @@ public:
// starting point for searching
// takes O(#alpha_shape) time
-
//------------------- GEOMETRIC PRIMITIVES ----------------------------
private:
NT squared_radius(const Cell_handle& s) const
- {
- return Compute_squared_radius_3()(*this)(
- this->point(s,0), this->point(s,1),
- this->point(s,2), this->point(s,3));
- }
+ {
+ return Compute_squared_radius_3()(*this)(point(s,0), point(s,1),
+ point(s,2), point(s,3));
+ }
NT squared_radius(const Cell_handle& s, const int& i) const
- {
- return Compute_squared_radius_3()(*this) (
- this->point(s,vertex_triple_index(i,0)),
- this->point(s,vertex_triple_index(i,1)),
- this->point(s,vertex_triple_index(i,2)) );
- }
+ {
+ return Compute_squared_radius_3()(*this)(point(s,vertex_triple_index(i,0)),
+ point(s,vertex_triple_index(i,1)),
+ point(s,vertex_triple_index(i,2)));
+ }
NT squared_radius(const Facet& f) const {
return squared_radius(f.first, f.second);
}
- NT squared_radius(const Cell_handle& s,
- const int& i, const int& j) const
- {
- return Compute_squared_radius_3()(*this)(
- this->point(s,i), this->point(s,j));
- }
+ NT squared_radius(const Cell_handle& s, const int& i, const int& j) const
+ {
+ return Compute_squared_radius_3()(*this)(point(s,i), point(s,j));
+ }
NT squared_radius(const Edge& e) const {
- return squared_radius(e.first,e.second,e.third);
+ return squared_radius(e.first,e.second,e.third);
}
NT squared_radius(const Vertex_handle& v) const {
return Compute_squared_radius_3()(*this)(v->point());
}
-
//---------------------------------------------------------------------
private:
@@ -792,6 +787,7 @@ private:
//---------------------------------------------------------------------
public:
#ifdef CGAL_USE_GEOMVIEW
+ void show_triangulation_edges(Geomview_stream &gv) const;
void show_alpha_shape_faces(Geomview_stream &gv) const;
#endif
diff --git a/Alpha_shapes_3/include/CGAL/IO/alpha_shape_geomview_ostream_3.h b/Alpha_shapes_3/include/CGAL/IO/alpha_shape_geomview_ostream_3.h
index 1013af9c64f..e7c5981177a 100644
--- a/Alpha_shapes_3/include/CGAL/IO/alpha_shape_geomview_ostream_3.h
+++ b/Alpha_shapes_3/include/CGAL/IO/alpha_shape_geomview_ostream_3.h
@@ -12,21 +12,20 @@
// This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
// WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
//
-// $URL$
-// $Id$
-//
-//
// Author(s) : Tran Kai Frank DA
+// Mael Rouxel-Labbé
#ifndef CGAL_IO_ALPHA_SHAPE_GEOMVIEW_OSTREAM_3_H
#define CGAL_IO_ALPHA_SHAPE_GEOMVIEW_OSTREAM_3_H
#include
-
#include
#include
+#include