Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGraphDelaunay.h
Go to the documentation of this file.
1// @(#)root/hist:$Id: TGraphDelaunay.h,v 1.00
2// Author: Olivier Couet, Luke Jones (Royal Holloway, University of London)
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
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#ifndef ROOT_TGraphDelaunay
13#define ROOT_TGraphDelaunay
14
15
16//////////////////////////////////////////////////////////////////////////
17// //
18// TGraphDelaunay //
19// //
20// This class uses the Delaunay triangles technique to interpolate and //
21// render the data set. //
22// //
23//////////////////////////////////////////////////////////////////////////
24
25#include "TNamed.h"
26
27class TGraph2D;
28class TView;
29
30class TGraphDelaunay : public TNamed {
31
32private:
33
36
37protected:
38
39 Int_t fNdt; ///<! Number of Delaunay triangles found
40 Int_t fNpoints; ///<! Number of data points in fGraph2D
41 Int_t fNhull; ///<! Number of points in the hull
42 Double_t *fX; ///<! Pointer to fGraph2D->fX
43 Double_t *fY; ///<! Pointer to fGraph2D->fY
44 Double_t *fZ; ///<! Pointer to fGraph2D->fZ
45 Double_t *fXN; ///<! fGraph2D vectors normalized of size fNpoints
46 Double_t *fYN; ///<! fGraph2D vectors normalized of size fNpoints
47 Double_t fXNmin; ///<! Minimum value of fXN
48 Double_t fXNmax; ///<! Maximum value of fXN
49 Double_t fYNmin; ///<! Minimum value of fYN
50 Double_t fYNmax; ///<! Maximum value of fYN
52 Double_t fYoffset; ///<! Parameters used to normalize user data
55 Double_t fZout; ///<! Histogram bin height for points lying outside the convex hull
56 Double_t *fDist; ///<! Array used to order mass points by distance
57 Int_t fMaxIter; ///<! Maximum number of iterations to find Delaunay triangles
58 Int_t fTriedSize; ///<! Real size of the fxTried arrays
59 Int_t *fPTried; ///<!
60 Int_t *fNTried; ///<! Delaunay triangles storage of size fNdt
61 Int_t *fMTried; ///<!
62 Int_t *fHullPoints; ///<! Hull points of size fNhull
63 Int_t *fOrder; ///<! Array used to order mass points by distance
64 Bool_t fAllTri; ///<! True if FindAllTriangles() has been performed on fGraph2D
65 Bool_t fInit; ///<! True if CreateTrianglesDataStructure() and FindHull() have been performed
66 TGraph2D *fGraph2D; ///<! 2D graph containing the user data
67
69 Bool_t Enclose(Int_t T1, Int_t T2, Int_t T3, Int_t Ex) const;
70 void FileIt(Int_t P, Int_t N, Int_t M);
71 void FindHull();
72 Bool_t InHull(Int_t E, Int_t X) const;
73 Double_t InterpolateOnPlane(Int_t TI1, Int_t TI2, Int_t TI3, Int_t E) const;
74
75public:
76
79
80 ~TGraphDelaunay() override;
81
83 void FindAllTriangles();
84 TGraph2D *GetGraph2D() const {return fGraph2D;}
86 Int_t GetNdt() const {return fNdt;}
87 Int_t *GetPTried() const {return fPTried;}
88 Int_t *GetNTried() const {return fNTried;}
89 Int_t *GetMTried() const {return fMTried;}
90 Double_t *GetXN() const {return fXN;}
91 Double_t *GetYN() const {return fYN;}
92 Double_t GetXNmin() const {return fXNmin;}
93 Double_t GetXNmax() const {return fXNmax;}
94 Double_t GetYNmin() const {return fYNmin;}
95 Double_t GetYNmax() const {return fYNmax;}
97 void SetMaxIter(Int_t n=100000);
99
100 ClassDefOverride(TGraphDelaunay,1) // Delaunay triangulation
101};
102
103#endif
#define g(i)
Definition RSha256.hxx:105
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
#define X(type, name)
#define N
Graphics object made of three arrays X, Y and Z with the same number of points each.
Definition TGraph2D.h:41
TGraphDelaunay generates a Delaunay triangulation of a TGraph2D.
Int_t fNpoints
! Number of data points in fGraph2D
Double_t ComputeZ(Double_t x, Double_t y)
Return the z value corresponding to the (x,y) point in fGraph2D.
Int_t * fNTried
! Delaunay triangles storage of size fNdt
Int_t GetNdt() const
Double_t * fY
! Pointer to fGraph2D->fY
Double_t fXNmax
! Maximum value of fXN
Double_t GetYNmax() const
Bool_t InHull(Int_t E, Int_t X) const
Is point e inside the hull defined by all points apart from x ?
Int_t * GetMTried() const
Double_t fXNmin
! Minimum value of fXN
Double_t InterpolateOnPlane(Int_t TI1, Int_t TI2, Int_t TI3, Int_t E) const
Finds the z-value at point e given that it lies on the plane defined by t1,t2,t3.
Double_t * fYN
! fGraph2D vectors normalized of size fNpoints
Double_t fYoffset
! Parameters used to normalize user data
Double_t GetXNmin() const
Double_t GetMarginBinsContent() const
void FindHull()
Finds those points which make up the convex hull of the set.
void SetMarginBinsContent(Double_t z=0.)
Sets the histogram bin height for points lying outside the convex hull ie: the bins in the margin.
Double_t fXoffset
!
Double_t * fDist
! Array used to order mass points by distance
TGraph2D * GetGraph2D() const
void FileIt(Int_t P, Int_t N, Int_t M)
Files the triangle defined by the 3 vertices p, n and m into the fxTried arrays.
Double_t GetXNmax() const
Double_t fYNmin
! Minimum value of fYN
Int_t fTriedSize
! Real size of the fxTried arrays
TGraphDelaunay & operator=(const TGraphDelaunay &)=delete
Double_t * fX
! Pointer to fGraph2D->fX
TGraph2D * fGraph2D
! 2D graph containing the user data
void FindAllTriangles()
Attempt to find all the Delaunay triangles of the point set.
TGraphDelaunay(const TGraphDelaunay &)=delete
Double_t * GetXN() const
Bool_t fAllTri
! True if FindAllTriangles() has been performed on fGraph2D
Int_t * GetPTried() const
void SetMaxIter(Int_t n=100000)
Defines the number of triangles tested for a Delaunay triangle (number of iterations) before abandoni...
Int_t fMaxIter
! Maximum number of iterations to find Delaunay triangles
Int_t * GetNTried() const
Bool_t fInit
! True if CreateTrianglesDataStructure() and FindHull() have been performed
Double_t GetYNmin() const
Int_t * fOrder
! Array used to order mass points by distance
Bool_t Enclose(Int_t T1, Int_t T2, Int_t T3, Int_t Ex) const
Is point e inside the triangle t1-t2-t3 ?
Double_t fXScaleFactor
!
Int_t fNdt
! Number of Delaunay triangles found
Double_t * fZ
! Pointer to fGraph2D->fZ
Double_t fYNmax
! Maximum value of fYN
~TGraphDelaunay() override
TGraphDelaunay destructor.
Double_t fZout
! Histogram bin height for points lying outside the convex hull
Double_t fYScaleFactor
!
void CreateTrianglesDataStructure()
Function used internally only.
Int_t fNhull
! Number of points in the hull
TGraphDelaunay()
TGraphDelaunay default constructor.
Int_t * fHullPoints
! Hull points of size fNhull
Double_t * GetYN() const
Double_t * fXN
! fGraph2D vectors normalized of size fNpoints
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
See TView3D.
Definition TView.h:25
#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
#define T2
Definition md5.inl:147
#define T3
Definition md5.inl:148
#define T1
Definition md5.inl:146