LArSoft  v06_85_00
Liquid Argon Software toolkit - http://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 <list>
14 #include <algorithm>
15 #include <queue>
16 
17 // LArSoft includes
21 
22 #include <boost/polygon/voronoi.hpp>
23 
24 // Algorithm includes
26 //------------------------------------------------------------------------------------------------------------------------------------------
27 
28 namespace voronoi2d
29 {
34 {
35 public:
36  using PointPair = std::pair<dcel2d::Point,dcel2d::Point>;
37  using MinMaxPointPair = std::pair<PointPair,PointPair>;
38 
45 
49  virtual ~VoronoiDiagram();
50 
57 
59 
63  const dcel2d::FaceList& getFaceList() const {return fFaceList;}
64 
68  const dcel2d::VertexList& getVertexList() const {return fVertexList;}
69 
73  double getVoronoiDiagramArea() const {return fVoronoiDiagramArea;}
74 
79 
83  PointPair findNearestEdge(const dcel2d::Point&, double&) const;
84 
88  double findNearestDistance(const dcel2d::Point&) const;
89 
90 private:
91 
92  using EventQueue = std::priority_queue<IEvent*, std::vector<IEvent*>, bool(*)(const IEvent*,const IEvent*)>;
93 
97  void handleSiteEvents(BeachLine&, EventQueue&, IEvent*);
98 
102  void handleCircleEvents(BeachLine&, EventQueue&, IEvent*);
103 
104  void makeLeftCircleEvent(EventQueue&, BSTNode*, double);
105  void makeRightCircleEvent(EventQueue&, BSTNode*, double);
106 
110  IEvent* makeCircleEvent(BSTNode*, BSTNode*, BSTNode*, double);
111 
112  bool computeCircleCenter(const dcel2d::Coords&, const dcel2d::Coords&, const dcel2d::Coords&, dcel2d::Coords&, double&, double&) const;
113 
114  bool computeCircleCenter2(const dcel2d::Coords&, const dcel2d::Coords&, const dcel2d::Coords&, dcel2d::Coords&, double&, double&) const;
115 
116  bool computeCircleCenter3(const dcel2d::Coords&, const dcel2d::Coords&, const dcel2d::Coords&, dcel2d::Coords&, double&, double&) const;
117 
121  void getConvexHull(const BSTNode*);
122 
126  void terminateInfiniteEdges(BeachLine&, double);
127 
131  bool isInsideConvexHull(const dcel2d::Vertex&) const;
132 
137 
141  void findBoundingBox(const dcel2d::VertexList&);
142 
146  using BoostEdgeToEdgeMap = std::map<const boost::polygon::voronoi_edge<double>*, dcel2d::HalfEdge*>;
147  using BoostVertexToVertexMap = std::map<const boost::polygon::voronoi_vertex<double>*, dcel2d::Vertex*>;
148  using BoostCellToFaceMap = std::map<const boost::polygon::voronoi_cell<double>*, dcel2d::Face*>;
149 
151  const boost::polygon::voronoi_edge<double>*,
152  const boost::polygon::voronoi_edge<double>*,
156 
161 
165  double ComputeFaceArea();
166 
170  double crossProduct(const dcel2d::Point& p0, const dcel2d::Point& p1, const dcel2d::Point& p2) const;
171 
175  double Area() const;
176 
180  bool isLeft(const dcel2d::Point& p0, const dcel2d::Point& p1, const dcel2d::Point& pCheck) const;
181 
185 
187  SiteEventList fSiteEventList; //< Container for site events
188  CircleEventList fCircleEventList; //< Container for circle events
189  BSTNodeList fCircleNodeList; //< Container for the circle "nodes"
190 
191  dcel2d::PointList fConvexHullList; //< Points representing the convex hull
192  dcel2d::Coords fConvexHullCenter; //< Center of the convex hull
193 
194  double fXMin; //< Bounding box min value x
195  double fXMax; //< Bounding box max value x
196  double fYMin; //< Bounding box min value y
197  double fYMax; //< Bounding box max value y
198  mutable int fNumBadCircles; //< Number bad circles
200 
201  EventUtilities fUtilities; //< Handy functions to operate on arcs
202 
203 };
204 
205 } // namespace lar_cluster3d
206 #endif
std::pair< PointPair, PointPair > MinMaxPointPair
Definition: Voronoi.h:37
double getVoronoiDiagramArea() const
recover the area of the convex hull
Definition: Voronoi.h:73
This defines the actual beach line. The idea is to implement this as a self balancing binary search t...
Definition: BeachLine.h:129
bool computeCircleCenter2(const dcel2d::Coords &, const dcel2d::Coords &, const dcel2d::Coords &, dcel2d::Coords &, double &, double &) const
Definition: Voronoi.cxx:764
std::list< SiteEvent > SiteEventList
Definition: SweepEvent.h:99
dcel2d::HalfEdgeList & fHalfEdgeList
Definition: Voronoi.h:182
std::list< HalfEdge > HalfEdgeList
Definition: DCEL.h:180
std::list< CircleEvent > CircleEventList
Definition: SweepEvent.h:100
Provides some basic functions operating in IEvent class objects.
BSTNode class definiton specifically for use in constructing Voronoi diagrams. We are trying to follo...
Definition: BeachLine.h:32
bool isInsideConvexHull(const dcel2d::Vertex &) const
Is a vertex inside the convex hull - meant to be a fast check.
Definition: Voronoi.cxx:1042
double Area() const
Computes the area of the enclosed convext hull.
Definition: Voronoi.cxx:106
dcel2d::FaceList & fFaceList
Definition: Voronoi.h:184
std::list< BSTNode > BSTNodeList
Definition: BeachLine.h:122
Internal class definitions to facilitate construction of diagram.
Represents the beachline implemented as a self balancing binary search tree.
virtual ~VoronoiDiagram()
Destructor.
Definition: Voronoi.cxx:78
void makeLeftCircleEvent(EventQueue &, BSTNode *, double)
Definition: Voronoi.cxx:559
bool computeCircleCenter(const dcel2d::Coords &, const dcel2d::Coords &, const dcel2d::Coords &, dcel2d::Coords &, double &, double &) const
Definition: Voronoi.cxx:692
BSTNodeList fCircleNodeList
Definition: Voronoi.h:189
IEvent * makeCircleEvent(BSTNode *, BSTNode *, BSTNode *, double)
There are two types of events in the queue, here we handle circle events.
Definition: Voronoi.cxx:655
void handleCircleEvents(BeachLine &, EventQueue &, IEvent *)
There are two types of events in the queue, here we handle circle events.
Definition: Voronoi.cxx:474
VoronoiDiagram class definiton.
Definition: Voronoi.h:33
intermediate_table::const_iterator const_iterator
std::list< Face > FaceList
Definition: DCEL.h:179
CircleEventList fCircleEventList
Definition: Voronoi.h:188
SiteEventList fSiteEventList
Definition: Voronoi.h:187
EventUtilities fUtilities
Definition: Voronoi.h:201
double ComputeFaceArea()
Compute the area of the faces.
Definition: Voronoi.cxx:1155
std::priority_queue< IEvent *, std::vector< IEvent * >, bool(*)(const IEvent *, const IEvent *)> EventQueue
Definition: Voronoi.h:92
std::tuple< double, double, const reco::ClusterHit3D * > Point
Definitions used by the VoronoiDiagram algorithm.
Definition: DCEL.h:34
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:1357
void findBoundingBox(const dcel2d::VertexList &)
Find the min/max values in x-y to use as a bounding box.
Definition: Voronoi.cxx:1279
dcel2d::VertexList & fVertexList
Definition: Voronoi.h:183
void mergeDegenerateVertices()
merge degenerate vertices (found by zero length edges)
Definition: Voronoi.cxx:1127
dcel2d::Coords fConvexHullCenter
Definition: Voronoi.h:192
void makeRightCircleEvent(EventQueue &, BSTNode *, double)
Definition: Voronoi.cxx:608
dcel2d::PointList fPointList
Definition: Voronoi.h:186
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:93
PointPair findNearestEdge(const dcel2d::Point &, double &) const
Given an input Point, find the nearest edge.
Definition: Voronoi.cxx:1296
dcel2d::PointList fConvexHullList
Definition: Voronoi.h:191
std::pair< dcel2d::Point, dcel2d::Point > PointPair
Definition: Voronoi.h:36
const dcel2d::PointList & getConvexHull() const
recover the point list representing the convex hull
Definition: Voronoi.h:78
void boostTranslation(const dcel2d::PointList &, const boost::polygon::voronoi_edge< double > *, const boost::polygon::voronoi_edge< double > *, BoostEdgeToEdgeMap &, BoostVertexToVertexMap &, BoostCellToFaceMap &)
Definition: Voronoi.cxx:323
std::map< const boost::polygon::voronoi_edge< double > *, dcel2d::HalfEdge * > BoostEdgeToEdgeMap
Translate boost to dcel.
Definition: Voronoi.h:146
std::map< const boost::polygon::voronoi_vertex< double > *, dcel2d::Vertex * > BoostVertexToVertexMap
Definition: Voronoi.h:147
bool computeCircleCenter3(const dcel2d::Coords &, const dcel2d::Coords &, const dcel2d::Coords &, dcel2d::Coords &, double &, double &) const
Definition: Voronoi.cxx:795
void terminateInfiniteEdges(BeachLine &, double)
this aims to process remaining elements in the beachline after the event queue has been depleted ...
Definition: Voronoi.cxx:831
std::list< Point > PointList
Definition: DCEL.h:35
void buildVoronoiDiagramBoost(const dcel2d::PointList &)
Definition: Voronoi.cxx:267
Eigen::Vector3f Coords
Definition: DCEL.h:36
const dcel2d::VertexList & getVertexList() const
Recover the list of vertices.
Definition: Voronoi.h:68
const dcel2d::FaceList & getFaceList() const
Recover the list of faces.
Definition: Voronoi.h:63
std::list< Vertex > VertexList
Definition: DCEL.h:178
std::map< const boost::polygon::voronoi_cell< double > *, dcel2d::Face * > BoostCellToFaceMap
Definition: Voronoi.h:148
VoronoiDiagram(dcel2d::HalfEdgeList &, dcel2d::VertexList &, dcel2d::FaceList &)
Constructor.
Definition: Voronoi.cxx:53
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:1071
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:84
void handleSiteEvents(BeachLine &, EventQueue &, IEvent *)
There are two types of events in the queue, here we handle site events.
Definition: Voronoi.cxx:417