summaryrefslogtreecommitdiff
path: root/src/common/include/gudhi
diff options
context:
space:
mode:
authorHind-M <hind.montassif@gmail.com>2021-09-30 11:05:17 +0200
committerHind-M <hind.montassif@gmail.com>2021-09-30 11:05:17 +0200
commit767d9fca5da2d3dd9698a5c27e9bedc159271f67 (patch)
treedefb85dddefbe6cd2064059baf9baeed0e72b87e /src/common/include/gudhi
parent44f754ee58aeee043891f4494892798b9807374b (diff)
Move minimal enclosing balls cache to Cech_complex.h instead of the blocker
Modify cech test to be more relevant regarding the current algorithm Do some cleaning
Diffstat (limited to 'src/common/include/gudhi')
-rw-r--r--src/common/include/gudhi/distance_functions.h20
-rw-r--r--src/common/include/gudhi/graph_simplicial_complex.h3
2 files changed, 6 insertions, 17 deletions
diff --git a/src/common/include/gudhi/distance_functions.h b/src/common/include/gudhi/distance_functions.h
index ae5168aa..a8ee4a75 100644
--- a/src/common/include/gudhi/distance_functions.h
+++ b/src/common/include/gudhi/distance_functions.h
@@ -65,19 +65,17 @@ class Euclidean_distance {
* The points are assumed to have the same dimension. */
class Minimal_enclosing_ball_radius {
public:
- /** \brief TODO
+ /** \brief Enclosing ball radius from two points using CGAL.
*
* @param[in] point_1
* @param[in] point_2
- * @return
- * \tparam Point
+ * @return Enclosing ball radius for the two points.
+ * \tparam Point must be a Kernel::Point_d from CGAL.
*
*/
- //typename FT = typename Kernel::FT,
template< typename Kernel = CGAL::Epeck_d<CGAL::Dynamic_dimension_tag>,
typename Point= typename Kernel::Point_d>
double operator()(const Point& point_1, const Point& point_2) const {
- std::clog << "Added template: distance betw points 1 and 2" << std::endl;
Kernel kernel_;
return std::sqrt(CGAL::to_double(kernel_.squared_distance_d_object()(point_1, point_2))) / 2.;
}
@@ -94,25 +92,21 @@ class Minimal_enclosing_ball_radius {
template< typename Point >
typename std::iterator_traits<typename boost::range_iterator<Point>::type>::value_type
operator()(const Point& point_1, const Point& point_2) const {
- std::clog << "Hind: Minimal_enclosing_ball_radius point1 et 2; Euclidean" << std::endl;
std::clog << "#" << *point_1.begin() << "##" << *point_2.begin() << std::endl;
return Euclidean_distance()(point_1, point_2) / 2.;
}
-
- /** \brief TODO
+ /** \brief Enclosing ball radius from a point cloud using CGAL.
*
* @param[in] point_cloud The points.
- * @return
- * \tparam Point_cloud
+ * @return Enclosing ball radius for the points.
+ * \tparam Point_cloud must be a range of Kernel::Point_d points from CGAL.
*
*/
- //typename FT = typename Kernel::FT,
template< typename Kernel = CGAL::Epeck_d<CGAL::Dynamic_dimension_tag>,
typename Point= typename Kernel::Point_d,
typename Point_cloud = std::vector<Point>>
double operator()(const Point_cloud& point_cloud) const {
- std::clog << "Added template: distance in point cloud" << std::endl;
Kernel kernel_;
return std::sqrt(CGAL::to_double(kernel_.compute_squared_radius_d_object()(point_cloud.begin(), point_cloud.end())));
}
@@ -133,8 +127,6 @@ class Minimal_enclosing_ball_radius {
typename Coordinate = typename std::iterator_traits<Coordinate_iterator>::value_type>
Coordinate
operator()(const Point_cloud& point_cloud) const {
- std::clog << "Hind: Minimal_enclosing_ball_radius point cloud; Miniball" << std::endl;
-
using Min_sphere = Miniball::Miniball<Miniball::CoordAccessor<Point_iterator, Coordinate_iterator>>;
Min_sphere ms(boost::size(*point_cloud.begin()), point_cloud.begin(), point_cloud.end());
diff --git a/src/common/include/gudhi/graph_simplicial_complex.h b/src/common/include/gudhi/graph_simplicial_complex.h
index 9190182c..da9dee7d 100644
--- a/src/common/include/gudhi/graph_simplicial_complex.h
+++ b/src/common/include/gudhi/graph_simplicial_complex.h
@@ -18,8 +18,6 @@
#include <map>
#include <tuple> // for std::tie
-#include <iostream>
-
namespace Gudhi {
/** @file
* @brief Graph simplicial complex methods
@@ -78,7 +76,6 @@ Proximity_graph<SimplicialComplexForProximityGraph> compute_proximity_graph(
for (auto it_u = points.begin(); it_u != points.end(); ++it_u) {
idx_v = idx_u + 1;
for (auto it_v = it_u + 1; it_v != points.end(); ++it_v, ++idx_v) {
- std::clog << "#idx_u" << idx_u << "#idx_v " << idx_v << std::endl;
fil = distance(*it_u, *it_v);
if (fil <= threshold) {
edges.emplace_back(idx_u, idx_v);