20 using namespace RooFit;
22 void rf310_sliceplot()
29 RooRealVar dt(
"dt",
"dt", -20, 20);
32 RooCategory mixState(
"mixState",
"B0/B0bar mixing state");
33 RooCategory tagFlav(
"tagFlav",
"Flavour of the tagged B0");
36 mixState.defineType(
"mixed", -1);
37 mixState.defineType(
"unmixed", 1);
38 tagFlav.defineType(
"B0", 1);
39 tagFlav.defineType(
"B0bar", -1);
42 RooRealVar dm(
"dm",
"delta m(B)", 0.472, 0., 1.0);
43 RooRealVar tau(
"tau",
"B0 decay time", 1.547, 1.0, 2.0);
44 RooRealVar w(
"w",
"Flavor Mistag rate", 0.03, 0.0, 1.0);
45 RooRealVar dw(
"dw",
"Flavor Mistag rate difference between B0 and B0bar", 0.01);
48 RooRealVar bias1(
"bias1",
"bias1", 0);
49 RooRealVar sigma1(
"sigma1",
"sigma1", 0.01);
50 RooGaussModel gm1(
"gm1",
"gauss model 1", dt, bias1, sigma1);
53 RooBMixDecay bmix_gm1(
"bmix",
"decay", dt, mixState, tagFlav, tau, dm, w, dw, gm1, RooBMixDecay::DoubleSided);
56 RooDataSet *data = bmix_gm1.generate(RooArgSet(dt, tagFlav, mixState), 20000);
62 RooPlot *frame = dt.frame(Title(
"Inclusive decay distribution"));
64 bmix_gm1.plotOn(frame);
70 RooPlot *frame2 = dt.frame(Title(
"Decay distribution of mixed events"));
71 data->plotOn(frame2, Cut(
"mixState==mixState::mixed"));
74 bmix_gm1.plotOn(frame2, Slice(mixState,
"mixed"));
77 RooPlot *frame3 = dt.frame(Title(
"Decay distribution of unmixed events"));
78 data->plotOn(frame3, Cut(
"mixState==mixState::unmixed"));
81 bmix_gm1.plotOn(frame3, Slice(mixState,
"unmixed"));
83 TCanvas *c =
new TCanvas(
"rf310_sliceplot",
"rf310_sliceplot", 1200, 400);
86 gPad->SetLeftMargin(0.15);
87 frame->GetYaxis()->SetTitleOffset(1.4);
91 gPad->SetLeftMargin(0.15);
92 frame2->GetYaxis()->SetTitleOffset(1.4);
96 gPad->SetLeftMargin(0.15);
97 frame3->GetYaxis()->SetTitleOffset(1.4);