ROOT logo
// Example of analysis class for the H1 data.
// =========================================
//
// This file uses 4 large data sets from the H1 collaboration at DESY Hamburg.
// One can access these data sets (277 MBytes) from the standard Root web site
// at:       ftp://root.cern.ch/root/h1analysis/
// The Physics plots below generated by this example cannot be produced when
// using smaller data sets.
//
// There are several ways to analyze data stored in a Root Tree
//   -Using TTree::Draw: This is very convenient and efficient for small tasks.
//     A TTree::Draw call produces one histogram at the time. The histogram
//     is automatically generated. The selection expression may be specified
//     in the command line.
//
//   -Using the TTreeViewer: This is a graphical interface to TTree::Draw
//      with the same functionality.
//
//   -Using the code generated by TTree::MakeClass: In this case, the user
//      creates an instance of the analysis class. He has the control over
//      the event loop and he can generate an unlimited number of histograms.
//
//   -Using the code generated by TTree::MakeSelector. Like for the code
//      generated by TTree::MakeClass, the user can do complex analysis.
//      However, he cannot control the event loop. The event loop is controlled
//      by TTree::Process called by the user. This solution is illustrated
//      by the current code. The advantage of this method is that it can be run
//      in a parallel environment using PROOF (the Parallel Root Facility).
//
// A chain of 4 files (originally converted from PAW ntuples) is used
// to illustrate the various ways to loop on Root data sets.
// Each data set contains a Root Tree named "h42"
// The class definition in h1analysis.h has been generated automatically
// by the Root utility TTree::MakeSelector using one of the files with the
// following statement:
//       h42->MakeSelector("h1analysis");
// This produces two files: h1analysis.h and h1analysis.C (skeleton of this file)
// The h1analysis class is derived from the Root class TSelector.
//
// The following members functions are called by the TTree::Process functions.
//    Begin():      called everytime a loop on the tree starts.
//                  a convenient place to create your histograms.
//    SlaveBegin():
//
//    Notify():     This function is called at the first entry of a new Tree
//                  in a chain.
//    Process():    called to analyze each entry.
//
//    SlaveTerminate():
//
//    Terminate():    called at the end of a loop on a TTree.
//                  a convenient place to draw/fit your histograms.
//
//   To use this file, try the following session
//
// Root > gROOT->Time(); //will show RT & CPU time per command
//
//==>   A-  create a TChain with the 4 H1 data files
// The chain can be created by executed the short macro h1chain.C below:
// {
//   TChain chain("h42");
//   chain.Add("$H1/dstarmb.root");  //  21330730 bytes  21920 events
//   chain.Add("$H1/dstarp1a.root"); //  71464503 bytes  73243 events
//   chain.Add("$H1/dstarp1b.root"); //  83827959 bytes  85597 events
//   chain.Add("$H1/dstarp2.root");  // 100675234 bytes 103053 events
//   //where $H1 is a system symbol pointing to the H1 data directory.
// }
//
// Root > .x h1chain.C
//
//==>   B- loop on all events
// Root > chain.Process("h1analysis.C")
//
//==>   C- same as B, but in addition fill the entry list with selected entries.
// The entry list is saved to a file "elist.root" by the Terminate function.
// To see the list of selected events, you can do elist->Print("all").
// The selection function has selected 7525 events out of the 283813 events
// in the chain of files. (2.65 per cent)
// Root > chain.Process("h1analysis.C","fillList")
//
//==>   D- Process only entries in the entry list
// The entry list is read from the file in elist.root generated by step C
// Root > chain.Process("h1analysis.C","useList")
//
//==>   E- the above steps have been executed via the interpreter.
//      You can repeat the steps B, C and D using the script compiler
//      by replacing "h1analysis.C" by "h1analysis.C+" or "h1analysis.C++"
//
// in a new session with ,eg:
//
//==>   F- Create the chain as in A, then execute
// Root > chain.Process("h1analysis.C+","useList")
//
// The commands executed with the 4 different methods B,C,D and E
// produce two canvases shown below:
// begin_html <a href="gif/h1analysis_dstar.gif" >the Dstar plot</a> end_html
// begin_html <a href="gif/h1analysis_tau.gif" >the Tau D0 plot</a> end_html
//
// The same analysis can be run on PROOF. For a quick try start a PROOF-Lite
// session
//
// Root > TProof *p = TProof::Open("")
//
// create (if mot already done) the chain by executing the 'h1chain.C' macro
// mentioned above, and then tell ROOT to use PROOF to process the chain:
//
// Root > chain.SetProof()
//
// You can then repeat step B above. Step C can also be executed in PROOF. However,
// step D cannot be executed in PROOF as in the local session (i.e. just passing
// option 'useList'): to use the entry list you have to
//
//==>   G- Load first in the session the list form the file
//
// Root > TFile f("elist.root")
// Root > TEntryList *elist = (TEntryList *) f.Get("elist")
//
// set it on the chain:
//
// Root > chain.SetEntryList(elist)
//
// call Process as in step B. Of course this works also for local processing.
//
//Author: Rene Brun

