ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
h1analysisProxy.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_tree
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: `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 ///
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 everytime 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 /// \macro_code
105 /// \author Philippe Canal from original h1analysis.C by Rene Brun
106 
107 TEntryList *elist;
108 Bool_t useList, fillList;
109 TH1F *hdmd;
110 TH2F *h2;
111 
112 
113 void h1analysisProxy_Begin(TTree *tree)
114 {
115 // function called before starting the event loop
116 // -it performs some cleanup
117 // -it creates histograms
118 // -it sets some initialisation for the event list
119 
120  //print the option specified in the Process function.
121  TString option = GetOption();
122  printf("Starting (begin) h1analysis with process option: %s\n",option.Data());
123 
124  //process cases with event list
125  fillList = kFALSE;
126  useList = kFALSE;
127  if (fChain) fChain->SetEntryList(0);
128  delete gDirectory->GetList()->FindObject("elist");
129 
130  // case when one creates/fills the event list
131  if (option.Contains("fillList")) {
132  fillList = kTRUE;
133  elist = new TEntryList("elist","H1 selection from Cut");
134  // Add to the input list for processing in PROOF, if needed
135  if (fInput) {
136  fInput->Add(new TNamed("fillList",""));
137  fInput->Add(elist);
138  }
139  } else elist = 0;
140 
141  // case when one uses the event list generated in a previous call
142  if (option.Contains("useList")) {
143  useList = kTRUE;
144  if (fInput) {
145  tree->SetEntryList(elist);
146  TFile f("elist.root");
147  elist = (TEntryList*)f.Get("elist");
148  if (elist) elist->SetDirectory(0); //otherwise the file destructor will delete elist
149  } else {
150  // Option "useList" not supported in PROOF directly
151  Warning("Begin", "option 'useList' not supported in PROOF - ignoring");
152  Warning("Begin", "the entry list must be set on the chain *before* calling Process");
153  }
154  }
155 }
156 
157 
158 void h1analysisProxy_SlaveBegin(TTree *tree)
159 {
160 // function called before starting the event loop
161 // -it performs some cleanup
162 // -it creates histograms
163 // -it sets some initialisation for the entry list
164 
165  //initialize the Tree branch addresses
166  Init(tree);
167 
168  //print the option specified in the Process function.
169  TString option = GetOption();
170  printf("Starting (slave) h1analysis with process option: %s\n",option.Data());
171 
172  //create histograms
173  hdmd = new TH1F("hdmd","dm_d",40,0.13,0.17);
174  h2 = new TH2F("h2","ptD0 vs dm_d",30,0.135,0.165,30,-3,6);
175 
176  fOutput->Add(hdmd);
177  fOutput->Add(h2);
178 
179  //process cases with entry list
180  fillList = kFALSE;
181  useList = kFALSE;
182 
183  // case when one creates/fills the entry list
184  if (option.Contains("fillList")) {
185  fillList = kTRUE;
186  // Get the list
187  if (fInput) {
188  if ((elist = (TEntryList *) fInput->FindObject("elist")))
189  // Need to clone to avoid problems when destroying the selector
190  elist = (TEntryList *) elist->Clone();
191  }
192  if (elist)
193  fOutput->Add(elist);
194  else
195  fillList = kFALSE;
196  } else elist = 0;
197 
198  // case when one uses the entry list generated in a previous call
199  if (option.Contains("useList")) {
200  useList = kTRUE;
201  TFile f("elist.root");
202  elist = (TEntryList*)f.Get("elist");
203  if (elist) elist->SetDirectory(0); //otherwise the file destructor will delete elist
204  if (tree) tree->SetEntryList(elist);
205  else {
206  // Option "useList" not supported in PROOF directly
207  Warning("Begin", "option 'useList' not supported in PROOF - ignoring");
208  Warning("Begin", "the entry list must be set on the chain *before* calling Process");
209  }
210  }
211 
212 }
213 
214 Double_t h1analysisProxy() {
215  return 0;
216 }
217 
218 
219 Bool_t h1analysisProxy_Process(Long64_t entry)
220 {
221 // entry is the entry number in the current Tree
222 // Selection function to select D* and D0.
223 
224  //in case one entry list is given in input, the selection has already been done.
225  if (!useList) {
226 
227  float f1 = md0_d;
228  float f2 = md0_d-1.8646;
229  bool test = TMath::Abs(md0_d-1.8646) >= 0.04;
230  if (gDebug>0) fprintf(stderr,"entry #%lld f1=%f f2=%f test=%d\n",
231  fChain->GetReadEntry(),f1,f2,test);
232 
233  if (TMath::Abs(md0_d-1.8646) >= 0.04) return kFALSE;
234  if (ptds_d <= 2.5) return kFALSE;
235  if (TMath::Abs(etads_d) >= 1.5) return kFALSE;
236 
237  int cik = ik-1; //original ik used f77 convention starting at 1
238  int cipi = ipi-1; //original ipi used f77 convention starting at 1
239 
240  f1 = nhitrp[cik];
241  f2 = nhitrp[cipi];
242  test = nhitrp[cik]*nhitrp[cipi] <= 1;
243  if (gDebug>0) fprintf(stderr,"entry #%lld f1=%f f2=%f test=%d\n",
244  fChain->GetReadEntry(),f1,f2,test);
245 
246  if (nhitrp[cik]*nhitrp[cipi] <= 1) return kFALSE;
247  if (rend[cik] -rstart[cik] <= 22) return kFALSE;
248  if (rend[cipi]-rstart[cipi] <= 22) return kFALSE;
249  if (nlhk[cik] <= 0.1) return kFALSE;
250  if (nlhpi[cipi] <= 0.1) return kFALSE;
251  // fix because read-only
252  if (nlhpi[ipis-1] <= 0.1) return kFALSE;
253  if (njets < 1) return kFALSE;
254 
255  }
256  // if option fillList, fill the event list
257  if (fillList) elist->Enter(entry);
258 
259  //fill some histograms
260  hdmd->Fill(dm_d);
261  h2->Fill(dm_d,rpd0_t/0.029979*1.8646/ptd0_d);
262 
263  return kTRUE;
264 }
265 
266 
267 
268 void h1analysisProxy_SlaveTerminate()
269 {
270  // nothing to be done
271  printf("Terminate (slave) h1analysis\n");
272 }
273 
274 
275 void h1analysisProxy_Terminate()
276 {
277  printf("Terminate (final) h1analysis\n");
278 
279  // function called at the end of the event loop
280 
281  hdmd = dynamic_cast<TH1F*>(fOutput->FindObject("hdmd"));
282  h2 = dynamic_cast<TH2F*>(fOutput->FindObject("h2"));
283 
284  if (hdmd == 0 || h2 == 0) {
285  Error("Terminate", "hdmd = %p , h2 = %p", hdmd, h2);
286  return;
287  }
288 
289  //create the canvas for the h1analysis fit
290  gStyle->SetOptFit();
291  TCanvas *c1 = new TCanvas("c1","h1analysis analysis",10,10,800,600);
292  c1->SetBottomMargin(0.15);
293  hdmd->GetXaxis()->SetTitle("m_{K#pi#pi} - m_{K#pi}[GeV/c^{2}]");
294  hdmd->GetXaxis()->SetTitleOffset(1.4);
295 
296  //fit histogram hdmd with function f5 using the loglikelihood option
297  TF1 *f5 = new TF1("f5",fdm5,0.139,0.17,5);
298  f5->SetParameters(1000000, .25, 2000, .1454, .001);
299  hdmd->Fit("f5","lr");
300 
301  //create the canvas for tau d0
302  gStyle->SetOptFit(0);
303  gStyle->SetOptStat(1100);
304  TCanvas *c2 = new TCanvas("c2","tauD0",100,100,800,600);
305  c2->SetGrid();
306  c2->SetBottomMargin(0.15);
307 
308  // Project slices of 2-d histogram h2 along X , then fit each slice
309  // with function f2 and make a histogram for each fit parameter
310  // Note that the generated histograms are added to the list of objects
311  // in the current directory.
312  TF1 *f2 = new TF1("f2",fdm2,0.139,0.17,2);
313  f2->SetParameters(10000, 10);
314  h2->FitSlicesX(f2,0,-1,1,"qln");
315  TH1D *h2_1 = (TH1D*)gDirectory->Get("h2_1");
316  h2_1->GetXaxis()->SetTitle("#tau[ps]");
317  h2_1->SetMarkerStyle(21);
318  h2_1->Draw();
319  c2->Update();
320  TLine *line = new TLine(0,0,0,c2->GetUymax());
321  line->Draw();
322 
323  // Have the number of entries on the first histogram (to cross check when running
324  // with entry lists)
325  TPaveStats *psdmd = (TPaveStats *)hdmd->GetListOfFunctions()->FindObject("stats");
326  psdmd->SetOptStat(1110);
327  c1->Modified();
328 
329  //save the entry list to a Root file if one was produced
330  if (fillList) {
331  elist = dynamic_cast<TEntryList*>(fOutput->FindObject("elist"));
332  if (elist) {
333  TFile efile("elist.root","recreate");
334  elist->Write();
335  } else {
336  Error("Terminate", "entry list requested but not found in output");
337  }
338  }
339 }
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
long long Long64_t
Definition: RtypesCore.h:69
R__EXTERN TStyle * gStyle
Definition: TStyle.h:423
tuple f2
Definition: surfaces.py:24
Double_t fdm2(Double_t *xx, Double_t *par)
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.
Basic string class.
Definition: TString.h:137
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
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
The histogram statistics painter class.
Definition: TPaveStats.h:28
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
TFile * f
virtual void SetDirectory(TDirectory *dir)
Add reference to directory dir. dir can be 0.
const char * Data() const
Definition: TString.h:349
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:33
virtual void SetGrid(Int_t valuex=1, Int_t valuey=1)
Definition: TPad.h:326
void Init(TClassEdit::TInterpreterLookupHelper *helper)
Definition: TClassEdit.cxx:118
TH2D * h2
Definition: fit2dHist.C:45
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition: TNamed.cxx:63
void Error(const char *location, const char *msgfmt,...)
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
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
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:256
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
A simple line.
Definition: TLine.h:41
virtual void SetMarkerStyle(Style_t mstyle=1)
Definition: TAttMarker.h:53
void Warning(const char *location, const char *msgfmt,...)
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:613
The Canvas class.
Definition: TCanvas.h:48
virtual void SetEntryList(TEntryList *list, Option_t *opt="")
Set an EntryList.
Definition: TTree.cxx:8141
double Double_t
Definition: RtypesCore.h:55
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
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)
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
1-Dim function class
Definition: TF1.h:149
TF1 * f1
Definition: legend1.C:11
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
R__EXTERN Int_t gDebug
Definition: Rtypes.h:128
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
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
void SetOptStat(Int_t stat=1)
Set the stat option.
Definition: TPaveStats.cxx:303