Loading [MathJax]/extensions/tex2jax.js
Logo ROOT  
Reference Guide
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
h1analysis.C File Reference

Detailed Description

View in nbviewer Open in SWAN 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. They have 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, they 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 every time a loop on the tree starts. A convenient place to create your histograms.
  • Notify(): This function is called at the first entry of a new Tree in a chain.
  • Process(): Called to analyze each entry.
  • 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 sessions

Root > gROOT->Time(); /// will show RT & CPU time per command
#define gROOT
Definition: TROOT.h:406
const char * Root
Definition: TXMLSetup.cxx:46

Case 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.
}
A chain is a collection of files containing TTree objects.
Definition: TChain.h:34

Case B: Loop on all events

Root > chain.Process("h1analysis.C")

Case 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")

Case 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")

Case 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 (see F).

Case F: Create the chain as in A, then execute

Root > chain.Process("h1analysis.C+","useList")

The same analysis can be run on PROOF. For a quick try start a PROOF-Lite session

This class controls a Parallel ROOT Facility, PROOF, cluster.
Definition: TProof.h:316
static TProof * Open(const char *url=0, const char *conffile=0, const char *confdir=0, Int_t loglevel=0)
Start a PROOF session on a specific cluster.
Definition: TProof.cxx:11573

create (if not 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

Case G: Load first in the session the list form the file

Root > TFile f("elist.root")
Root > TEntryList *elist = (TEntryList *) f.Get("elist")
#define f(i)
Definition: RSha256.hxx:104
A List of entry numbers in a TTree or TChain.
Definition: TEntryList.h:26
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition: TFile.h:53

set it on the chain:

Root > chain.SetEntryList(elist)

call Process as in step B. Of course this works also for local processing.

#include "h1analysis.h"
#include "TH2.h"
#include "TF1.h"
#include "TStyle.h"
#include "TBranch.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 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 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
// This is needed when re-processing the object
Reset();
//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
delete gDirectory->GetList()->FindObject("elist");
// case when one creates/fills the entry list
if (option.Contains("fillList")) {
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",""));
// We send a clone to avoid double deletes when importing the result
// This is needed to avoid warnings from output-to-members mapping
elist = 0;
}
Info("Begin", "creating an entry-list");
}
// case when one uses the entry list generated in a previous call
if (option.Contains("useList")) {
if (fInput) {
// In PROOF option "useList" is processed in SlaveBegin and we do not need
// to do anything here
} else {
TFile f("elist.root");
elist = (TEntryList*)f.Get("elist");
if (elist) elist->SetDirectory(0); //otherwise the file destructor will delete elist
}
}
}
{
// 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
//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);
// Entry list stuff (re-parse option because on PROOF only SlaveBegin is called)
if (option.Contains("fillList")) {
// Get the list
if (fInput) {
if ((elist = (TEntryList *) fInput->FindObject("elist")))
// Need to clone to avoid problems when destroying the selector
if (elist)
else
}
}
if (fillList) Info("SlaveBegin", "creating an entry-list");
if (option.Contains("useList")) useList = kTRUE;
}
{
// 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_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
h2->Fill(dm_d,rpd0_t/0.029979*1.8646/ptd0_d);
// Count the number of selected events
return kTRUE;
}
{
// nothing to be done
}
{
// 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
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}]");
//fit histogram hdmd with function f5 using the log-likelihood 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
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)
psdmd->SetOptStat(1110);
c1->Modified();
//save the entry list to a Root file if one was produced
if (fillList) {
if (!elist)
elist = dynamic_cast<TEntryList*>(fOutput->FindObject("elist"));
if (elist) {
Printf("Entry list 'elist' created:");
TFile efile("elist.root","recreate");
} else {
Error("Terminate", "entry list requested but not found in output");
}
}
// Notify the amount of processed events
if (!fInput) Info("Terminate", "processed %lld events", fProcessed);
}
const Bool_t kFALSE
Definition: RtypesCore.h:90
bool Bool_t
Definition: RtypesCore.h:61
double Double_t
Definition: RtypesCore.h:57
long long Long64_t
Definition: RtypesCore.h:71
const Bool_t kTRUE
Definition: RtypesCore.h:89
#define gDirectory
Definition: TDirectory.h:229
void Printf(const char *fmt,...)
R__EXTERN TStyle * gStyle
Definition: TStyle.h:410
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:294
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition: TAttMarker.h:40
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:1591
The Canvas class.
Definition: TCanvas.h:27
virtual void SetDirectory(TDirectory *dir)
Add reference to directory dir. dir can be 0.
virtual Bool_t Enter(Long64_t entry, TTree *tree=0)
Add entry #entry to the list.
Definition: TEntryList.cxx:560
virtual void Print(const Option_t *option="") const
Print this list.
Definition: TEntryList.cxx:999
1-Dim function class
Definition: TF1.h:210
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:638
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:614
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:571
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition: TH1.h:316
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition: TH1.cxx:3808
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3275
TList * GetListOfFunctions() const
Definition: TH1.h:239
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2998
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:251
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:294
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:897
TObject * FindObject(const char *name) const
Find object using its name.
Definition: THashList.cxx:262
A simple line.
Definition: TLine.h:23
virtual void Add(TObject *obj)
Definition: TList.h:87
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition: TList.cxx:577
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:164
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition: TNamed.cxx:74
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:796
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:891
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition: TObject.cxx:195
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition: TObject.cxx:865
The histogram statistics painter class.
Definition: TPaveStats.h:18
void SetOptStat(Int_t stat=1)
Set the stat option.
Definition: TPaveStats.cxx:303
TList * fInput
List of objects available during processing.
Definition: TSelector.h:43
TSelectorList * fOutput
! List of objects created during processing
Definition: TSelector.h:44
Long64_t fStatus
Selector status.
Definition: TSelector.h:39
virtual const char * GetOption() const
Definition: TSelector.h:59
Basic string class.
Definition: TString.h:131
const char * Data() const
Definition: TString.h:364
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:619
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:1590
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:1542
A TTree represents a columnar dataset.
Definition: TTree.h:78
virtual void SetEntryList(TEntryList *list, Option_t *opt="")
Set an EntryList.
Definition: TTree.cxx:8893
TBranch * b_md0_d
Definition: h1analysis.h:265
TBranch * b_dm_d
Definition: h1analysis.h:254
Float_t rpd0_t
Definition: h1analysis.h:131
TBranch * b_nlhpi
Definition: h1analysis.h:320
TBranch * b_ipis
Definition: h1analysis.h:261
void Terminate()
Bool_t Process(Long64_t entry)
TBranch * b_ipi
Definition: h1analysis.h:260
Float_t nlhpi[200]
Definition: h1analysis.h:166
Bool_t useList
Definition: h1analysis.h:30
Int_t nhitrp[200]
Definition: h1analysis.h:157
TBranch * b_ik
Definition: h1analysis.h:259
TBranch * b_rstart
Definition: h1analysis.h:315
TBranch * b_rend
Definition: h1analysis.h:316
TBranch * b_rpd0_t
Definition: h1analysis.h:285
Bool_t fillList
Definition: h1analysis.h:31
Float_t md0_d
Definition: h1analysis.h:111
void Init(TTree *tree)
Definition: h1analysis.h:390
TBranch * b_etads_d
Definition: h1analysis.h:253
TTree * fChain
Definition: h1analysis.h:35
TBranch * b_nhitrp
Definition: h1analysis.h:311
void Begin(TTree *tree)
Long64_t fProcessed
Definition: h1analysis.h:33
Int_t ik
Definition: h1analysis.h:105
TH1F * hdmd
Definition: h1analysis.h:27
Float_t dm_d
Definition: h1analysis.h:100
void Reset()
Definition: h1analysis.h:376
Float_t rend[200]
Definition: h1analysis.h:162
TBranch * b_nlhk
Definition: h1analysis.h:319
TBranch * b_ptds_d
Definition: h1analysis.h:252
TBranch * b_ntracks
Definition: h1analysis.h:303
Float_t rstart[200]
Definition: h1analysis.h:161
TBranch * b_njets
Definition: h1analysis.h:328
void SlaveTerminate()
Int_t ipis
Definition: h1analysis.h:107
Float_t ptd0_d
Definition: h1analysis.h:109
Int_t njets
Definition: h1analysis.h:174
Float_t nlhk[200]
Definition: h1analysis.h:165
void SlaveBegin(TTree *tree)
TEntryList * elist
Definition: h1analysis.h:32
Float_t etads_d
Definition: h1analysis.h:99
TH2F * h2
Definition: h1analysis.h:28
TBranch * b_ptd0_d
Definition: h1analysis.h:263
Float_t ptds_d
Definition: h1analysis.h:98
Int_t ipi
Definition: h1analysis.h:106
TLine * line
Double_t fdm5(Double_t *xx, Double_t *par)
const Double_t sigma
const Double_t dxbin
Double_t fdm2(Double_t *xx, Double_t *par)
return c1
Definition: legend1.C:41
Double_t x[n]
Definition: legend1.C:17
return c2
Definition: legend2.C:14
Double_t Exp(Double_t x)
Definition: TMath.h:717
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:725
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
Definition: tree.py:1
Author
Rene Brun

Definition in file h1analysis.C.