Logo ROOT   6.18/05
Reference Guide
h1analysisTreeReader.C File Reference

Detailed Description

View in nbviewer Open in SWAN H1 analysis example expressed in terms of TTreeReader

#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 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;
}
//_____________________________________________________________________
// 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) {
// 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
h2->Fill(*fDm_d,*fRpd0_t/0.029979*1.8646/ *fPtd0_d);
return kTRUE;
}
// 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")) {
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
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);
// 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;
}
// 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 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
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);
}
}
// called when loading a new file
// get branch pointers
Info("Notify","processing file: %s",myTreeReader.GetTree()->GetCurrentFile()->GetName());
if (fillList) {
} else if (useList) {
}
}
return kTRUE;
}
#define f(i)
Definition: RSha256.hxx:104
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
long long Long64_t
Definition: RtypesCore.h:69
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define gDirectory
Definition: TDirectory.h:218
#define gROOT
Definition: TROOT.h:414
void Printf(const char *fmt,...)
R__EXTERN TStyle * gStyle
Definition: TStyle.h:406
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title Offset is a correction factor with respect to the "s...
Definition: TAttAxis.cxx:294
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition: TAttMarker.h:40
The Canvas class.
Definition: TCanvas.h:31
A List of entry numbers in a TTree or TChain.
Definition: TEntryList.h:26
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,...
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:558
virtual void Print(const Option_t *option="") const
Print this list.
Definition: TEntryList.cxx:997
1-Dim function class
Definition: TF1.h:211
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:638
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition: TFile.h:48
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:3791
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3258
TList * GetListOfFunctions() const
Definition: TH1.h:239
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2981
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:248
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:292
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:851
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:575
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 const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
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
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
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:854
The histogram statistics painter class.
Definition: TPaveStats.h:18
void SetOptStat(Int_t stat=1)
Set the stat option.
Definition: TPaveStats.cxx:302
TList * fInput
List of objects available during processing.
Definition: TSelector.h:43
TSelectorList * fOutput
! List of objects created during processing
Definition: TSelector.h:44
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:1444
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:1396
T & At(std::size_t idx)
TTree * GetTree() const
Definition: TTreeReader.h:174
EEntryStatus SetLocalEntry(Long64_t entry)
Set the next local tree entry.
Definition: TTreeReader.h:202
A TTree represents a columnar dataset.
Definition: TTree.h:71
TFile * GetCurrentFile() const
Return pointer to the current file.
Definition: TTree.cxx:5263
virtual void SetEntryList(TEntryList *list, Option_t *opt="")
Set an EntryList.
Definition: TTree.cxx:8702
TTreeReaderValue< Int_t > fIpi
void Init(TTree *myTree)
TTreeReaderArray< Int_t > fNhitrp
Bool_t Process(Long64_t entry)
TTreeReaderValue< Float_t > fRpd0_t
TTreeReaderArray< Float_t > fNlhpi
TTreeReaderValue< Float_t > fMd0_d
void SlaveBegin(TTree *)
TTreeReaderValue< Float_t > fPtd0_d
TTreeReaderValue< Int_t > fIk
Bool_t Notify()
This method must be overridden to handle object notification.
TTreeReaderValue< Int_t > fIpis
TTreeReaderValue< Int_t > fNjets
TTreeReaderArray< Float_t > fRend
TTreeReaderValue< Float_t > fEtads_d
void Begin(TTree *)
TTreeReaderArray< Float_t > fRstart
TTreeReaderValue< Float_t > fDm_d
TTreeReaderArray< Float_t > fNlhk
TTreeReaderValue< Float_t > fPtds_d
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:715
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:723
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
Author
Anders Eie, 2013

Definition in file h1analysisTreeReader.C.