rf313_paramranges.C

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #313
00004 // 
00005 // Working with parameterized ranges to define non-rectangular regions
00006 // for fitting and integration
00007 //
00008 //
00009 //
00010 // 07/2008 - Wouter Verkerke 
00011 // 
00012 /////////////////////////////////////////////////////////////////////////
00013 
00014 #ifndef __CINT__
00015 #include "RooGlobalFunc.h"
00016 #endif
00017 #include "RooRealVar.h"
00018 #include "RooDataSet.h"
00019 #include "RooGaussian.h"
00020 #include "RooConstVar.h"
00021 #include "RooPolynomial.h"
00022 #include "RooProdPdf.h"
00023 #include "TCanvas.h"
00024 #include "TAxis.h"
00025 #include "RooPlot.h"
00026 using namespace RooFit ;
00027 
00028 
00029 void rf313_paramranges()
00030 {
00031 
00032   // C r e a t e   3 D   p d f 
00033   // -------------------------
00034 
00035   // Define observable (x,y,z)
00036   RooRealVar x("x","x",0,10) ;
00037   RooRealVar y("y","y",0,10) ;
00038   RooRealVar z("z","z",0,10) ;
00039 
00040   // Define 3 dimensional pdf
00041   RooRealVar z0("z0","z0",-0.1,1) ;
00042   RooPolynomial px("px","px",x,RooConst(0)) ;
00043   RooPolynomial py("py","py",y,RooConst(0)) ;
00044   RooPolynomial pz("pz","pz",z,z0) ;
00045   RooProdPdf pxyz("pxyz","pxyz",RooArgSet(px,py,pz)) ;
00046 
00047 
00048 
00049   // D e f i n e d   n o n - r e c t a n g u l a r   r e g i o n   R   i n   ( x , y , z ) 
00050   // -------------------------------------------------------------------------------------
00051 
00052   //
00053   // R = Z[0 - 0.1*Y^2] * Y[0.1*X - 0.9*X] * X[0 - 10]
00054   //
00055 
00056   // Construct range parameterized in "R" in y [ 0.1*x, 0.9*x ]
00057   RooFormulaVar ylo("ylo","0.1*x",x) ;
00058   RooFormulaVar yhi("yhi","0.9*x",x) ;
00059   y.setRange("R",ylo,yhi) ;
00060 
00061   // Construct parameterized ranged "R" in z [ 0, 0.1*y^2 ]
00062   RooFormulaVar zlo("zlo","0.0*y",y) ;
00063   RooFormulaVar zhi("zhi","0.1*y*y",y) ;
00064   z.setRange("R",zlo,zhi) ;
00065 
00066 
00067 
00068   // C a l c u l a t e   i n t e g r a l   o f   n o r m a l i z e d   p d f   i n   R 
00069   // ----------------------------------------------------------------------------------
00070 
00071   // Create integral over normalized pdf model over x,y,z in "R" region
00072   RooAbsReal* intPdf = pxyz.createIntegral(RooArgSet(x,y,z),RooArgSet(x,y,z),"R") ;
00073 
00074   // Plot value of integral as function of pdf parameter z0
00075   RooPlot* frame = z0.frame(Title("Integral of pxyz over x,y,z in region R")) ;
00076   intPdf->plotOn(frame) ;
00077 
00078 
00079 
00080   new TCanvas("rf313_paramranges","rf313_paramranges",600,600) ;
00081   gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;
00082 
00083   return ;
00084 }

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