Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
h1analysisProxy.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_tree_analysis
3
4/// Example of analysis class for the H1 data using code generated by MakeProxy.
5///
6/// This file uses 4 large data sets from the H1 collaboration at DESY Hamburg.
7/// One can access these data sets (277 MBytes) from the standard Root web site
8/// at: [https://root.cern/files/h1/](https://root.cern/files/h1/)
9/// The Physics plots below generated by this example cannot be produced when
10/// using smaller data sets.
11///
12/// There are several ways to analyze data stored in a Root Tree
13/// - Using TTree::Draw: This is very convenient and efficient for small tasks.
14/// A TTree::Draw call produces one histogram at the time. The histogram
15/// is automatically generated. The selection expression may be specified
16/// in the command line.
17///
18/// - Using the TTreeViewer: This is a graphical interface to TTree::Draw
19/// with the same functionality.
20///
21/// - Using the code generated by TTree::MakeClass: In this case, the user
22/// creates an instance of the analysis class. They have the control over
23/// the event loop and he can generate an unlimited number of histograms.
24///
25/// - Using the code generated by TTree::MakeSelector. Like for the code
26/// generated by TTree::MakeClass, the user can do complex analysis.
27/// However, they cannot control the event loop. The event loop is controlled
28/// by TTree::Process called by the user. This solution is illustrated
29/// by the current code. The advantage of this method is that it can be run
30/// in a parallel environment using PROOF (the Parallel Root Facility).
31///
32/// A chain of 4 files (originally converted from PAW ntuples) is used
33/// to illustrate the various ways to loop on Root data sets.
34/// Each data set contains a Root Tree named "h42"
35///
36/// h1analysProxy.C can be used either via TTree::Draw:
37/// ~~~{.cpp}
38/// h42->Draw("h1analysisProxy.C");
39/// ~~~
40/// or it can be used directly with TTree::MakeProxy, for example to generate a
41/// shared library. TTree::MakeProxy will generate a TSelector skeleton that
42/// include h1analysProxy.C:
43/// ~~~{.cpp}
44/// h42->MakeProxy("h1sel","h1analysisProxy.C");
45/// ~~~
46/// This produces one file: h1sel.h which does a `#include "h1analysProxy.C"`
47/// The h1sel class is derived from the Root class TSelector and can then
48/// be used as:
49/// ~~~{.cpp}
50/// h42->Process("h1sel.h+");
51/// ~~~
52///
53/// The following members functions are called by the TTree::Process functions.
54/// - **h1analysProxy_Begin()**: Called every time a loop on the tree starts.
55/// a convenient place to create your histograms.
56///
57/// - **h1analysProxy_Notify()**: This function is called at the first entry of a new Tree
58/// in a chain.
59/// - **h1analysProxy_Process()**: called to analyze each entry.
60///
61/// - **h1analysProxy_Terminate()**: called at the end of a loop on a TTree.
62/// a convenient place to draw/fit your histograms.
63///
64/// To use this file, try the following session
65/// ~~~{.cpp}
66/// Root > gROOT->Time(); ///will show RT & CPU time per command
67/// ~~~
68/// ### Case A: Create a TChain with the 4 H1 data files
69/// The chain can be created by executed the short macro h1chain.C below:
70/// ~~~{.cpp}
71/// {
72/// TChain chain("h42");
73/// chain.Add("$H1/dstarmb.root"); /// 21330730 bytes 21920 events
74/// chain.Add("$H1/dstarp1a.root"); /// 71464503 bytes 73243 events
75/// chain.Add("$H1/dstarp1b.root"); /// 83827959 bytes 85597 events
76/// chain.Add("$H1/dstarp2.root"); /// 100675234 bytes 103053 events
77/// ///where $H1 is a system symbol pointing to the H1 data directory.
78/// }
79/// ~~~
80///
81/// ### Case B: Loop on all events
82/// ~~~{.cpp}
83/// Root > chain.Draw("h1analysisProxy.C")
84/// ~~~
85///
86/// ### Case C: Same as B, but in addition fill the event list with selected entries.
87/// The event list is saved to a file "elist.root" by the Terminate function.
88/// To see the list of selected events, you can do elist->Print("all").
89/// The selection function has selected 7525 events out of the 283813 events
90/// in the chain of files. (2.65 per cent)
91/// ~~~{.cpp}
92/// Root > chain.Draw("h1analysisProxy.C","","fillList")
93/// ~~~
94/// ### Case D: Process only entries in the event list
95/// The event list is read from the file in elist.root generated by step C
96/// ~~~{.cpp}
97/// Root > chain.Draw("h1analysisProxy.C","","useList")
98/// ~~~
99///
100/// The commands executed with the 3 different methods B,C and D
101/// produce two canvases shown below:
102/// begin_html <a href="gif/h1analysis_dstar.gif" >the Dstar plot</a> end_html
103/// begin_html <a href="gif/h1analysis_tau.gif" >the Tau D0 plot</a> end_html
104///
105/// \macro_code
106///
107/// \author Philippe Canal from original h1analysis.C by Rene Brun
108
109TEntryList *elist;
110Bool_t useList, fillList;
111TH1F *hdmd;
112h2;
113
114
116{
117// function called before starting the event loop
118// -it performs some cleanup
119// -it creates histograms
120// -it sets some initialisation for the event list
121
122 //print the option specified in the Process function.
123 TString option = GetOption();
124 printf("Starting (begin) h1analysis with process option: %s\n",option.Data());
125
126 //process cases with event list
127 fillList = kFALSE;
128 useList = kFALSE;
129 if (fChain) fChain->SetEntryList(0);
130 delete gDirectory->GetList()->FindObject("elist");
131
132 // case when one creates/fills the event list
133 if (option.Contains("fillList")) {
134 fillList = kTRUE;
135 elist = new TEntryList("elist","H1 selection from Cut");
136 // Add to the input list for processing in PROOF, if needed
137 if (fInput) {
138 fInput->Add(new TNamed("fillList",""));
139 fInput->Add(elist);
140 }
141 } else elist = 0;
142
143 // case when one uses the event list generated in a previous call
144 if (option.Contains("useList")) {
145 useList = kTRUE;
146 if (fInput) {
147 tree->SetEntryList(elist);
148 TFile f("elist.root");
149 elist = (TEntryList*)f.Get("elist");
150 if (elist) elist->SetDirectory(0); //otherwise the file destructor will delete elist
151 } else {
152 // Option "useList" not supported in PROOF directly
153 Warning("Begin", "option 'useList' not supported in PROOF - ignoring");
154 Warning("Begin", "the entry list must be set on the chain *before* calling Process");
155 }
156 }
157}
158
159
161{
162// function called before starting the event loop
163// -it performs some cleanup
164// -it creates histograms
165// -it sets some initialisation for the entry list
166
167 //initialize the Tree branch addresses
168 Init(tree);
169
170 //print the option specified in the Process function.
171 TString option = GetOption();
172 printf("Starting (slave) h1analysis with process option: %s\n",option.Data());
173
174 //create histograms
175 hdmd = new TH1F("hdmd","dm_d",40,0.13,0.17);
176 TH2F("h2","ptD0 vs dm_d",30,0.135,0.165,30,-3,6);
177
178 fOutput->Add(hdmd);
179 fOutput->Add(h2);
180
181 //process cases with entry list
182 fillList = kFALSE;
183 useList = kFALSE;
184
185 // case when one creates/fills the entry list
186 if (option.Contains("fillList")) {
187 fillList = kTRUE;
188 // Get the list
189 if (fInput) {
190 if ((elist = (TEntryList *) fInput->FindObject("elist")))
191 // Need to clone to avoid problems when destroying the selector
192 elist = (TEntryList *) elist->Clone();
193 }
194 if (elist)
195 fOutput->Add(elist);
196 else
197 fillList = kFALSE;
198 } else elist = 0;
199
200 // case when one uses the entry list generated in a previous call
201 if (option.Contains("useList")) {
202 useList = kTRUE;
203 TFile f("elist.root");
204 elist = (TEntryList*)f.Get("elist");
205 if (elist) elist->SetDirectory(0); //otherwise the file destructor will delete elist
206 if (tree) tree->SetEntryList(elist);
207 else {
208 // Option "useList" not supported in PROOF directly
209 Warning("Begin", "option 'useList' not supported in PROOF - ignoring");
210 Warning("Begin", "the entry list must be set on the chain *before* calling Process");
211 }
212 }
213
214}
215
217 return 0;
218}
219
220
222{
223// entry is the entry number in the current Tree
224// Selection function to select D* and D0.
225
226 //in case one entry list is given in input, the selection has already been done.
227 if (!useList) {
228
229 float f1 = md0_d;
230 float f2 = md0_d-1.8646;
231 bool test = TMath::Abs(md0_d-1.8646) >= 0.04;
232 if (gDebug>0) fprintf(stderr,"entry #%lld f1=%f f2=%f test=%d\n",
233 fChain->GetReadEntry(),f1,f2,test);
234
235 if (TMath::Abs(md0_d-1.8646) >= 0.04) return kFALSE;
236 if (ptds_d <= 2.5) return kFALSE;
237 if (TMath::Abs(etads_d) >= 1.5) return kFALSE;
238
239 int cik = ik-1; //original ik used f77 convention starting at 1
240 int cipi = ipi-1; //original ipi used f77 convention starting at 1
241
242 f1 = nhitrp[cik];
243 f2 = nhitrp[cipi];
244 test = nhitrp[cik]*nhitrp[cipi] <= 1;
245 if (gDebug>0) fprintf(stderr,"entry #%lld f1=%f f2=%f test=%d\n",
246 fChain->GetReadEntry(),f1,f2,test);
247
248 if (nhitrp[cik]*nhitrp[cipi] <= 1) return kFALSE;
249 if (rend[cik] -rstart[cik] <= 22) return kFALSE;
250 if (rend[cipi]-rstart[cipi] <= 22) return kFALSE;
251 if (nlhk[cik] <= 0.1) return kFALSE;
252 if (nlhpi[cipi] <= 0.1) return kFALSE;
253 // fix because read-only
254 if (nlhpi[ipis-1] <= 0.1) return kFALSE;
255 if (njets < 1) return kFALSE;
256
257 }
258 // if option fillList, fill the event list
259 if (fillList) elist->Enter(entry);
260
261 //fill some histograms
262 hdmd->Fill(dm_d);
263 h2->Fill(dm_d,rpd0_t/0.029979*1.8646/ptd0_d);
264
265 return kTRUE;
266}
267
268
269
271{
272 // nothing to be done
273 printf("Terminate (slave) h1analysis\n");
274}
275
276
278{
279 printf("Terminate (final) h1analysis\n");
280
281 // function called at the end of the event loop
282
283 hdmd = dynamic_cast<TH1F*>(fOutput->FindObject("hdmd"));
284 TH2F*>(fOutput->FindObject("h2"));
285
286 if (hdmd == 0 || h2 == 0) {
287 Error("Terminate", "hdmd = %p , h2 = %p", hdmd, h2);
288 return;
289 }
290
291 //create the canvas for the h1analysis fit
292 gStyle->SetOptFit();
293 TCanvas *c1 = new TCanvas("c1","h1analysis analysis",10,10,800,600);
294 c1->SetBottomMargin(0.15);
295 hdmd->GetXaxis()->SetTitle("m_{K#pi#pi} - m_{K#pi}[GeV/c^{2}]");
296 hdmd->GetXaxis()->SetTitleOffset(1.4);
297
298 //fit histogram hdmd with function f5 using the log-likelihood option
299 TF1 *f5 = new TF1("f5",fdm5,0.139,0.17,5);
300 f5->SetParameters(1000000, .25, 2000, .1454, .001);
301 hdmd->Fit("f5","lr");
302
303 //create the canvas for tau d0
304 gStyle->SetOptFit(0);
305 gStyle->SetOptStat(1100);
306 TCanvas *c2 = new TCanvas("c2","tauD0",100,100,800,600);
307 c2->SetGrid();
308 c2->SetBottomMargin(0.15);
309
310 // Project slices of 2-d histogram h2 along X , then fit each slice
311 // with function f2 and make a histogram for each fit parameter
312 // Note that the generated histograms are added to the list of objects
313 // in the current directory.
314 TF1 *f2 = new TF1("f2",fdm2,0.139,0.17,2);
315 f2->SetParameters(10000, 10);
316 h2->FitSlicesX(f2,0,-1,1,"qln");
317 TH1D *h2_1 = (TH1D*)gDirectory->Get("h2_1");
318 h2_1->GetXaxis()->SetTitle("#tau[ps]");
319 h2_1->SetMarkerStyle(21);
320 h2_1->Draw();
321 c2->Update();
322 TLine *line = new TLine(0,0,0,c2->GetUymax());
323 line->Draw();
324
325 // Have the number of entries on the first histogram (to cross check when running
326 // with entry lists)
328 psdmd->SetOptStat(1110);
329 c1->Modified();
330
331 //save the entry list to a Root file if one was produced
332 if (fillList) {
333 elist = dynamic_cast<TEntryList*>(fOutput->FindObject("elist"));
334 if (elist) {
335 TFile efile("elist.root","recreate");
336 elist->Write();
337 } else {
338 Error("Terminate", "entry list requested but not found in output");
339 }
340 }
341}
#define f(i)
Definition RSha256.hxx:104
bool Bool_t
Definition RtypesCore.h:63
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
double Double_t
Definition RtypesCore.h:59
long long Long64_t
Definition RtypesCore.h:69
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define gDirectory
Definition TDirectory.h:384
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:185
void Warning(const char *location, const char *msgfmt,...)
Use this function in warning situations.
Definition TError.cxx:229
Option_t Option_t option
Int_t Init(const char *masterurl, const char *conffile, const char *confdir, Int_t loglevel, const char *alias=0)
Int_t gDebug
Definition TROOT.cxx:622
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:280
The Canvas class.
Definition TCanvas.h:23
A List of entry numbers in a TTree or TChain.
Definition TEntryList.h:26
virtual bool Enter(Long64_t entry, TTree *tree=nullptr)
Add entry #entry to the list.
virtual void SetDirectory(TDirectory *dir)
Add reference to directory dir. dir can be 0.
1-Dim function class
Definition TF1.h:234
virtual void SetParameters(const Double_t *params)
Definition TF1.h:685
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:131
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:925
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:877
TAxis * GetXaxis()
Definition TH1.h:571
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition TH1.cxx:3876
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3316
TList * GetListOfFunctions() const
Definition TH1.h:488
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:307
Use the TLine constructor to create a simple line.
Definition TLine.h:22
TObject * FindObject(const char *name) const override
Find an object in this list using its name.
Definition TList.cxx:576
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
TObject * Clone(const char *newname="") const override
Make a clone of an object using the Streamer facility.
Definition TNamed.cxx:74
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:174
virtual Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition TObject.cxx:964
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:293
The histogram statistics painter class.
Definition TPaveStats.h:18
Basic string class.
Definition TString.h:139
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:1642
void SetOptFit(Int_t fit=1)
The type of information about fit parameters printed in the histogram statistics box can be selected ...
Definition TStyle.cxx:1595
A TTree represents a columnar dataset.
Definition TTree.h:84
TLine * line
Double_t fdm5(Double_t *xx, Double_t *par)
Double_t fdm2(Double_t *xx, Double_t *par)
return c1
Definition legend1.C:41
TF1 * f1
Definition legend1.C:11
return c2
Definition legend2.C:14
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123