LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
lar_cluster3d::ConvexHull Class Reference

ConvexHull class definiton. More...

#include "ConvexHull.h"

Public Types

using Point = std::tuple< float, float, const reco::ClusterHit3D * >
 Definitions used by the ConvexHull algorithm. More...
 
using PointList = std::list< Point >
 
using PointPair = std::pair< Point, Point >
 
using MinMaxPointPair = std::pair< PointPair, PointPair >
 

Public Member Functions

 ConvexHull (const PointList &)
 Constructor. More...
 
virtual ~ConvexHull ()
 Destructor. More...
 
const PointListgetPointsList ()
 recover the list of points used to build convex hull More...
 
const PointListgetConvexHull () const
 recover the list of convex hull vertices More...
 
const MinMaxPointPairgetMinMaxPointPair () const
 find the ends of the convex hull (along its x axis) More...
 
float getConvexHullArea () const
 recover the area of the convex hull More...
 
PointPair findNearestEdge (const Point &, float &) const
 Given an input Point, find the nearest edge. More...
 
float findNearestDistance (const Point &) const
 Given an input point, find the distance to the nearest edge/point. More...
 

Private Member Functions

void getConvexHull (const PointList &)
 Given an input set of 2D points build a convex hull around them. More...
 
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. More...
 
float Area () const
 Computes the area of the enclosed convext hull. More...
 
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. More...
 

Private Attributes

const PointListfPoints
 
PointList fConvexHull
 
MinMaxPointPair fMinMaxPointPair
 
float fConvexHullArea
 

Detailed Description

ConvexHull class definiton.

Definition at line 25 of file ConvexHull.h.

Member Typedef Documentation

Definition at line 34 of file ConvexHull.h.

using lar_cluster3d::ConvexHull::Point = std::tuple<float,float,const reco::ClusterHit3D*>

Definitions used by the ConvexHull algorithm.

Definition at line 31 of file ConvexHull.h.

Definition at line 32 of file ConvexHull.h.

Definition at line 33 of file ConvexHull.h.

Constructor & Destructor Documentation

lar_cluster3d::ConvexHull::ConvexHull ( const PointList pointListIn)

Constructor.

Parameters
pset

Definition at line 28 of file ConvexHull.cxx.

References Area(), fConvexHull, fConvexHullArea, fPoints, and getConvexHull().

28  : fPoints(pointListIn), fConvexHullArea(0.)
29 {
30  fConvexHull.clear();
31 
32  // Build out the convex hull around the input points
34 
35  // And the area
37 }
float Area() const
Computes the area of the enclosed convext hull.
Definition: ConvexHull.cxx:69
const PointList & getConvexHull() const
recover the list of convex hull vertices
Definition: ConvexHull.h:56
const PointList & fPoints
Definition: ConvexHull.h:102
lar_cluster3d::ConvexHull::~ConvexHull ( )
virtual

Destructor.

Definition at line 41 of file ConvexHull.cxx.

42 {
43 }

Member Function Documentation

float lar_cluster3d::ConvexHull::Area ( ) const
private

Computes the area of the enclosed convext hull.

Definition at line 69 of file ConvexHull.cxx.

References crossProduct(), and fConvexHull.

Referenced by ConvexHull(), and getConvexHullArea().

70 {
71  float area(0.);
72 
73  // Compute the area by taking advantage of
74  // 1) the ability to decompose a convex hull into triangles,
75  // 2) the ability to use the cross product to calculate the area
76  // So, the technique is to pick a point (for technical reasons we use 0,0)
77  // and then sum the signed area of triangles formed from this point to two adjecent
78  // vertices on the convex hull.
79  Point center(0.,0.,0);
80  Point lastPoint = fConvexHull.front();
81 
82  for(const auto& point : fConvexHull)
83  {
84  if (point != lastPoint) area += 0.5 * crossProduct(center,lastPoint,point);
85 
86  lastPoint = point;
87  }
88 
89  return area;
90 }
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:56
std::tuple< float, float, const reco::ClusterHit3D * > Point
Definitions used by the ConvexHull algorithm.
Definition: ConvexHull.h:31
float lar_cluster3d::ConvexHull::crossProduct ( const Point p0,
const Point p1,
const Point p2 
) const
private

Gets the cross product of line from p0 to p1 and p0 to p2.

Definition at line 56 of file ConvexHull.cxx.

Referenced by Area(), getConvexHullArea(), and isLeft().

57 {
58  // Define a quick 2D cross product here since it will used quite a bit!
59  float deltaX = std::get<0>(p1) - std::get<0>(p0);
60  float deltaY = std::get<1>(p1) - std::get<1>(p0);
61  float dCheckX = std::get<0>(p2) - std::get<0>(p0);
62  float dCheckY = std::get<1>(p2) - std::get<1>(p0);
63 
64  return ((deltaX * dCheckY) - (deltaY * dCheckX));
65 }
float lar_cluster3d::ConvexHull::findNearestDistance ( const Point point) const

Given an input point, find the distance to the nearest edge/point.

Definition at line 255 of file ConvexHull.cxx.

References findNearestEdge().

