23 using namespace RooFit;
25 void rf707_kernelestimation()
31 RooRealVar x(
"x",
"x", 0, 20);
32 RooPolynomial p(
"p",
"p", x, RooArgList(RooConst(0.01), RooConst(-0.01), RooConst(0.0004)));
35 RooDataSet *data1 = p.generate(x, 200);
43 RooKeysPdf kest1(
"kest1",
"kest1", x, *data1, RooKeysPdf::MirrorBoth);
47 RooKeysPdf kest2(
"kest2",
"kest2", x, *data1, RooKeysPdf::NoMirror);
51 RooKeysPdf kest3(
"kest1",
"kest1", x, *data1, RooKeysPdf::MirrorBoth, 2);
54 RooPlot *frame = x.frame(Title(
"Adaptive kernel estimation pdf with and w/o mirroring"), Bins(20));
57 kest2.plotOn(frame, LineStyle(kDashed), LineColor(kRed));
60 RooPlot *frame2 = x.frame(Title(
"Adaptive kernel estimation pdf with regular, increased bandwidth"));
62 kest3.plotOn(frame2, LineColor(kMagenta));
68 RooRealVar y(
"y",
"y", 0, 20);
69 RooPolynomial py(
"py",
"py", y, RooArgList(RooConst(0.01), RooConst(0.01), RooConst(-0.0004)));
70 RooProdPdf pxy(
"pxy",
"pxy", RooArgSet(p, py));
71 RooDataSet *data2 = pxy.generate(RooArgSet(x, y), 1000);
77 RooNDKeysPdf kest4(
"kest4",
"kest4", RooArgSet(x, y), *data2,
"am");
80 RooNDKeysPdf kest5(
"kest5",
"kest5", RooArgSet(x, y), *data2,
"am", 2);
83 TH1 *hh_data = data2->createHistogram(
"hh_data", x, Binning(10), YVar(y, Binning(10)));
86 TH1 *hh_pdf = kest4.createHistogram(
"hh_pdf", x, Binning(25), YVar(y, Binning(25)));
87 TH1 *hh_pdf2 = kest5.createHistogram(
"hh_pdf2", x, Binning(25), YVar(y, Binning(25)));
88 hh_pdf->SetLineColor(kBlue);
89 hh_pdf2->SetLineColor(kMagenta);
91 TCanvas *c =
new TCanvas(
"rf707_kernelestimation",
"rf707_kernelestimation", 800, 800);
94 gPad->SetLeftMargin(0.15);
95 frame->GetYaxis()->SetTitleOffset(1.4);
98 gPad->SetLeftMargin(0.15);
99 frame2->GetYaxis()->SetTitleOffset(1.8);
102 gPad->SetLeftMargin(0.15);
103 hh_data->GetZaxis()->SetTitleOffset(1.4);
104 hh_data->Draw(
"lego");
106 gPad->SetLeftMargin(0.20);
107 hh_pdf->GetZaxis()->SetTitleOffset(2.4);
108 hh_pdf->Draw(
"surf");
109 hh_pdf2->Draw(
"surfsame");