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: https://root.cern/files/h1/ 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
Case A: Create a TChain with the 4 H1 data files
The chain can be created by executed the short macro h1chain.C below:
{
chain.Add("$H1/dstarmb.root");
chain.Add("$H1/dstarp1a.root");
chain.Add("$H1/dstarp1b.root");
chain.Add("$H1/dstarp2.root");
}
A chain is a collection of files containing TTree objects.
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
Root > TProof *p = TProof::Open("")
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:
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
A file, usually with extension .root, that stores data and code in the form of serialized objects in ...
set it on the chain:
Root > chain.SetEntryList(elist)
call Process as in step B. Of course this works also for local processing.
{
if (
x <= 0.13957)
return 0;
+ par[2] / 2.5066/par[4]*
TMath::Exp(-xp3/2/par[4]/par[4]));
return res;
}
{
if (
x <= 0.13957)
return 0;
return res;
}
{
Info(
"Begin",
"starting h1analysis with process option: %s", option.
Data());
delete gDirectory->GetList()->FindObject(
"elist");
elist =
new TEntryList(
"elist",
"H1 selection from Cut");
fInput->Add(
new TNamed(
"fillList",
""));
}
Info(
"Begin",
"creating an entry-list");
}
} else {
elist = (TEntryList*)
f.Get(
"elist");
}
}
}
{
"starting h1analysis with process option: %s (tree: %p)", option.
Data(), tree);
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);
if ((
elist = (TEntryList *)
fInput->FindObject(
"elist")))
else
}
}
}
{
}
}
{
}
{
hdmd =
dynamic_cast<TH1F*
>(
fOutput->FindObject(
"hdmd"));
h2 =
dynamic_cast<TH2F*
>(
fOutput->FindObject(
"h2"));
if (
hdmd ==
nullptr ||
h2 ==
nullptr) {
return;
}
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);
if (
gROOT->GetListOfFunctions()->FindObject(
"f5"))
delete gROOT->GetFunction(
"f5");
TF1 *f5 =
new TF1(
"f5",
fdm5,0.139,0.17,5);
TCanvas *
c2 =
new TCanvas(
"c2",
"tauD0",100,100,800,600);
c2->SetBottomMargin(0.15);
if (
gROOT->GetListOfFunctions()->FindObject(
"f2"))
delete gROOT->GetFunction(
"f2");
TF1 *f2 =
new TF1(
"f2",
fdm2,0.139,0.17,2);
h2->FitSlicesX(f2,0,-1,1,
"qln");
TLine *
line =
new TLine(0,0,0,
c2->GetUymax());
TPaveStats *psdmd = (TPaveStats *)
hdmd->GetListOfFunctions()->FindObject(
"stats");
elist =
dynamic_cast<TEntryList*
>(
fOutput->FindObject(
"elist"));
Printf(
"Entry list 'elist' created:");
TFile efile("elist.root","recreate");
} else {
Error(
"Terminate",
"entry list requested but not found in output");
}
}
}
bool Bool_t
Boolean (0=false, 1=true) (bool).
double Double_t
Double 8 bytes.
long long Long64_t
Portable signed long integer 8 bytes.
void Printf(const char *fmt,...)
Formats a string in a circular formatting buffer and prints the string.
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
virtual void SetParameters(const Double_t *params)
void Draw(Option_t *option="") override
Draw this histogram with options.
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
void SetOptStat(Int_t stat=1)
TList * fInput
List of objects available during processing.
TSelectorList * fOutput
! List of objects created during processing
Long64_t fStatus
Selector status.
const char * GetOption() const override
const char * Data() const
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
A TTree represents a columnar dataset.
void SlaveTerminate() override
void Init(TTree *tree) override
TTree * fChain
! //pointer to the analyzed TTree or TChain
Bool_t Process(Long64_t entry) override
The Process() function is called for each entry in the tree to be processed.
void Begin(TTree *tree) override
void SlaveBegin(TTree *tree) override
void Terminate() override
Double_t fdm5(Double_t *xx, Double_t *par)
Double_t fdm2(Double_t *xx, Double_t *par)
Double_t Exp(Double_t x)
Returns the base-e exponential function of x, which is e raised to the power x.
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
- Author
- Rene Brun
Definition in file h1analysis.C.