summaryrefslogtreecommitdiff
path: root/src/Alpha_complex/include/gudhi/Alpha_complex.h
diff options
context:
space:
mode:
authorROUVREAU Vincent <vincent.rouvreau@inria.fr>2019-12-03 16:51:09 +0100
committerROUVREAU Vincent <vincent.rouvreau@inria.fr>2019-12-03 16:51:09 +0100
commitf87a1026d5d1accbb6a558ef5d50f5c4c035e0d8 (patch)
tree6e1246c10f520557e43ad4b5b97260bc49068fe5 /src/Alpha_complex/include/gudhi/Alpha_complex.h
parent7df0113ddf2892b0cccf025a836214021b301072 (diff)
Fix #134
Diffstat (limited to 'src/Alpha_complex/include/gudhi/Alpha_complex.h')
-rw-r--r--src/Alpha_complex/include/gudhi/Alpha_complex.h10
1 files changed, 6 insertions, 4 deletions
diff --git a/src/Alpha_complex/include/gudhi/Alpha_complex.h b/src/Alpha_complex/include/gudhi/Alpha_complex.h
index 6b4d8463..6f19cb6c 100644
--- a/src/Alpha_complex/include/gudhi/Alpha_complex.h
+++ b/src/Alpha_complex/include/gudhi/Alpha_complex.h
@@ -121,8 +121,8 @@ class Alpha_complex {
// size_type type from CGAL.
typedef typename Delaunay_triangulation::size_type size_type;
- // Map type to switch from simplex tree vertex handle to CGAL vertex iterator.
- typedef typename std::map< std::size_t, CGAL_vertex_iterator > Vector_vertex_iterator;
+ // Structure to switch from simplex tree vertex handle to CGAL vertex iterator.
+ typedef typename std::vector< CGAL_vertex_iterator > Vector_vertex_iterator;
private:
/** \brief Vertex iterator vector to switch from simplex tree vertex handle to CGAL vertex iterator.
@@ -238,14 +238,16 @@ class Alpha_complex {
hint = pos->full_cell();
}
// --------------------------------------------------------------------------------------------
- // double map to retrieve simplex tree vertex handles from CGAL vertex iterator and vice versa
+ // 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());
// Loop on triangulation vertices list
for (CGAL_vertex_iterator vit = triangulation_->vertices_begin(); vit != triangulation_->vertices_end(); ++vit) {
if (!triangulation_->is_infinite(*vit)) {
#ifdef DEBUG_TRACES
std::cout << "Vertex insertion - " << vit->data() << " -> " << vit->point() << std::endl;
#endif // DEBUG_TRACES
- vertex_handle_to_iterator_.emplace(vit->data(), vit);
+ vertex_handle_to_iterator_[vit->data()] = vit;
}
}
// --------------------------------------------------------------------------------------------