Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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//#include "RtypesCore.h"
25
26#include <map>
27#include <vector>
28#include <set>
29#include <functional>
30
31#ifdef HAS_CGAL
32 /* CGAL uses the name PTR as member name in its Handle class
33 * but its a macro defined in mmalloc.h of ROOT
34 * Safe it, disable it and then re-enable it later on*/
35 #pragma push_macro("PTR")
36 #undef PTR
37
38 #include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
39 #include <CGAL/Delaunay_triangulation_2.h>
40 #include <CGAL/Triangulation_vertex_base_with_info_2.h>
41 #include <CGAL/Interpolation_traits_2.h>
42 #include <CGAL/natural_neighbor_coordinates_2.h>
43 #include <CGAL/interpolation_functions.h>
44
45 #pragma pop_macro("PTR")
46#endif
47
48#ifdef THREAD_SAFE
49 #include<atomic> //atomic operations for thread safety
50#endif
51
52
53namespace ROOT {
54
55
56
57 namespace Math {
58
59/**
60
61 Class to generate a Delaunay triangulation of a 2D set of points.
62 Algorithm based on **Triangle**, a two-dimensional quality mesh generator and
63 Delaunay triangulator from Jonathan Richard Shewchuk.
64
65 See [http://www.cs.cmu.edu/~quake/triangle.html]
66
67 After having found the triangles using the above library, barycentric coordinates are used
68 to test whether a point is inside a triangle (inTriangle test) and for interpolation.
69 All this below is implemented in the DoInterpolateNormalized function.
70
71 Given triangle ABC and point P, P can be expressed by
72
73 P.x = la * A.x + lb * B.x + lc * C.x
74 P.y = la * A.y + lb * B.y + lc * C.y
75
76 with lc = 1 - la - lb
77
78 P.x = la * A.x + lb * B.x + (1-la-lb) * C.x
79 P.y = la * A.y + lb * B.y + (1-la-lb) * C.y
80
81 Rearranging yields
82
83 la * (A.x - C.x) + lb * (B.x - C.x) = P.x - C.x
84 la * (A.y - C.y) + lb * (B.y - C.y) = P.y - C.y
85
86 Thus
87
88 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) )
89 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) )
90 lc = 1 - la - lb
91
92 We save the inverse denominator to speedup computation
93
94 invDenom = 1 / ( (B.y - C.y)*(A.x - C.x) + (C.x - B.x)*(A.y - C.y) )
95
96 P is in triangle (including edges if
97
98 0 <= [la, lb, lc] <= 1
99
100 The interpolation of P.z is
101
102 P.z = la * A.z + lb * B.z + lc * C.z
103
104 To speed up localisation of points (to see to which triangle belong) a grid is laid over the internal coordinate space.
105 A reference to triangle ABC is added to _all_ grid cells that include ABC's bounding box.
106 The size of the grid is defined to be 25x25
107
108 Optionally (if the compiler macro `HAS_GCAL` is defined ) the triangle findings and interpolation can be computed
109 using the GCAL library. This is however not supported when using the class within ROOT
110
111 \ingroup MathCore
112 */
113
114
116
117public:
118
119 struct Triangle {
120 double x[3]; // x of triangle vertices
121 double y[3]; // y of triangle vertices
122 unsigned int idx[3]; // point corresponding to vertices
123
124 #ifndef HAS_CGAL
125 double invDenom; // cached inv denominator for computing barycentric coordinates (see above)
126 #endif
127 };
128
129 typedef std::vector<Triangle> Triangles;
130
131public:
132
133
134 Delaunay2D(int n, const double *x, const double * y, const double * z, double xmin=0, double xmax=0, double ymin=0, double ymax=0);
135
136 /// set the input points for building the graph
137 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);
138
139 /// Return the Interpolated z value corresponding to the given (x,y) point.
140 /// Note that in case no Delaunay triangles are found, for example when the
141 /// points are aligned, then a default value of zero is always return.
142 /// See the class documentation for how the interpolation is computed.
143 double Interpolate(double x, double y);
144
145 /// Find all triangles
146 void FindAllTriangles();
147
148 /// return the number of triangles
149 int NumberOfTriangles() const {return fNdt;}
150
151 double XMin() const {return fXNmin;}
152 double XMax() const {return fXNmax;}
153 double YMin() const {return fYNmin;}
154 double YMax() const {return fYNmax;}
155
156 /// set z value to be returned for points outside the region
157 void SetZOuterValue(double z=0.) { fZout = z; }
158
159 /// return the user defined Z-outer value
160 double ZOuterValue() const { return fZout; }
161
162 // iterators on the found triangles
163 Triangles::const_iterator begin() const { return fTriangles.begin(); }
164 Triangles::const_iterator end() const { return fTriangles.end(); }
165
166
167private:
168
169 // internal methods
170
171
172 inline double Linear_transform(double x, double offset, double factor){
173 return (x+offset)*factor;
174 }
175
176 /// internal function to normalize the points.
177 /// See the class documentation for the details on how it is computed.
178 void DoNormalizePoints();
179
180 /// internal function to find the triangle
181 /// use Triangle or CGAL if flag is set
182 void DoFindTriangles();
183
184 /// internal method to compute the interpolation
185 double DoInterpolateNormalized(double x, double y);
186
187
188
189private:
190 // class is not copyable
191 Delaunay2D(const Delaunay2D&); // Not implemented
192 Delaunay2D& operator=(const Delaunay2D&); // Not implemented
193
194protected:
195
196 int fNdt; ///<! Number of Delaunay triangles found
197 int fNpoints; ///<! Number of data points
198
199 const double *fX; ///<! Pointer to X array (managed externally)
200 const double *fY; ///<! Pointer to Y array
201 const double *fZ; ///<! Pointer to Z array
202
203 double fXNmin; ///<! Minimum value of fXN
204 double fXNmax; ///<! Maximum value of fXN
205 double fYNmin; ///<! Minimum value of fYN
206 double fYNmax; ///<! Maximum value of fYN
207
208 double fOffsetX; ///<! Normalization offset X
209 double fOffsetY; ///<! Normalization offset Y
210
211 double fScaleFactorX; ///<! Normalization factor X
212 double fScaleFactorY; ///<! Normalization factor Y
213
214 double fZout; ///<! Height for points lying outside the convex hull
215
216 bool fInit; ///<! True if FindAllTriangles() has been performed
217
218
219 Triangles fTriangles; ///<! Triangles of Triangulation
220
221#ifndef HAS_CGAL
222
223 //using triangle library
224
225 std::vector<double> fXN; ///<! normalized X
226 std::vector<double> fYN; ///<! normalized Y
227
228 static const int fNCells = 25; ///<! number of cells to divide the normalized space
229 double fXCellStep; ///<! inverse denominator to calculate X cell = fNCells / (fXNmax - fXNmin)
230 double fYCellStep; ///<! inverse denominator to calculate X cell = fNCells / (fYNmax - fYNmin)
231 std::set<unsigned int> fCells[(fNCells+1)*(fNCells+1)]; ///<! grid cells with containing triangles
232
233 inline unsigned int Cell(unsigned int x, unsigned int y) const {
234 return x*(fNCells+1) + y;
235 }
236
237 inline int CellX(double x) const {
238 return (x - fXNmin) * fXCellStep;
239 }
240
241 inline int CellY(double y) const {
242 return (y - fYNmin) * fYCellStep;
243 }
244
245#else // HAS_CGAL
246 // case of using GCAL
247 //Functor class for accessing the function values/gradients
248 template< class PointWithInfoMap, typename ValueType >
249 struct Data_access : public std::unary_function< typename PointWithInfoMap::key_type,
250 std::pair<ValueType, bool> >
251 {
252
253 Data_access(const PointWithInfoMap& points, const ValueType * values)
254 : _points(points), _values(values){};
255
256 std::pair< ValueType, bool>
257 operator()(const typename PointWithInfoMap::key_type& p) const {
258 typename PointWithInfoMap::const_iterator mit = _points.find(p);
259 if(mit!= _points.end())
260 return std::make_pair(_values[mit->second], true);
261 return std::make_pair(ValueType(), false);
262 };
263
264 const PointWithInfoMap& _points;
265 const ValueType * _values;
266 };
267
268 typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
269 typedef CGAL::Triangulation_vertex_base_with_info_2<uint, K> Vb;
270 typedef CGAL::Triangulation_data_structure_2<Vb> Tds;
271 typedef CGAL::Delaunay_triangulation_2<K, Tds> Delaunay;
272 typedef CGAL::Interpolation_traits_2<K> Traits;
273 typedef K::FT Coord_type;
274 typedef K::Point_2 Point;
275 typedef std::map<Point, Vb::Info, K::Less_xy_2> PointWithInfoMap;
276 typedef Data_access< PointWithInfoMap, double > Value_access;
277
278 Delaunay fCGALdelaunay; //! CGAL delaunay triangulation object
279 PointWithInfoMap fNormalizedPoints; //! Normalized function values
280
281#endif //HAS_CGAL
282
283
284};
285
286
287} // namespace Math
288} // namespace ROOT
289
290
291#endif
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h offset
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
float xmin
float ymin
float xmax
float ymax
TRObject operator()(const T1 &t1) const
Class to generate a Delaunay triangulation of a 2D set of points.
Definition Delaunay2D.h:115
Triangles::const_iterator end() const
Definition Delaunay2D.h:164
int fNpoints
! Number of data points
Definition Delaunay2D.h:197
double fXNmin
! Minimum value of fXN
Definition Delaunay2D.h:203
std::set< unsigned int > fCells[(fNCells+1) *(fNCells+1)]
! grid cells with containing triangles
Definition Delaunay2D.h:231
int fNdt
! Number of Delaunay triangles found
Definition Delaunay2D.h:196
double DoInterpolateNormalized(double x, double y)
internal method to compute the interpolation
int CellY(double y) const
Definition Delaunay2D.h:241
Triangles::const_iterator begin() const
Definition Delaunay2D.h:163
Delaunay2D(const Delaunay2D &)
const double * fZ
! Pointer to Z array
Definition Delaunay2D.h:201
double fXNmax
! Maximum value of fXN
Definition Delaunay2D.h:204
void SetZOuterValue(double z=0.)
set z value to be returned for points outside the region
Definition Delaunay2D.h:157
double fZout
! Height for points lying outside the convex hull
Definition Delaunay2D.h:214
double YMin() const
Definition Delaunay2D.h:153
double fOffsetY
! Normalization offset Y
Definition Delaunay2D.h:209
Delaunay2D & operator=(const Delaunay2D &)
double XMax() const
Definition Delaunay2D.h:152
double ZOuterValue() const
return the user defined Z-outer value
Definition Delaunay2D.h:160
Triangles fTriangles
! Triangles of Triangulation
Definition Delaunay2D.h:219
double Linear_transform(double x, double offset, double factor)
Definition Delaunay2D.h:172
double fYCellStep
! inverse denominator to calculate X cell = fNCells / (fYNmax - fYNmin)
Definition Delaunay2D.h:230
std::vector< double > fYN
! normalized Y
Definition Delaunay2D.h:226
void FindAllTriangles()
Find all triangles.
static const int fNCells
! number of cells to divide the normalized space
Definition Delaunay2D.h:228
double YMax() const
Definition Delaunay2D.h:154
bool fInit
! True if FindAllTriangles() has been performed
Definition Delaunay2D.h:216
unsigned int Cell(unsigned int x, unsigned int y) const
Definition Delaunay2D.h:233
std::vector< double > fXN
! normalized X
Definition Delaunay2D.h:225
double fXCellStep
! inverse denominator to calculate X cell = fNCells / (fXNmax - fXNmin)
Definition Delaunay2D.h:229
double fYNmin
! Minimum value of fYN
Definition Delaunay2D.h:205
double fOffsetX
! Normalization offset X
Definition Delaunay2D.h:208
double XMin() const
Definition Delaunay2D.h:151
std::vector< Triangle > Triangles
Definition Delaunay2D.h:129
const double * fY
! Pointer to Y array
Definition Delaunay2D.h:200
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
const double * fX
! Pointer to X array (managed externally)
Definition Delaunay2D.h:199
void DoNormalizePoints()
internal function to normalize the points.
double fYNmax
! Maximum value of fYN
Definition Delaunay2D.h:206
int NumberOfTriangles() const
return the number of triangles
Definition Delaunay2D.h:149
double fScaleFactorY
! Normalization factor Y
Definition Delaunay2D.h:212
double fScaleFactorX
! Normalization factor X
Definition Delaunay2D.h:211
void DoFindTriangles()
internal function to find the triangle use Triangle or CGAL if flag is set
int CellX(double x) const
Definition Delaunay2D.h:237
#define Interpolate(a, x, b, y)
Definition geom.c:179
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
Namespace for new Math classes and functions.
This file contains a specialised ROOT message handler to test for diagnostic in unit tests.
constexpr Double_t K()
Boltzmann's constant in : .
Definition TMath.h:247