Skip to Content.
Sympa Menu

cgal-discuss - Re: [cgal-discuss] convex hull, indices of extreme points

Subject: CGAL users discussion list

List archive

Re: [cgal-discuss] convex hull, indices of extreme points


Chronological Thread 
  • From: Andreas Fabri <>
  • To:
  • Subject: Re: [cgal-discuss] convex hull, indices of extreme points
  • Date: Sat, 29 Oct 2011 23:16:57 +0200
  • Organization: GeometryFactory

On 29/10/2011 22:20,

wrote:
Dear CGAL community,

following
http://www.cgal.org/Manual/latest/doc_html/cgal_manual/Convex_hull_2/Chapter_main.html#Section_15.3
to perform a convec hull computation was quite simple and straight forward -
thank you very much.

Starting from there: Is it possible to receive the indices of the resulting
extreme points in the original set of points?

Thanks a lot,
Thomas


The traits classes are all about this kind of plumbing.
Instead of the transform iterator you could have first
created a container with point/index pairs.

See the attached code.

best regards,

andreas

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/convex_hull_2.h>

#include <vector>
#include <boost/iterator/transform_iterator.hpp>

typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
typedef K::Point_2 Point_2;

typedef std::pair<Point_2,int> PointWithIndex;

typedef std::vector<Point_2> Points;


class pwi {
private:
  const Point_2* base;
public:
  typedef PointWithIndex result_type;

  pwi(): base(0) {}

  pwi(const Point_2* base): base(base) {}

  PointWithIndex operator()(const Point_2& p) const
  {
    return PointWithIndex(p, &p-base);
  }
};

typedef boost::transform_iterator<pwi, Points::iterator> Iterator;


struct My_traits {
  typedef PointWithIndex Point_2;

  struct Less_xy_2 {
    bool operator()(const Point_2& p, const Point_2& q) const
    {
      K::Less_xy_2 less_xy_2;
      return less_xy_2(p.first, q.first);
    }
  };

  struct Left_turn_2 {
    bool operator()(const Point_2& p, const Point_2& q, const Point_2& r) const
    {
      K::Left_turn_2 left_turn_2;
      return left_turn_2(p.first, q.first, r.first);
    }
  };

  struct Equal_2 {
    bool operator()(const Point_2& p, const Point_2& q) const
    {
      K::Equal_2 equal_2;
      return equal_2(p.first, q.first);
    }
  };

  Equal_2 equal_2_object() const { return Equal_2();}
  Less_xy_2 less_xy_2_object() const { return Less_xy_2();}
  Left_turn_2 left_turn_2_object() const { return Left_turn_2();}
};

int main()
{
  Points points;
  std::vector<PointWithIndex> result;
  points.push_back(Point_2(0,0));
  points.push_back(Point_2(10,0));
  points.push_back(Point_2(6,5));
  points.push_back(Point_2(4,6));
  points.push_back(Point_2(10,10));
  points.push_back(Point_2(4,1));

  Iterator it(points.begin(), pwi(&points[0])), end(points.end(), pwi());

  CGAL::ch_graham_andrew( it, end, std::back_inserter(result),My_traits() );
  std::cout << result.size() << " points on the convex hull" << std::endl;

  for(std::vector<PointWithIndex>::iterator it = result.begin(); it != result.end(); ++it){
    std::cout << it->first << " |  " << it->second << std::endl; 
  }
  return 0;
}



Archive powered by MHonArc 2.6.16.

Top of Page