LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
ConvexHull.h
Go to the documentation of this file.
1 
9 #ifndef ConvexHull_h
10 #define ConvexHull_h
11 
12 // Algorithm includes
14 
15 // std includes
16 #include <list>
17 #include <tuple>
18 #include <utility>
19 
20 //------------------------------------------------------------------------------------------------------------------------------------------
21 
22 namespace lar_cluster3d {
26  class ConvexHull {
27  public:
31  using Point =
32  std::tuple<float, float, const reco::ClusterHit3D*>;
33  using PointList = std::list<Point>;
34  using PointPair = std::pair<Point, Point>;
35  using MinMaxPointPair = std::pair<PointPair, PointPair>;
36 
42  ConvexHull(const PointList&, float = 0.85, float = 0.35);
43 
47  ~ConvexHull();
48 
52  const PointList& getPointsList() { return fPoints; }
53 
57  const PointList& getConvexHull() const { return fConvexHull; }
58 
63 
67  const PointList& getExtremePoints();
68 
73 
77  float getConvexHullArea() const { return fConvexHullArea; }
78 
82  PointPair findNearestEdge(const Point&, float&) const;
83 
87  float findNearestDistance(const Point&) const;
88 
89  private:
95  void getConvexHull(const PointList&);
96 
100  float crossProduct(const Point& p0, const Point& p1, const Point& p2) const;
101 
105  float Area() const;
106 
110  bool isLeft(const Point& p0, const Point& p1, const Point& pCheck) const;
111 
112  float fKinkAngle;
114 
121  };
122 
123 } // namespace lar_cluster3d
124 #endif
const PointList & getPointsList()
recover the list of points used to build convex hull
Definition: ConvexHull.h:52
std::tuple< float, float, const reco::ClusterHit3D * > Point
projected x,y position and 3D hit
Definition: ConvexHull.h:32
MinMaxPointPair fMinMaxPointPair
Definition: ConvexHull.h:117
std::list< Point > PointList
The list of the projected points.
Definition: ConvexHull.h:33
float crossProduct(const Point &p0, const Point &p1, const Point &p2) const
Gets the cross product of line from p0 to p1 and p0 to p2.
Definition: ConvexHull.cxx:59
std::pair< Point, Point > PointPair
Definition: ConvexHull.h:34
const MinMaxPointPair & getMinMaxPointPair() const
find the ends of the convex hull (along its x axis)
Definition: ConvexHull.h:62
std::pair< PointPair, PointPair > MinMaxPointPair
Definition: ConvexHull.h:35
float Area() const
Computes the area of the enclosed convext hull.
Definition: ConvexHull.cxx:72
PointPair findNearestEdge(const Point &, float &) const
Given an input Point, find the nearest edge.
Definition: ConvexHull.cxx:329
float findNearestDistance(const Point &) const
Given an input point, find the distance to the nearest edge/point.
Definition: ConvexHull.cxx:389
const reco::ConvexHullKinkTupleList & getKinkPoints()
Find the points with the largest angles.
Definition: ConvexHull.cxx:274
std::list< ConvexHullKinkTuple > ConvexHullKinkTupleList
Definition: Cluster3D.h:349
float getConvexHullArea() const
recover the area of the convex hull
Definition: ConvexHull.h:77
const PointList & getConvexHull() const
recover the list of convex hull vertices
Definition: ConvexHull.h:57
ConvexHull class definiton.
Definition: ConvexHull.h:26
reco::ConvexHullKinkTupleList fKinkPoints
Definition: ConvexHull.h:120
const PointList & getExtremePoints()
Find the two points on the hull which are furthest apart.
Definition: ConvexHull.cxx:206
bool isLeft(const Point &p0, const Point &p1, const Point &pCheck) const
Determines whether a point is to the left, right or on line specifiec by p0 and p1.
Definition: ConvexHull.cxx:50
ConvexHull(const PointList &, float=0.85, float=0.35)
Constructor.
Definition: ConvexHull.cxx:29
const PointList & fPoints
Definition: ConvexHull.h:115
~ConvexHull()
Destructor.
Definition: ConvexHull.cxx:46