rf312_multirangefit.C

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #312
00004 // 
00005 // Performing fits in multiple (disjoint) ranges in one or more dimensions
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 "RooProdPdf.h"
00021 #include "RooAddPdf.h"
00022 #include "RooPolynomial.h"
00023 #include "TCanvas.h"
00024 #include "TAxis.h"
00025 #include "RooPlot.h"
00026 #include "RooFitResult.h"
00027 using namespace RooFit ;
00028 
00029 
00030 void rf312_multirangefit()
00031 {
00032 
00033   // C r e a t e   2 D   p d f   a n d   d a t a 
00034   // -------------------------------------------
00035 
00036   // Define observables x,y
00037   RooRealVar x("x","x",-10,10) ;
00038   RooRealVar y("y","y",-10,10) ;
00039 
00040   // Construct the signal pdf gauss(x)*gauss(y)
00041   RooRealVar mx("mx","mx",1,-10,10) ;
00042   RooRealVar my("my","my",1,-10,10) ;
00043 
00044   RooGaussian gx("gx","gx",x,mx,RooConst(1)) ;
00045   RooGaussian gy("gy","gy",y,my,RooConst(1)) ;
00046 
00047   RooProdPdf sig("sig","sig",gx,gy) ;
00048 
00049   // Construct the background pdf (flat in x,y)
00050   RooPolynomial px("px","px",x) ;
00051   RooPolynomial py("py","py",y) ;
00052   RooProdPdf bkg("bkg","bkg",px,py) ;
00053 
00054   // Construct the composite model sig+bkg
00055   RooRealVar f("f","f",0.,1.) ;
00056   RooAddPdf model("model","model",RooArgList(sig,bkg),f) ;
00057 
00058   // Sample 10000 events in (x,y) from the model
00059   RooDataSet* modelData = model.generate(RooArgSet(x,y),10000) ;
00060 
00061 
00062 
00063   // D e f i n e   s i g n a l   a n d   s i d e b a n d   r e g i o n s
00064   // -------------------------------------------------------------------
00065 
00066   // Construct the SideBand1,SideBand2,Signal regions
00067   //
00068   //                    |
00069   //      +-------------+-----------+                 
00070   //      |             |           |             
00071   //      |    Side     |   Sig     |        
00072   //      |    Band1    |   nal     |             
00073   //      |             |           |            
00074   //    --+-------------+-----------+--   
00075   //      |                         |       
00076   //      |           Side          |        
00077   //      |           Band2         |            
00078   //      |                         |          
00079   //      +-------------+-----------+            
00080   //                    |                       
00081 
00082   x.setRange("SB1",-10,+10) ;
00083   y.setRange("SB1",-10,0) ;
00084 
00085   x.setRange("SB2",-10,0) ;
00086   y.setRange("SB2",0,+10) ;
00087 
00088   x.setRange("SIG",0,+10) ;
00089   y.setRange("SIG",0,+10) ;
00090 
00091   x.setRange("FULL",-10,+10) ;
00092   y.setRange("FULL",-10,+10) ;
00093 
00094 
00095   // P e r f o r m   f i t s   i n   i n d i v i d u a l   s i d e b a n d   r e g i o n s
00096   // -------------------------------------------------------------------------------------
00097 
00098   // Perform fit in SideBand1 region (RooAddPdf coefficients will be interpreted in full range)
00099   RooFitResult* r_sb1 = model.fitTo(*modelData,Range("SB1"),Save()) ;
00100 
00101   // Perform fit in SideBand2 region (RooAddPdf coefficients will be interpreted in full range)
00102   RooFitResult* r_sb2 = model.fitTo(*modelData,Range("SB2"),Save()) ;
00103 
00104 
00105 
00106   // P e r f o r m   f i t s   i n   j o i n t    s i d e b a n d   r e g i o n s
00107   // -----------------------------------------------------------------------------
00108 
00109   // Now perform fit to joint 'L-shaped' sideband region 'SB1|SB2' 
00110   // (RooAddPdf coefficients will be interpreted in full range)
00111   RooFitResult* r_sb12 = model.fitTo(*modelData,Range("SB1,SB2"),Save()) ;
00112 
00113 
00114   // Print results for comparison
00115   r_sb1->Print() ;
00116   r_sb2->Print() ;
00117   r_sb12->Print() ;
00118   
00119 
00120 }

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