12 auto Select = [](ROOT::RDataFrame &dataFrame) {
13 using Farray_t = ROOT::VecOps::RVec<float>;
14 using Iarray_t = ROOT::VecOps::RVec<int>;
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");
33 const Double_t dxbin = (0.17 - 0.13) / 40;
35 Double_t fdm5(Double_t *xx, Double_t *par)
40 Double_t xp3 = (x - par[3]) * (x - par[3]);
42 dxbin * (par[0] * pow(x - 0.13957, par[1]) + par[2] / 2.5066 / par[4] * exp(-xp3 / 2 / par[4] / par[4]));
46 Double_t fdm2(Double_t *xx, Double_t *par)
48 static const Double_t sigma = 0.0012;
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));
57 void FitAndPlotHdmd(TH1 &hdmd)
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);
66 auto f5 =
new TF1(
"f5", fdm5, 0.139, 0.17, 5);
67 f5->SetParameters(1000000, .25, 2000, .1454, .001);
73 void FitAndPlotH2(TH2 &h2)
76 auto c2 =
new TCanvas(
"c2",
"tauD0", 100, 100, 800, 600);
79 c2->SetBottomMargin(0.15);
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");
90 auto h2_1 = (TH1D *)gDirectory->Get(
"h2_1");
91 h2_1->GetXaxis()->SetTitle(
"#tau [ps]");
92 h2_1->SetMarkerStyle(21);
96 auto line =
new TLine(0, 0, 0, c2->GetUymax());
100 void df101_h1Analysis()
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");
108 ROOT::EnableImplicitMT(4);
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");
116 FitAndPlotHdmd(*hdmdARP);
117 FitAndPlotH2(*h2ARP);