#include "h1analysis.h"
#include "TH2.h"
#include "TF1.h"
#include "TStyle.h"
#include "TCanvas.h"
#include "TPaveStats.h"
#include "TLine.h"
#include "TMath.h"

const Double_t dxbin = (0.17-0.13)/40;   // Bin-width
const Double_t sigma = 0.0012;

//_____________________________________________________________________
Double_t fdm5(Double_t *xx, Double_t *par)
{
   Double_t x = xx[0];
   if (x <= 0.13957) return 0;
   Double_t xp3 = (x-par[3])*(x-par[3]);
   Double_t res = dxbin*(par[0]*TMath::Power(x-0.13957, par[1])
       + par[2] / 2.5066/par[4]*TMath::Exp(-xp3/2/par[4]/par[4]));
   return res;
}

//_____________________________________________________________________
Double_t fdm2(Double_t *xx, Double_t *par)
{
   Double_t x = xx[0];
   if (x <= 0.13957) return 0;
   Double_t xp3 = (x-0.1454)*(x-0.1454);
   Double_t res = dxbin*(par[0]*TMath::Power(x-0.13957, 0.25)
       + par[1] / 2.5066/sigma*TMath::Exp(-xp3/2/sigma/sigma));
   return res;
}

//_____________________________________________________________________
void h1analysis::Begin(TTree * /*tree*/)
{
// function called before starting the event loop
//  -it performs some cleanup
//  -it creates histograms
//  -it sets some initialisation for the entry list

   //print the option specified in the Process function.
   TString option = GetOption();
   Info("Begin", "starting h1analysis with process option: %s", option.Data());

   //process cases with entry list
   if (fChain) fChain->SetEntryList(0);
   delete gDirectory->GetList()->FindObject("elist");

   // case when one creates/fills the entry list
   if (option.Contains("fillList")) {
      fillList = kTRUE;
      elist = new TEntryList("elist", "H1 selection from Cut");
      // Add to the input list for processing in PROOF, if needed
      if (fInput) {
         fInput->Add(new TNamed("fillList",""));
         fInput->Add(elist);
      }
   }
   // case when one uses the entry list generated in a previous call
   if (option.Contains("useList")) {
      useList  = kTRUE;
      if (fInput) {
         // Option "useList" not supported in PROOF directly
         Warning("Begin", "option 'useList' not supported in PROOF - ignoring");
         Warning("Begin", "the entry list must be set on the chain *before* calling Process");
      } else {
         TFile f("elist.root");
         elist = (TEntryList*)f.Get("elist");
         if (elist) elist->SetDirectory(0); //otherwise the file destructor will delete elist
      }
   }
}