Referenced by getConvexHullArea().

256 {
257  float closestDistance;
258 
259  findNearestEdge(point,closestDistance);
260 
261  return closestDistance;
262 }
PointPair findNearestEdge(const Point &, float &) const
Given an input Point, find the nearest edge.
Definition: ConvexHull.cxx:194
ConvexHull::PointPair lar_cluster3d::ConvexHull::findNearestEdge ( const Point point,
float &  closestDistance 
) const

Given an input Point, find the nearest edge.

Definition at line 194 of file ConvexHull.cxx.

References fConvexHull, isLeft(), and max.

Referenced by findNearestDistance(), and getConvexHullArea().

195 {
196  // The idea is to find the nearest edge of the convex hull, defined by
197  // two adjacent vertices of the hull, to the input point.
198  // As near as I can tell, the best way to do this is to apply brute force...
199  // Idea will be to iterate over pairs of points
200  PointList::const_iterator curPointItr = fConvexHull.begin();
201  Point prevPoint = *curPointItr++;
202  Point curPoint = *curPointItr;
203 
204  // Set up the winner
205  PointPair closestEdge(prevPoint,curPoint);
206 
207  closestDistance = std::numeric_limits<float>::max();
208 
209  // curPointItr is meant to point to the second point
210  while(curPointItr != fConvexHull.end())
211  {
212  if (curPoint != prevPoint)
213  {
214  // Dereference some stuff
215  float xPrevToPoint = (std::get<0>(point) - std::get<0>(prevPoint));
216  float yPrevToPoint = (std::get<1>(point) - std::get<1>(prevPoint));
217  float xPrevToCur = (std::get<0>(curPoint) - std::get<0>(prevPoint));
218  float yPrevToCur = (std::get<1>(curPoint) - std::get<1>(prevPoint));
219  float edgeLength = std::sqrt(xPrevToCur * xPrevToCur + yPrevToCur * yPrevToCur);
220 
221  // Find projection onto convex hull edge
222  float projection = ((xPrevToPoint * xPrevToCur) + (yPrevToPoint * yPrevToCur)) / edgeLength;
223 
224  // DOCA point
225  Point docaPoint(std::get<0>(prevPoint) + projection * xPrevToCur / edgeLength,
226  std::get<1>(prevPoint) + projection * yPrevToCur / edgeLength, 0);
227 
228  if (projection > edgeLength) docaPoint = curPoint;
229  if (projection < 0) docaPoint = prevPoint;
230 
231  float xDocaDist = std::get<0>(point) - std::get<0>(docaPoint);
232  float yDocaDist = std::get<1>(point) - std::get<1>(docaPoint);
233  float docaDist = xDocaDist * xDocaDist + yDocaDist * yDocaDist;
234 
235  if (docaDist < closestDistance)
236  {
237  closestEdge = PointPair(prevPoint,curPoint);
238  closestDistance = docaDist;
239  }
240  }
241 
242  prevPoint = curPoint;
243  curPoint = *curPointItr++;
244  }
245 
246  closestDistance = std::sqrt(closestDistance);
247 
248  // Convention is convex hull vertices sorted in counter clockwise fashion so if the point
249  // lays to the left of the nearest edge then it must be an interior point
250  if (isLeft(closestEdge.first,closestEdge.second,point)) closestDistance = -closestDistance;
251 
252  return closestEdge;
253 }
std::pair< Point, Point > PointPair
Definition: ConvexHull.h:33
Int_t max
Definition: plot.C:27
std::tuple< float, float, const reco::ClusterHit3D * > Point
Definitions used by the ConvexHull algorithm.
Definition: ConvexHull.h:31
intermediate_table::const_iterator const_iterator
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:47
const PointList& lar_cluster3d::ConvexHull::getConvexHull ( ) const
inline
void lar_cluster3d::ConvexHull::getConvexHull ( const PointList pointList)
private

Given an input set of 2D points build a convex hull around them.

Parameters
PointListThe input list of 2D points to build hull around

Definition at line 94 of file ConvexHull.cxx.

References fConvexHull, fMinMaxPointPair, and isLeft().

