rf306_condpereventerrors.C

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #306
00004 // 
00005 // Complete example with use of conditional 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 "RooConstVar.h"
00020 #include "RooGaussModel.h"
00021 #include "RooDecay.h"
00022 #include "RooLandau.h"
00023 #include "RooPlot.h"
00024 #include "TCanvas.h"
00025 #include "TAxis.h"
00026 #include "TH2D.h"
00027 using namespace RooFit ;
00028 
00029 
00030 
00031 void rf306_condpereventerrors()
00032 {
00033   // 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
00034   // ----------------------------------------------------------------------------------------------
00035 
00036   // Observables
00037   RooRealVar dt("dt","dt",-10,10) ;
00038   RooRealVar dterr("dterr","per-event error on dt",0.01,10) ;
00039 
00040   // Build a gaussian resolution model scaled by the per-event error = gauss(dt,bias,sigma*dterr)
00041   RooRealVar bias("bias","bias",0,-10,10) ;
00042   RooRealVar sigma("sigma","per-event error scale factor",1,0.1,10) ;
00043   RooGaussModel gm("gm1","gauss model scaled bt per-event error",dt,bias,sigma,dterr) ;
00044 
00045   // Construct decay(dt) (x) gauss1(dt|dterr)
00046   RooRealVar tau("tau","tau",1.548) ;
00047   RooDecay decay_gm("decay_gm","decay",dt,tau,gm,RooDecay::DoubleSided) ;
00048 
00049 
00050 
00051   // C o n s t r u c t   f a k e   ' e x t e r n a l '   d a t a    w i t h   p e r - e v e n t   e r r o r
00052   // ------------------------------------------------------------------------------------------------------
00053 
00054   // Use landau p.d.f to get somewhat realistic distribution with long tail
00055   RooLandau pdfDtErr("pdfDtErr","pdfDtErr",dterr,RooConst(1),RooConst(0.25)) ;
00056   RooDataSet* expDataDterr = pdfDtErr.generate(dterr,10000) ;
00057 
00058 
00059 
00060   // S a m p l e   d a t a   f r o m   c o n d i t i o n a l   d e c a y _ g m ( d t | d t e r r )
00061   // ---------------------------------------------------------------------------------------------
00062 
00063   // Specify external dataset with dterr values to use decay_dm as conditional p.d.f.
00064   RooDataSet* data = decay_gm.generate(dt,ProtoData(*expDataDterr)) ;
00065 
00066   
00067 
00068   // 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 )
00069   // ---------------------------------------------------------------------
00070 
00071   // Specify dterr as conditional observable
00072   decay_gm.fitTo(*data,ConditionalObservables(dterr)) ;
00073 
00074 
00075   
00076   // 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 )
00077   // ---------------------------------------------------------------------
00078 
00079 
00080   // Make two-dimensional plot of conditional p.d.f in (dt,dterr)
00081   TH1* hh_decay = decay_gm.createHistogram("hh_decay",dt,Binning(50),YVar(dterr,Binning(50))) ;
00082   hh_decay->SetLineColor(kBlue) ;
00083 
00084 
00085   // Plot decay_gm(dt|dterr) at various values of dterr
00086   RooPlot* frame = dt.frame(Title("Slices of decay(dt|dterr) at various dterr")) ;
00087   for (Int_t ibin=0 ; ibin<100 ; ibin+=20) {
00088     dterr.setBin(ibin) ;
00089     decay_gm.plotOn(frame,Normalization(5.)) ;
00090   }
00091 
00092 
00093   // Make projection of data an dt
00094   RooPlot* frame2 = dt.frame(Title("Projection of decay(dt|dterr) on dt")) ;
00095   data->plotOn(frame2) ;
00096 
00097   // Make projection of decay(dt|dterr) on dt. 
00098   //
00099   // Instead of integrating out dterr, make a weighted average of curves
00100   // at values dterr_i as given in the external dataset. 
00101   // (The kTRUE argument bins the data before projection to speed up the process)
00102   decay_gm.plotOn(frame2,ProjWData(*expDataDterr,kTRUE)) ;
00103 
00104 
00105 
00106   // Draw all frames on canvas
00107   TCanvas* c = new TCanvas("rf306_condpereventerrors","rf306_condperventerrors",1200, 400);
00108   c->Divide(3) ;
00109   c->cd(1) ; gPad->SetLeftMargin(0.20) ; hh_decay->GetZaxis()->SetTitleOffset(2.5) ; hh_decay->Draw("surf") ;
00110   c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;
00111   c->cd(3) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.6) ; frame2->Draw() ;
00112 
00113 
00114 }
00115 

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