Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
h1analysisTreeReader.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_tree
3 /// \notebook
4 /// H1 analysis example expressed in terms of TTreeReader
5 ///
6 /// \macro_code
7 ///
8 /// \author Anders Eie, 2013
9 
10 #include "h1analysisTreeReader.h"
11 #include "TStyle.h"
12 #include "TCanvas.h"
13 #include "TPaveStats.h"
14 #include "TLine.h"
15 #include "TMath.h"
16 #include "TFile.h"
17 #include "TROOT.h"
18 
19 
20 const Double_t dxbin = (0.17-0.13)/40; // Bin-width
21 const Double_t sigma = 0.0012;
22 
23 //_____________________________________________________________________
24 Double_t fdm5(Double_t *xx, Double_t *par)
25 {
26  Double_t x = xx[0];
27  if (x <= 0.13957) return 0;
28  Double_t xp3 = (x-par[3])*(x-par[3]);
29  Double_t res = dxbin*(par[0]*TMath::Power(x-0.13957, par[1])
30  + par[2] / 2.5066/par[4]*TMath::Exp(-xp3/2/par[4]/par[4]));
31  return res;
32 }
33 
34 //_____________________________________________________________________
35 Double_t fdm2(Double_t *xx, Double_t *par)
36 {
37  Double_t x = xx[0];
38  if (x <= 0.13957) return 0;
39  Double_t xp3 = (x-0.1454)*(x-0.1454);
40  Double_t res = dxbin*(par[0]*TMath::Power(x-0.13957, 0.25)
41  + par[1] / 2.5066/sigma*TMath::Exp(-xp3/2/sigma/sigma));
42  return res;
43 }
44 //_____________________________________________________________________
45 Bool_t h1analysisTreeReader::Process(Long64_t entry){
46 // entry is the entry number in the current Tree
47 // Selection function to select D* and D0.
48  myTreeReader.SetLocalEntry(entry);
49  fProcessed++;
50  //in case one entry list is given in input, the selection has already been done.
51  if (!useList) {
52  // Return as soon as a bad entry is detected
53  if (TMath::Abs(*fMd0_d-1.8646) >= 0.04) return kFALSE;
54  if (*fPtds_d <= 2.5) return kFALSE;
55  if (TMath::Abs(*fEtads_d) >= 1.5) return kFALSE;
56  (*fIk)--; //original fIk used f77 convention starting at 1
57  (*fIpi)--;
58 
59 
60  if (fNhitrp.At(*fIk)* fNhitrp.At(*fIpi) <= 1) return kFALSE;
61 
62 
63  if (fRend.At(*fIk) -fRstart.At(*fIk) <= 22) return kFALSE;
64  if (fRend.At(*fIpi)-fRstart.At(*fIpi) <= 22) return kFALSE;
65  if (fNlhk.At(*fIk) <= 0.1) return kFALSE;
66  if (fNlhpi.At(*fIpi) <= 0.1) return kFALSE;
67  (*fIpis)--; if (fNlhpi.At(*fIpis) <= 0.1) return kFALSE;
68  if (*fNjets < 1) return kFALSE;
69  }
70  // if option fillList, fill the entry list
71  if (fillList) elist->Enter(entry);
72 
73  //fill some histograms
74  hdmd->Fill(*fDm_d);
75  h2->Fill(*fDm_d,*fRpd0_t/0.029979*1.8646/ *fPtd0_d);
76 
77  return kTRUE;
78 }
79 
80 void h1analysisTreeReader::Begin(TTree* /*myTree*/) {
81 // function called before starting the event loop
82 // -it performs some cleanup
83 // -it creates histograms
84 // -it sets some initialisation for the entry list
85 
86  Reset();
87 
88  //print the option specified in the Process function.
89  TString option = GetOption();
90  Info("Begin", "starting h1analysis with process option: %s", option.Data());
91 
92  delete gDirectory->GetList()->FindObject("elist");
93 
94  // case when one creates/fills the entry list
95  if (option.Contains("fillList")) {
96  fillList = kTRUE;
97  elist = new TEntryList("elist", "H1 selection from Cut");
98  // Add to the input list for processing in PROOF, if needed
99  if (fInput) {
100  fInput->Add(new TNamed("fillList",""));
101  // We send a clone to avoid double deletes when importing the result
102  fInput->Add(elist);
103  // This is needed to avoid warnings from output-to-members mapping
104  elist = 0;
105  }
106  Info("Begin", "creating an entry-list");
107  }
108  // case when one uses the entry list generated in a previous call
109  if (option.Contains("useList")) {
110  useList = kTRUE;
111  if (fInput) {
112  // In PROOF option "useList" is processed in SlaveBegin and we do not need
113  // to do anything here
114  } else {
115  TFile f("elist.root");
116  elist = (TEntryList*)f.Get("elist");
117  if (elist) elist->SetDirectory(0); //otherwise the file destructor will delete elist
118  }
119  }
120 }
121 
122 void h1analysisTreeReader::SlaveBegin(TTree *myTree){
123 
124 // function called before starting the event loop
125 // -it performs some cleanup
126 // -it creates histograms
127 // -it sets some initialisation for the entry list
128 
129  Init(myTree);
130 
131  //print the option specified in the Process function.
132  TString option = GetOption();
133  Info("SlaveBegin",
134  "starting h1analysis with process option: %s (tree: %p)", option.Data(), myTree);
135 
136  //create histograms
137  hdmd = new TH1F("hdmd","Dm_d",40,0.13,0.17);
138  h2 = new TH2F("h2","ptD0 vs Dm_d",30,0.135,0.165,30,-3,6);
139 
140  fOutput->Add(hdmd);
141  fOutput->Add(h2);
142 
143  // Entry list stuff (re-parse option because on PROOF only SlaveBegin is called)
144  if (option.Contains("fillList")) {
145  fillList = kTRUE;
146  // Get the list
147  if (fInput) {
148  if ((elist = (TEntryList *) fInput->FindObject("elist")))
149  // Need to clone to avoid problems when destroying the selector
150  elist = (TEntryList *) elist->Clone();
151  if (elist)
152  fOutput->Add(elist);
153  else
154  fillList = kFALSE;
155  }
156  }
157  if (fillList) Info("SlaveBegin", "creating an entry-list");
158  if (option.Contains("useList")) useList = kTRUE;
159 }
160 
161 void h1analysisTreeReader::Terminate() {
162  // function called at the end of the event loop
163 
164  hdmd = dynamic_cast<TH1F*>(fOutput->FindObject("hdmd"));
165  h2 = dynamic_cast<TH2F*>(fOutput->FindObject("h2"));
166 
167  if (hdmd == 0 || h2 == 0) {
168  Error("Terminate", "hdmd = %p , h2 = %p", hdmd, h2);
169  return;
170  }
171 
172  //create the canvas for the h1analysis fit
173  gStyle->SetOptFit();
174  TCanvas *c1 = new TCanvas("c1","h1analysis analysis",10,10,800,600);
175  c1->SetBottomMargin(0.15);
176  hdmd->GetXaxis()->SetTitle("m_{K#pi#pi} - m_{K#pi}[GeV/c^{2}]");
177  hdmd->GetXaxis()->SetTitleOffset(1.4);
178 
179  //fit histogram hdmd with function f5 using the loglfIkelihood option
180  if (gROOT->GetListOfFunctions()->FindObject("f5"))
181  delete gROOT->GetFunction("f5");
182  TF1 *f5 = new TF1("f5",fdm5,0.139,0.17,5);
183  f5->SetParameters(1000000, .25, 2000, .1454, .001);
184  hdmd->Fit("f5","lr");
185 
186  //create the canvas for tau d0
187  gStyle->SetOptFit(0);
188  gStyle->SetOptStat(1100);
189  TCanvas *c2 = new TCanvas("c2","tauD0",100,100,800,600);
190  c2->SetGrid();
191  c2->SetBottomMargin(0.15);
192 
193  // Project slices of 2-d histogram h2 along X , then fit each slice
194  // with function f2 and make a histogram for each fit parameter
195  // Note that the generated histograms are added to the list of objects
196  // in the current directory.
197  if (gROOT->GetListOfFunctions()->FindObject("f2"))
198  delete gROOT->GetFunction("f2");
199  TF1 *f2 = new TF1("f2",fdm2,0.139,0.17,2);
200  f2->SetParameters(10000, 10);
201  h2->FitSlicesX(f2,0,-1,1,"qln");
202  TH1D *h2_1 = (TH1D*)gDirectory->Get("h2_1");
203  h2_1->GetXaxis()->SetTitle("#tau[ps]");
204  h2_1->SetMarkerStyle(21);
205  h2_1->Draw();
206  c2->Update();
207  TLine *line = new TLine(0,0,0,c2->GetUymax());
208  line->Draw();
209 
210  // Have the number of entries on the first histogram (to cross check when running
211  // with entry lists)
212  TPaveStats *psdmd = (TPaveStats *)hdmd->GetListOfFunctions()->FindObject("stats");
213  psdmd->SetOptStat(1110);
214  c1->Modified();
215 
216  //save the entry list to a Root file if one was produced
217  if (fillList) {
218  if (!elist)
219  elist = dynamic_cast<TEntryList*>(fOutput->FindObject("elist"));
220  if (elist) {
221  Printf("Entry list 'elist' created:");
222  elist->Print();
223  TFile efile("elist.root","recreate");
224  elist->Write();
225  } else {
226  Error("Terminate", "entry list requested but not found in output");
227  }
228  }
229  // Notify the amount of processed events
230  if (!fInput) Info("Terminate", "processed %lld events", fProcessed);
231 }
232 
233 void h1analysisTreeReader::SlaveTerminate(){
234 
235 }
236 
237 Bool_t h1analysisTreeReader::Notify() {
238 // called when loading a new file
239 // get branch pointers
240 
241  Info("Notify","processing file: %s",myTreeReader.GetTree()->GetCurrentFile()->GetName());
242 
243  if (elist && myTreeReader.GetTree()) {
244  if (fillList) {
245  elist->SetTree(myTreeReader.GetTree());
246  } else if (useList) {
247  myTreeReader.GetTree()->SetEntryList(elist);
248  }
249  }
250  return kTRUE;
251 }