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