rf312_multirangefit.cxx

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 "RooProdPdf.h"
00020 #include "RooAddPdf.h"
00021 #include "RooPolynomial.h"
00022 #include "TCanvas.h"
00023 #include "RooPlot.h"
00024 #include "RooFitResult.h"
00025 using namespace RooFit ;
00026 
00027 
00028 class TestBasic312 : public RooFitTestUnit
00029 {
00030 public: 
00031   TestBasic312(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Fit in multiple rectangular ranges",refFile,writeRef,verbose) {} ;
00032   Bool_t testCode() {
00033 
00034   // C r e a t e   2 D   p d f   a n d   d a t a 
00035   // -------------------------------------------
00036 
00037   // Define observables x,y
00038   RooRealVar x("x","x",-10,10) ;
00039   RooRealVar y("y","y",-10,10) ;
00040 
00041   // Construct the signal pdf gauss(x)*gauss(y)
00042   RooRealVar mx("mx","mx",1,-10,10) ;
00043   RooRealVar my("my","my",1,-10,10) ;
00044 
00045   RooGaussian gx("gx","gx",x,mx,RooConst(1)) ;
00046   RooGaussian gy("gy","gy",y,my,RooConst(1)) ;
00047 
00048   RooProdPdf sig("sig","sig",gx,gy) ;
00049 
00050   // Construct the background pdf (flat in x,y)
00051   RooPolynomial px("px","px",x) ;
00052   RooPolynomial py("py","py",y) ;
00053   RooProdPdf bkg("bkg","bkg",px,py) ;
00054 
00055   // Construct the composite model sig+bkg
00056   RooRealVar f("f","f",0.,1.) ;
00057   RooAddPdf model("model","model",RooArgList(sig,bkg),f) ;
00058 
00059   // Sample 10000 events in (x,y) from the model
00060   RooDataSet* modelData = model.generate(RooArgSet(x,y),10000) ;
00061 
00062 
00063 
00064   // 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
00065   // -------------------------------------------------------------------
00066 
00067   // Construct the SideBand1,SideBand2,Signal regions
00068   //
00069   //                    |
00070   //      +-------------+-----------+                 
00071   //      |             |           |             
00072   //      |    Side     |   Sig     |        
00073   //      |    Band1    |   nal     |             
00074   //      |             |           |            
00075   //    --+-------------+-----------+--   
00076   //      |                         |       
00077   //      |           Side          |        
00078   //      |           Band2         |            
00079   //      |                         |          
00080   //      +-------------+-----------+            
00081   //                    |                       
00082 
00083   x.setRange("SB1",-10,+10) ;
00084   y.setRange("SB1",-10,0) ;
00085 
00086   x.setRange("SB2",-10,0) ;
00087   y.setRange("SB2",0,+10) ;
00088 
00089   x.setRange("SIG",0,+10) ;
00090   y.setRange("SIG",0,+10) ;
00091 
00092   x.setRange("FULL",-10,+10) ;
00093   y.setRange("FULL",-10,+10) ;
00094 
00095 
00096   // 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
00097   // -------------------------------------------------------------------------------------
00098 
00099   // Perform fit in SideBand1 region (RooAddPdf coefficients will be interpreted in full range)
00100   RooFitResult* r_sb1 = model.fitTo(*modelData,Range("SB1"),Save()) ;
00101 
00102   // Perform fit in SideBand2 region (RooAddPdf coefficients will be interpreted in full range)
00103   RooFitResult* r_sb2 = model.fitTo(*modelData,Range("SB2"),Save()) ;
00104 
00105 
00106 
00107   // 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
00108   // -----------------------------------------------------------------------------
00109 
00110   // Now perform fit to joint 'L-shaped' sideband region 'SB1|SB2' 
00111   // (RooAddPdf coefficients will be interpreted in full range)
00112   RooFitResult* r_sb12 = model.fitTo(*modelData,Range("SB1,SB2"),Save()) ;
00113 
00114 
00115   regResult(r_sb1,"rf312_fit_sb1") ;
00116   regResult(r_sb2,"rf312_fit_sb2") ;
00117   regResult(r_sb12,"rf312_fit_sb12") ;
00118 
00119   delete modelData ;
00120 
00121   return kTRUE ;
00122 
00123   }
00124 } ;

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