ROOT  6.06/09
Reference Guide
Delaunay2D.h
Go to the documentation of this file.
1 // @(#)root/mathcore:$Id: Delaunay2D.h,v 1.00
2 // Author: Daniel Funke, Lorenzo Moneta
3 
4 /*************************************************************************
5  * Copyright (C) 2015 ROOT Math Team *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 // Header file for class Delaunay2D
13 
14 #ifndef ROOT_Math_Delaunay2D
15 #define ROOT_Math_Delaunay2D
16 
17 //for testing purposes HAS_CGAL can be [un]defined here
18 //#define HAS_CGAL
19 
20 //for testing purposes THREAD_SAFE can [un]defined here
21 //#define THREAD_SAFE
22 
23 
24 #ifndef ROOT_TNamed
25 #include "TNamed.h"
26 #endif
27 
28 #include <map>
29 #include <vector>
30 #include <set>
31 #include <functional>
32 
33 #ifdef HAS_CGAL
34  /* CGAL uses the name PTR as member name in its Handle class
35  * but its a macro defined in mmalloc.h of ROOT
36  * Safe it, disable it and then re-enable it later on*/
37  #pragma push_macro("PTR")
38  #undef PTR
39 
40  #include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
41  #include <CGAL/Delaunay_triangulation_2.h>
42  #include <CGAL/Triangulation_vertex_base_with_info_2.h>
43  #include <CGAL/Interpolation_traits_2.h>
44  #include <CGAL/natural_neighbor_coordinates_2.h>
45  #include <CGAL/interpolation_functions.h>
46 
47  #pragma pop_macro("PTR")
48 #endif
49 
50 #ifdef THREAD_SAFE
51  #include<atomic> //atomic operations for thread safety
52 #endif
53 
54 
55 namespace ROOT {
56 
57 
58 
59  namespace Math {
60 
61 /**
62 
63  Class to generate a Delaunay triangulation of a 2D set of points.
64  Algorithm based on **Triangle**, a two-dimensional quality mesh generator and
65  Delaunay triangulator from Jonathan Richard Shewchuk.
66 
67  See [http://www.cs.cmu.edu/~quake/triangle.html]
68 
69  \ingroup MathCore
70  */
71 
72 
73 class Delaunay2D {
74 
75 public:
76 
77  struct Triangle {
78  double x[3];
79  double y[3];
80  UInt_t idx[3];
81 
82  #ifndef HAS_CGAL
83  double invDenom; //see comment below in CGAL fall back section
84  #endif
85  };
86 
87  typedef std::vector<Triangle> Triangles;
88 
89 public:
90 
91 
92  Delaunay2D(int n, const double *x, const double * y, const double * z, double xmin=0, double xmax=0, double ymin=0, double ymax=0);
93 
94  /// set the input points for building the graph
95  void SetInputPoints(int n, const double *x, const double * y, const double * z, double xmin=0, double xmax=0, double ymin=0, double ymax=0);
96 
97  /// Return the Interpolated z value corresponding to the (x,y) point
98  double Interpolate(double x, double y);
99 
100  /// Find all triangles
101  void FindAllTriangles();
102 
103  /// return the number of triangles
104  Int_t NumberOfTriangles() const {return fNdt;}
105 
106  double XMin() const {return fXNmin;}
107  double XMax() const {return fXNmax;}
108  double YMin() const {return fYNmin;}
109  double YMax() const {return fYNmax;}
110 
111  /// set z value to be returned for points outside the region
112  void SetZOuterValue(double z=0.) { fZout = z; }
113 
114  /// return the user defined Z-outer value
115  double ZOuterValue() const { return fZout; }
116 
117  // iterators on the found triangles
118  Triangles::const_iterator begin() const { return fTriangles.begin(); }
119  Triangles::const_iterator end() const { return fTriangles.end(); }
120 
121 
122 private:
123 
124  // internal methods
125 
126 
127  inline double Linear_transform(double x, double offset, double factor){
128  return (x+offset)*factor;
129  }
130 
131  /// internal function to normalize the points
132  void DoNormalizePoints();
133 
134  /// internal function to find the triangle
135  /// use Triangle or CGAL if flag is set
136  void DoFindTriangles();
137 
138  /// internal method to compute the interpolation
139  double DoInterpolateNormalized(double x, double y);
140 
141 
142 
143 private:
144  // class is not copyable
145  Delaunay2D(const Delaunay2D&); // Not implemented
146  Delaunay2D& operator=(const Delaunay2D&); // Not implemented
147 
148 protected:
149 
150  //typedef std::function<double(double)> Transformer;
151 
152  Int_t fNdt; //!Number of Delaunay triangles found
153  Int_t fNpoints; //!Number of data points
154 
155  const double *fX; //!Pointer to X array (managed externally)
156  const double *fY; //!Pointer to Y array
157  const double *fZ; //!Pointer to Z array
158 
159  double fXNmin; //!Minimum value of fXN
160  double fXNmax; //!Maximum value of fXN
161  double fYNmin; //!Minimum value of fYN
162  double fYNmax; //!Maximum value of fYN
163 
164  //Transformer xTransformer; //!transform x values to mapped space
165  //Transformer yTransformer; //!transform y values to mapped space
166 
167  double fOffsetX; //!Normalization offset X
168  double fOffsetY; //!Normalization offset Y
169 
170  double fScaleFactorX; //!Normalization factor X
171  double fScaleFactorY; //!Normalization factor Y
172 
173  double fZout; //!Height for points lying outside the convex hull
174 
175 #ifdef THREAD_SAFE
176 
177  enum class Initialization : char {UNINITIALIZED, INITIALIZING, INITIALIZED};
178  std::atomic<Initialization> fInit; //!Indicate initialization state
179 
180 #else
181  Bool_t fInit; //!True if FindAllTriangels() has been performed
182 #endif
183 
184 
185  Triangles fTriangles; //!Triangles of Triangulation
186 
187 #ifdef HAS_CGAL
188 
189  //Functor class for accessing the function values/gradients
190  template< class PointWithInfoMap, typename ValueType >
191  struct Data_access : public std::unary_function< typename PointWithInfoMap::key_type,
192  std::pair<ValueType, bool> >
193  {
194 
195  Data_access(const PointWithInfoMap& points, const ValueType * values)
196  : _points(points), _values(values){};
197 
198  std::pair< ValueType, bool>
199  operator()(const typename PointWithInfoMap::key_type& p) const {
200  typename PointWithInfoMap::const_iterator mit = _points.find(p);
201  if(mit!= _points.end())
202  return std::make_pair(_values[mit->second], true);
203  return std::make_pair(ValueType(), false);
204  };
205 
206  const PointWithInfoMap& _points;
207  const ValueType * _values;
208  };
209 
210  typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
211  typedef CGAL::Triangulation_vertex_base_with_info_2<uint, K> Vb;
212  typedef CGAL::Triangulation_data_structure_2<Vb> Tds;
213  typedef CGAL::Delaunay_triangulation_2<K, Tds> Delaunay;
214  typedef CGAL::Interpolation_traits_2<K> Traits;
215  typedef K::FT Coord_type;
216  typedef K::Point_2 Point;
217  typedef std::map<Point, Vb::Info, K::Less_xy_2> PointWithInfoMap;
218  typedef Data_access< PointWithInfoMap, double > Value_access;
219 
220  Delaunay fCGALdelaunay; //! CGAL delaunay triangulation object
221  PointWithInfoMap fNormalizedPoints; //! Normalized function values
222 
223 #else // HAS_CGAL
224  //fallback to triangle library
225 
226  /* Using barycentric coordinates for inTriangle test and interpolation
227  *
228  * Given triangle ABC and point P, P can be expressed by
229  *
230  * P.x = la * A.x + lb * B.x + lc * C.x
231  * P.y = la * A.y + lb * B.y + lc * C.y
232  *
233  * with lc = 1 - la - lb
234  *
235  * P.x = la * A.x + lb * B.x + (1-la-lb) * C.x
236  * P.y = la * A.y + lb * B.y + (1-la-lb) * C.y
237  *
238  * Rearranging yields
239  *
240  * la * (A.x - C.x) + lb * (B.x - C.x) = P.x - C.x
241  * la * (A.y - C.y) + lb * (B.y - C.y) = P.y - C.y
242  *
243  * Thus
244  *
245  * la = ( (B.y - C.y)*(P.x - C.x) + (C.x - B.x)*(P.y - C.y) ) / ( (B.y - C.y)*(A.x - C.x) + (C.x - B.x)*(A.y - C.y) )
246  * lb = ( (C.y - A.y)*(P.x - C.x) + (A.x - C.x)*(P.y - C.y) ) / ( (B.y - C.y)*(A.x - C.x) + (C.x - B.x)*(A.y - C.y) )
247  * lc = 1 - la - lb
248  *
249  * We save the inverse denominator to speedup computation
250  *
251  * invDenom = 1 / ( (B.y - C.y)*(A.x - C.x) + (C.x - B.x)*(A.y - C.y) )
252  *
253  * P is in triangle (including edges if
254  *
255  * 0 <= [la, lb, lc] <= 1
256  *
257  * The interpolation of P.z is
258  *
259  * P.z = la * A.z + lb * B.z + lc * C.z
260  *
261  */
262 
263  std::vector<double> fXN; //! normalized X
264  std::vector<double> fYN; //! normalized Y
265 
266  /* To speed up localisation of points a grid is layed over normalized space
267  *
268  * A reference to triangle ABC is added to _all_ grid cells that include ABC's bounding box
269  */
270 
271  static const int fNCells = 25; //! number of cells to divide the normalized space
272  double fXCellStep; //! inverse denominator to calculate X cell = fNCells / (fXNmax - fXNmin)
273  double fYCellStep; //! inverse denominator to calculate X cell = fNCells / (fYNmax - fYNmin)
274  std::set<UInt_t> fCells[(fNCells+1)*(fNCells+1)]; //! grid cells with containing triangles
275 
276  inline unsigned int Cell(UInt_t x, UInt_t y) const {
277  return x*(fNCells+1) + y;
278  }
279 
280  inline int CellX(double x) const {
281  return (x - fXNmin) * fXCellStep;
282  }
283 
284  inline int CellY(double y) const {
285  return (y - fYNmin) * fYCellStep;
286  }
287 
288 #endif //HAS_CGAL
289 
290 
291 };
292 
293 
294 } // namespace Math
295 } // namespace ROOT
296 
297 
298 #endif
unsigned int Cell(UInt_t x, UInt_t y) const
grid cells with containing triangles
Definition: Delaunay2D.h:276
double fZout
Normalization factor Y.
Definition: Delaunay2D.h:173
void DoNormalizePoints()
internal function to normalize the points
Definition: Delaunay2D.cxx:241
float xmin
Definition: THbookFile.cxx:93
std::vector< Triangle > Triangles
Definition: Delaunay2D.h:87
Namespace for new ROOT classes and functions.
Definition: ROOT.py:1
float ymin
Definition: THbookFile.cxx:93
Bool_t fInit
Height for points lying outside the convex hull.
Definition: Delaunay2D.h:181
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
double XMax() const
Definition: Delaunay2D.h:107
void FindAllTriangles()
Find all triangles.
Definition: Delaunay2D.cxx:125
Triangles::const_iterator end() const
Definition: Delaunay2D.h:119
void SetZOuterValue(double z=0.)
set z value to be returned for points outside the region
Definition: Delaunay2D.h:112
double fScaleFactorX
Normalization offset Y.
Definition: Delaunay2D.h:170
double DoInterpolateNormalized(double x, double y)
internal method to compute the interpolation
Definition: Delaunay2D.cxx:376
Double_t x[n]
Definition: legend1.C:17
double YMin() const
Definition: Delaunay2D.h:108
double fOffsetY
Normalization offset X.
Definition: Delaunay2D.h:168
double fScaleFactorY
Normalization factor X.
Definition: Delaunay2D.h:171
double fOffsetX
Maximum value of fYN.
Definition: Delaunay2D.h:167
void SetInputPoints(int n, const double *x, const double *y, const double *z, double xmin=0, double xmax=0, double ymin=0, double ymax=0)
set the input points for building the graph
Definition: Delaunay2D.cxx:56
std::set< UInt_t > fCells[(fNCells+1)*(fNCells+1)]
inverse denominator to calculate X cell = fNCells / (fYNmax - fYNmin)
Definition: Delaunay2D.h:274
double fXNmin
Pointer to Z array.
Definition: Delaunay2D.h:159
point * points
Definition: X3DBuffer.c:20
float ymax
Definition: THbookFile.cxx:93
const double * fY
Pointer to X array (managed externally)
Definition: Delaunay2D.h:156
int CellX(double x) const
Definition: Delaunay2D.h:280
double fXNmax
Minimum value of fXN.
Definition: Delaunay2D.h:160
unsigned int UInt_t
Definition: RtypesCore.h:42
Double_t K()
Definition: TMath.h:95
Delaunay2D & operator=(const Delaunay2D &)
std::vector< double > fYN
normalized X
Definition: Delaunay2D.h:264
float xmax
Definition: THbookFile.cxx:93
double fYNmax
Minimum value of fYN.
Definition: Delaunay2D.h:162
double fXCellStep
number of cells to divide the normalized space
Definition: Delaunay2D.h:272
TRObject operator()(const T1 &t1) const
Int_t fNpoints
Number of Delaunay triangles found.
Definition: Delaunay2D.h:153
double fYCellStep
inverse denominator to calculate X cell = fNCells / (fXNmax - fXNmin)
Definition: Delaunay2D.h:273
Double_t y[n]
Definition: legend1.C:17
Triangles::const_iterator begin() const
Definition: Delaunay2D.h:118
Namespace for new Math classes and functions.
Int_t NumberOfTriangles() const
return the number of triangles
Definition: Delaunay2D.h:104
double YMax() const
Definition: Delaunay2D.h:109
static const Int_t UNINITIALIZED
static const int fNCells
normalized Y
Definition: Delaunay2D.h:271
void DoFindTriangles()
internal function to find the triangle use Triangle or CGAL if flag is set
Definition: Delaunay2D.cxx:253
std::vector< double > fXN
Triangles of Triangulation.
Definition: Delaunay2D.h:263
const double * fX
Number of data points.
Definition: Delaunay2D.h:155
Class to generate a Delaunay triangulation of a 2D set of points.
Definition: Delaunay2D.h:73
double fYNmin
Maximum value of fXN.
Definition: Delaunay2D.h:161
int CellY(double y) const
Definition: Delaunay2D.h:284
double Interpolate(double x, double y)
Return the Interpolated z value corresponding to the (x,y) point.
Definition: Delaunay2D.cxx:101
const Int_t n
Definition: legend1.C:16
double ZOuterValue() const
return the user defined Z-outer value
Definition: Delaunay2D.h:115
Delaunay2D(int n, const double *x, const double *y, const double *z, double xmin=0, double xmax=0, double ymin=0, double ymax=0)
class constructor from array of data points
Definition: Delaunay2D.cxx:32
double XMin() const
Definition: Delaunay2D.h:106
double Linear_transform(double x, double offset, double factor)
Definition: Delaunay2D.h:127
const double * fZ
Pointer to Y array.
Definition: Delaunay2D.h:157
Triangles fTriangles
True if FindAllTriangels() has been performed.
Definition: Delaunay2D.h:185