//_____________________________________________________________________
void h1analysis::SlaveBegin(TTree *tree)
{
// function called before starting the event loop
//  -it performs some cleanup
//  -it creates histograms
//  -it sets some initialisation for the entry list

   //initialize the Tree branch addresses
   Init(tree);

   //print the option specified in the Process function.
   TString option = GetOption();
   Info("SlaveBegin",
        "starting h1analysis with process option: %s (tree: %p)", option.Data(), tree);

   //create histograms
   hdmd = new TH1F("hdmd","dm_d",40,0.13,0.17);
   h2   = new TH2F("h2","ptD0 vs dm_d",30,0.135,0.165,30,-3,6);

   fOutput->Add(hdmd);
   fOutput->Add(h2);

   // Entry list stuff (re-parse option because on PROOF only SlaveBegin is called)
   if (option.Contains("fillList")) {
      fillList = kTRUE;
      if (fInput) {
         elist = (TEntryList *) fInput->FindObject("elist");
         fInput->Remove(elist);
      }
      if (elist)
         fOutput->Add(elist);
      else
         fillList = kFALSE;
   }
}

//_____________________________________________________________________
Bool_t h1analysis::Process(Long64_t entry)
{
// entry is the entry number in the current Tree
// Selection function to select D* and D0.

   //in case one entry list is given in input, the selection has already been done.
   if (!useList) {
      // Read only the necessary branches to select entries.
      // return as soon as a bad entry is detected
      // to read complete event, call fChain->GetTree()->GetEntry(entry)
      b_md0_d->GetEntry(entry);   if (TMath::Abs(md0_d-1.8646) >= 0.04) return kFALSE;
      b_ptds_d->GetEntry(entry);  if (ptds_d <= 2.5) return kFALSE;
      b_etads_d->GetEntry(entry); if (TMath::Abs(etads_d) >= 1.5) return kFALSE;
      b_ik->GetEntry(entry);  ik--; //original ik used f77 convention starting at 1
      b_ipi->GetEntry(entry); ipi--;
      b_ntracks->GetEntry(entry);
      b_nhitrp->GetEntry(entry);
      if (nhitrp[ik]*nhitrp[ipi] <= 1) return kFALSE;
      b_rend->GetEntry(entry);
      b_rstart->GetEntry(entry);
      if (rend[ik] -rstart[ik]  <= 22) return kFALSE;
      if (rend[ipi]-rstart[ipi] <= 22) return kFALSE;
      b_nlhk->GetEntry(entry);         if (nlhk[ik] <= 0.1)    return kFALSE;
      b_nlhpi->GetEntry(entry);        if (nlhpi[ipi] <= 0.1)  return kFALSE;
      b_ipis->GetEntry(entry); ipis--; if (nlhpi[ipis] <= 0.1) return kFALSE;
      b_njets->GetEntry(entry);        if (njets < 1)          return kFALSE;
   }
   // if option fillList, fill the entry list
   if (fillList) elist->Enter(entry);

   // to read complete event, call fChain->GetTree()->GetEntry(entry)
   // read branches not processed in ProcessCut
   b_dm_d->GetEntry(entry);         //read branch holding dm_d
   b_rpd0_t->GetEntry(entry);       //read branch holding rpd0_t
   b_ptd0_d->GetEntry(entry);       //read branch holding ptd0_d

   //fill some histograms
   hdmd->Fill(dm_d);
   h2->Fill(dm_d,rpd0_t/0.029979*1.8646/ptd0_d);

   return kTRUE;
}


//_____________________________________________________________________
void h1analysis::SlaveTerminate()
{
   // nothing to be done
}

