00001 // Example of analysis class for the H1 data. 00002 // ========================================= 00003 // 00004 // This file uses 4 large data sets from the H1 collaboration at DESY Hamburg. 00005 // One can access these data sets (277 MBytes) from the standard Root web site 00006 // at: ftp://root.cern.ch/root/h1analysis/ 00007 // The Physics plots below generated by this example cannot be produced when 00008 // using smaller data sets. 00009 // 00010 // There are several ways to analyze data stored in a Root Tree 00011 // -Using TTree::Draw: This is very convenient and efficient for small tasks. 00012 // A TTree::Draw call produces one histogram at the time. The histogram 00013 // is automatically generated. The selection expression may be specified 00014 // in the command line. 00015 // 00016 // -Using the TTreeViewer: This is a graphical interface to TTree::Draw 00017 // with the same functionality. 00018 // 00019 // -Using the code generated by TTree::MakeClass: In this case, the user 00020 // creates an instance of the analysis class. He has the control over 00021 // the event loop and he can generate an unlimited number of histograms. 00022 // 00023 // -Using the code generated by TTree::MakeSelector. Like for the code 00024 // generated by TTree::MakeClass, the user can do complex analysis. 00025 // However, he cannot control the event loop. The event loop is controlled 00026 // by TTree::Process called by the user. This solution is illustrated 00027 // by the current code. The advantage of this method is that it can be run 00028 // in a parallel environment using PROOF (the Parallel Root Facility). 00029 // 00030 // A chain of 4 files (originally converted from PAW ntuples) is used 00031 // to illustrate the various ways to loop on Root data sets. 00032 // Each data set contains a Root Tree named "h42" 00033 // The class definition in h1analysis.h has been generated automatically 00034 // by the Root utility TTree::MakeSelector using one of the files with the 00035 // following statement: 00036 // h42->MakeSelector("h1analysis"); 00037 // This produces two files: h1analysis.h and h1analysis.C (skeleton of this file) 00038 // The h1analysis class is derived from the Root class TSelector. 00039 // 00040 // The following members functions are called by the TTree::Process functions. 00041 // Begin(): called everytime a loop on the tree starts. 00042 // a convenient place to create your histograms. 00043 // SlaveBegin(): 00044 // 00045 // Notify(): This function is called at the first entry of a new Tree 00046 // in a chain. 00047 // Process(): called to analyze each entry. 00048 // 00049 // SlaveTerminate(): 00050 // 00051 // Terminate(): called at the end of a loop on a TTree. 00052 // a convenient place to draw/fit your histograms. 00053 // 00054 // To use this file, try the following session 00055 // 00056 // Root > gROOT->Time(); //will show RT & CPU time per command 00057 // 00058 //==> A- create a TChain with the 4 H1 data files 00059 // The chain can be created by executed the short macro h1chain.C below: 00060 // { 00061 // TChain chain("h42"); 00062 // chain.Add("$H1/dstarmb.root"); // 21330730 bytes 21920 events 00063 // chain.Add("$H1/dstarp1a.root"); // 71464503 bytes 73243 events 00064 // chain.Add("$H1/dstarp1b.root"); // 83827959 bytes 85597 events 00065 // chain.Add("$H1/dstarp2.root"); // 100675234 bytes 103053 events 00066 // //where $H1 is a system symbol pointing to the H1 data directory. 00067 // } 00068 // 00069 // Root > .x h1chain.C 00070 // 00071 //==> B- loop on all events 00072 // Root > chain.Process("h1analysis.C") 00073 // 00074 //==> C- same as B, but in addition fill the entry list with selected entries. 00075 // The entry list is saved to a file "elist.root" by the Terminate function. 00076 // To see the list of selected events, you can do elist->Print("all"). 00077 // The selection function has selected 7525 events out of the 283813 events 00078 // in the chain of files. (2.65 per cent) 00079 // Root > chain.Process("h1analysis.C","fillList") 00080 // 00081 //==> D- Process only entries in the entry list 00082 // The entry list is read from the file in elist.root generated by step C 00083 // Root > chain.Process("h1analysis.C","useList") 00084 // 00085 //==> E- the above steps have been executed via the interpreter. 00086 // You can repeat the steps B, C and D using the script compiler 00087 // by replacing "h1analysis.C" by "h1analysis.C+" or "h1analysis.C++" 00088 // 00089 // in a new session with ,eg: 00090 // 00091 //==> F- Create the chain as in A, then execute 00092 // Root > chain.Process("h1analysis.C+","useList") 00093 // 00094 // The commands executed with the 4 different methods B,C,D and E 00095 // produce two canvases shown below: 00096 // begin_html <a href="gif/h1analysis_dstar.gif" >the Dstar plot</a> end_html 00097 // begin_html <a href="gif/h1analysis_tau.gif" >the Tau D0 plot</a> end_html 00098 // 00099 // The same analysis can be run on PROOF. For a quick try start a PROOF-Lite 00100 // session 00101 // 00102 // Root > TProof *p = TProof::Open("") 00103 // 00104 // create (if mot already done) the chain by executing the 'h1chain.C' macro 00105 // mentioned above, and then tell ROOT to use PROOF to process the chain: 00106 // 00107 // Root > chain.SetProof() 00108 // 00109 // You can then repeat step B above. Step C can also be executed in PROOF. However, 00110 // step D cannot be executed in PROOF as in the local session (i.e. just passing 00111 // option 'useList'): to use the entry list you have to 00112 // 00113 //==> G- Load first in the session the list form the file 00114 // 00115 // Root > TFile f("elist.root") 00116 // Root > TEntryList *elist = (TEntryList *) f.Get("elist") 00117 // 00118 // set it on the chain: 00119 // 00120 // Root > chain.SetEntryList(elist) 00121 // 00122 // call Process as in step B. Of course this works also for local processing. 00123 // 00124 //Author: Rene Brun 00125 00126 #include "h1analysis.h" 00127 #include "TH2.h" 00128 #include "TF1.h" 00129 #include "TStyle.h" 00130 #include "TCanvas.h" 00131 #include "TPaveStats.h" 00132 #include "TLine.h" 00133 #include "TMath.h" 00134 00135 const Double_t dxbin = (0.17-0.13)/40; // Bin-width 00136 const Double_t sigma = 0.0012; 00137 00138 //_____________________________________________________________________ 00139 Double_t fdm5(Double_t *xx, Double_t *par) 00140 { 00141 Double_t x = xx[0]; 00142 if (x <= 0.13957) return 0; 00143 Double_t xp3 = (x-par[3])*(x-par[3]); 00144 Double_t res = dxbin*(par[0]*TMath::Power(x-0.13957, par[1]) 00145 + par[2] / 2.5066/par[4]*TMath::Exp(-xp3/2/par[4]/par[4])); 00146 return res; 00147 } 00148 00149 //_____________________________________________________________________ 00150 Double_t fdm2(Double_t *xx, Double_t *par) 00151 { 00152 Double_t x = xx[0]; 00153 if (x <= 0.13957) return 0; 00154 Double_t xp3 = (x-0.1454)*(x-0.1454); 00155 Double_t res = dxbin*(par[0]*TMath::Power(x-0.13957, 0.25) 00156 + par[1] / 2.5066/sigma*TMath::Exp(-xp3/2/sigma/sigma)); 00157 return res; 00158 } 00159 00160 //_____________________________________________________________________ 00161 void h1analysis::Begin(TTree * /*tree*/) 00162 { 00163 // function called before starting the event loop 00164 // -it performs some cleanup 00165 // -it creates histograms 00166 // -it sets some initialisation for the entry list 00167 00168 //print the option specified in the Process function. 00169 TString option = GetOption(); 00170 Info("Begin", "starting h1analysis with process option: %s", option.Data()); 00171 00172 //process cases with entry list 00173 if (fChain) fChain->SetEntryList(0); 00174 delete gDirectory->GetList()->FindObject("elist"); 00175 00176 // case when one creates/fills the entry list 00177 if (option.Contains("fillList")) { 00178 fillList = kTRUE; 00179 // Add to the input list for processing in PROOF, if needed 00180 if (fInput) { 00181 fInput->Add(new TNamed("fillList","")); 00182 // We send a clone to avoid double deletes when importing the result 00183 fInput->Add(new TEntryList("elist", "H1 selection from Cut")); 00184 } 00185 } 00186 // case when one uses the entry list generated in a previous call 00187 if (option.Contains("useList")) { 00188 useList = kTRUE; 00189 if (fInput) { 00190 // Option "useList" not supported in PROOF directly 00191 Warning("Begin", "option 'useList' not supported in PROOF - ignoring"); 00192 Warning("Begin", "the entry list must be set on the chain *before* calling Process"); 00193 } else { 00194 TFile f("elist.root"); 00195 elist = (TEntryList*)f.Get("elist"); 00196 if (elist) elist->SetDirectory(0); //otherwise the file destructor will delete elist 00197 } 00198 } 00199 } 00200 00201 //_____________________________________________________________________ 00202 void h1analysis::SlaveBegin(TTree *tree) 00203 { 00204 // function called before starting the event loop 00205 // -it performs some cleanup 00206 // -it creates histograms 00207 // -it sets some initialisation for the entry list 00208 00209 //initialize the Tree branch addresses 00210 Init(tree); 00211 00212 //print the option specified in the Process function. 00213 TString option = GetOption(); 00214 Info("SlaveBegin", 00215 "starting h1analysis with process option: %s (tree: %p)", option.Data(), tree); 00216 00217 //create histograms 00218 hdmd = new TH1F("hdmd","dm_d",40,0.13,0.17); 00219 h2 = new TH2F("h2","ptD0 vs dm_d",30,0.135,0.165,30,-3,6); 00220 00221 fOutput->Add(hdmd); 00222 fOutput->Add(h2); 00223 00224 // Entry list stuff (re-parse option because on PROOF only SlaveBegin is called) 00225 if (option.Contains("fillList")) { 00226 fillList = kTRUE; 00227 // Get the list 00228 if (fInput) { 00229 if ((elist = (TEntryList *) fInput->FindObject("elist"))) 00230 // Need to clone to avoid problems when destroying the selector 00231 elist = (TEntryList *) elist->Clone(); 00232 } 00233 if (elist) 00234 fOutput->Add(elist); 00235 else 00236 fillList = kFALSE; 00237 } 00238 } 00239 00240 //_____________________________________________________________________ 00241 Bool_t h1analysis::Process(Long64_t entry) 00242 { 00243 // entry is the entry number in the current Tree 00244 // Selection function to select D* and D0. 00245 00246 //in case one entry list is given in input, the selection has already been done. 00247 if (!useList) { 00248 // Read only the necessary branches to select entries. 00249 // return as soon as a bad entry is detected 00250 // to read complete event, call fChain->GetTree()->GetEntry(entry) 00251 b_md0_d->GetEntry(entry); if (TMath::Abs(md0_d-1.8646) >= 0.04) return kFALSE; 00252 b_ptds_d->GetEntry(entry); if (ptds_d <= 2.5) return kFALSE; 00253 b_etads_d->GetEntry(entry); if (TMath::Abs(etads_d) >= 1.5) return kFALSE; 00254 b_ik->GetEntry(entry); ik--; //original ik used f77 convention starting at 1 00255 b_ipi->GetEntry(entry); ipi--; 00256 b_ntracks->GetEntry(entry); 00257 b_nhitrp->GetEntry(entry); 00258 if (nhitrp[ik]*nhitrp[ipi] <= 1) return kFALSE; 00259 b_rend->GetEntry(entry); 00260 b_rstart->GetEntry(entry); 00261 if (rend[ik] -rstart[ik] <= 22) return kFALSE; 00262 if (rend[ipi]-rstart[ipi] <= 22) return kFALSE; 00263 b_nlhk->GetEntry(entry); if (nlhk[ik] <= 0.1) return kFALSE; 00264 b_nlhpi->GetEntry(entry); if (nlhpi[ipi] <= 0.1) return kFALSE; 00265 b_ipis->GetEntry(entry); ipis--; if (nlhpi[ipis] <= 0.1) return kFALSE; 00266 b_njets->GetEntry(entry); if (njets < 1) return kFALSE; 00267 } 00268 // if option fillList, fill the entry list 00269 if (fillList) elist->Enter(entry); 00270 00271 // to read complete event, call fChain->GetTree()->GetEntry(entry) 00272 // read branches not processed in ProcessCut 00273 b_dm_d->GetEntry(entry); //read branch holding dm_d 00274 b_rpd0_t->GetEntry(entry); //read branch holding rpd0_t 00275 b_ptd0_d->GetEntry(entry); //read branch holding ptd0_d 00276 00277 //fill some histograms 00278 hdmd->Fill(dm_d); 00279 h2->Fill(dm_d,rpd0_t/0.029979*1.8646/ptd0_d); 00280 00281 return kTRUE; 00282 } 00283 00284 00285 //_____________________________________________________________________ 00286 void h1analysis::SlaveTerminate() 00287 { 00288 // nothing to be done 00289 } 00290 00291 //_____________________________________________________________________ 00292 void h1analysis::Terminate() 00293 { 00294 // function called at the end of the event loop 00295 00296 hdmd = dynamic_cast<TH1F*>(fOutput->FindObject("hdmd")); 00297 h2 = dynamic_cast<TH2F*>(fOutput->FindObject("h2")); 00298 00299 if (hdmd == 0 || h2 == 0) { 00300 Error("Terminate", "hdmd = %p , h2 = %p", hdmd, h2); 00301 return; 00302 } 00303 00304 //create the canvas for the h1analysis fit 00305 gStyle->SetOptFit(); 00306 TCanvas *c1 = new TCanvas("c1","h1analysis analysis",10,10,800,600); 00307 c1->SetBottomMargin(0.15); 00308 hdmd->GetXaxis()->SetTitle("m_{K#pi#pi} - m_{K#pi}[GeV/c^{2}]"); 00309 hdmd->GetXaxis()->SetTitleOffset(1.4); 00310 00311 //fit histogram hdmd with function f5 using the loglikelihood option 00312 if (gROOT->GetListOfFunctions()->FindObject("f5")) 00313 delete gROOT->GetFunction("f5"); 00314 TF1 *f5 = new TF1("f5",fdm5,0.139,0.17,5); 00315 f5->SetParameters(1000000, .25, 2000, .1454, .001); 00316 hdmd->Fit("f5","lr"); 00317 00318 //create the canvas for tau d0 00319 gStyle->SetOptFit(0); 00320 gStyle->SetOptStat(1100); 00321 TCanvas *c2 = new TCanvas("c2","tauD0",100,100,800,600); 00322 c2->SetGrid(); 00323 c2->SetBottomMargin(0.15); 00324 00325 // Project slices of 2-d histogram h2 along X , then fit each slice 00326 // with function f2 and make a histogram for each fit parameter 00327 // Note that the generated histograms are added to the list of objects 00328 // in the current directory. 00329 if (gROOT->GetListOfFunctions()->FindObject("f2")) 00330 delete gROOT->GetFunction("f2"); 00331 TF1 *f2 = new TF1("f2",fdm2,0.139,0.17,2); 00332 f2->SetParameters(10000, 10); 00333 h2->FitSlicesX(f2,0,-1,1,"qln"); 00334 TH1D *h2_1 = (TH1D*)gDirectory->Get("h2_1"); 00335 h2_1->GetXaxis()->SetTitle("#tau[ps]"); 00336 h2_1->SetMarkerStyle(21); 00337 h2_1->Draw(); 00338 c2->Update(); 00339 TLine *line = new TLine(0,0,0,c2->GetUymax()); 00340 line->Draw(); 00341 00342 // Have the number of entries on the first histogram (to cross check when running 00343 // with entry lists) 00344 TPaveStats *psdmd = (TPaveStats *)hdmd->GetListOfFunctions()->FindObject("stats"); 00345 psdmd->SetOptStat(1110); 00346 c1->Modified(); 00347 00348 //save the entry list to a Root file if one was produced 00349 if (fillList) { 00350 elist = dynamic_cast<TEntryList*>(fOutput->FindObject("elist")); 00351 if (elist) { 00352 TFile efile("elist.root","recreate"); 00353 elist->Write(); 00354 } else { 00355 Error("Terminate", "entry list requested but not found in output"); 00356 } 00357 } 00358 }