rf703_effpdfprod.C

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'SPECIAL PDFS' RooFit tutorial macro #703
00004 // 
00005 // Using a product of an (acceptance) efficiency and a p.d.f as p.d.f.
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 "RooExponential.h"
00021 #include "RooEffProd.h"
00022 #include "RooFormulaVar.h"
00023 #include "TCanvas.h"
00024 #include "TAxis.h"
00025 #include "RooPlot.h"
00026 using namespace RooFit ;
00027 
00028 
00029 void rf703_effpdfprod()
00030 {
00031   // D e f i n e   o b s e r v a b l e s   a n d   d e c a y   p d f
00032   // ---------------------------------------------------------------
00033 
00034   // Declare observables
00035   RooRealVar t("t","t",0,5) ;
00036 
00037   // Make pdf
00038   RooRealVar tau("tau","tau",-1.54,-4,-0.1) ;
00039   RooExponential model("model","model",t,tau) ;
00040 
00041 
00042 
00043   // D e f i n e   e f f i c i e n c y   f u n c t i o n
00044   // ---------------------------------------------------
00045   
00046   // Use error function to simulate turn-on slope
00047   RooFormulaVar eff("eff","0.5*(TMath::Erf((t-1)/0.5)+1)",t) ;
00048 
00049 
00050 
00051   // D e f i n e   d e c a y   p d f   w i t h   e f f i c i e n c y 
00052   // ---------------------------------------------------------------
00053 
00054   // Multiply pdf(t) with efficiency in t
00055   RooEffProd modelEff("modelEff","model with efficiency",model,eff) ;
00056 
00057   
00058 
00059   // P l o t   e f f i c i e n c y ,   p d f  
00060   // ----------------------------------------
00061 
00062   RooPlot* frame1 = t.frame(Title("Efficiency")) ;
00063   eff.plotOn(frame1,LineColor(kRed)) ;
00064 
00065   RooPlot* frame2 = t.frame(Title("Pdf with and without efficiency")) ;
00066 
00067   model.plotOn(frame2,LineStyle(kDashed)) ;
00068   modelEff.plotOn(frame2) ;
00069 
00070 
00071 
00072   // G e n e r a t e   t o y   d a t a ,   f i t   m o d e l E f f   t o   d a t a
00073   // ------------------------------------------------------------------------------
00074 
00075   // Generate events. If the input pdf has an internal generator, the internal generator
00076   // is used and an accept/reject sampling on the efficiency is applied. 
00077   RooDataSet* data = modelEff.generate(t,10000) ;
00078 
00079   // Fit pdf. The normalization integral is calculated numerically. 
00080   modelEff.fitTo(*data) ;
00081 
00082   // Plot generated data and overlay fitted pdf
00083   RooPlot* frame3 = t.frame(Title("Fitted pdf with efficiency")) ;
00084   data->plotOn(frame3) ;
00085   modelEff.plotOn(frame3) ;
00086 
00087   
00088   TCanvas* c = new TCanvas("rf703_effpdfprod","rf703_effpdfprod",1200,400) ;
00089   c->Divide(3) ;
00090   c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.4) ; frame1->Draw() ;
00091   c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.6) ; frame2->Draw() ;
00092   c->cd(3) ; gPad->SetLeftMargin(0.15) ; frame3->GetYaxis()->SetTitleOffset(1.6) ; frame3->Draw() ;
00093 
00094 }

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