rf307_fullpereventerrors.C

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #307
00004 // 
00005 // Complete example with use of full p.d.f. with per-event errors
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 "RooGaussModel.h"
00020 #include "RooConstVar.h"
00021 #include "RooDecay.h"
00022 #include "RooLandau.h"
00023 #include "RooProdPdf.h"
00024 #include "RooHistPdf.h"
00025 #include "RooPlot.h"
00026 #include "TCanvas.h"
00027 #include "TAxis.h"
00028 #include "TH1.h"
00029 using namespace RooFit ;
00030 
00031 
00032 void rf307_fullpereventerrors()
00033 {
00034   // B - p h y s i c s   p d f   w i t h   p e r - e v e n t  G a u s s i a n   r e s o l u t i o n
00035   // ----------------------------------------------------------------------------------------------
00036 
00037   // Observables
00038   RooRealVar dt("dt","dt",-10,10) ;
00039   RooRealVar dterr("dterr","per-event error on dt",0.01,10) ;
00040 
00041   // Build a gaussian resolution model scaled by the per-event error = gauss(dt,bias,sigma*dterr)
00042   RooRealVar bias("bias","bias",0,-10,10) ;
00043   RooRealVar sigma("sigma","per-event error scale factor",1,0.1,10) ;
00044   RooGaussModel gm("gm1","gauss model scaled bt per-event error",dt,bias,sigma,dterr) ;
00045 
00046   // Construct decay(dt) (x) gauss1(dt|dterr)
00047   RooRealVar tau("tau","tau",1.548) ;
00048   RooDecay decay_gm("decay_gm","decay",dt,tau,gm,RooDecay::DoubleSided) ;
00049 
00050 
00051 
00052   // C o n s t r u c t   e m p i r i c a l   p d f   f o r   p e r - e v e n t   e r r o r
00053   // -----------------------------------------------------------------
00054 
00055   // Use landau p.d.f to get empirical distribution with long tail
00056   RooLandau pdfDtErr("pdfDtErr","pdfDtErr",dterr,RooConst(1),RooConst(0.25)) ;
00057   RooDataSet* expDataDterr = pdfDtErr.generate(dterr,10000) ;
00058 
00059   // Construct a histogram pdf to describe the shape of the dtErr distribution
00060   RooDataHist* expHistDterr = expDataDterr->binnedClone() ;
00061   RooHistPdf pdfErr("pdfErr","pdfErr",dterr,*expHistDterr) ;
00062 
00063 
00064   // C o n s t r u c t   c o n d i t i o n a l   p r o d u c t   d e c a y _ d m ( d t | d t e r r ) * p d f ( d t e r r )
00065   // ----------------------------------------------------------------------------------------------------------------------
00066 
00067   // Construct production of conditional decay_dm(dt|dterr) with empirical pdfErr(dterr)
00068   RooProdPdf model("model","model",pdfErr,Conditional(decay_gm,dt)) ;
00069 
00070   // (Alternatively you could also use the landau shape pdfDtErr)
00071   //RooProdPdf model("model","model",pdfDtErr,Conditional(decay_gm,dt)) ;
00072 
00073   
00074 
00075   // S a m p l e,   f i t   a n d   p l o t   p r o d u c t   m o d e l 
00076   // ------------------------------------------------------------------
00077 
00078   // Specify external dataset with dterr values to use model_dm as conditional p.d.f.
00079   RooDataSet* data = model.generate(RooArgSet(dt,dterr),10000) ;
00080 
00081   
00082 
00083   // F i t   c o n d i t i o n a l   d e c a y _ d m ( d t | d t e r r )
00084   // ---------------------------------------------------------------------
00085 
00086   // Specify dterr as conditional observable
00087   model.fitTo(*data) ;
00088 
00089 
00090   
00091   // P l o t   c o n d i t i o n a l   d e c a y _ d m ( d t | d t e r r )
00092   // ---------------------------------------------------------------------
00093 
00094 
00095   // Make two-dimensional plot of conditional p.d.f in (dt,dterr)
00096   TH1* hh_model = model.createHistogram("hh_model",dt,Binning(50),YVar(dterr,Binning(50))) ;
00097   hh_model->SetLineColor(kBlue) ;
00098 
00099 
00100   // Make projection of data an dt
00101   RooPlot* frame = dt.frame(Title("Projection of model(dt|dterr) on dt")) ;
00102   data->plotOn(frame) ;
00103   model.plotOn(frame) ;
00104 
00105 
00106   // Draw all frames on canvas
00107   TCanvas* c = new TCanvas("rf307_fullpereventerrors","rf307_fullperventerrors",800, 400);
00108   c->Divide(2) ;
00109   c->cd(1) ; gPad->SetLeftMargin(0.20) ; hh_model->GetZaxis()->SetTitleOffset(2.5) ; hh_model->Draw("surf") ;
00110   c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;
00111 
00112 
00113 
00114 }

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