LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
Voronoi.h
Go to the documentation of this file.
1 
9 #ifndef VoronoiDiagram_h
10 #define VoronoiDiagram_h
11 
12 // std includes
13 #include <queue>
14 
15 // LArSoft includes
19 
20 #include <boost/polygon/voronoi.hpp>
21 
22 // Algorithm includes
24 //------------------------------------------------------------------------------------------------------------------------------------------
25 
26 namespace voronoi2d {
31  public:
32  using PointPair = std::pair<dcel2d::Point, dcel2d::Point>;
33  using MinMaxPointPair = std::pair<PointPair, PointPair>;
34 
41 
46 
53 
55 
59  const dcel2d::FaceList& getFaceList() const { return fFaceList; }
60 
64  const dcel2d::VertexList& getVertexList() const { return fVertexList; }
65 
69  double getVoronoiDiagramArea() const { return fVoronoiDiagramArea; }
70 
74  const dcel2d::PointList& getConvexHull() const { return fConvexHullList; }
75 
80 
84  PointPair findNearestEdge(const dcel2d::Point&, double&) const;
85 
89  double findNearestDistance(const dcel2d::Point&) const;
90 
91  private:
92  using EventQueue =
93  std::priority_queue<IEvent*, std::vector<IEvent*>, bool (*)(const IEvent*, const IEvent*)>;
94 
98  void handleSiteEvents(BeachLine&, EventQueue&, IEvent*);
99 
103  void handleCircleEvents(BeachLine&, EventQueue&, IEvent*);
104 
105  void makeLeftCircleEvent(EventQueue&, BSTNode*, double);
106  void makeRightCircleEvent(EventQueue&, BSTNode*, double);
107 
111  IEvent* makeCircleEvent(BSTNode*, BSTNode*, BSTNode*, double);
112 
114  const dcel2d::Coords&,
115  const dcel2d::Coords&,
117  double&,
118  double&) const;
119 
121  const dcel2d::Coords&,
122  const dcel2d::Coords&,
124  double&) const;
125 
127  const dcel2d::Coords&,
128  const dcel2d::Coords&,
130  double&) const;
131 
135  void getConvexHull(const BSTNode*);
136 
140  void terminateInfiniteEdges(BeachLine&, double);
141 
145  bool isInsideConvexHull(const dcel2d::Vertex&) const;
146 
153  double&) const;
154 
158  void findBoundingBox(const dcel2d::VertexList&);
159 
163  using BoostEdgeToEdgeMap =
164  std::map<const boost::polygon::voronoi_edge<double>*, dcel2d::HalfEdge*>;
165  using BoostVertexToVertexMap =
166  std::map<const boost::polygon::voronoi_vertex<double>*, dcel2d::Vertex*>;
167  using BoostCellToFaceMap = std::map<const boost::polygon::voronoi_cell<double>*, dcel2d::Face*>;
168 
170  const boost::polygon::voronoi_edge<double>*,
171  const boost::polygon::voronoi_edge<double>*,
175 
180 
184  double ComputeFaceArea();
185 
189  double crossProduct(const dcel2d::Point& p0,
190  const dcel2d::Point& p1,
191  const dcel2d::Point& p2) const;
192 
196  double Area() const;
197 
201  bool isLeft(const dcel2d::Point& p0,
202  const dcel2d::Point& p1,
203  const dcel2d::Point& pCheck) const;
204 
208 
210  SiteEventList fSiteEventList; //< Container for site events
211  CircleEventList fCircleEventList; //< Container for circle events
212  BSTNodeList fCircleNodeList; //< Container for the circle "nodes"
213 
214  dcel2d::PointList fConvexHullList; //< Points representing the convex hull
215  dcel2d::Coords fConvexHullCenter; //< Center of the convex hull
216 
217  double fXMin; //< Bounding box min value x
218  double fXMax; //< Bounding box max value x
219  double fYMin; //< Bounding box min value y
220  double fYMax; //< Bounding box max value y
221  mutable int fNumBadCircles; //< Number bad circles
223 
224  EventUtilities fUtilities; //< Handy functions to operate on arcs
225  };
226 
227 } // namespace lar_cluster3d
228 #endif
double getVoronoiDiagramArea() const
recover the area of the convex hull
Definition: Voronoi.h:69
This defines the actual beach line. The idea is to implement this as a self balancing binary search t...
Definition: BeachLine.h:125
std::list< SiteEvent > SiteEventList
Definition: SweepEvent.h:99
dcel2d::HalfEdgeList & fHalfEdgeList
Definition: Voronoi.h:205
std::list< HalfEdge > HalfEdgeList
Definition: DCEL.h:171
std::list< CircleEvent > CircleEventList
Definition: SweepEvent.h:100
Provides some basic functions operating in IEvent class objects.
std::pair< dcel2d::Point, dcel2d::Point > PointPair
Definition: Voronoi.h:32
BSTNode class definiton specifically for use in constructing Voronoi diagrams. We are trying to follo...
Definition: BeachLine.h:34
std::pair< PointPair, PointPair > MinMaxPointPair
Definition: Voronoi.h:33
bool isInsideConvexHull(const dcel2d::Vertex &) const
Is a vertex inside the convex hull - meant to be a fast check.
Definition: Voronoi.cxx:1158
double Area() const
Computes the area of the enclosed convext hull.
Definition: Voronoi.cxx:104
bool computeCircleCenter2(const dcel2d::Coords &, const dcel2d::Coords &, const dcel2d::Coords &, dcel2d::Coords &, double &) const
Definition: Voronoi.cxx:785
intermediate_table::const_iterator const_iterator
dcel2d::FaceList & fFaceList
Definition: Voronoi.h:207
std::list< BSTNode > BSTNodeList
Definition: BeachLine.h:118
Internal class definitions to facilitate construction of diagram.
Represents the beachline implemented as a self balancing binary search tree.
~VoronoiDiagram()
Destructor.
Definition: Voronoi.cxx:74
void makeLeftCircleEvent(EventQueue &, BSTNode *, double)
Definition: Voronoi.cxx:591
bool computeCircleCenter(const dcel2d::Coords &, const dcel2d::Coords &, const dcel2d::Coords &, dcel2d::Coords &, double &, double &) const
Definition: Voronoi.cxx:714
BSTNodeList fCircleNodeList
Definition: Voronoi.h:212
IEvent * makeCircleEvent(BSTNode *, BSTNode *, BSTNode *, double)
There are two types of events in the queue, here we handle circle events.
Definition: Voronoi.cxx:673
void handleCircleEvents(BeachLine &, EventQueue &, IEvent *)
There are two types of events in the queue, here we handle circle events.
Definition: Voronoi.cxx:503
VoronoiDiagram class definiton.
Definition: Voronoi.h:30
std::list< Face > FaceList
Definition: DCEL.h:170
CircleEventList fCircleEventList
Definition: Voronoi.h:211
SiteEventList fSiteEventList
Definition: Voronoi.h:210
EventUtilities fUtilities
Definition: Voronoi.h:224
double ComputeFaceArea()
Compute the area of the faces.
Definition: Voronoi.cxx:1270
void buildVoronoiDiagram(const dcel2d::PointList &)
Given an input set of 2D points construct a 2D voronoi diagram.
Definition: Voronoi.cxx:135
double findNearestDistance(const dcel2d::Point &) const
Given an input point, find the distance to the nearest edge/point.
Definition: Voronoi.cxx:1494
void findBoundingBox(const dcel2d::VertexList &)
Find the min/max values in x-y to use as a bounding box.
Definition: Voronoi.cxx:1408
dcel2d::VertexList & fVertexList
Definition: Voronoi.h:206
std::priority_queue< IEvent *, std::vector< IEvent * >, bool(*)(const IEvent *, const IEvent *)> EventQueue
Definition: Voronoi.h:93
void mergeDegenerateVertices()
merge degenerate vertices (found by zero length edges)
Definition: Voronoi.cxx:1240
dcel2d::Coords fConvexHullCenter
Definition: Voronoi.h:215
void makeRightCircleEvent(EventQueue &, BSTNode *, double)
Definition: Voronoi.cxx:633
dcel2d::PointList fPointList
Definition: Voronoi.h:209
double crossProduct(const dcel2d::Point &p0, const dcel2d::Point &p1, const dcel2d::Point &p2) const
Gets the cross product of line from p0 to p1 and p0 to p2.
Definition: Voronoi.cxx:89
PointPair findNearestEdge(const dcel2d::Point &, double &) const
Given an input Point, find the nearest edge.
Definition: Voronoi.cxx:1433
dcel2d::PointList fConvexHullList
Definition: Voronoi.h:214
std::tuple< double, double, const reco::ClusterHit3D * > Point
Definitions used by the VoronoiDiagram algorithm.
Definition: DCEL.h:42
const dcel2d::PointList & getConvexHull() const
recover the point list representing the convex hull
Definition: Voronoi.h:74
void boostTranslation(const dcel2d::PointList &, const boost::polygon::voronoi_edge< double > *, const boost::polygon::voronoi_edge< double > *, BoostEdgeToEdgeMap &, BoostVertexToVertexMap &, BoostCellToFaceMap &)
Definition: Voronoi.cxx:355
PointPair getExtremePoints() const
Given an input Point, find the nearest edge.
Definition: Voronoi.cxx:1097
std::map< const boost::polygon::voronoi_edge< double > *, dcel2d::HalfEdge * > BoostEdgeToEdgeMap
Translate boost to dcel.
Definition: Voronoi.h:164
std::map< const boost::polygon::voronoi_vertex< double > *, dcel2d::Vertex * > BoostVertexToVertexMap
Definition: Voronoi.h:166
void terminateInfiniteEdges(BeachLine &, double)
this aims to process remaining elements in the beachline after the event queue has been depleted ...
Definition: Voronoi.cxx:854
bool computeCircleCenter3(const dcel2d::Coords &, const dcel2d::Coords &, const dcel2d::Coords &, dcel2d::Coords &, double &) const
Definition: Voronoi.cxx:818
std::list< Point > PointList
Definition: DCEL.h:43
void buildVoronoiDiagramBoost(const dcel2d::PointList &)
Definition: Voronoi.cxx:282
Eigen::Vector3f Coords
Definition: DCEL.h:44
const dcel2d::VertexList & getVertexList() const
Recover the list of vertices.
Definition: Voronoi.h:64
const dcel2d::FaceList & getFaceList() const
Recover the list of faces.
Definition: Voronoi.h:59
std::list< Vertex > VertexList
Definition: DCEL.h:169
std::map< const boost::polygon::voronoi_cell< double > *, dcel2d::Face * > BoostCellToFaceMap
Definition: Voronoi.h:167
VoronoiDiagram(dcel2d::HalfEdgeList &, dcel2d::VertexList &, dcel2d::FaceList &)
Constructor.
Definition: Voronoi.cxx:47
bool isOutsideConvexHull(const dcel2d::Vertex &, dcel2d::PointList::const_iterator, dcel2d::Coords &, double &) const
is vertex outside the convex hull and if so return some useful information
Definition: Voronoi.cxx:1185
bool isLeft(const dcel2d::Point &p0, const dcel2d::Point &p1, const dcel2d::Point &pCheck) const
Determines whether a point is to the left, right or on line specifiec by p0 and p1.
Definition: Voronoi.cxx:78
void handleSiteEvents(BeachLine &, EventQueue &, IEvent *)
There are two types of events in the queue, here we handle site events.
Definition: Voronoi.cxx:445