Class to generate a Delaunay triangulation of a 2D set of points.
Algorithm based on Triangle, a two-dimensional quality mesh generator and Delaunay triangulator from Jonathan Richard Shewchuk.
See [http://www.cs.cmu.edu/~quake/triangle.html]
After having found the triangles using the above library, barycentric coordinates are used to test whether a point is inside a triangle (inTriangle test) and for interpolation. All this below is implemented in the DoInterpolateNormalized function.
Given triangle ABC and point P, P can be expressed by
P.x = la * A.x + lb * B.x + lc * C.x P.y = la * A.y + lb * B.y + lc * C.y
with lc = 1 - la - lb
P.x = la * A.x + lb * B.x + (1-la-lb) * C.x P.y = la * A.y + lb * B.y + (1-la-lb) * C.y
Rearranging yields
la * (A.x - C.x) + lb * (B.x - C.x) = P.x - C.x la * (A.y - C.y) + lb * (B.y - C.y) = P.y - C.y
Thus
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) ) 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) ) lc = 1 - la - lb
We save the inverse denominator to speedup computation
invDenom = 1 / ( (B.y - C.y)*(A.x - C.x) + (C.x - B.x)*(A.y - C.y) )
P is in triangle (including edges if
0 <= [la, lb, lc] <= 1
The interpolation of P.z is
P.z = la * A.z + lb * B.z + lc * C.z
To speed up localisation of points (to see to which triangle belong) a grid is laid over the internal coordinate space. A reference to triangle ABC is added to all grid cells that include ABC's bounding box. The size of the grid is defined to be 25x25
Optionally (if the compiler macro HAS_GCAL
is defined ) the triangle findings and interpolation can be computed using the GCAL library. This is however not supported when using the class within ROOT
Definition at line 115 of file Delaunay2D.h.
Classes | |
struct | Triangle |
Public Types | |
typedef std::vector< Triangle > | Triangles |
Public Member Functions | |
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 | |
Triangles::const_iterator | begin () const |
Triangles::const_iterator | end () const |
void | FindAllTriangles () |
Find all triangles. | |
double | Interpolate (double x, double y) |
Return the Interpolated z value corresponding to the given (x,y) point. | |
int | NumberOfTriangles () const |
return the number of triangles | |
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 | |
void | SetZOuterValue (double z=0.) |
set z value to be returned for points outside the region | |
double | XMax () const |
double | XMin () const |
double | YMax () const |
double | YMin () const |
double | ZOuterValue () const |
return the user defined Z-outer value | |
Protected Member Functions | |
unsigned int | Cell (unsigned int x, unsigned int y) const |
int | CellX (double x) const |
int | CellY (double y) const |
Protected Attributes | |
std::set< unsigned int > | fCells [(fNCells+1) *(fNCells+1)] |
! grid cells with containing triangles | |
bool | fInit |
! True if FindAllTriangles() has been performed | |
int | fNdt |
! Number of Delaunay triangles found | |
int | fNpoints |
! Number of data points | |
double | fOffsetX |
! Normalization offset X | |
double | fOffsetY |
! Normalization offset Y | |
double | fScaleFactorX |
! Normalization factor X | |
double | fScaleFactorY |
! Normalization factor Y | |
Triangles | fTriangles |
! Triangles of Triangulation | |
const double * | fX |
! Pointer to X array (managed externally) | |
double | fXCellStep |
! inverse denominator to calculate X cell = fNCells / (fXNmax - fXNmin) | |
std::vector< double > | fXN |
! normalized X | |
double | fXNmax |
! Maximum value of fXN | |
double | fXNmin |
! Minimum value of fXN | |
const double * | fY |
! Pointer to Y array | |
double | fYCellStep |
! inverse denominator to calculate X cell = fNCells / (fYNmax - fYNmin) | |
std::vector< double > | fYN |
! normalized Y | |
double | fYNmax |
! Maximum value of fYN | |
double | fYNmin |
! Minimum value of fYN | |
const double * | fZ |
! Pointer to Z array | |
double | fZout |
! Height for points lying outside the convex hull | |
Static Protected Attributes | |
static const int | fNCells = 25 |
! number of cells to divide the normalized space | |
Private Member Functions | |
Delaunay2D (const Delaunay2D &) | |
void | DoFindTriangles () |
internal function to find the triangle use Triangle or CGAL if flag is set | |
double | DoInterpolateNormalized (double x, double y) |
internal method to compute the interpolation | |
void | DoNormalizePoints () |
internal function to normalize the points. | |
double | Linear_transform (double x, double offset, double factor) |
Delaunay2D & | operator= (const Delaunay2D &) |
#include <Math/Delaunay2D.h>
typedef std::vector<Triangle> ROOT::Math::Delaunay2D::Triangles |
Definition at line 129 of file Delaunay2D.h.
ROOT::Math::Delaunay2D::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 at line 36 of file Delaunay2D.cxx.
|
private |
|
inline |
Definition at line 163 of file Delaunay2D.h.
Definition at line 233 of file Delaunay2D.h.
Definition at line 237 of file Delaunay2D.h.
Definition at line 241 of file Delaunay2D.h.
|
private |
internal function to find the triangle use Triangle or CGAL if flag is set
Triangle implementation for finding all the triangles.
Definition at line 177 of file Delaunay2D.cxx.
internal method to compute the interpolation
Triangle implementation for interpolation Finds the Delaunay triangle that the point (xi,yi) sits in (if any) and calculate a z-value for it by linearly interpolating the z-values that make up that triangle.
Relay that all the triangles have been found before see comment in class description (in Delaunay2D.h) for implementation details: finding barycentric coordinates and computing the interpolation
Definition at line 303 of file Delaunay2D.cxx.
|
private |
internal function to normalize the points.
Triangle implementation for points normalization.
See the class documentation for the details on how it is computed.
Definition at line 165 of file Delaunay2D.cxx.
|
inline |
Definition at line 164 of file Delaunay2D.h.
void ROOT::Math::Delaunay2D::FindAllTriangles | ( | ) |
Find all triangles.
Definition at line 136 of file Delaunay2D.cxx.
Return the Interpolated z value corresponding to the given (x,y) point.
Note that in case no Delaunay triangles are found, for example when the points are aligned, then a default value of zero is always return. See the class documentation for how the interpolation is computed.
Definition at line 106 of file Delaunay2D.cxx.
|
inlineprivate |
Definition at line 172 of file Delaunay2D.h.
|
inline |
return the number of triangles
Definition at line 149 of file Delaunay2D.h.
|
private |
void ROOT::Math::Delaunay2D::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
set the input points
Definition at line 61 of file Delaunay2D.cxx.
|
inline |
set z value to be returned for points outside the region
Definition at line 157 of file Delaunay2D.h.
|
inline |
Definition at line 152 of file Delaunay2D.h.
|
inline |
Definition at line 151 of file Delaunay2D.h.
|
inline |
Definition at line 154 of file Delaunay2D.h.
|
inline |
Definition at line 153 of file Delaunay2D.h.
|
inline |
return the user defined Z-outer value
Definition at line 160 of file Delaunay2D.h.
! grid cells with containing triangles
Definition at line 231 of file Delaunay2D.h.
|
protected |
! True if FindAllTriangles() has been performed
Definition at line 216 of file Delaunay2D.h.
|
staticprotected |
! number of cells to divide the normalized space
Definition at line 228 of file Delaunay2D.h.
|
protected |
! Number of Delaunay triangles found
Definition at line 196 of file Delaunay2D.h.
|
protected |
! Number of data points
Definition at line 197 of file Delaunay2D.h.
|
protected |
! Normalization offset X
Definition at line 208 of file Delaunay2D.h.
|
protected |
! Normalization offset Y
Definition at line 209 of file Delaunay2D.h.
|
protected |
! Normalization factor X
Definition at line 211 of file Delaunay2D.h.
|
protected |
! Normalization factor Y
Definition at line 212 of file Delaunay2D.h.
|
protected |
! Triangles of Triangulation
Definition at line 219 of file Delaunay2D.h.
|
protected |
! Pointer to X array (managed externally)
Definition at line 199 of file Delaunay2D.h.
|
protected |
! inverse denominator to calculate X cell = fNCells / (fXNmax - fXNmin)
Definition at line 229 of file Delaunay2D.h.
|
protected |
! normalized X
Definition at line 225 of file Delaunay2D.h.
|
protected |
! Maximum value of fXN
Definition at line 204 of file Delaunay2D.h.
|
protected |
! Minimum value of fXN
Definition at line 203 of file Delaunay2D.h.
|
protected |
! Pointer to Y array
Definition at line 200 of file Delaunay2D.h.
|
protected |
! inverse denominator to calculate X cell = fNCells / (fYNmax - fYNmin)
Definition at line 230 of file Delaunay2D.h.
|
protected |
! normalized Y
Definition at line 226 of file Delaunay2D.h.
|
protected |
! Maximum value of fYN
Definition at line 206 of file Delaunay2D.h.
|
protected |
! Minimum value of fYN
Definition at line 205 of file Delaunay2D.h.
|
protected |
! Pointer to Z array
Definition at line 201 of file Delaunay2D.h.
|
protected |
! Height for points lying outside the convex hull
Definition at line 214 of file Delaunay2D.h.