rf203_ranges.C

Go to the documentation of this file.
00001 /////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'ADDITION AND CONVOLUTION' RooFit tutorial macro #203
00004 // 
00005 // Fitting and plotting in sub ranges
00006 //
00007 //
00008 // 07/2008 - Wouter Verkerke 
00009 //
00010 /////////////////////////////////////////////////////////////////////////
00011 
00012 #ifndef __CINT__
00013 #include "RooGlobalFunc.h"
00014 #endif
00015 #include "RooRealVar.h"
00016 #include "RooDataSet.h"
00017 #include "RooGaussian.h"
00018 #include "RooConstVar.h"
00019 #include "RooPolynomial.h"
00020 #include "RooAddPdf.h"
00021 #include "RooFitResult.h"
00022 #include "RooPlot.h"
00023 #include "TCanvas.h"
00024 #include "TAxis.h"
00025 #include "TH1.h"
00026 using namespace RooFit ;
00027 
00028 
00029 void rf203_ranges()
00030 {
00031   // S e t u p   m o d e l 
00032   // ---------------------
00033 
00034   // Construct observables x
00035   RooRealVar x("x","x",-10,10) ;
00036   
00037   // Construct gaussx(x,mx,1)
00038   RooRealVar mx("mx","mx",0,-10,10) ;
00039   RooGaussian gx("gx","gx",x,mx,RooConst(1)) ;
00040 
00041   // Construct px = 1 (flat in x)
00042   RooPolynomial px("px","px",x) ;
00043 
00044   // Construct model = f*gx + (1-f)px
00045   RooRealVar f("f","f",0.,1.) ;
00046   RooAddPdf model("model","model",RooArgList(gx,px),f) ;
00047 
00048   // Generated 10000 events in (x,y) from p.d.f. model
00049   RooDataSet* modelData = model.generate(x,10000) ;
00050 
00051   // F i t   f u l l   r a n g e 
00052   // ---------------------------
00053 
00054   // Fit p.d.f to all data
00055   RooFitResult* r_full = model.fitTo(*modelData,Save(kTRUE)) ;
00056 
00057 
00058   // F i t   p a r t i a l   r a n g e 
00059   // ----------------------------------
00060 
00061   // Define "signal" range in x as [-3,3]
00062   x.setRange("signal",-3,3) ;  
00063 
00064   // Fit p.d.f only to data in "signal" range
00065   RooFitResult* r_sig = model.fitTo(*modelData,Save(kTRUE),Range("signal")) ;
00066 
00067 
00068   // P l o t   /   p r i n t   r e s u l t s 
00069   // ---------------------------------------
00070 
00071   // Make plot frame in x and add data and fitted model
00072   RooPlot* frame = x.frame(Title("Fitting a sub range")) ;
00073   modelData->plotOn(frame) ;
00074   model.plotOn(frame,Range("Full"),LineStyle(kDashed),LineColor(kRed)) ; // Add shape in full ranged dashed
00075   model.plotOn(frame) ; // By default only fitted range is shown
00076   
00077   // Print fit results 
00078   cout << "result of fit on all data " << endl ;
00079   r_full->Print() ;  
00080   cout << "result of fit in in signal region (note increased error on signal fraction)" << endl ;
00081   r_sig->Print() ;
00082 
00083   // Draw frame on canvas
00084   new TCanvas("rf203_ranges","rf203_ranges",600,600) ;
00085   gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
00086 
00087   return ;
00088   
00089 }

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