//_____________________________________________________________________
void h1analysis::Terminate()
{
// function called at the end of the event loop

   hdmd = dynamic_cast<TH1F*>(fOutput->FindObject("hdmd"));
   h2 = dynamic_cast<TH2F*>(fOutput->FindObject("h2"));

   if (hdmd == 0 || h2 == 0) {
      Error("Terminate", "hdmd = %p , h2 = %p", hdmd, h2);
      return;
   }

   //create the canvas for the h1analysis fit
   gStyle->SetOptFit();
   TCanvas *c1 = new TCanvas("c1","h1analysis analysis",10,10,800,600);
   c1->SetBottomMargin(0.15);
   hdmd->GetXaxis()->SetTitle("m_{K#pi#pi} - m_{K#pi}[GeV/c^{2}]");
   hdmd->GetXaxis()->SetTitleOffset(1.4);

   //fit histogram hdmd with function f5 using the loglikelihood option
   if (gROOT->GetListOfFunctions()->FindObject("f5"))
      delete gROOT->GetFunction("f5");
   TF1 *f5 = new TF1("f5",fdm5,0.139,0.17,5);
   f5->SetParameters(1000000, .25, 2000, .1454, .001);
   hdmd->Fit("f5","lr");

   //create the canvas for tau d0
   gStyle->SetOptFit(0);
   gStyle->SetOptStat(1100);
   TCanvas *c2 = new TCanvas("c2","tauD0",100,100,800,600);
   c2->SetGrid();
   c2->SetBottomMargin(0.15);

   // Project slices of 2-d histogram h2 along X , then fit each slice
   // with function f2 and make a histogram for each fit parameter
   // Note that the generated histograms are added to the list of objects
   // in the current directory.
   if (gROOT->GetListOfFunctions()->FindObject("f2"))
      delete gROOT->GetFunction("f2");
   TF1 *f2 = new TF1("f2",fdm2,0.139,0.17,2);
   f2->SetParameters(10000, 10);
   h2->FitSlicesX(f2,0,-1,1,"qln");
   TH1D *h2_1 = (TH1D*)gDirectory->Get("h2_1");
   h2_1->GetXaxis()->SetTitle("#tau[ps]");
   h2_1->SetMarkerStyle(21);
   h2_1->Draw();
   c2->Update();
   TLine *line = new TLine(0,0,0,c2->GetUymax());
   line->Draw();

   // Have the number of entries on the first histogram (to cross check when running
   // with entry lists)
   TPaveStats *psdmd = (TPaveStats *)hdmd->GetListOfFunctions()->FindObject("stats");
   psdmd->SetOptStat(1110);
   c1->Modified();

   //save the entry list to a Root file if one was produced
   if (fillList) {
      elist = dynamic_cast<TEntryList*>(fOutput->FindObject("elist"));
      if (elist) {
         TFile efile("elist.root","recreate");
         elist->Write();
      } else {
         Error("Terminate", "entry list requested but not found in output");
      }
   }
}
 h1analysis.C:1
 h1analysis.C:2
 h1analysis.C:3
 h1analysis.C:4
 h1analysis.C:5
 h1analysis.C:6
 h1analysis.C:7
 h1analysis.C:8
 h1analysis.C:9
 h1analysis.C:10
 h1analysis.C:11
 h1analysis.C:12
 h1analysis.C:13
 h1analysis.C:14
 h1analysis.C:15
 h1analysis.C:16
 h1analysis.C:17
 h1analysis.C:18
 h1analysis.C:19
 h1analysis.C:20
 h1analysis.C:21
 h1analysis.C:22
 h1analysis.C:23
 h1analysis.C:24
 h1analysis.C:25
 h1analysis.C:26
 h1analysis.C:27
 h1analysis.C:28
 h1analysis.C:29
 h1analysis.C:30
 h1analysis.C:31
 h1analysis.C:32
 h1analysis.C:33
 h1analysis.C:34
 h1analysis.C:35
 h1analysis.C:36
 h1analysis.C:37
 h1analysis.C:38
 h1analysis.C:39
 h1analysis.C:40
 h1analysis.C:41
 h1analysis.C:42
 h1analysis.C:43
 h1analysis.C:44
 h1analysis.C:45
 h1analysis.C:46
 h1analysis.C:47
 h1analysis.C:48
 h1analysis.C:49
 h1analysis.C:50
 h1analysis.C:51
 h1analysis.C:52
 h1analysis.C:53
 h1analysis.C:54
 h1analysis.C:55
 h1analysis.C:56
 h1analysis.C:57
 h1analysis.C:58
 h1analysis.C:59
 h1analysis.C:60
 h1analysis.C:61
 h1analysis.C:62
 h1analysis.C:63
 h1analysis.C:64
 h1analysis.C:65
 h1analysis.C:66
 h1analysis.C:67
 h1analysis.C:68
 h1analysis.C:69
 h1analysis.C:70
 h1analysis.C:71
 h1analysis.C:72
 h1analysis.C:73
 h1analysis.C:74
 h1analysis.C:75
 h1analysis.C:76
 h1analysis.C:77
 h1analysis.C:78
 h1analysis.C:79
 h1analysis.C:80
 h1analysis.C:81
 h1analysis.C:82
 h1analysis.C:83
 h1analysis.C:84
 h1analysis.C:85
 h1analysis.C:86
 h1analysis.C:87
 h1analysis.C:88
 h1analysis.C:89
 h1analysis.C:90
 h1analysis.C:91
 h1analysis.C:92
 h1analysis.C:93
 h1analysis.C:94
 h1analysis.C:95
 h1analysis.C:96
 h1analysis.C:97
 h1analysis.C:98
 h1analysis.C:99
 h1analysis.C:100
 h1analysis.C:101
 h1analysis.C:102
 h1analysis.C:103
 h1analysis.C:104
 h1analysis.C:105
 h1analysis.C:106
 h1analysis.C:107
 h1analysis.C:108
 h1analysis.C:109
 h1analysis.C:110
 h1analysis.C:111
 h1analysis.C:112
 h1analysis.C:113
 h1analysis.C:114
 h1analysis.C:115
 h1analysis.C:116
 h1analysis.C:117
 h1analysis.C:118
 h1analysis.C:119
 h1analysis.C:120
 h1analysis.C:121
 h1analysis.C:122
 h1analysis.C:123
 h1analysis.C:124
 h1analysis.C:125
 h1analysis.C:126
 h1analysis.C:127
 h1analysis.C:128
 h1analysis.C:129
 h1analysis.C:130
 h1analysis.C:131
 h1analysis.C:132
 h1analysis.C:133
 h1analysis.C:134
 h1analysis.C:135
 h1analysis.C:136
 h1analysis.C:137
 h1analysis.C:138
 h1analysis.C:139
 h1analysis.C:140
 h1analysis.C:141
 h1analysis.C:142
 h1analysis.C:143
 h1analysis.C:144
 h1analysis.C:145
 h1analysis.C:146
 h1analysis.C:147
 h1analysis.C:148
 h1analysis.C:149
 h1analysis.C:150
 h1analysis.C:151
 h1analysis.C:152
 h1analysis.C:153
 h1analysis.C:154
 h1analysis.C:155
 h1analysis.C:156
 h1analysis.C:157
 h1analysis.C:158
 h1analysis.C:159
 h1analysis.C:160
 h1analysis.C:161
 h1analysis.C:162
 h1analysis.C:163
 h1analysis.C:164
 h1analysis.C:165
 h1analysis.C:166
 h1analysis.C:167
 h1analysis.C:168
 h1analysis.C:169
 h1analysis.C:170
 h1analysis.C:171
 h1analysis.C:172
 h1analysis.C:173
 h1analysis.C:174
 h1analysis.C:175
 h1analysis.C:176
 h1analysis.C:177
 h1analysis.C:178
 h1analysis.C:179
 h1analysis.C:180
 h1analysis.C:181
 h1analysis.C:182
 h1analysis.C:183
 h1analysis.C:184
 h1analysis.C:185
 h1analysis.C:186
 h1analysis.C:187
 h1analysis.C:188
 h1analysis.C:189
 h1analysis.C:190
 h1analysis.C:191
 h1analysis.C:192
 h1analysis.C:193
 h1analysis.C:194
 h1analysis.C:195
 h1analysis.C:196
 h1analysis.C:197
 h1analysis.C:198
 h1analysis.C:199
 h1analysis.C:200
 h1analysis.C:201
 h1analysis.C:202
 h1analysis.C:203
 h1analysis.C:204
 h1analysis.C:205
 h1analysis.C:206
 h1analysis.C:207
 h1analysis.C:208
 h1analysis.C:209
 h1analysis.C:210
 h1analysis.C:211
 h1analysis.C:212
 h1analysis.C:213
 h1analysis.C:214
 h1analysis.C:215
 h1analysis.C:216
 h1analysis.C:217
 h1analysis.C:218
 h1analysis.C:219
 h1analysis.C:220
 h1analysis.C:221
 h1analysis.C:222
 h1analysis.C:223
 h1analysis.C:224
 h1analysis.C:225
 h1analysis.C:226
 h1analysis.C:227
 h1analysis.C:228
 h1analysis.C:229
 h1analysis.C:230
 h1analysis.C:231
 h1analysis.C:232
 h1analysis.C:233
 h1analysis.C:234
 h1analysis.C:235
 h1analysis.C:236
 h1analysis.C:237
 h1analysis.C:238
 h1analysis.C:239
 h1analysis.C:240
 h1analysis.C:241
 h1analysis.C:242
 h1analysis.C:243
 h1analysis.C:244
 h1analysis.C:245
 h1analysis.C:246
 h1analysis.C:247
 h1analysis.C:248
 h1analysis.C:249
 h1analysis.C:250
 h1analysis.C:251
 h1analysis.C:252
 h1analysis.C:253
 h1analysis.C:254
 h1analysis.C:255
 h1analysis.C:256
 h1analysis.C:257
 h1analysis.C:258
 h1analysis.C:259
 h1analysis.C:260
 h1analysis.C:261
 h1analysis.C:262
 h1analysis.C:263
 h1analysis.C:264
 h1analysis.C:265
 h1analysis.C:266
 h1analysis.C:267
 h1analysis.C:268
 h1analysis.C:269
 h1analysis.C:270
 h1analysis.C:271
 h1analysis.C:272
 h1analysis.C:273
 h1analysis.C:274
 h1analysis.C:275
 h1analysis.C:276
 h1analysis.C:277
 h1analysis.C:278
 h1analysis.C:279
 h1analysis.C:280
 h1analysis.C:281
 h1analysis.C:282
 h1analysis.C:283
 h1analysis.C:284
 h1analysis.C:285
 h1analysis.C:286
 h1analysis.C:287
 h1analysis.C:288
 h1analysis.C:289
 h1analysis.C:290
 h1analysis.C:291
 h1analysis.C:292
 h1analysis.C:293
 h1analysis.C:294
 h1analysis.C:295
 h1analysis.C:296
 h1analysis.C:297
 h1analysis.C:298
 h1analysis.C:299
 h1analysis.C:300
 h1analysis.C:301
 h1analysis.C:302
 h1analysis.C:303
 h1analysis.C:304
 h1analysis.C:305
 h1analysis.C:306
 h1analysis.C:307
 h1analysis.C:308
 h1analysis.C:309
 h1analysis.C:310
 h1analysis.C:311
 h1analysis.C:312
 h1analysis.C:313
 h1analysis.C:314
 h1analysis.C:315
 h1analysis.C:316
 h1analysis.C:317
 h1analysis.C:318
 h1analysis.C:319
 h1analysis.C:320
 h1analysis.C:321
 h1analysis.C:322
 h1analysis.C:323
 h1analysis.C:324
 h1analysis.C:325
 h1analysis.C:326
 h1analysis.C:327
 h1analysis.C:328
 h1analysis.C:329
 h1analysis.C:330
 h1analysis.C:331
 h1analysis.C:332
 h1analysis.C:333
 h1analysis.C:334
 h1analysis.C:335
 h1analysis.C:336
 h1analysis.C:337
 h1analysis.C:338
 h1analysis.C:339
 h1analysis.C:340
 h1analysis.C:341
 h1analysis.C:342
 h1analysis.C:343
 h1analysis.C:344
 h1analysis.C:345
 h1analysis.C:346
 h1analysis.C:347
 h1analysis.C:348
 h1analysis.C:349
 h1analysis.C:350
 h1analysis.C:351
 h1analysis.C:352
 h1analysis.C:353
 h1analysis.C:354
 h1analysis.C:355
 h1analysis.C:356
 h1analysis.C:357