ROOT  6.06/09
Reference Guide
rf313_paramranges.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////
2 //
3 // 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #313
4 //
5 // Working with parameterized ranges to define non-rectangular regions
6 // for fitting and integration
7 //
8 //
9 //
10 // 07/2008 - Wouter Verkerke
11 //
12 /////////////////////////////////////////////////////////////////////////
13 
14 #ifndef __CINT__
15 #include "RooGlobalFunc.h"
16 #endif
17 #include "RooRealVar.h"
18 #include "RooDataSet.h"
19 #include "RooGaussian.h"
20 #include "RooPolynomial.h"
21 #include "RooProdPdf.h"
22 #include "TCanvas.h"
23 #include "RooPlot.h"
24 using namespace RooFit ;
25 
26 
27 class TestBasic313 : public RooFitTestUnit
28 {
29 public:
30  TestBasic313(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Integration over non-rectangular regions",refFile,writeRef,verbose) {} ;
31  Bool_t testCode() {
32 
33  // C r e a t e 3 D p d f
34  // -------------------------
35 
36  // Define observable (x,y,z)
37  RooRealVar x("x","x",0,10) ;
38  RooRealVar y("y","y",0,10) ;
39  RooRealVar z("z","z",0,10) ;
40 
41  // Define 3 dimensional pdf
42  RooRealVar z0("z0","z0",-0.1,1) ;
43  RooPolynomial px("px","px",x,RooConst(0)) ;
44  RooPolynomial py("py","py",y,RooConst(0)) ;
45  RooPolynomial pz("pz","pz",z,z0) ;
46  RooProdPdf pxyz("pxyz","pxyz",RooArgSet(px,py,pz)) ;
47 
48 
49 
50  // 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 )
51  // -------------------------------------------------------------------------------------
52 
53  //
54  // R = Z[0 - 0.1*Y^2] * Y[0.1*X - 0.9*X] * X[0 - 10]
55  //
56 
57  // Construct range parameterized in "R" in y [ 0.1*x, 0.9*x ]
58  RooFormulaVar ylo("ylo","0.1*x",x) ;
59  RooFormulaVar yhi("yhi","0.9*x",x) ;
60  y.setRange("R",ylo,yhi) ;
61 
62  // Construct parameterized ranged "R" in z [ 0, 0.1*y^2 ]
63  RooFormulaVar zlo("zlo","0.0*y",y) ;
64  RooFormulaVar zhi("zhi","0.1*y*y",y) ;
65  z.setRange("R",zlo,zhi) ;
66 
67 
68 
69  // 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
70  // ----------------------------------------------------------------------------------
71 
72  // Create integral over normalized pdf model over x,y,z in "R" region
73  RooAbsReal* intPdf = pxyz.createIntegral(RooArgSet(x,y,z),RooArgSet(x,y,z),"R") ;
74 
75  // Plot value of integral as function of pdf parameter z0
76  RooPlot* frame = z0.frame(Title("Integral of pxyz over x,y,z in region R")) ;
77  intPdf->plotOn(frame) ;
78 
79 
80  regPlot(frame,"rf313_plot1") ;
81 
82  delete intPdf ;
83 
84  return kTRUE;
85  }
86 } ;
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg(), const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg(), const RooCmdArg &arg9=RooCmdArg(), const RooCmdArg &arg10=RooCmdArg()) const
Plot (project) PDF on specified frame.
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:45
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
RooCmdArg Title(const char *name)
Double_t x[n]
Definition: legend1.C:17
bool verbose
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
Double_t y[n]
Definition: legend1.C:17
RooConstVar & RooConst(Double_t val)
const Bool_t kTRUE
Definition: Rtypes.h:91
RooAbsReal * createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables listed in ...
Definition: RooAbsReal.cxx:503