00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013 #ifndef __CINT__
00014 #include "RooGlobalFunc.h"
00015 #endif
00016 #include "RooRealVar.h"
00017 #include "RooDataSet.h"
00018 #include "RooGaussian.h"
00019 #include "RooCategory.h"
00020 #include "RooBMixDecay.h"
00021 #include "RooBCPEffDecay.h"
00022 #include "RooBDecay.h"
00023 #include "RooFormulaVar.h"
00024 #include "RooTruthModel.h"
00025 #include "TCanvas.h"
00026 #include "RooPlot.h"
00027 using namespace RooFit ;
00028
00029
00030 class TestBasic708 : public RooFitTestUnit
00031 {
00032 public:
00033 TestBasic708(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("B Physics p.d.f.s",refFile,writeRef,verbose) {} ;
00034 Bool_t testCode() {
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044 RooRealVar dt("dt","dt",-10,10) ;
00045 dt.setBins(40) ;
00046
00047
00048 RooRealVar dm("dm","delta m(B0)",0.472) ;
00049 RooRealVar tau("tau","tau (B0)",1.547) ;
00050 RooRealVar w("w","flavour mistag rate",0.1) ;
00051 RooRealVar dw("dw","delta mistag rate for B0/B0bar",0.1) ;
00052
00053 RooCategory mixState("mixState","B0/B0bar mixing state") ;
00054 mixState.defineType("mixed",-1) ;
00055 mixState.defineType("unmixed",1) ;
00056
00057 RooCategory tagFlav("tagFlav","Flavour of the tagged B0") ;
00058 tagFlav.defineType("B0",1) ;
00059 tagFlav.defineType("B0bar",-1) ;
00060
00061
00062 RooTruthModel tm("tm","truth model",dt) ;
00063
00064
00065 RooBMixDecay bmix("bmix","decay",dt,mixState,tagFlav,tau,dm,w,dw,tm,RooBMixDecay::DoubleSided) ;
00066
00067
00068
00069
00070
00071
00072
00073 RooDataSet* data = bmix.generate(RooArgSet(dt,mixState,tagFlav),10000) ;
00074
00075
00076
00077
00078 RooPlot* frame1 = dt.frame(Title("B decay distribution with mixing (B0/B0bar)")) ;
00079
00080 data->plotOn(frame1,Cut("tagFlav==tagFlav::B0")) ;
00081 bmix.plotOn(frame1,Slice(tagFlav,"B0")) ;
00082
00083 data->plotOn(frame1,Cut("tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
00084 bmix.plotOn(frame1,Slice(tagFlav,"B0bar"),LineColor(kCyan),Name("alt")) ;
00085
00086
00087
00088 RooPlot* frame2 = dt.frame(Title("B decay distribution of mixed events (B0/B0bar)")) ;
00089
00090 data->plotOn(frame2,Cut("mixState==mixState::mixed&&tagFlav==tagFlav::B0")) ;
00091 bmix.plotOn(frame2,Slice(tagFlav,"B0"),Slice(mixState,"mixed")) ;
00092
00093 data->plotOn(frame2,Cut("mixState==mixState::mixed&&tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
00094 bmix.plotOn(frame2,Slice(tagFlav,"B0bar"),Slice(mixState,"mixed"),LineColor(kCyan),Name("alt")) ;
00095
00096
00097
00098 RooPlot* frame3 = dt.frame(Title("B decay distribution of unmixed events (B0/B0bar)")) ;
00099
00100 data->plotOn(frame3,Cut("mixState==mixState::unmixed&&tagFlav==tagFlav::B0")) ;
00101 bmix.plotOn(frame3,Slice(tagFlav,"B0"),Slice(mixState,"unmixed")) ;
00102
00103 data->plotOn(frame3,Cut("mixState==mixState::unmixed&&tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
00104 bmix.plotOn(frame3,Slice(tagFlav,"B0bar"),Slice(mixState,"unmixed"),LineColor(kCyan),Name("alt")) ;
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118 RooRealVar CPeigen("CPeigen","CP eigen value",-1) ;
00119 RooRealVar absLambda("absLambda","|lambda|",1,0,2) ;
00120 RooRealVar argLambda("absLambda","|lambda|",0.7,-1,1) ;
00121 RooRealVar effR("effR","B0/B0bar reco efficiency ratio",1) ;
00122
00123
00124 RooBCPEffDecay bcp("bcp","bcp", dt, tagFlav, tau, dm, w, CPeigen, absLambda, argLambda, effR, dw, tm, RooBCPEffDecay::DoubleSided) ;
00125
00126
00127
00128
00129
00130
00131
00132 RooDataSet* data2 = bcp.generate(RooArgSet(dt,tagFlav),10000) ;
00133
00134
00135 RooPlot* frame4 = dt.frame(Title("B decay distribution with CPV(|l|=1,Im(l)=0.7) (B0/B0bar)")) ;
00136
00137 data2->plotOn(frame4,Cut("tagFlav==tagFlav::B0")) ;
00138 bcp.plotOn(frame4,Slice(tagFlav,"B0")) ;
00139
00140 data2->plotOn(frame4,Cut("tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
00141 bcp.plotOn(frame4,Slice(tagFlav,"B0bar"),LineColor(kCyan),Name("alt")) ;
00142
00143
00144
00145
00146
00147
00148 absLambda=0.7 ;
00149
00150
00151 RooDataSet* data3 = bcp.generate(RooArgSet(dt,tagFlav),10000) ;
00152
00153
00154 RooPlot* frame5 = dt.frame(Title("B decay distribution with CPV(|l|=0.7,Im(l)=0.7) (B0/B0bar)")) ;
00155
00156 data3->plotOn(frame5,Cut("tagFlav==tagFlav::B0")) ;
00157 bcp.plotOn(frame5,Slice(tagFlav,"B0")) ;
00158
00159 data3->plotOn(frame5,Cut("tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
00160 bcp.plotOn(frame5,Slice(tagFlav,"B0bar"),LineColor(kCyan),Name("alt")) ;
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172 RooRealVar DGbG("DGbG","DGamma/GammaAvg",0.5,-1,1);
00173 RooRealVar Adir("Adir","-[1-abs(l)**2]/[1+abs(l)**2]",0);
00174 RooRealVar Amix("Amix","2Im(l)/[1+abs(l)**2]",0.7);
00175 RooRealVar Adel("Adel","2Re(l)/[1+abs(l)**2]",0.7);
00176
00177
00178 RooFormulaVar DG("DG","Delta Gamma","@1/@0",RooArgList(tau,DGbG));
00179
00180
00181 RooFormulaVar fsin("fsin","fsin","@0*@1*(1-2*@2)",RooArgList(Amix,tagFlav,w));
00182 RooFormulaVar fcos("fcos","fcos","@0*@1*(1-2*@2)",RooArgList(Adir,tagFlav,w));
00183 RooFormulaVar fsinh("fsinh","fsinh","@0",RooArgList(Adel));
00184
00185
00186 RooBDecay bcpg("bcpg","bcpg",dt,tau,DG,RooConst(1),fsinh,fcos,fsin,dm,tm,RooBDecay::DoubleSided);
00187
00188
00189
00190
00191
00192
00193
00194 RooDataSet* data4 = bcpg.generate(RooArgSet(dt,tagFlav),10000) ;
00195
00196
00197 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)")) ;
00198
00199 data4->plotOn(frame6,Cut("tagFlav==tagFlav::B0")) ;
00200 bcpg.plotOn(frame6,Slice(tagFlav,"B0")) ;
00201
00202 data4->plotOn(frame6,Cut("tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
00203 bcpg.plotOn(frame6,Slice(tagFlav,"B0bar"),LineColor(kCyan),Name("alt")) ;
00204
00205
00206 regPlot(frame1,"rf708_plot1") ;
00207 regPlot(frame2,"rf708_plot2") ;
00208 regPlot(frame3,"rf708_plot3") ;
00209 regPlot(frame4,"rf708_plot4") ;
00210 regPlot(frame5,"rf708_plot5") ;
00211 regPlot(frame6,"rf708_plot6") ;
00212
00213 delete data ;
00214 delete data2 ;
00215 delete data3 ;
00216 delete data4 ;
00217
00218 return kTRUE ;
00219 }
00220 } ;