95 {
96  // Start by identifying the min/max points
97  Point pMinMin = pointList.front();
98 
99  // Find the point with maximum y sharing the same x value...
100  PointList::const_iterator pMinMaxItr = std::find_if(pointList.begin(),pointList.end(),[pMinMin](const auto& elem){return std::get<0>(elem) != std::get<0>(pMinMin);});
101 
102  Point pMinMax = *(--pMinMaxItr); // This could be / probably is the same point
103 
104  fMinMaxPointPair.first.first = pMinMin;
105  fMinMaxPointPair.first.second = pMinMax;
106 
107  Point pMaxMax = pointList.back(); // get the last point
108 
109  // Find the point with minimum y sharing the same x value...
110  PointList::const_reverse_iterator pMaxMinItr = std::find_if(pointList.rbegin(),pointList.rend(),[pMaxMax](const auto& elem){return std::get<0>(elem) != std::get<0>(pMaxMax);});
111 
112  Point pMaxMin = *(--pMaxMinItr); // This could be / probably is the same point
113 
114  fMinMaxPointPair.second.first = pMaxMin;
115  fMinMaxPointPair.second.second = pMaxMax;
116 
117  // Get the lower convex hull
118  PointList lowerHullList;
119 
120  lowerHullList.push_back(pMinMin);
121 
122  // loop over points in the set
123  for(auto curPoint : pointList)
124  {
125  // First check that we even want to consider this point
126  if (isLeft(pMinMin,pMaxMin,curPoint)) continue;
127 
128  // In words: treat "lowerHullVec" as a stack. If there is only one entry then
129  // push the current point onto the stack. If more than one entry then check if
130  // the current point is to the left of the line from the last two points in the
131  // stack. If it is then add the current point to the stack, if it is not then
132  // pop the last point off the stack and try again.
133  while(lowerHullList.size() > 1)
134  {
135  Point lastPoint = lowerHullList.back();
136 
137  lowerHullList.pop_back();
138 
139  Point nextToLastPoint = lowerHullList.back();
140 
141  if (isLeft(nextToLastPoint,lastPoint,curPoint))
142  {
143  lowerHullList.push_back(lastPoint);
144  break;
145  }
146  }
147 
148  lowerHullList.push_back(curPoint);
149  }
150 
151  // Now get the upper hull
152  PointList upperHullList;
153 
154  upperHullList.push_back(pMaxMax);
155 
156  for(auto curPoint : boost::adaptors::reverse(pointList))
157  {
158  // First check that we even want to consider this point
159  // Remember that we are going "backwards" so still want
160  // curPoint to lie to the "left"
161  if (isLeft(pMaxMax,pMinMax,curPoint)) continue;
162 
163  // Replicate the above but going the other direction...
164  while(upperHullList.size() > 1)
165  {
166  Point lastPoint = upperHullList.back();
167 
168  upperHullList.pop_back();
169 
170  Point nextToLastPoint = upperHullList.back();
171 
172  if (isLeft(nextToLastPoint,lastPoint,curPoint))
173  {
174  upperHullList.push_back(lastPoint);
175  break;
176  }
177  }
178 
179  upperHullList.push_back(curPoint);
180  }
181 
182  // Now we merge the two lists into the output list
183  std::copy(lowerHullList.begin(),lowerHullList.end(),std::back_inserter(fConvexHull));
184 
185  if (pMaxMin == pMaxMax) upperHullList.pop_front();
186 
187  std::copy(upperHullList.begin(),upperHullList.end(),std::back_inserter(fConvexHull));
188 
189  if (pMinMin != pMinMax) fConvexHull.push_back(pMinMin);
190 
191  return;
192 }
MinMaxPointPair fMinMaxPointPair
Definition: ConvexHull.h:104
std::list< Point > PointList
Definition: ConvexHull.h:32
std::tuple< float, float, const reco::ClusterHit3D * > Point
Definitions used by the ConvexHull algorithm.
Definition: ConvexHull.h:31
intermediate_table::const_iterator const_iterator
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:47
float lar_cluster3d::ConvexHull::getConvexHullArea ( ) const
inline
const MinMaxPointPair& lar_cluster3d::ConvexHull::getMinMaxPointPair ( ) const
inline

find the ends of the convex hull (along its x axis)

Definition at line 61 of file ConvexHull.h.

References fMinMaxPointPair.

61 {return fMinMaxPointPair;}
MinMaxPointPair fMinMaxPointPair
Definition: ConvexHull.h:104
const PointList& lar_cluster3d::ConvexHull::getPointsList ( )
inline

recover the list of points used to build convex hull

Definition at line 51 of file ConvexHull.h.

References fPoints.

51 {return fPoints;}
const PointList & fPoints
Definition: ConvexHull.h:102
bool lar_cluster3d::ConvexHull::isLeft ( const Point p0,
const Point p1,
const Point pCheck 
) const
private

Determines whether a point is to the left, right or on line specifiec by p0 and p1.

Definition at line 47 of file ConvexHull.cxx.

References crossProduct().

Referenced by findNearestEdge(), getConvexHull(), and getConvexHullArea().

48 {
49  // Use the cross product to determine if the check point lies to the left, on or right
50  // of the line defined by points p0 and p1
51  return crossProduct(p0, p1, pCheck) > 0;
52 }
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:56

Member Data Documentation

PointList lar_cluster3d::ConvexHull::fConvexHull
private

Definition at line 103 of file ConvexHull.h.

Referenced by Area(), ConvexHull(), findNearestEdge(), and getConvexHull().

float lar_cluster3d::ConvexHull::fConvexHullArea
private

Definition at line 105 of file ConvexHull.h.

Referenced by ConvexHull(), and getConvexHullArea().

MinMaxPointPair lar_cluster3d::ConvexHull::fMinMaxPointPair
private

Definition at line 104 of file ConvexHull.h.

Referenced by getConvexHull(), and getMinMaxPointPair().

const PointList& lar_cluster3d::ConvexHull::fPoints
private

Definition at line 102 of file ConvexHull.h.

Referenced by ConvexHull(), and getPointsList().


The documentation for this class was generated from the following files: