ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
h1analysis.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_tree
3 ///
4 /// Example of analysis class for the H1 data.
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: `ftp:/// root.cern.ch/root/h1analysis`
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 /// The class definition in h1analysis.h has been generated automatically
36 /// by the Root utility TTree::MakeSelector using one of the files with the
37 /// following statement:
38 /// ~~~{.cpp}
39 /// h42->MakeSelector("h1analysis");
40 /// ~~~
41 /// This produces two files: h1analysis.h and h1analysis.C (skeleton of this file)
42 /// The h1analysis class is derived from the Root class TSelector.
43 ///
44 /// The following members functions are called by the TTree::Process functions.
45 /// - **Begin()**: Called everytime a loop on the tree starts.
46 /// A convenient place to create your histograms.
47 /// - **Notify()**: This function is called at the first entry of a new Tree
48 /// in a chain.
49 /// - **Process()**: Called to analyze each entry.
50 ///
51 /// - **Terminate()**: Called at the end of a loop on a TTree.
52 /// A convenient place to draw/fit your histograms.
53 ///
54 /// To use this file, try the following sessions
55 /// ~~~{.cpp}
56 /// Root > gROOT->Time(); /// will show RT & CPU time per command
57 /// ~~~
58 /// ## Case A: Create a TChain with the 4 H1 data files
59 /// The chain can be created by executed the short macro h1chain.C below:
60 /// ~~~{.cpp}
61 /// {
62 /// TChain chain("h42");
63 /// chain.Add("$H1/dstarmb.root"); /// 21330730 bytes 21920 events
64 /// chain.Add("$H1/dstarp1a.root"); /// 71464503 bytes 73243 events
65 /// chain.Add("$H1/dstarp1b.root"); /// 83827959 bytes 85597 events
66 /// chain.Add("$H1/dstarp2.root"); /// 100675234 bytes 103053 events
67 /// /// where $H1 is a system symbol pointing to the H1 data directory.
68 /// }
69 /// ~~~
70 /// ## Case B: Loop on all events
71 /// ~~~{.cpp}
72 /// Root > chain.Process("h1analysis.C")
73 /// ~~~
74 /// ## Case C: Same as B, but in addition fill the entry list with selected entries.
75 /// The entry list is saved to a file "elist.root" by the Terminate function.
76 /// To see the list of selected events, you can do `elist->Print("all")`.
77 /// The selection function has selected 7525 events out of the 283813 events
78 /// in the chain of files. (2.65 per cent)
79 /// ~~~{.cpp}
80 /// Root > chain.Process("h1analysis.C","fillList")
81 /// ~~~
82 /// ## Case D: Process only entries in the entry list
83 /// The entry list is read from the file in elist.root generated by step C
84 /// ~~~{.cpp}
85 /// Root > chain.Process("h1analysis.C","useList")
86 /// ~~~
87 /// ## Case E: The above steps have been executed via the interpreter.
88 /// You can repeat the steps B, C and D using the script compiler
89 /// by replacing "h1analysis.C" by "h1analysis.C+" or "h1analysis.C++"
90 /// in a new session (see F).
91 ///
92 /// ## Case F: Create the chain as in A, then execute
93 /// ~~~{.cpp}
94 /// Root > chain.Process("h1analysis.C+","useList")
95 /// ~~~
96 ///
97 /// The same analysis can be run on PROOF. For a quick try start a PROOF-Lite
98 /// session
99 /// ~~~{.cpp}
100 /// Root > TProof *p = TProof::Open("")
101 /// ~~~
102 /// create (if mot already done) the chain by executing the 'h1chain.C' macro
103 /// mentioned above, and then tell ROOT to use PROOF to process the chain:
104 /// ~~~{.cpp}
105 /// Root > chain.SetProof()
106 /// ~~~
107 /// You can then repeat step B above. Step C can also be executed in PROOF. However,
108 /// step D cannot be executed in PROOF as in the local session (i.e. just passing
109 /// option 'useList'): to use the entry list you have to
110 ///
111 /// ## Case G: Load first in the session the list form the file
112 /// ~~~{.cpp}
113 /// Root > TFile f("elist.root")
114 /// Root > TEntryList *elist = (TEntryList *) f.Get("elist")
115 /// ~~~
116 /// set it on the chain:
117 /// ~~~{.cpp}
118 /// Root > chain.SetEntryList(elist)
119 /// ~~~
120 /// call Process as in step B. Of course this works also for local processing.
121 /// \macro_code
122 /// \author Rene Brun
123 
124 #include "h1analysis.h"
125 #include "TH2.h"
126 #include "TF1.h"
127 #include "TStyle.h"
128 #include "TCanvas.h"
129 #include "TPaveStats.h"
130 #include "TLine.h"
131 #include "TMath.h"
132 
133 const Double_t dxbin = (0.17-0.13)/40; // Bin-width
134 const Double_t sigma = 0.0012;
135 
136 
137 Double_t fdm5(Double_t *xx, Double_t *par)
138 {
139  Double_t x = xx[0];
140  if (x <= 0.13957) return 0;
141  Double_t xp3 = (x-par[3])*(x-par[3]);
142  Double_t res = dxbin*(par[0]*TMath::Power(x-0.13957, par[1])
143  + par[2] / 2.5066/par[4]*TMath::Exp(-xp3/2/par[4]/par[4]));
144  return res;
145 }
146 
147 
148 Double_t fdm2(Double_t *xx, Double_t *par)
149 {
150  Double_t x = xx[0];
151  if (x <= 0.13957) return 0;
152  Double_t xp3 = (x-0.1454)*(x-0.1454);
153  Double_t res = dxbin*(par[0]*TMath::Power(x-0.13957, 0.25)
154  + par[1] / 2.5066/sigma*TMath::Exp(-xp3/2/sigma/sigma));
155  return res;
156 }
157 
158 
159 void h1analysis::Begin(TTree * /*tree*/)
160 {
161 // function called before starting the event loop
162 // -it performs some cleanup
163 // -it creates histograms
164 // -it sets some initialisation for the entry list
165 
166  // This is needed when re-processing the object
167  Reset();
168 
169  //print the option specified in the Process function.
170  TString option = GetOption();
171  Info("Begin", "starting h1analysis with process option: %s", option.Data());
172 
173  //process cases with entry list
174  if (fChain) fChain->SetEntryList(0);
175  delete gDirectory->GetList()->FindObject("elist");
176 
177  // case when one creates/fills the entry list
178  if (option.Contains("fillList")) {
179  fillList = kTRUE;
180  elist = new TEntryList("elist", "H1 selection from Cut");
181  // Add to the input list for processing in PROOF, if needed
182  if (fInput) {
183  fInput->Add(new TNamed("fillList",""));
184  // We send a clone to avoid double deletes when importing the result
185  fInput->Add(elist);
186  // This is needed to avoid warnings from output-to-members mapping
187  elist = 0;
188  }
189  }
190  if (fillList) Info("Begin", "creating an entry-list");
191  // case when one uses the entry list generated in a previous call
192  if (option.Contains("useList")) {
193  useList = kTRUE;
194  if (fInput) {
195  // Option "useList" not supported in PROOF directly
196  Warning("Begin", "option 'useList' not supported in PROOF - ignoring");
197  Warning("Begin", "the entry list must be set on the chain *before* calling Process");
198  } else {
199  TFile f("elist.root");
200  elist = (TEntryList*)f.Get("elist");
201  if (elist) elist->SetDirectory(0); //otherwise the file destructor will delete elist
202  }
203  }
204 }
205 
206 
207 void h1analysis::SlaveBegin(TTree *tree)
208 {
209 // function called before starting the event loop
210 // -it performs some cleanup
211 // -it creates histograms
212 // -it sets some initialisation for the entry list
213 
214  //initialize the Tree branch addresses
215  Init(tree);
216 
217  //print the option specified in the Process function.
218  TString option = GetOption();
219  Info("SlaveBegin",
220  "starting h1analysis with process option: %s (tree: %p)", option.Data(), tree);
221 
222  //create histograms
223  hdmd = new TH1F("hdmd","dm_d",40,0.13,0.17);
224  h2 = new TH2F("h2","ptD0 vs dm_d",30,0.135,0.165,30,-3,6);
225 
226  fOutput->Add(hdmd);
227  fOutput->Add(h2);
228 
229  // Entry list stuff (re-parse option because on PROOF only SlaveBegin is called)
230  if (option.Contains("fillList")) {
231  fillList = kTRUE;
232  // Get the list
233  if (fInput) {
234  if ((elist = (TEntryList *) fInput->FindObject("elist")))
235  // Need to clone to avoid problems when destroying the selector
236  elist = (TEntryList *) elist->Clone();
237  if (elist)
238  fOutput->Add(elist);
239  else
240  fillList = kFALSE;
241  }
242  }
243  if (fillList) Info("SlaveBegin", "creating an entry-list");
244 }
245 
246 
248 {
249 // entry is the entry number in the current Tree
250 // Selection function to select D* and D0.
251 
252  fProcessed++;
253  //in case one entry list is given in input, the selection has already been done.
254  if (!useList) {
255  // Read only the necessary branches to select entries.
256  // return as soon as a bad entry is detected
257  // to read complete event, call fChain->GetTree()->GetEntry(entry)
258  b_md0_d->GetEntry(entry); if (TMath::Abs(md0_d-1.8646) >= 0.04) return kFALSE;
259  b_ptds_d->GetEntry(entry); if (ptds_d <= 2.5) return kFALSE;
260  b_etads_d->GetEntry(entry); if (TMath::Abs(etads_d) >= 1.5) return kFALSE;
261  b_ik->GetEntry(entry); ik--; //original ik used f77 convention starting at 1
262  b_ipi->GetEntry(entry); ipi--;
263  b_ntracks->GetEntry(entry);
264  b_nhitrp->GetEntry(entry);
265  if (nhitrp[ik]*nhitrp[ipi] <= 1) return kFALSE;
266  b_rend->GetEntry(entry);
267  b_rstart->GetEntry(entry);
268  if (rend[ik] -rstart[ik] <= 22) return kFALSE;
269  if (rend[ipi]-rstart[ipi] <= 22) return kFALSE;
270  b_nlhk->GetEntry(entry); if (nlhk[ik] <= 0.1) return kFALSE;
271  b_nlhpi->GetEntry(entry); if (nlhpi[ipi] <= 0.1) return kFALSE;
272  b_ipis->GetEntry(entry); ipis--; if (nlhpi[ipis] <= 0.1) return kFALSE;
273  b_njets->GetEntry(entry); if (njets < 1) return kFALSE;
274  }
275  // if option fillList, fill the entry list
276  if (fillList) elist->Enter(entry);
277 
278  // to read complete event, call fChain->GetTree()->GetEntry(entry)
279  // read branches not processed in ProcessCut
280  b_dm_d->GetEntry(entry); //read branch holding dm_d
281  b_rpd0_t->GetEntry(entry); //read branch holding rpd0_t
282  b_ptd0_d->GetEntry(entry); //read branch holding ptd0_d
283 
284  //fill some histograms
285  hdmd->Fill(dm_d);
286  h2->Fill(dm_d,rpd0_t/0.029979*1.8646/ptd0_d);
287 
288  // Count the number of selected events
289  fStatus++;
290 
291  return kTRUE;
292 }
293 
294 
295 
297 {
298  // nothing to be done
299 }
300 
301 
303 {
304 // function called at the end of the event loop
305 
306  hdmd = dynamic_cast<TH1F*>(fOutput->FindObject("hdmd"));
307  h2 = dynamic_cast<TH2F*>(fOutput->FindObject("h2"));
308 
309  if (hdmd == 0 || h2 == 0) {
310  Error("Terminate", "hdmd = %p , h2 = %p", hdmd, h2);
311  return;
312  }
313 
314  //create the canvas for the h1analysis fit
315  gStyle->SetOptFit();
316  TCanvas *c1 = new TCanvas("c1","h1analysis analysis",10,10,800,600);
317  c1->SetBottomMargin(0.15);
318  hdmd->GetXaxis()->SetTitle("m_{K#pi#pi} - m_{K#pi}[GeV/c^{2}]");
319  hdmd->GetXaxis()->SetTitleOffset(1.4);
320 
321  //fit histogram hdmd with function f5 using the loglikelihood option
322  if (gROOT->GetListOfFunctions()->FindObject("f5"))
323  delete gROOT->GetFunction("f5");
324  TF1 *f5 = new TF1("f5",fdm5,0.139,0.17,5);
325  f5->SetParameters(1000000, .25, 2000, .1454, .001);
326  hdmd->Fit("f5","lr");
327 
328  //create the canvas for tau d0
329  gStyle->SetOptFit(0);
330  gStyle->SetOptStat(1100);
331  TCanvas *c2 = new TCanvas("c2","tauD0",100,100,800,600);
332  c2->SetGrid();
333  c2->SetBottomMargin(0.15);
334 
335  // Project slices of 2-d histogram h2 along X , then fit each slice
336  // with function f2 and make a histogram for each fit parameter
337  // Note that the generated histograms are added to the list of objects
338  // in the current directory.
339  if (gROOT->GetListOfFunctions()->FindObject("f2"))
340  delete gROOT->GetFunction("f2");
341  TF1 *f2 = new TF1("f2",fdm2,0.139,0.17,2);
342  f2->SetParameters(10000, 10);
343  // Restrict to three bins in this example
344  Info("Fit Slices","Restricting fit to two bins only in this example...");
345  h2->FitSlicesX(f2,10,20,10,"g5 l");
346  TH1D *h2_1 = (TH1D*)gDirectory->Get("h2_1");
347  h2_1->GetXaxis()->SetTitle("#tau[ps]");
348  h2_1->SetMarkerStyle(21);
349  h2_1->Draw();
350  c2->Update();
351  TLine *line = new TLine(0,0,0,c2->GetUymax());
352  line->Draw();
353 
354  // Have the number of entries on the first histogram (to cross check when running
355  // with entry lists)
356  TPaveStats *psdmd = (TPaveStats *)hdmd->GetListOfFunctions()->FindObject("stats");
357  psdmd->SetOptStat(1110);
358  c1->Modified();
359 
360  //save the entry list to a Root file if one was produced
361  if (fillList) {
362  if (!elist)
363  elist = dynamic_cast<TEntryList*>(fOutput->FindObject("elist"));
364  if (elist) {
365  Printf("Entry list 'elist' created:");
366  elist->Print();
367  TFile efile("elist.root","recreate");
368  elist->Write();
369  } else {
370  Error("Terminate", "entry list requested but not found in output");
371  }
372  }
373  // Notify the amount of processed events
374  if (!fInput) Info("Terminate", "processed %lld events", fProcessed);
375 }
virtual const char * GetOption() const
Definition: TSelector.h:65
TBranch * b_md0_d
Definition: h1analysis.h:265
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition: TObject.cxx:823
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:432
TSelectorList * fOutput
Definition: TSelector.h:50
long long Long64_t
Definition: RtypesCore.h:69
Bool_t Process(Long64_t entry)
TBranch * b_rstart
Definition: h1analysis.h:315
TObject * FindObject(const char *name) const
Find object using its name.
Definition: THashList.cxx:213
Float_t rstart[200]
Definition: h1analysis.h:161
R__EXTERN TStyle * gStyle
Definition: TStyle.h:423
Int_t njets
Definition: h1analysis.h:174
Double_t fdm2(Double_t *xx, Double_t *par)
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition: TObject.cxx:892
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:45
Double_t GetUymax() const
Definition: TPad.h:229
virtual TObject * Get(const char *namecycle)
Return pointer to object identified by namecycle.
TBranch * b_ptds_d
Definition: h1analysis.h:252
#define gROOT
Definition: TROOT.h:344
Float_t etads_d
Definition: h1analysis.h:99
Basic string class.
Definition: TString.h:137
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
TBranch * b_rpd0_t
Definition: h1analysis.h:285
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition: TObject.cxx:254
void Terminate()
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition: TList.cxx:497
The histogram statistics painter class.
Definition: TPaveStats.h:28
TBranch * b_ptd0_d
Definition: h1analysis.h:263
void Reset()
Definition: h1analysis.h:376
TBranch * b_ipi
Definition: h1analysis.h:260
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
TFile * f
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
TBranch * b_rend
Definition: h1analysis.h:316
TEntryList * elist
Definition: h1analysis.h:32
void SlaveBegin(TTree *tree)
virtual void SetDirectory(TDirectory *dir)
Add reference to directory dir. dir can be 0.
const char * Data() const
Definition: TString.h:349
Float_t nlhpi[200]
Definition: h1analysis.h:166
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:33
TH2F * h2
Definition: h1analysis.h:28
Float_t ptd0_d
Definition: h1analysis.h:109
Float_t md0_d
Definition: h1analysis.h:111
virtual void SetGrid(Int_t valuex=1, Int_t valuey=1)
Definition: TPad.h:326
TBranch * b_etads_d
Definition: h1analysis.h:253
const Double_t sigma
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition: TNamed.cxx:63
TBranch * b_njets
Definition: h1analysis.h:328
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
Int_t nhitrp[200]
Definition: h1analysis.h:157
Long64_t fProcessed
Definition: h1analysis.h:33
void SlaveTerminate()
TBranch * b_ik
Definition: h1analysis.h:259
virtual void SetBottomMargin(Float_t bottommargin)
Set Pad bottom margin in fraction of the pad height.
Definition: TAttPad.cxx:97
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
Long64_t fStatus
Definition: TSelector.h:45
void Init(TTree *tree)
Definition: h1analysis.h:390
virtual void FitSlicesX(TF1 *f1=0, Int_t firstybin=0, Int_t lastybin=-1, Int_t cut=0, Option_t *option="QNR", TObjArray *arr=0)
Project slices along X in case of a 2-D histogram, then fit each slice with function f1 and make a hi...
Definition: TH2.cxx:891
void Begin(TTree *tree)
TBranch * b_nlhk
Definition: h1analysis.h:319
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:256
Float_t dm_d
Definition: h1analysis.h:100
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:1204
const Double_t dxbin
A simple line.
Definition: TLine.h:41
virtual Int_t GetEntry(Long64_t entry=0, Int_t getall=0)
Read all leaves of entry and return total number of bytes read.
Definition: TBranch.cxx:1199
virtual void SetMarkerStyle(Style_t mstyle=1)
Definition: TAttMarker.h:53
TTree * fChain
Definition: h1analysis.h:35
Int_t ipi
Definition: h1analysis.h:106
Bool_t useList
Definition: h1analysis.h:30
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:613
Float_t rend[200]
Definition: h1analysis.h:162
#define Printf
Definition: TGeoToOCC.h:18
The Canvas class.
Definition: TCanvas.h:48
Double_t Exp(Double_t x)
Definition: TMath.h:495
tuple tree
Definition: tree.py:24
Int_t ik
Definition: h1analysis.h:105
Float_t nlhk[200]
Definition: h1analysis.h:165
virtual void SetEntryList(TEntryList *list, Option_t *opt="")
Set an EntryList.
Definition: TTree.cxx:8141
TBranch * b_ipis
Definition: h1analysis.h:261
double Double_t
Definition: RtypesCore.h:55
TBranch * b_dm_d
Definition: h1analysis.h:254
TBranch * b_ntracks
Definition: h1analysis.h:303
Float_t rpd0_t
Definition: h1analysis.h:131
Int_t ipis
Definition: h1analysis.h:107
virtual Bool_t Enter(Long64_t entry, TTree *tree=0)
Add entry entry to the list.
Definition: TEntryList.cxx:558
Double_t fdm5(Double_t *xx, Double_t *par)
TList * fInput
Current object if processing object (vs. TTree)
Definition: TSelector.h:49
TBranch * b_nlhpi
Definition: h1analysis.h:320
virtual void Add(TObject *obj)
Definition: TList.h:81
Float_t ptds_d
Definition: h1analysis.h:98
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
1-Dim function class
Definition: TF1.h:149
TBranch * b_nhitrp
Definition: h1analysis.h:311
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:1252
A TTree object has a header with a name and a title.
Definition: TTree.h:98
#define gDirectory
Definition: TDirectory.h:221
virtual void Update()
Update canvas pad buffers.
Definition: TCanvas.cxx:2179
virtual void Print(const Option_t *option="") const
Print this list.
Definition: TEntryList.cxx:989
const Bool_t kTRUE
Definition: Rtypes.h:91
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:287
virtual void SetTitle(const char *title="")
Change (i.e. set) the title of the TNamed.
Definition: TNamed.cxx:152
virtual TObject * FindObject(const char *name) const
Must be redefined in derived classes.
Definition: TObject.cxx:379
A List of entry numbers in a TTree or TChain.
Definition: TEntryList.h:27
void Modified(Bool_t flag=1)
Definition: TPad.h:407
TAxis * GetXaxis()
Definition: TH1.h:319
Bool_t fillList
Definition: h1analysis.h:31
void SetOptStat(Int_t stat=1)
Set the stat option.
Definition: TPaveStats.cxx:303
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:904