Logo ROOT  
Reference Guide
df101_h1Analysis.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_dataframe
3/// \notebook
4/// This tutorial illustrates how to express the H1 analysis with a RDataFrame.
5///
6/// \macro_code
7/// \macro_image
8///
9/// \date December 2016
10/// \author Axel Naumann, Danilo Piparo (CERN)
11
12auto Select = [](ROOT::RDataFrame &dataFrame) {
13 using Farray_t = ROOT::VecOps::RVec<float>;
14 using Iarray_t = ROOT::VecOps::RVec<int>;
15
16 auto ret = dataFrame.Filter("TMath::Abs(md0_d - 1.8646) < 0.04")
17 .Filter("ptds_d > 2.5")
18 .Filter("TMath::Abs(etads_d) < 1.5")
19 .Filter([](int ik, int ipi, Iarray_t& nhitrp) { return nhitrp[ik - 1] * nhitrp[ipi - 1] > 1; },
20 {"ik", "ipi", "nhitrp"})
21 .Filter([](int ik, Farray_t& rstart, Farray_t& rend) { return rend[ik - 1] - rstart[ik - 1] > 22; },
22 {"ik", "rstart", "rend"})
23 .Filter([](int ipi, Farray_t& rstart, Farray_t& rend) { return rend[ipi - 1] - rstart[ipi - 1] > 22; },
24 {"ipi", "rstart", "rend"})
25 .Filter([](int ik, Farray_t& nlhk) { return nlhk[ik - 1] > 0.1; }, {"ik", "nlhk"})
26 .Filter([](int ipi, Farray_t& nlhpi) { return nlhpi[ipi - 1] > 0.1; }, {"ipi", "nlhpi"})
27 .Filter([](int ipis, Farray_t& nlhpi) { return nlhpi[ipis - 1] > 0.1; }, {"ipis", "nlhpi"})
28 .Filter("njets >= 1");
29
30 return ret;
31};
32
33const Double_t dxbin = (0.17 - 0.13) / 40; // Bin-width
34
36{
37 Double_t x = xx[0];
38 if (x <= 0.13957)
39 return 0;
40 Double_t xp3 = (x - par[3]) * (x - par[3]);
41 Double_t res =
42 dxbin * (par[0] * pow(x - 0.13957, par[1]) + par[2] / 2.5066 / par[4] * exp(-xp3 / 2 / par[4] / par[4]));
43 return res;
44}
45
47{
48 static const Double_t sigma = 0.0012;
49 Double_t x = xx[0];
50 if (x <= 0.13957)
51 return 0;
52 Double_t xp3 = (x - 0.1454) * (x - 0.1454);
53 Double_t res = dxbin * (par[0] * pow(x - 0.13957, 0.25) + par[1] / 2.5066 / sigma * exp(-xp3 / 2 / sigma / sigma));
54 return res;
55}
56
57void FitAndPlotHdmd(TH1 &hdmd)
58{
59 // create the canvas for the h1analysis fit
61 auto c1 = new TCanvas("c1", "h1analysis analysis", 10, 10, 800, 600);
62 hdmd.GetXaxis()->SetTitle("m_{K#pi#pi} - m_{K#pi}[GeV/c^{2}]");
63 hdmd.GetXaxis()->SetTitleOffset(1.4);
64
65 // fit histogram hdmd with function f5 using the loglikelihood option
66 auto f5 = new TF1("f5", fdm5, 0.139, 0.17, 5);
67 f5->SetParameters(1000000, .25, 2000, .1454, .001);
68 hdmd.Fit("f5", "lr");
69
70 hdmd.DrawClone();
71}
72
73void FitAndPlotH2(TH2 &h2)
74{
75 // create the canvas for tau d0
76 auto c2 = new TCanvas("c2", "tauD0", 100, 100, 800, 600);
77
78 c2->SetGrid();
79 c2->SetBottomMargin(0.15);
80
81 // Project slices of 2-d histogram h2 along X , then fit each slice
82 // with function f2 and make a histogram for each fit parameter
83 // Note that the generated histograms are added to the list of objects
84 // in the current directory.
85 auto f2 = new TF1("f2", fdm2, 0.139, 0.17, 2);
86 f2->SetParameters(10000, 10);
87 h2.FitSlicesX(f2, 0, -1, 1, "qln");
88
89 // See TH2::FitSlicesX documentation
90 auto h2_1 = (TH1D *)gDirectory->Get("h2_1");
91 h2_1->GetXaxis()->SetTitle("#tau [ps]");
92 h2_1->SetMarkerStyle(21);
93 h2_1->DrawClone();
94 c2->Update();
95
96 auto line = new TLine(0, 0, 0, c2->GetUymax());
97 line->Draw();
98}
99
100void df101_h1Analysis()
101{
102 TChain chain("h42");
103 chain.Add("root://eospublic.cern.ch//eos/root-eos/h1/dstarmb.root");
104 chain.Add("root://eospublic.cern.ch//eos/root-eos/h1/dstarp1a.root");
105 chain.Add("root://eospublic.cern.ch//eos/root-eos/h1/dstarp1b.root");
106 chain.Add("root://eospublic.cern.ch//eos/root-eos/h1/dstarp2.root");
107
109
110 ROOT::RDataFrame dataFrame(chain);
111 auto selected = Select(dataFrame);
112 auto hdmdARP = selected.Histo1D({"hdmd", "Dm_d", 40, 0.13, 0.17}, "dm_d");
113 auto selectedAddedBranch = selected.Define("h2_y", "rpd0_t / 0.029979f * 1.8646f / ptd0_d");
114 auto h2ARP = selectedAddedBranch.Histo2D({"h2", "ptD0 vs Dm_d", 30, 0.135, 0.165, 30, -3, 6}, "dm_d", "h2_y");
115
116 FitAndPlotHdmd(*hdmdARP);
117 FitAndPlotH2(*h2ARP);
118}
double Double_t
Definition: RtypesCore.h:55
#define gDirectory
Definition: TDirectory.h:223
double pow(double, double)
double exp(double)
R__EXTERN TStyle * gStyle
Definition: TStyle.h:407
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTrees,...
Definition: RDataFrame.hxx:42
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:294
The Canvas class.
Definition: TCanvas.h:31
A chain is a collection of files containing TTree objects.
Definition: TChain.h:34
1-Dim function class
Definition: TF1.h:211
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:614
The TH1 histogram class.
Definition: TH1.h:56
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
Service class for 2-Dim histogram classes.
Definition: TH2.h:30
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
A simple line.
Definition: TLine.h:23
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:164
virtual TObject * DrawClone(Option_t *option="") const
Draw a clone of this object in the current selected pad for instance with: gROOT->SetSelectedPad(gPad...
Definition: TObject.cxx:219
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition: TObject.cxx:195
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
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
RVec< T > Filter(const RVec< T > &v, F &&f)
Create a new collection with the elements passing the filter expressed by the predicate.
Definition: RVec.hxx:938
void EnableImplicitMT(UInt_t numthreads=0)
Enable ROOT's implicit multi-threading for all objects and methods that provide an internal paralleli...
Definition: TROOT.cxx:580