rf307_fullpereventerrors.cxx

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 "RooDecay.h"
00021 #include "RooLandau.h"
00022 #include "RooProdPdf.h"
00023 #include "RooHistPdf.h"
00024 #include "RooPlot.h"
00025 #include "TCanvas.h"
00026 #include "TH1.h"
00027 using namespace RooFit ;
00028 
00029 
00030 class TestBasic307 : public RooFitTestUnit
00031 {
00032 public: 
00033   TestBasic307(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Full per-event error p.d.f. F(t|dt)G(dt)",refFile,writeRef,verbose) {} ;
00034   Bool_t testCode() {
00035 
00036   // 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
00037   // ----------------------------------------------------------------------------------------------
00038 
00039   // Observables
00040   RooRealVar dt("dt","dt",-10,10) ;
00041   RooRealVar dterr("dterr","per-event error on dt",0.1,10) ;
00042 
00043   // Build a gaussian resolution model scaled by the per-event error = gauss(dt,bias,sigma*dterr)
00044   RooRealVar bias("bias","bias",0,-10,10) ;
00045   RooRealVar sigma("sigma","per-event error scale factor",1,0.1,10) ;
00046   RooGaussModel gm("gm1","gauss model scaled bt per-event error",dt,bias,sigma,dterr) ;
00047 
00048   // Construct decay(dt) (x) gauss1(dt|dterr)
00049   RooRealVar tau("tau","tau",1.548) ;
00050   RooDecay decay_gm("decay_gm","decay",dt,tau,gm,RooDecay::DoubleSided) ;
00051 
00052 
00053 
00054   // 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
00055   // -----------------------------------------------------------------
00056 
00057   // Use landau p.d.f to get empirical distribution with long tail
00058   RooLandau pdfDtErr("pdfDtErr","pdfDtErr",dterr,RooConst(1),RooConst(0.25)) ;
00059   RooDataSet* expDataDterr = pdfDtErr.generate(dterr,10000) ;
00060 
00061   // Construct a histogram pdf to describe the shape of the dtErr distribution
00062   RooDataHist* expHistDterr = expDataDterr->binnedClone() ;
00063   RooHistPdf pdfErr("pdfErr","pdfErr",dterr,*expHistDterr) ;
00064 
00065 
00066   // 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 )
00067   // ----------------------------------------------------------------------------------------------------------------------
00068 
00069   // Construct production of conditional decay_dm(dt|dterr) with empirical pdfErr(dterr)
00070   RooProdPdf model("model","model",pdfErr,Conditional(decay_gm,dt)) ;
00071 
00072   // (Alternatively you could also use the landau shape pdfDtErr)
00073   //RooProdPdf model("model","model",pdfDtErr,Conditional(decay_gm,dt)) ;
00074 
00075   
00076 
00077   // 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 
00078   // ------------------------------------------------------------------
00079 
00080   // Specify external dataset with dterr values to use model_dm as conditional p.d.f.
00081   RooDataSet* data = model.generate(RooArgSet(dt,dterr),10000) ;
00082 
00083   
00084 
00085   // 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 )
00086   // ---------------------------------------------------------------------
00087 
00088   // Specify dterr as conditional observable
00089   model.fitTo(*data) ;
00090 
00091 
00092   
00093   // 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 )
00094   // ---------------------------------------------------------------------
00095 
00096 
00097   // Make projection of data an dt
00098   RooPlot* frame = dt.frame(Title("Projection of model(dt|dterr) on dt")) ;
00099   data->plotOn(frame) ;
00100   model.plotOn(frame) ;
00101 
00102 
00103   regPlot(frame,"rf307_plot1") ;
00104 
00105   delete expDataDterr ;
00106   delete expHistDterr ;
00107   delete data ;
00108   
00109   return kTRUE ;
00110 
00111   }
00112 } ;

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