rf201_composite.C

Go to the documentation of this file.
00001 /////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'ADDITION AND CONVOLUTION' RooFit tutorial macro #201
00004 // 
00005 // Composite p.d.f with signal and background component
00006 //
00007 // pdf = f_bkg * bkg(x,a0,a1) + (1-fbkg) * (f_sig1 * sig1(x,m,s1 + (1-f_sig1) * sig2(x,m,s2)))
00008 //
00009 //
00010 // 07/2008 - Wouter Verkerke 
00011 //
00012 /////////////////////////////////////////////////////////////////////////
00013 
00014 #ifndef __CINT__
00015 #include "RooGlobalFunc.h"
00016 #endif
00017 #include "RooRealVar.h"
00018 #include "RooDataSet.h"
00019 #include "RooGaussian.h"
00020 #include "RooChebychev.h"
00021 #include "RooAddPdf.h"
00022 #include "TCanvas.h"
00023 #include "TAxis.h"
00024 #include "RooPlot.h"
00025 using namespace RooFit ;
00026 
00027 
00028 void rf201_composite()
00029 {
00030   // S e t u p   c o m p o n e n t   p d f s 
00031   // ---------------------------------------
00032 
00033   // Declare observable x
00034   RooRealVar x("x","x",0,10) ;
00035 
00036   // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their paramaters
00037   RooRealVar mean("mean","mean of gaussians",5) ;
00038   RooRealVar sigma1("sigma1","width of gaussians",0.5) ;
00039   RooRealVar sigma2("sigma2","width of gaussians",1) ;
00040 
00041   RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;  
00042   RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;  
00043   
00044   // Build Chebychev polynomial p.d.f.  
00045   RooRealVar a0("a0","a0",0.5,0.,1.) ;
00046   RooRealVar a1("a1","a1",-0.2,0.,1.) ;
00047   RooChebychev bkg("bkg","Background",x,RooArgSet(a0,a1)) ;
00048 
00049 
00050   ////////////////////////////////////////////////////
00051   // M E T H O D   1 - T w o   R o o A d d P d f s  //
00052   ////////////////////////////////////////////////////
00053 
00054 
00055   // A d d   s i g n a l   c o m p o n e n t s 
00056   // ------------------------------------------
00057 
00058   // Sum the signal components into a composite signal p.d.f.
00059   RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
00060   RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;
00061 
00062 
00063   // A d d  s i g n a l   a n d   b a c k g r o u n d
00064   // ------------------------------------------------
00065 
00066   // Sum the composite signal and background 
00067   RooRealVar bkgfrac("bkgfrac","fraction of background",0.5,0.,1.) ;
00068   RooAddPdf  model("model","g1+g2+a",RooArgList(bkg,sig),bkgfrac) ;
00069 
00070 
00071   // S a m p l e ,   f i t   a n d   p l o t   m o d e l
00072   // ---------------------------------------------------
00073 
00074   // Generate a data sample of 1000 events in x from model
00075   RooDataSet *data = model.generate(x,1000) ;
00076 
00077   // Fit model to data
00078   model.fitTo(*data) ;
00079   
00080   // Plot data and PDF overlaid
00081   RooPlot* xframe = x.frame(Title("Example of composite pdf=(sig1+sig2)+bkg")) ;
00082   data->plotOn(xframe) ;
00083   model.plotOn(xframe) ;
00084 
00085   // Overlay the background component of model with a dashed line
00086   model.plotOn(xframe,Components(bkg),LineStyle(kDashed)) ;
00087 
00088   // Overlay the background+sig2 components of model with a dotted line
00089   model.plotOn(xframe,Components(RooArgSet(bkg,sig2)),LineStyle(kDotted)) ;
00090 
00091   // Print structure of composite p.d.f.
00092   model.Print("t") ;
00093 
00094 
00095   ////////////////////////////////////////////////////////////////////////////////////////////////////
00096   // M E T H O D   2 - O n e   R o o A d d P d f   w i t h   r e c u r s i v e   f r a c t i o n s  //
00097   ////////////////////////////////////////////////////////////////////////////////////////////////////
00098 
00099   // Construct sum of models on one go using recursive fraction interpretations
00100   //
00101   //   model2 = bkg + (sig1 + sig2)
00102   //
00103   RooAddPdf  model2("model","g1+g2+a",RooArgList(bkg,sig1,sig2),RooArgList(bkgfrac,sig1frac),kTRUE) ;    
00104 
00105   // NB: Each coefficient is interpreted as the fraction of the 
00106   // left-hand component of the i-th recursive sum, i.e.
00107   //
00108   //   sum4 = A + ( B + ( C + D)  with fraction fA, fB and fC expands to
00109   //
00110   //   sum4 = fA*A + (1-fA)*(fB*B + (1-fB)*(fC*C + (1-fC)*D))
00111 
00112 
00113   // P l o t   r e c u r s i v e   a d d i t i o n   m o d e l
00114   // ---------------------------------------------------------
00115   model2.plotOn(xframe,LineColor(kRed),LineStyle(kDashed)) ;
00116   model2.plotOn(xframe,Components(RooArgSet(bkg,sig2)),LineColor(kRed),LineStyle(kDashed)) ;
00117   model2.Print("t") ;
00118 
00119   
00120   // Draw the frame on the canvas
00121   new TCanvas("rf201_composite","rf201_composite",600,600) ;
00122   gPad->SetLeftMargin(0.15) ; xframe->GetYaxis()->SetTitleOffset(1.4) ; xframe->Draw() ;
00123 
00124 
00125 }
00126 

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