Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ContourList.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_hist
3/// \notebook
4/// Getting Contours From TH2D.
5///
6/// #### Image produced by `.x ContourList.C`
7/// The contours values are drawn next to each contour.
8/// \macro_image
9///
10/// #### Output produced by `.x ContourList.C`
11/// It shows that 6 contours and 12 graphs were found.
12/// \macro_output
13///
14/// #### `ContourList.C`
15/// \macro_code
16///
17/// \authors Josh de Bever (CSI Medical Physics Group, The University of Western Ontario, London, Ontario, Canada), Olivier Couet
18
19Double_t SawTooth(Double_t x, Double_t WaveLen);
20
21TCanvas *ContourList(){
22
23 const Double_t PI = TMath::Pi();
24
25 TCanvas* c = new TCanvas("c","Contour List",0,0,600,600);
26 c->SetRightMargin(0.15);
27 c->SetTopMargin(0.15);
28
29 Int_t i, j;
30
31 Int_t nZsamples = 80;
32 Int_t nPhiSamples = 80;
33
34 Double_t HofZwavelength = 4.0; // 4 meters
35 Double_t dZ = HofZwavelength/(Double_t)(nZsamples - 1);
36 Double_t dPhi = 2*PI/(Double_t)(nPhiSamples - 1);
37
38 TArrayD z(nZsamples);
39 TArrayD HofZ(nZsamples);
40 TArrayD phi(nPhiSamples);
41 TArrayD FofPhi(nPhiSamples);
42
43 // Discretized Z and Phi Values
44 for ( i = 0; i < nZsamples; i++) {
45 z[i] = (i)*dZ - HofZwavelength/2.0;
46 HofZ[i] = SawTooth(z[i], HofZwavelength);
47 }
48
49 for(Int_t i=0; i < nPhiSamples; i++){
50 phi[i] = (i)*dPhi;
51 FofPhi[i] = sin(phi[i]);
52 }
53
54 // Create Histogram
55 TH2D *HistStreamFn = new TH2D("HstreamFn",
56 "#splitline{Histogram with negative and positive contents. Six contours are defined.}{It is plotted with options CONT LIST to retrieve the contours points in TGraphs}",
57 nZsamples, z[0], z[nZsamples-1], nPhiSamples, phi[0], phi[nPhiSamples-1]);
58
59 // Load Histogram Data
60 for (Int_t i = 0; i < nZsamples; i++) {
61 for(Int_t j = 0; j < nPhiSamples; j++){
62 HistStreamFn->SetBinContent(i,j, HofZ[i]*FofPhi[j]);
63 }
64 }
65
67 gStyle->SetTitleW(0.99);
68 gStyle->SetTitleH(0.08);
69
70 Double_t contours[6];
71 contours[0] = -0.7;
72 contours[1] = -0.5;
73 contours[2] = -0.1;
74 contours[3] = 0.1;
75 contours[4] = 0.4;
76 contours[5] = 0.8;
77
78 HistStreamFn->SetContour(6, contours);
79
80 // Draw contours as filled regions, and Save points
81 HistStreamFn->Draw("CONT Z LIST");
82 c->Update(); // Needed to force the plotting and retrieve the contours in TGraphs
83
84 // Get Contours
85 TObjArray *conts = (TObjArray*)gROOT->GetListOfSpecials()->FindObject("contours");
86
87 if (!conts){
88 printf("*** No Contours Were Extracted!\n");
89 return nullptr;
90 }
91
92 TList* contLevel = nullptr;
93 TGraph* curv = nullptr;
94 TGraph* gc = nullptr;
95
96 Int_t nGraphs = 0;
97 Int_t TotalConts = conts->GetSize();
98
99 printf("TotalConts = %d\n", TotalConts);
100
101 for(i = 0; i < TotalConts; i++){
102 contLevel = (TList*)conts->At(i);
103 printf("Contour %d has %d Graphs\n", i, contLevel->GetSize());
104 nGraphs += contLevel->GetSize();
105 }
106
107 nGraphs = 0;
108
109 TCanvas* c1 = new TCanvas("c1","Contour List",610,0,600,600);
110 c1->SetTopMargin(0.15);
111 TH2F *hr = new TH2F("hr",
112 "#splitline{Negative contours are returned first (highest to lowest). Positive contours are returned from}{lowest to highest. On this plot Negative contours are drawn in red and positive contours in blue.}",
113 2, -2, 2, 2, 0, 6.5);
114
115 hr->Draw();
116 Double_t xval0, yval0, zval0;
117 TLatex l;
118 l.SetTextSize(0.03);
119 char val[20];
120
121 for(i = 0; i < TotalConts; i++){
122 contLevel = (TList*)conts->At(i);
123 if (i<3) zval0 = contours[2-i];
124 else zval0 = contours[i];
125 printf("Z-Level Passed in as: Z = %f\n", zval0);
126
127 // Get first graph from list on curves on this level
128 curv = (TGraph*)contLevel->First();
129 for(j = 0; j < contLevel->GetSize(); j++){
130 curv->GetPoint(0, xval0, yval0);
131 if (zval0<0) curv->SetLineColor(kRed);
132 if (zval0>0) curv->SetLineColor(kBlue);
133 nGraphs ++;
134 printf("\tGraph: %d -- %d Elements\n", nGraphs,curv->GetN());
135
136 // Draw clones of the graphs to avoid deletions in case the 1st
137 // pad is redrawn.
138 gc = (TGraph*)curv->Clone();
139 gc->Draw("C");
140
141 sprintf(val,"%g",zval0);
142 l.DrawLatex(xval0,yval0,val);
143 curv = (TGraph*)contLevel->After(curv); // Get Next graph
144 }
145 }
146 c1->Update();
147 printf("\n\n\tExtracted %d Contours and %d Graphs \n", TotalConts, nGraphs );
148 gStyle->SetTitleW(0.);
149 gStyle->SetTitleH(0.);
150 return c1;
151}
152
153
154Double_t SawTooth(Double_t x, Double_t WaveLen){
155
156// This function is specific to a sawtooth function with period
157// WaveLen, symmetric about x = 0, and with amplitude = 1. Each segment
158// is 1/4 of the wavelength.
159//
160// |
161// /\ |
162// / \ |
163// / \ |
164// / \
165// /--------\--------/------------
166// |\ /
167// | \ /
168// | \ /
169// | \/
170//
171
172 Double_t y;
173 if ( (x < -WaveLen/2) || (x > WaveLen/2)) y = -99999999; // Error X out of bounds
174 if (x <= -WaveLen/4) {
175 y = x + 2.0;
176 } else if ((x > -WaveLen/4) && (x <= WaveLen/4)) {
177 y = -x ;
178 } else if (( x > WaveLen/4) && (x <= WaveLen/2)) {
179 y = x - 2.0;
180 }
181 return y;
182}
#define c(i)
Definition RSha256.hxx:101
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
@ kRed
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
#define PI
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void gc
#define gROOT
Definition TROOT.h:406
R__EXTERN TStyle * gStyle
Definition TStyle.h:433
Array of doubles (64 bits per element).
Definition TArrayD.h:27
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
The Canvas class.
Definition TCanvas.h:23
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
Int_t GetN() const
Definition TGraph.h:131
virtual Int_t GetPoint(Int_t i, Double_t &x, Double_t &y) const
Get x and y values for point number i.
Definition TGraph.cxx:1516
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3066
virtual void SetContour(Int_t nlevels, const Double_t *levels=nullptr)
Set the number and values of contour levels.
Definition TH1.cxx:8451
2-D histogram with a double per channel (see TH1 documentation)
Definition TH2.h:357
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:307
void SetBinContent(Int_t bin, Double_t content) override
Set bin content.
Definition TH2.cxx:2616
To draw Mathematical Formula.
Definition TLatex.h:18
A doubly linked list.
Definition TList.h:38
TObject * After(const TObject *obj) const override
Returns the object after object obj.
Definition TList.cxx:328
TObject * First() const override
Return the first object in the list. Returns 0 when list is empty.
Definition TList.cxx:657
TObject * Clone(const char *newname="") const override
Make a clone of an object using the Streamer facility.
Definition TNamed.cxx:74
An array of TObjects.
Definition TObjArray.h:31
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
TObject * FindObject(const char *name) const override
Find an object in this collection using its name.
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition TStyle.cxx:1636
void SetTitleW(Float_t w=0)
Definition TStyle.h:412
void SetTitleH(Float_t h=0)
Definition TStyle.h:413
RVec< PromoteType< T > > sin(const RVec< T > &v)
Definition RVec.hxx:1851
Double_t y[n]
Definition legend1.C:17
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
constexpr Double_t Pi()
Definition TMath.h:37
TLine l
Definition textangle.C:4