Multidimensional models: projecting pdf and data slices in discrete observables
{
RooCategory mixState(
"mixState",
"B0/B0bar mixing state");
RooCategory tagFlav(
"tagFlav",
"Flavour of the tagged B0");
mixState.defineType("mixed", -1);
mixState.defineType("unmixed", 1);
tagFlav.defineType("B0", 1);
tagFlav.defineType("B0bar", -1);
RooRealVar dm(
"dm",
"delta m(B)", 0.472, 0., 1.0);
RooRealVar tau(
"tau",
"B0 decay time", 1.547, 1.0, 2.0);
RooRealVar w(
"w",
"Flavor Mistag rate", 0.03, 0.0, 1.0);
RooRealVar dw(
"dw",
"Flavor Mistag rate difference between B0 and B0bar", 0.01);
RooBMixDecay bmix_gm1(
"bmix",
"decay", dt, mixState, tagFlav, tau, dm,
w, dw, gm1,
RooBMixDecay::DoubleSided);
std::unique_ptr<RooDataSet>
data{bmix_gm1.generate({dt, tagFlav, mixState}, 20000)};
bmix_gm1.plotOn(frame);
data->plotOn(frame2,
Cut(
"mixState==mixState::mixed"));
bmix_gm1.plotOn(frame2,
Slice(mixState,
"mixed"));
data->plotOn(frame3,
Cut(
"mixState==mixState::unmixed"));
bmix_gm1.plotOn(frame3,
Slice(mixState,
"unmixed"));
gPad->SetLeftMargin(0.15);
gPad->SetLeftMargin(0.15);
gPad->SetLeftMargin(0.15);
}
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Class RooBMixDecay is a RooAbsAnaConvPdf implementation that describes the decay of B mesons with the...
Object to represent discrete states.
Class RooGaussModel implements a RooResolutionModel that models a Gaussian distribution.
A RooPlot is a plot frame and a container for graphics objects within that frame.
static RooPlot * frame(const RooAbsRealLValue &var, double xmin, double xmax, Int_t nBins)
Create a new frame for a given variable in x.
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
RooRealVar represents a variable that can be changed from the outside.
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
RooCmdArg Slice(const RooArgSet &sliceSet)
RooCmdArg Cut(const char *cutSpec)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
[#1] INFO:Plotting -- RooAbsReal::plotOn(bmix) plot on dt integrates over variables (tagFlav,mixState)
[#1] INFO:Plotting -- RooTreeData::plotOn: plotting 3787 events out of 20000 total events
[#1] INFO:Plotting -- RooAbsReal::plotOn(bmix) plot on dt represents a slice in (mixState)
[#1] INFO:Plotting -- RooAbsReal::plotOn(bmix) plot on dt integrates over variables (tagFlav)
[#1] INFO:Plotting -- RooTreeData::plotOn: plotting 16213 events out of 20000 total events
[#1] INFO:Plotting -- RooAbsReal::plotOn(bmix) plot on dt represents a slice in (mixState)
[#1] INFO:Plotting -- RooAbsReal::plotOn(bmix) plot on dt integrates over variables (tagFlav)
- Date
- July 2008
- Author
- Wouter Verkerke
Definition in file rf310_sliceplot.C.