30 class TestBasic708 :
public RooFitTestUnit
33 TestBasic708(TFile* refFile,
Bool_t writeRef,
Int_t verbose) : RooFitTestUnit(
"B Physics p.d.f.s",refFile,writeRef,verbose) {} ;
51 RooRealVar dw(
"dw",
"delta mistag rate for B0/B0bar",0.1) ;
53 RooCategory mixState(
"mixState",
"B0/B0bar mixing state") ;
54 mixState.defineType(
"mixed",-1) ;
55 mixState.defineType(
"unmixed",1) ;
57 RooCategory tagFlav(
"tagFlav",
"Flavour of the tagged B0") ;
58 tagFlav.defineType(
"B0",1) ;
59 tagFlav.defineType(
"B0bar",-1) ;
65 RooBMixDecay bmix(
"bmix",
"decay",dt,mixState,tagFlav,
tau,dm,w,dw,tm,
RooBMixDecay::DoubleSided) ;
78 RooPlot* frame1 = dt.frame(
Title(
"B decay distribution with mixing (B0/B0bar)")) ;
80 data->
plotOn(frame1,
Cut(
"tagFlav==tagFlav::B0")) ;
81 bmix.plotOn(frame1,
Slice(tagFlav,
"B0")) ;
88 RooPlot* frame2 = dt.frame(
Title(
"B decay distribution of mixed events (B0/B0bar)")) ;
90 data->
plotOn(frame2,
Cut(
"mixState==mixState::mixed&&tagFlav==tagFlav::B0")) ;
91 bmix.plotOn(frame2,
Slice(tagFlav,
"B0"),
Slice(mixState,
"mixed")) ;
98 RooPlot* frame3 = dt.frame(
Title(
"B decay distribution of unmixed events (B0/B0bar)")) ;
100 data->
plotOn(frame3,
Cut(
"mixState==mixState::unmixed&&tagFlav==tagFlav::B0")) ;
101 bmix.plotOn(frame3,
Slice(tagFlav,
"B0"),
Slice(mixState,
"unmixed")) ;
118 RooRealVar CPeigen(
"CPeigen",
"CP eigen value",-1) ;
119 RooRealVar absLambda(
"absLambda",
"|lambda|",1,0,2) ;
120 RooRealVar argLambda(
"absLambda",
"|lambda|",0.7,-1,1) ;
121 RooRealVar effR(
"effR",
"B0/B0bar reco efficiency ratio",1) ;
124 RooBCPEffDecay bcp(
"bcp",
"bcp", dt, tagFlav,
tau, dm, w, CPeigen, absLambda, argLambda, effR, dw, tm,
RooBCPEffDecay::DoubleSided) ;
135 RooPlot* frame4 = dt.frame(
Title(
"B decay distribution with CPV(|l|=1,Im(l)=0.7) (B0/B0bar)")) ;
137 data2->
plotOn(frame4,
Cut(
"tagFlav==tagFlav::B0")) ;
138 bcp.plotOn(frame4,
Slice(tagFlav,
"B0")) ;
154 RooPlot* frame5 = dt.frame(
Title(
"B decay distribution with CPV(|l|=0.7,Im(l)=0.7) (B0/B0bar)")) ;
156 data3->
plotOn(frame5,
Cut(
"tagFlav==tagFlav::B0")) ;
157 bcp.plotOn(frame5,
Slice(tagFlav,
"B0")) ;
172 RooRealVar DGbG(
"DGbG",
"DGamma/GammaAvg",0.5,-1,1);
173 RooRealVar Adir(
"Adir",
"-[1-abs(l)**2]/[1+abs(l)**2]",0);
174 RooRealVar Amix(
"Amix",
"2Im(l)/[1+abs(l)**2]",0.7);
175 RooRealVar Adel(
"Adel",
"2Re(l)/[1+abs(l)**2]",0.7);
186 RooBDecay bcpg(
"bcpg",
"bcpg",dt,
tau,DG,
RooConst(1),fsinh,fcos,fsin,dm,tm,
RooBDecay::DoubleSided);
197 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)")) ;
199 data4->
plotOn(frame6,
Cut(
"tagFlav==tagFlav::B0")) ;
200 bcpg.plotOn(frame6,
Slice(tagFlav,
"B0")) ;
206 regPlot(frame1,
"rf708_plot1") ;
207 regPlot(frame2,
"rf708_plot2") ;
208 regPlot(frame3,
"rf708_plot3") ;
209 regPlot(frame4,
"rf708_plot4") ;
210 regPlot(frame5,
"rf708_plot5") ;
211 regPlot(frame6,
"rf708_plot6") ;
RooTruthModel is an implementation of RooResolution model that provides a delta-function resolution m...
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Plot dataset on specified frame.
RooCmdArg Cut(const char *cutSpec)
RooCmdArg LineColor(Color_t color)
PDF describing decay time distribution of B meson including effects of standard model CP violation...
RooCmdArg Title(const char *name)
Most general description of B decay time distribution with effects of CP violation, mixing and life time differences.
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t tau
Class RooBMixDecay is a RooAbsAnaConvPdf implementation that describes the decay of B mesons with the...
RooRealVar represents a fundamental (non-derived) real valued object.
RooCmdArg Name(const char *name)
RooDataSet is a container class to hold unbinned data.
RooCategory represents a fundamental (non-derived) discrete value object.
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
RooConstVar & RooConst(Double_t val)
RooCmdArg MarkerColor(Color_t color)
RooCmdArg Slice(const RooArgSet &sliceSet)