rf708_bphysics.C

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'SPECIAL PDFS' RooFit tutorial macro #708
00004 // 
00005 // Special decay pdf for B physics with mixing and/or CP violation
00006 //
00007 //
00008 //
00009 // 07/2008 - Wouter Verkerke 
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 "RooConstVar.h"
00020 #include "RooCategory.h"
00021 #include "RooBMixDecay.h"
00022 #include "RooBCPEffDecay.h"
00023 #include "RooBDecay.h"
00024 #include "RooFormulaVar.h"
00025 #include "RooTruthModel.h"
00026 #include "TCanvas.h"
00027 #include "TAxis.h"
00028 #include "RooPlot.h"
00029 using namespace RooFit ;
00030 
00031 void rf708_bphysics()
00032 {
00033   ////////////////////////////////////////////////////
00034   // B - D e c a y   w i t h   m i x i n g          //
00035   ////////////////////////////////////////////////////
00036 
00037   // C o n s t r u c t   p d f 
00038   // -------------------------
00039   
00040   // Observable
00041   RooRealVar dt("dt","dt",-10,10) ;
00042   dt.setBins(40) ;
00043 
00044   // Parameters
00045   RooRealVar dm("dm","delta m(B0)",0.472) ;
00046   RooRealVar tau("tau","tau (B0)",1.547) ;
00047   RooRealVar w("w","flavour mistag rate",0.1) ;
00048   RooRealVar dw("dw","delta mistag rate for B0/B0bar",0.1) ;
00049 
00050   RooCategory mixState("mixState","B0/B0bar mixing state") ;
00051   mixState.defineType("mixed",-1) ;
00052   mixState.defineType("unmixed",1) ;
00053 
00054   RooCategory tagFlav("tagFlav","Flavour of the tagged B0") ;
00055   tagFlav.defineType("B0",1) ;
00056   tagFlav.defineType("B0bar",-1) ;
00057 
00058   // Use delta function resolution model
00059   RooTruthModel tm("tm","truth model",dt) ;
00060 
00061   // Construct Bdecay with mixing
00062   RooBMixDecay bmix("bmix","decay",dt,mixState,tagFlav,tau,dm,w,dw,tm,RooBMixDecay::DoubleSided) ;
00063 
00064 
00065 
00066   // P l o t   p d f   i n   v a r i o u s   s l i c e s
00067   // ---------------------------------------------------
00068 
00069   // Generate some data
00070   RooDataSet* data = bmix.generate(RooArgSet(dt,mixState,tagFlav),10000) ;
00071 
00072   // Plot B0 and B0bar tagged data separately
00073   // For all plots below B0 and B0 tagged data will look somewhat differently
00074   // if the flavor tagging mistag rate for B0 and B0 is different (i.e. dw!=0)
00075   RooPlot* frame1 = dt.frame(Title("B decay distribution with mixing (B0/B0bar)")) ;
00076 
00077   data->plotOn(frame1,Cut("tagFlav==tagFlav::B0")) ;
00078   bmix.plotOn(frame1,Slice(tagFlav,"B0")) ;
00079 
00080   data->plotOn(frame1,Cut("tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
00081   bmix.plotOn(frame1,Slice(tagFlav,"B0bar"),LineColor(kCyan)) ;
00082 
00083 
00084   // Plot mixed slice for B0 and B0bar tagged data separately
00085   RooPlot* frame2 = dt.frame(Title("B decay distribution of mixed events (B0/B0bar)")) ;
00086 
00087   data->plotOn(frame2,Cut("mixState==mixState::mixed&&tagFlav==tagFlav::B0")) ;
00088   bmix.plotOn(frame2,Slice(tagFlav,"B0"),Slice(mixState,"mixed")) ;
00089 
00090   data->plotOn(frame2,Cut("mixState==mixState::mixed&&tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
00091   bmix.plotOn(frame2,Slice(tagFlav,"B0bar"),Slice(mixState,"mixed"),LineColor(kCyan)) ;
00092 
00093 
00094   // Plot unmixed slice for B0 and B0bar tagged data separately
00095   RooPlot* frame3 = dt.frame(Title("B decay distribution of unmixed events (B0/B0bar)")) ;
00096 
00097   data->plotOn(frame3,Cut("mixState==mixState::unmixed&&tagFlav==tagFlav::B0")) ;
00098   bmix.plotOn(frame3,Slice(tagFlav,"B0"),Slice(mixState,"unmixed")) ;
00099 
00100   data->plotOn(frame3,Cut("mixState==mixState::unmixed&&tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
00101   bmix.plotOn(frame3,Slice(tagFlav,"B0bar"),Slice(mixState,"unmixed"),LineColor(kCyan)) ;
00102 
00103 
00104 
00105 
00106 
00107   ///////////////////////////////////////////////////////
00108   // B - D e c a y   w i t h   C P   v i o l a t i o n //
00109   ///////////////////////////////////////////////////////
00110 
00111   // C o n s t r u c t   p d f 
00112   // -------------------------
00113   
00114   // Additional parameters needed for B decay with CPV
00115   RooRealVar CPeigen("CPeigen","CP eigen value",-1) ;
00116   RooRealVar absLambda("absLambda","|lambda|",1,0,2) ;
00117   RooRealVar argLambda("absLambda","|lambda|",0.7,-1,1) ;
00118   RooRealVar effR("effR","B0/B0bar reco efficiency ratio",1) ;
00119 
00120   // Construct Bdecay with CP violation
00121   RooBCPEffDecay bcp("bcp","bcp", dt, tagFlav, tau, dm, w, CPeigen, absLambda, argLambda, effR, dw, tm, RooBCPEffDecay::DoubleSided) ;
00122       
00123 
00124 
00125   // P l o t   s c e n a r i o   1   -   s i n ( 2 b )   =   0 . 7 ,   | l | = 1 
00126   // ---------------------------------------------------------------------------
00127 
00128   // Generate some data
00129   RooDataSet* data2 = bcp.generate(RooArgSet(dt,tagFlav),10000) ;
00130 
00131   // Plot B0 and B0bar tagged data separately
00132   RooPlot* frame4 = dt.frame(Title("B decay distribution with CPV(|l|=1,Im(l)=0.7) (B0/B0bar)")) ;
00133 
00134   data2->plotOn(frame4,Cut("tagFlav==tagFlav::B0")) ;
00135   bcp.plotOn(frame4,Slice(tagFlav,"B0")) ;
00136 
00137   data2->plotOn(frame4,Cut("tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
00138   bcp.plotOn(frame4,Slice(tagFlav,"B0bar"),LineColor(kCyan)) ;
00139 
00140 
00141 
00142   // P l o t   s c e n a r i o   2   -   s i n ( 2 b )   =   0 . 7 ,   | l | = 0 . 7 
00143   // -------------------------------------------------------------------------------
00144 
00145   absLambda=0.7 ;
00146 
00147   // Generate some data
00148   RooDataSet* data3 = bcp.generate(RooArgSet(dt,tagFlav),10000) ;
00149 
00150   // Plot B0 and B0bar tagged data separately (sin2b = 0.7 plus direct CPV |l|=0.5)
00151   RooPlot* frame5 = dt.frame(Title("B decay distribution with CPV(|l|=0.7,Im(l)=0.7) (B0/B0bar)")) ;  
00152 
00153   data3->plotOn(frame5,Cut("tagFlav==tagFlav::B0")) ;
00154   bcp.plotOn(frame5,Slice(tagFlav,"B0")) ;
00155 
00156   data3->plotOn(frame5,Cut("tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
00157   bcp.plotOn(frame5,Slice(tagFlav,"B0bar"),LineColor(kCyan)) ;
00158 
00159 
00160 
00161   //////////////////////////////////////////////////////////////////////////////////
00162   // 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  //
00163   //////////////////////////////////////////////////////////////////////////////////
00164 
00165   // C o n s t r u c t   p d f 
00166   // -------------------------
00167   
00168   // Model parameters
00169   RooRealVar DGbG("DGbG","DGamma/GammaAvg",0.5,-1,1);
00170   RooRealVar Adir("Adir","-[1-abs(l)**2]/[1+abs(l)**2]",0);
00171   RooRealVar Amix("Amix","2Im(l)/[1+abs(l)**2]",0.7);
00172   RooRealVar Adel("Adel","2Re(l)/[1+abs(l)**2]",0.7);
00173   
00174   // Derived input parameters for pdf
00175   RooFormulaVar DG("DG","Delta Gamma","@1/@0",RooArgList(tau,DGbG));
00176   
00177   // Construct coefficient functions for sin,cos,sinh modulations of decay distribution
00178   RooFormulaVar fsin("fsin","fsin","@0*@1*(1-2*@2)",RooArgList(Amix,tagFlav,w));
00179   RooFormulaVar fcos("fcos","fcos","@0*@1*(1-2*@2)",RooArgList(Adir,tagFlav,w));
00180   RooFormulaVar fsinh("fsinh","fsinh","@0",RooArgList(Adel));
00181   
00182   // Construct generic B decay pdf using above user coefficients
00183   RooBDecay bcpg("bcpg","bcpg",dt,tau,DG,RooConst(1),fsinh,fcos,fsin,dm,tm,RooBDecay::DoubleSided);
00184   
00185   
00186   
00187   // P l o t   -   I m ( l ) = 0 . 7 ,   R e ( l ) = 0 . 7   | l | = 1,   d G / G = 0 . 5 
00188   // -------------------------------------------------------------------------------------
00189   
00190   // Generate some data
00191   RooDataSet* data4 = bcpg.generate(RooArgSet(dt,tagFlav),10000) ;
00192   
00193   // Plot B0 and B0bar tagged data separately 
00194   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)")) ;  
00195   
00196   data4->plotOn(frame6,Cut("tagFlav==tagFlav::B0")) ;
00197   bcpg.plotOn(frame6,Slice(tagFlav,"B0")) ;
00198   
00199   data4->plotOn(frame6,Cut("tagFlav==tagFlav::B0bar"),MarkerColor(kCyan)) ;
00200   bcpg.plotOn(frame6,Slice(tagFlav,"B0bar"),LineColor(kCyan)) ;
00201   
00202  
00203  
00204 
00205   TCanvas* c = new TCanvas("rf708_bphysics","rf708_bphysics",1200,800) ;
00206   c->Divide(3,2) ;
00207   c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.6) ; frame1->Draw() ;
00208   c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.6) ; frame2->Draw() ;
00209   c->cd(3) ; gPad->SetLeftMargin(0.15) ; frame3->GetYaxis()->SetTitleOffset(1.6) ; frame3->Draw() ;
00210   c->cd(4) ; gPad->SetLeftMargin(0.15) ; frame4->GetYaxis()->SetTitleOffset(1.6) ; frame4->Draw() ;
00211   c->cd(5) ; gPad->SetLeftMargin(0.15) ; frame5->GetYaxis()->SetTitleOffset(1.6) ; frame5->Draw() ;
00212   c->cd(6) ; gPad->SetLeftMargin(0.15) ; frame6->GetYaxis()->SetTitleOffset(1.6) ; frame6->Draw() ;
00213   
00214 }

Generated on Tue Jul 5 15:45:09 2011 for ROOT_528-00b_version by  doxygen 1.5.1