Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf310_sliceplot.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook
4/// Multidimensional models: projecting pdf and data slices in discrete observables
5///
6/// \macro_image
7/// \macro_code
8/// \macro_output
9///
10/// \date July 2008
11/// \author Wouter Verkerke
12
13#include "RooRealVar.h"
14#include "RooDataSet.h"
15#include "RooGaussModel.h"
16#include "RooDecay.h"
17#include "RooBMixDecay.h"
18#include "RooCategory.h"
19#include "TCanvas.h"
20#include "TAxis.h"
21#include "RooPlot.h"
22using namespace RooFit;
23
24void rf310_sliceplot()
25{
26
27 // C r e a t e B d e c a y p d f w it h m i x i n g
28 // ----------------------------------------------------------
29
30 // Decay time observables
31 RooRealVar dt("dt", "dt", -20, 20);
32
33 // Discrete observables mixState (B0tag==B0reco?) and tagFlav (B0tag==B0(bar)?)
34 RooCategory mixState("mixState", "B0/B0bar mixing state");
35 RooCategory tagFlav("tagFlav", "Flavour of the tagged B0");
36
37 // Define state labels of discrete observables
38 mixState.defineType("mixed", -1);
39 mixState.defineType("unmixed", 1);
40 tagFlav.defineType("B0", 1);
41 tagFlav.defineType("B0bar", -1);
42
43 // Model parameters
44 RooRealVar dm("dm", "delta m(B)", 0.472, 0., 1.0);
45 RooRealVar tau("tau", "B0 decay time", 1.547, 1.0, 2.0);
46 RooRealVar w("w", "Flavor Mistag rate", 0.03, 0.0, 1.0);
47 RooRealVar dw("dw", "Flavor Mistag rate difference between B0 and B0bar", 0.01);
48
49 // Build a gaussian resolution model
50 RooRealVar bias1("bias1", "bias1", 0);
51 RooRealVar sigma1("sigma1", "sigma1", 0.01);
52 RooGaussModel gm1("gm1", "gauss model 1", dt, bias1, sigma1);
53
54 // Construct a decay pdf, smeared with single gaussian resolution model
55 RooBMixDecay bmix_gm1("bmix", "decay", dt, mixState, tagFlav, tau, dm, w, dw, gm1, RooBMixDecay::DoubleSided);
56
57 // Generate BMixing data with above set of event errors
58 std::unique_ptr<RooDataSet> data{bmix_gm1.generate({dt, tagFlav, mixState}, 20000)};
59
60 // P l o t f u l l d e c a y d i s t r i b u t i o n
61 // ----------------------------------------------------------
62
63 // Create frame, plot data and pdf projection (integrated over tagFlav and mixState)
64 RooPlot *frame = dt.frame(Title("Inclusive decay distribution"));
65 data->plotOn(frame);
66 bmix_gm1.plotOn(frame);
67
68 // P l o t d e c a y d i s t r . f o r m i x e d a n d u n m i x e d s l i c e o f m i x S t a t e
69 // ------------------------------------------------------------------------------------------------------------------
70
71 // Create frame, plot data (mixed only)
72 RooPlot *frame2 = dt.frame(Title("Decay distribution of mixed events"));
73 data->plotOn(frame2, Cut("mixState==mixState::mixed"));
74
75 // Position slice in mixState at "mixed" and plot slice of pdf in mixstate over data (integrated over tagFlav)
76 bmix_gm1.plotOn(frame2, Slice(mixState, "mixed"));
77
78 // Create frame, plot data (unmixed only)
79 RooPlot *frame3 = dt.frame(Title("Decay distribution of unmixed events"));
80 data->plotOn(frame3, Cut("mixState==mixState::unmixed"));
81
82 // Position slice in mixState at "unmixed" and plot slice of pdf in mixstate over data (integrated over tagFlav)
83 bmix_gm1.plotOn(frame3, Slice(mixState, "unmixed"));
84
85 TCanvas *c = new TCanvas("rf310_sliceplot", "rf310_sliceplot", 1200, 400);
86 c->Divide(3);
87 c->cd(1);
88 gPad->SetLeftMargin(0.15);
89 frame->GetYaxis()->SetTitleOffset(1.4);
90 gPad->SetLogy();
91 frame->Draw();
92 c->cd(2);
93 gPad->SetLeftMargin(0.15);
94 frame2->GetYaxis()->SetTitleOffset(1.4);
95 gPad->SetLogy();
96 frame2->Draw();
97 c->cd(3);
98 gPad->SetLeftMargin(0.15);
99 frame3->GetYaxis()->SetTitleOffset(1.4);
100 gPad->SetLogy();
101 frame3->Draw();
102}
#define c(i)
Definition RSha256.hxx:101
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
#define gPad
Class RooBMixDecay is a RooAbsAnaConvPdf implementation that describes the decay of B mesons with the...
Object to represent discrete states.
Definition RooCategory.h:28
Class RooGaussModel implements a RooResolutionModel that models a Gaussian distribution.
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
static RooPlot * frame(const RooAbsRealLValue &var, double xmin, double xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition RooPlot.cxx:185
TAxis * GetYaxis() const
Definition RooPlot.cxx:1224
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:597
Variable that can be changed from the outside.
Definition RooRealVar.h:37
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:298
The Canvas class.
Definition TCanvas.h:23
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...
Definition CodegenImpl.h:64
const char * Title
Definition TXMLSetup.cxx:68