Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf708_bphysics.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// Special pdf's: special decay pdf for B physics with mixing and/or CP violation
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 "RooGaussian.h"
16#include "RooConstVar.h"
17#include "RooCategory.h"
18#include "RooBMixDecay.h"
19#include "RooBCPEffDecay.h"
20#include "RooBDecay.h"
21#include "RooFormulaVar.h"
22#include "RooTruthModel.h"
23#include "TCanvas.h"
24#include "TAxis.h"
25#include "RooPlot.h"
26using namespace RooFit;
27
28void rf708_bphysics()
29{
30 // -------------------------------------
31 // B - D e c a y w i t h m i x i n g
32 // =====================================
33
34 // C o n s t r u c t p d f
35 // -------------------------
36
37 // Observable
38 RooRealVar dt("dt", "dt", -10, 10);
39 dt.setBins(40);
40
41 // Parameters
42 RooRealVar dm("dm", "delta m(B0)", 0.472);
43 RooRealVar tau("tau", "tau (B0)", 1.547);
44 RooRealVar w("w", "flavour mistag rate", 0.1);
45 RooRealVar dw("dw", "delta mistag rate for B0/B0bar", 0.1);
46
47 RooCategory mixState("mixState", "B0/B0bar mixing state");
48 mixState.defineType("mixed", -1);
49 mixState.defineType("unmixed", 1);
50
51 RooCategory tagFlav("tagFlav", "Flavour of the tagged B0");
52 tagFlav.defineType("B0", 1);
53 tagFlav.defineType("B0bar", -1);
54
55 // Use delta function resolution model
56 RooTruthModel truthModel("tm", "truth model", dt);
57
58 // Construct Bdecay with mixing
59 RooBMixDecay bmix("bmix", "decay", dt, mixState, tagFlav, tau, dm, w, dw, truthModel, RooBMixDecay::DoubleSided);
60
61 // P l o t p d f i n v a r i o u s s l i c e s
62 // ---------------------------------------------------
63
64 // Generate some data
65 std::unique_ptr<RooDataSet> data{bmix.generate({dt, mixState, tagFlav}, 10000)};
66
67 // Plot B0 and B0bar tagged data separately
68 // For all plots below B0 and B0 tagged data will look somewhat differently
69 // if the flavor tagging mistag rate for B0 and B0 is different (i.e. dw!=0)
70 RooPlot *frame1 = dt.frame(Title("B decay distribution with mixing (B0/B0bar)"));
71
72 data->plotOn(frame1, Cut("tagFlav==tagFlav::B0"));
73 bmix.plotOn(frame1, Slice(tagFlav, "B0"));
74
75 data->plotOn(frame1, Cut("tagFlav==tagFlav::B0bar"), MarkerColor(kCyan));
76 bmix.plotOn(frame1, Slice(tagFlav, "B0bar"), LineColor(kCyan));
77
78 // Plot mixed slice for B0 and B0bar tagged data separately
79 RooPlot *frame2 = dt.frame(Title("B decay distribution of mixed events (B0/B0bar)"));
80
81 data->plotOn(frame2, Cut("mixState==mixState::mixed&&tagFlav==tagFlav::B0"));
82 bmix.plotOn(frame2, Slice({{&tagFlav, "B0"}, {&mixState, "mixed"}}));
83
84 data->plotOn(frame2, Cut("mixState==mixState::mixed&&tagFlav==tagFlav::B0bar"), MarkerColor(kCyan));
85 bmix.plotOn(frame2, Slice({{&tagFlav, "B0bar"}, {&mixState, "mixed"}}), LineColor(kCyan));
86
87 // Plot unmixed slice for B0 and B0bar tagged data separately
88 RooPlot *frame3 = dt.frame(Title("B decay distribution of unmixed events (B0/B0bar)"));
89
90 data->plotOn(frame3, Cut("mixState==mixState::unmixed&&tagFlav==tagFlav::B0"));
91 bmix.plotOn(frame3, Slice({{&tagFlav, "B0"}, {&mixState, "unmixed"}}));
92
93 data->plotOn(frame3, Cut("mixState==mixState::unmixed&&tagFlav==tagFlav::B0bar"), MarkerColor(kCyan));
94 bmix.plotOn(frame3, Slice({{&tagFlav, "B0bar"}, {&mixState, "unmixed"}}), LineColor(kCyan));
95
96 // -------------------------------------------------
97 // B - D e c a y w i t h C P v i o l a t i o n
98 // =================================================
99
100 // C o n s t r u c t p d f
101 // -------------------------
102
103 // Additional parameters needed for B decay with CPV
104 RooRealVar CPeigen("CPeigen", "CP eigen value", -1);
105 RooRealVar absLambda("absLambda", "|lambda|", 1, 0, 2);
106 RooRealVar argLambda("absLambda", "|lambda|", 0.7, -1, 1);
107 RooRealVar effR("effR", "B0/B0bar reco efficiency ratio", 1);
108
109 // Construct Bdecay with CP violation
110 RooBCPEffDecay bcp("bcp", "bcp", dt, tagFlav, tau, dm, w, CPeigen, absLambda, argLambda, effR, dw, truthModel,
112
113 // P l o t s c e n a r i o 1 - s i n ( 2 b ) = 0 . 7 , | l | = 1
114 // ---------------------------------------------------------------------------
115
116 // Generate some data
117 std::unique_ptr<RooDataSet> data2{bcp.generate({dt, tagFlav}, 10000)};
118
119 // Plot B0 and B0bar tagged data separately
120 RooPlot *frame4 = dt.frame(Title("B decay distribution with CPV(|l|=1,Im(l)=0.7) (B0/B0bar)"));
121
122 data2->plotOn(frame4, Cut("tagFlav==tagFlav::B0"));
123 bcp.plotOn(frame4, Slice(tagFlav, "B0"));
124
125 data2->plotOn(frame4, Cut("tagFlav==tagFlav::B0bar"), MarkerColor(kCyan));
126 bcp.plotOn(frame4, Slice(tagFlav, "B0bar"), LineColor(kCyan));
127
128 // P l o t s c e n a r i o 2 - s i n ( 2 b ) = 0 . 7 , | l | = 0 . 7
129 // -------------------------------------------------------------------------------
130
131 absLambda = 0.7;
132
133 // Generate some data
134 std::unique_ptr<RooDataSet> data3{bcp.generate({dt, tagFlav}, 10000)};
135
136 // Plot B0 and B0bar tagged data separately (sin2b = 0.7 plus direct CPV |l|=0.5)
137 RooPlot *frame5 = dt.frame(Title("B decay distribution with CPV(|l|=0.7,Im(l)=0.7) (B0/B0bar)"));
138
139 data3->plotOn(frame5, Cut("tagFlav==tagFlav::B0"));
140 bcp.plotOn(frame5, Slice(tagFlav, "B0"));
141
142 data3->plotOn(frame5, Cut("tagFlav==tagFlav::B0bar"), MarkerColor(kCyan));
143 bcp.plotOn(frame5, Slice(tagFlav, "B0bar"), LineColor(kCyan));
144
145 // ---------------------------------------------------------------------------
146 // G e n e r i c B d e c a y w i t h u s e r c o e f f i c i e n t s
147 // ===========================================================================
148
149 // C o n s t r u c t p d f
150 // -------------------------
151
152 // Model parameters
153 RooRealVar DGbG("DGbG", "DGamma/GammaAvg", 0.5, -1, 1);
154 RooRealVar Adir("Adir", "-[1-abs(l)**2]/[1+abs(l)**2]", 0);
155 RooRealVar Amix("Amix", "2Im(l)/[1+abs(l)**2]", 0.7);
156 RooRealVar Adel("Adel", "2Re(l)/[1+abs(l)**2]", 0.7);
157
158 // Derived input parameters for pdf
159 RooFormulaVar DG("DG", "Delta Gamma", "@1/@0", RooArgList(tau, DGbG));
160
161 // Construct coefficient functions for sin,cos,sinh modulations of decay distribution
162 RooFormulaVar fsin("fsin", "fsin", "@0*@1*(1-2*@2)", RooArgList(Amix, tagFlav, w));
163 RooFormulaVar fcos("fcos", "fcos", "@0*@1*(1-2*@2)", RooArgList(Adir, tagFlav, w));
164 RooFormulaVar fsinh("fsinh", "fsinh", "@0", RooArgList(Adel));
165
166 // Construct generic B decay pdf using above user coefficients
167 RooBDecay bcpg("bcpg", "bcpg", dt, tau, DG, RooConst(1), fsinh, fcos, fsin, dm, truthModel, RooBDecay::DoubleSided);
168
169 // P l o t - I m ( l ) = 0 . 7 , R e ( l ) = 0 . 7 | l | = 1, d G / G = 0 . 5
170 // -------------------------------------------------------------------------------------
171
172 // Generate some data
173 std::unique_ptr<RooDataSet> data4{bcpg.generate({dt, tagFlav}, 10000)};
174
175 // Plot B0 and B0bar tagged data separately
176 RooPlot *frame6 = dt.frame(Title("B decay distribution with CPV(Im(l)=0.7,Re(l)=0.7,|l|=1,dG/G=0.5) (B0/B0bar)"));
177
178 data4->plotOn(frame6, Cut("tagFlav==tagFlav::B0"));
179 bcpg.plotOn(frame6, Slice(tagFlav, "B0"));
180
181 data4->plotOn(frame6, Cut("tagFlav==tagFlav::B0bar"), MarkerColor(kCyan));
182 bcpg.plotOn(frame6, Slice(tagFlav, "B0bar"), LineColor(kCyan));
183
184 TCanvas *c = new TCanvas("rf708_bphysics", "rf708_bphysics", 1200, 800);
185 c->Divide(3, 2);
186 c->cd(1);
187 gPad->SetLeftMargin(0.15);
188 frame1->GetYaxis()->SetTitleOffset(1.6);
189 frame1->Draw();
190 c->cd(2);
191 gPad->SetLeftMargin(0.15);
192 frame2->GetYaxis()->SetTitleOffset(1.6);
193 frame2->Draw();
194 c->cd(3);
195 gPad->SetLeftMargin(0.15);
196 frame3->GetYaxis()->SetTitleOffset(1.6);
197 frame3->Draw();
198 c->cd(4);
199 gPad->SetLeftMargin(0.15);
200 frame4->GetYaxis()->SetTitleOffset(1.6);
201 frame4->Draw();
202 c->cd(5);
203 gPad->SetLeftMargin(0.15);
204 frame5->GetYaxis()->SetTitleOffset(1.6);
205 frame5->Draw();
206 c->cd(6);
207 gPad->SetLeftMargin(0.15);
208 frame6->GetYaxis()->SetTitleOffset(1.6);
209 frame6->Draw();
210}
#define c(i)
Definition RSha256.hxx:101
@ kCyan
Definition Rtypes.h:66
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
#define gPad
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
PDF describing decay time distribution of B meson including effects of standard model CP violation.
Most general description of B decay time distribution with effects of CP violation,...
Definition RooBDecay.h:25
@ DoubleSided
Definition RooBDecay.h:29
Class RooBMixDecay is a RooAbsAnaConvPdf implementation that describes the decay of B mesons with the...
Object to represent discrete states.
Definition RooCategory.h:28
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:45
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:225
TAxis * GetYaxis() const
Definition RooPlot.cxx:1264
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:637
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Implements a RooResolution model that corresponds to a delta function.
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
RooConstVar & RooConst(double val)
RooCmdArg MarkerColor(Color_t color)
RooCmdArg Slice(const RooArgSet &sliceSet)
RooCmdArg Cut(const char *cutSpec)
RooCmdArg LineColor(Color_t color)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
const char * Title
Definition TXMLSetup.cxx:68