rf314_paramfitrange.C

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #314
00004 // 
00005 // Working with parameterized ranges in a fit. This an example of a
00006 // fit with an acceptance that changes per-event 
00007 //
00008 //  pdf = exp(-t/tau) with t[tmin,5]
00009 //
00010 //  where t and tmin are both observables in the dataset
00011 //
00012 // 07/2008 - Wouter Verkerke 
00013 // 
00014 /////////////////////////////////////////////////////////////////////////
00015 
00016 #ifndef __CINT__
00017 #include "RooGlobalFunc.h"
00018 #endif
00019 #include "RooRealVar.h"
00020 #include "RooDataSet.h"
00021 #include "RooGaussian.h"
00022 #include "RooConstVar.h"
00023 #include "RooExponential.h"
00024 #include "TCanvas.h"
00025 #include "TAxis.h"
00026 #include "RooPlot.h"
00027 #include "RooFitResult.h"
00028 
00029 using namespace RooFit ;
00030 
00031 
00032 void rf314_paramfitrange()
00033 {
00034 
00035   // 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 
00036   // ---------------------------------------------------------------
00037 
00038   // Declare observables
00039   RooRealVar t("t","t",0,5) ;
00040   RooRealVar tmin("tmin","tmin",0,0,5) ;
00041 
00042   // Make parameterized range in t : [tmin,5]
00043   t.setRange(tmin,RooConst(t.getMax())) ;
00044 
00045   // Make pdf
00046   RooRealVar tau("tau","tau",-1.54,-10,-0.1) ;
00047   RooExponential model("model","model",t,tau) ;
00048 
00049 
00050 
00051   // C r e a t e   i n p u t   d a t a 
00052   // ------------------------------------
00053 
00054   // Generate complete dataset without acceptance cuts (for reference)
00055   RooDataSet* dall = model.generate(t,10000) ;
00056 
00057   // Generate a (fake) prototype dataset for acceptance limit values
00058   RooDataSet* tmp = RooGaussian("gmin","gmin",tmin,RooConst(0),RooConst(0.5)).generate(tmin,5000) ;
00059 
00060   // Generate dataset with t values that observe (t>tmin)
00061   RooDataSet* dacc = model.generate(t,ProtoData(*tmp)) ;
00062 
00063 
00064 
00065   // F i t   p d f   t o   d a t a   i n   a c c e p t a n c e   r e g i o n
00066   // -----------------------------------------------------------------------
00067 
00068   RooFitResult* r = model.fitTo(*dacc,Save()) ;
00069 
00070 
00071  
00072   // P l o t   f i t t e d   p d f   o n   f u l l   a n d   a c c e p t e d   d a t a 
00073   // ---------------------------------------------------------------------------------
00074 
00075   // Make plot frame, add datasets and overlay model
00076   RooPlot* frame = t.frame(Title("Fit to data with per-event acceptance")) ;
00077   dall->plotOn(frame,MarkerColor(kRed),LineColor(kRed)) ;
00078   model.plotOn(frame) ;
00079   dacc->plotOn(frame) ;
00080 
00081   // Print fit results to demonstrate absence of bias
00082   r->Print("v") ;
00083 
00084 
00085   new TCanvas("rf314_paramranges","rf314_paramranges",600,600) ;
00086   gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;
00087 
00088   return ;
00089 }

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