ROOT
v6-20
Reference Guide
h1analysisTreeReader.C File Reference
Tutorials
»
Tree tutorials
Detailed Description
H1 analysis example expressed in terms of
TTreeReader
#include "
h1analysisTreeReader.h
"
#include "
TStyle.h
"
#include "
TCanvas.h
"
#include "
TPaveStats.h
"
#include "
TLine.h
"
#include "
TMath.h
"
#include "
TFile.h
"
#include "
TROOT.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;
}
//_____________________________________________________________________
Bool_t
h1analysisTreeReader::Process
(
Long64_t
entry){
// entry is the entry number in the current Tree
// Selection function to select D* and D0.
myTreeReader
.
SetLocalEntry
(entry);
fProcessed
++;
//in case one entry list is given in input, the selection has already been done.
if
(!
useList
) {
// Return as soon as a bad entry is detected
if
(
TMath::Abs
(*
fMd0_d
-1.8646) >= 0.04)
return
kFALSE
;
if
(*
fPtds_d
<= 2.5)
return
kFALSE
;
if
(
TMath::Abs
(*
fEtads_d
) >= 1.5)
return
kFALSE
;
(*fIk)--;
//original fIk used f77 convention starting at 1
(*fIpi)--;
if
(
fNhitrp
.
At
(*
fIk
)*
fNhitrp
.
At
(*
fIpi
) <= 1)
return
kFALSE
;
if
(
fRend
.
At
(*
fIk
) -
fRstart
.
At
(*
fIk
) <= 22)
return
kFALSE
;
if
(
fRend
.
At
(*
fIpi
)-
fRstart
.
At
(*
fIpi
) <= 22)
return
kFALSE
;
if
(
fNlhk
.
At
(*
fIk
) <= 0.1)
return
kFALSE
;
if
(
fNlhpi
.
At
(*
fIpi
) <= 0.1)
return
kFALSE
;
(*fIpis)--;
if
(
fNlhpi
.
At
(*
fIpis
) <= 0.1)
return
kFALSE
;
if
(*
fNjets
< 1)
return
kFALSE
;
}
// if option fillList, fill the entry list
if
(
fillList
)
elist
->
Enter
(entry);
//fill some histograms
hdmd
->
Fill
(*
fDm_d
);
h2
->
Fill
(*
fDm_d
,*
fRpd0_t
/0.029979*1.8646/ *
fPtd0_d
);
return
kTRUE
;
}
void
h1analysisTreeReader::Begin
(
TTree
*
/*myTree*/
) {
// function called before starting the event loop
// -it performs some cleanup
// -it creates histograms
// -it sets some initialisation for the entry list
Reset
();
//print the option specified in the Process function.
TString
option =
GetOption
();
Info
(
"Begin"
,
"starting h1analysis with process option: %s"
, option.
Data
());
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"
,
""
));
// We send a clone to avoid double deletes when importing the result
fInput
->
Add
(
elist
);
// 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"
)) {
useList
=
kTRUE
;
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
}
}
}
void
h1analysisTreeReader::SlaveBegin
(
TTree
*myTree){
// function called before starting the event loop
// -it performs some cleanup
// -it creates histograms
// -it sets some initialisation for the entry list
Init
(myTree);
//print the option specified in the Process function.
TString
option =
GetOption
();
Info
(
"SlaveBegin"
,
"starting h1analysis with process option: %s (tree: %p)"
, option.
Data
(), myTree);
//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
;
// Get the list
if
(
fInput
) {
if
((
elist
= (
TEntryList
*)
fInput
->
FindObject
(
"elist"
)))
// Need to clone to avoid problems when destroying the selector
elist
= (
TEntryList
*)
elist
->
Clone
();
if
(
elist
)
fOutput
->
Add
(
elist
);
else
fillList
=
kFALSE
;
}
}
if
(
fillList
)
Info
(
"SlaveBegin"
,
"creating an entry-list"
);
if
(option.
Contains
(
"useList"
))
useList
=
kTRUE
;
}
void
h1analysisTreeReader::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 loglfIkelihood 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
) {
if
(!
elist
)
elist
=
dynamic_cast<
TEntryList
*
>
(
fOutput
->
FindObject
(
"elist"
));
if
(
elist
) {
Printf
(
"Entry list 'elist' created:"
);
elist
->
Print
();
TFile
efile(
"elist.root"
,
"recreate"
);
elist
->
Write
();
}
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
);
}
void
h1analysisTreeReader::SlaveTerminate
(){
}
Bool_t
h1analysisTreeReader::Notify
() {
// called when loading a new file
// get branch pointers
Info
(
"Notify"
,
"processing file: %s"
,
myTreeReader
.
GetTree
()->
GetCurrentFile
()->
GetName
());
if
(
elist
&&
myTreeReader
.
GetTree
()) {
if
(
fillList
) {
elist
->
SetTree
(
myTreeReader
.
GetTree
());
}
else
if
(
useList
) {
myTreeReader
.
GetTree
()->
SetEntryList
(
elist
);
}
}
return
kTRUE
;
}
f
#define f(i)
Definition:
RSha256.hxx:104
kFALSE
const Bool_t kFALSE
Definition:
RtypesCore.h:88
Bool_t
bool Bool_t
Definition:
RtypesCore.h:59
Double_t
double Double_t
Definition:
RtypesCore.h:55
Long64_t
long long Long64_t
Definition:
RtypesCore.h:69
kTRUE
const Bool_t kTRUE
Definition:
RtypesCore.h:87
TCanvas.h
gDirectory
#define gDirectory
Definition:
TDirectory.h:223
TFile.h
TLine.h
TMath.h
TPaveStats.h
TROOT.h
gROOT
#define gROOT
Definition:
TROOT.h:415
Printf
void Printf(const char *fmt,...)
TStyle.h
gStyle
R__EXTERN TStyle * gStyle
Definition:
TStyle.h:407
TAttAxis::SetTitleOffset
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition:
TAttAxis.cxx:294
TAttMarker::SetMarkerStyle
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition:
TAttMarker.h:40
TCanvas
The Canvas class.
Definition:
TCanvas.h:31
TEntryList
A List of entry numbers in a TTree or TChain.
Definition:
TEntryList.h:26
TEntryList::SetTree
virtual void SetTree(const TTree *tree)
If a list for a tree with such name and filename exists, sets it as the current sublist If not,...
Definition:
TEntryList.cxx:1192
TEntryList::SetDirectory
virtual void SetDirectory(TDirectory *dir)
Add reference to directory dir. dir can be 0.
Definition:
TEntryList.cxx:1066
TEntryList::Enter
virtual Bool_t Enter(Long64_t entry, TTree *tree=0)
Add entry #entry to the list.
Definition:
TEntryList.cxx:558
TEntryList::Print
virtual void Print(const Option_t *option="") const
Print this list.
Definition:
TEntryList.cxx:997
TF1
1-Dim function class
Definition:
TF1.h:211
TF1::SetParameters
virtual void SetParameters(const Double_t *params)
Definition:
TF1.h:638
TFile
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition:
TFile.h:48
TH1D
1-D histogram with a double per channel (see TH1 documentation)}
Definition:
TH1.h:614
TH1F
1-D histogram with a float per channel (see TH1 documentation)}
Definition:
TH1.h:571
TH1::GetXaxis
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition:
TH1.h:316
TH1::Fit
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
TH1::Fill
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition:
TH1.cxx:3275
TH1::GetListOfFunctions
TList * GetListOfFunctions() const
Definition:
TH1.h:239
TH1::Draw
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition:
TH1.cxx:2998
TH2F
2-D histogram with a float per channel (see TH1 documentation)}
Definition:
TH2.h:251
TH2::Fill
Int_t Fill(Double_t)
Invalid Fill method.
Definition:
TH2.cxx:292
TH2::FitSlicesX
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:858
THashList::FindObject
TObject * FindObject(const char *name) const
Find object using its name.
Definition:
THashList.cxx:262
TLine
A simple line.
Definition:
TLine.h:23
TList::Add
virtual void Add(TObject *obj)
Definition:
TList.h:87
TList::FindObject
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition:
TList.cxx:575
TNamed
The TNamed class is the base class for all named ROOT classes.
Definition:
TNamed.h:29
TNamed::SetTitle
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition:
TNamed.cxx:164
TNamed::Clone
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition:
TNamed.cxx:74
TNamed::GetName
virtual const char * GetName() const
Returns name of object.
Definition:
TNamed.h:47
TObject::Write
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:785
TObject::Error
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition:
TObject.cxx:880
TObject::Draw
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition:
TObject.cxx:195
TObject::Info
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition:
TObject.cxx:854
TPaveStats
The histogram statistics painter class.
Definition:
TPaveStats.h:18
TPaveStats::SetOptStat
void SetOptStat(Int_t stat=1)
Set the stat option.
Definition:
TPaveStats.cxx:302
TSelector::fInput
TList * fInput
List of objects available during processing.
Definition:
TSelector.h:43
TSelector::fOutput
TSelectorList * fOutput
! List of objects created during processing
Definition:
TSelector.h:44
TSelector::GetOption
virtual const char * GetOption() const
Definition:
TSelector.h:59
TString
Basic string class.
Definition:
TString.h:131
TString::Data
const char * Data() const
Definition:
TString.h:364
TString::Contains
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition:
TString.h:619
TStyle::SetOptStat
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:1450
TStyle::SetOptFit
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:1402
TTreeReaderArray::At
T & At(std::size_t idx)
Definition:
TTreeReaderArray.h:205
TTreeReader::GetTree
TTree * GetTree() const
Definition:
TTreeReader.h:173
TTreeReader::SetLocalEntry
EEntryStatus SetLocalEntry(Long64_t entry)
Set the next local tree entry.
Definition:
TTreeReader.h:201
TTree
A TTree represents a columnar dataset.
Definition:
TTree.h:72
TTree::GetCurrentFile
TFile * GetCurrentFile() const
Return pointer to the current file.
Definition:
TTree.cxx:5338
TTree::SetEntryList
virtual void SetEntryList(TEntryList *list, Option_t *opt="")
Set an EntryList.
Definition:
TTree.cxx:8783
h1analysisTreeReader::myTreeReader
TTreeReader myTreeReader
Definition:
h1analysisTreeReader.h:14
h1analysisTreeReader::fIpi
TTreeReaderValue< Int_t > fIpi
Definition:
h1analysisTreeReader.h:20
h1analysisTreeReader::Init
void Init(TTree *myTree)
Definition:
h1analysisTreeReader.h:67
h1analysisTreeReader::fNhitrp
TTreeReaderArray< Int_t > fNhitrp
Definition:
h1analysisTreeReader.h:25
h1analysisTreeReader::Terminate
void Terminate()
h1analysisTreeReader::Process
Bool_t Process(Long64_t entry)
h1analysisTreeReader::fRpd0_t
TTreeReaderValue< Float_t > fRpd0_t
Definition:
h1analysisTreeReader.h:24
h1analysisTreeReader::fNlhpi
TTreeReaderArray< Float_t > fNlhpi
Definition:
h1analysisTreeReader.h:29
h1analysisTreeReader::fMd0_d
TTreeReaderValue< Float_t > fMd0_d
Definition:
h1analysisTreeReader.h:23
h1analysisTreeReader::SlaveBegin
void SlaveBegin(TTree *)
h1analysisTreeReader::SlaveTerminate
void SlaveTerminate()
h1analysisTreeReader::useList
Bool_t useList
Definition:
h1analysisTreeReader.h:35
h1analysisTreeReader::elist
TEntryList * elist
Definition:
h1analysisTreeReader.h:37
h1analysisTreeReader::fPtd0_d
TTreeReaderValue< Float_t > fPtd0_d
Definition:
h1analysisTreeReader.h:22
h1analysisTreeReader::fProcessed
Long64_t fProcessed
Definition:
h1analysisTreeReader.h:38
h1analysisTreeReader::fillList
Bool_t fillList
Definition:
h1analysisTreeReader.h:36
h1analysisTreeReader::fIk
TTreeReaderValue< Int_t > fIk
Definition:
h1analysisTreeReader.h:19
h1analysisTreeReader::Reset
void Reset()
Definition:
h1analysisTreeReader.h:81
h1analysisTreeReader::Notify
Bool_t Notify()
This method must be overridden to handle object notification.
h1analysisTreeReader::fIpis
TTreeReaderValue< Int_t > fIpis
Definition:
h1analysisTreeReader.h:21
h1analysisTreeReader::fNjets
TTreeReaderValue< Int_t > fNjets
Definition:
h1analysisTreeReader.h:30
h1analysisTreeReader::h2
TH2F * h2
Definition:
h1analysisTreeReader.h:33
h1analysisTreeReader::fRend
TTreeReaderArray< Float_t > fRend
Definition:
h1analysisTreeReader.h:27
h1analysisTreeReader::fEtads_d
TTreeReaderValue< Float_t > fEtads_d
Definition:
h1analysisTreeReader.h:17
h1analysisTreeReader::hdmd
TH1F * hdmd
Definition:
h1analysisTreeReader.h:32
h1analysisTreeReader::Begin
void Begin(TTree *)
h1analysisTreeReader::fRstart
TTreeReaderArray< Float_t > fRstart
Definition:
h1analysisTreeReader.h:26
h1analysisTreeReader::fDm_d
TTreeReaderValue< Float_t > fDm_d
Definition:
h1analysisTreeReader.h:18
h1analysisTreeReader::fNlhk
TTreeReaderArray< Float_t > fNlhk
Definition:
h1analysisTreeReader.h:28
h1analysisTreeReader::fPtds_d
TTreeReaderValue< Float_t > fPtds_d
Definition:
h1analysisTreeReader.h:16
line
TLine * line
Definition:
entrylistblock_figure1.C:235
fdm5
Double_t fdm5(Double_t *xx, Double_t *par)
Definition:
h1analysisProxy.h:14
sigma
const Double_t sigma
Definition:
h1analysisProxy.h:11
dxbin
const Double_t dxbin
Definition:
h1analysisProxy.h:10
fdm2
Double_t fdm2(Double_t *xx, Double_t *par)
Definition:
h1analysisProxy.h:25
h1analysisTreeReader.h
c1
return c1
Definition:
legend1.C:41
x
Double_t x[n]
Definition:
legend1.C:17
c2
return c2
Definition:
legend2.C:14
TMath::Exp
Double_t Exp(Double_t x)
Definition:
TMath.h:717
TMath::Power
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition:
TMath.h:725
TMath::Abs
Short_t Abs(Short_t d)
Definition:
TMathBase.h:120
Author
Anders Eie, 2013
Definition in file
h1analysisTreeReader.C
.
tutorials
tree
h1analysisTreeReader.C
ROOT v6-20 - Reference Guide Generated on Fri Apr 1 2022 00:24:48 (GVA Time) using Doxygen 1.9.4