Logo ROOT   6.10/09
Reference Guide
rf314_paramfitrange.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////
2 //
3 // 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #314
4 //
5 // Working with parameterized ranges in a fit. This an example of a
6 // fit with an acceptance that changes per-event
7 //
8 // pdf = exp(-t/tau) with t[tmin,5]
9 //
10 // where t and tmin are both observables in the dataset
11 //
12 // 07/2008 - Wouter Verkerke
13 //
14 /////////////////////////////////////////////////////////////////////////
15 
16 #ifndef __CINT__
17 #include "RooGlobalFunc.h"
18 #endif
19 #include "RooRealVar.h"
20 #include "RooDataSet.h"
21 #include "RooGaussian.h"
22 #include "RooExponential.h"
23 #include "TCanvas.h"
24 #include "RooPlot.h"
25 #include "RooFitResult.h"
26 
27 using namespace RooFit ;
28 
29 
30 class TestBasic314 : public RooFitTestUnit
31 {
32 public:
33  TestBasic314(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Fit with non-rectangular observable boundaries",refFile,writeRef,verbose) {} ;
34  Bool_t testCode() {
35 
36  // D e f i n e o b s e r v a b l e s a n d d e c a y p d f
37  // ---------------------------------------------------------------
38 
39  // Declare observables
40  RooRealVar t("t","t",0,5) ;
41  RooRealVar tmin("tmin","tmin",0,0,5) ;
42 
43  // Make parameterized range in t : [tmin,5]
44  t.setRange(tmin,RooConst(t.getMax())) ;
45 
46  // Make pdf
47  RooRealVar tau("tau","tau",-1.54,-10,-0.1) ;
48  RooExponential model("model","model",t,tau) ;
49 
50 
51 
52  // C r e a t e i n p u t d a t a
53  // ------------------------------------
54 
55  // Generate complete dataset without acceptance cuts (for reference)
56  RooDataSet* dall = model.generate(t,10000) ;
57 
58  // Generate a (fake) prototype dataset for acceptance limit values
59  RooDataSet* tmp = RooGaussian("gmin","gmin",tmin,RooConst(0),RooConst(0.5)).generate(tmin,5000) ;
60 
61  // Generate dataset with t values that observe (t>tmin)
62  RooDataSet* dacc = model.generate(t,ProtoData(*tmp)) ;
63 
64 
65 
66  // F i t p d f t o d a t a i n a c c e p t a n c e r e g i o n
67  // -----------------------------------------------------------------------
68 
69  RooFitResult* r = model.fitTo(*dacc,Save()) ;
70 
71 
72 
73  // P l o t f i t t e d p d f o n f u l l a n d a c c e p t e d d a t a
74  // ---------------------------------------------------------------------------------
75 
76  // Make plot frame, add datasets and overlay model
77  RooPlot* frame = t.frame(Title("Fit to data with per-event acceptance")) ;
78  dall->plotOn(frame,MarkerColor(kRed),LineColor(kRed)) ;
79  model.plotOn(frame) ;
80  dacc->plotOn(frame,Name("dacc")) ;
81 
82  // Print fit results to demonstrate absence of bias
83  regResult(r,"rf314_fit") ;
84  regPlot(frame,"rf314_plot1") ;
85 
86  delete tmp ;
87  delete dacc ;
88  delete dall ;
89 
90  return kTRUE;
91  }
92 } ;
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), 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
Plot dataset on specified frame.
Definition: RooAbsData.cxx:552
RooCmdArg LineColor(Color_t color)
Definition: Rtypes.h:56
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
RooCmdArg Title(const char *name)
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t tau
Definition: TRolke.cxx:630
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
Exponential p.d.f.
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
RooCmdArg Name(const char *name)
TRandom2 r(17)
bool verbose
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
RooCmdArg ProtoData(const RooDataSet &protoData, Bool_t randomizeOrder=kFALSE, Bool_t resample=kFALSE)
RooCmdArg Save(Bool_t flag=kTRUE)
RooConstVar & RooConst(Double_t val)
RooDataSet * generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none())
Generate a new dataset containing the specified variables with events sampled from our distribution...
Definition: RooAbsPdf.cxx:1702
RooCmdArg MarkerColor(Color_t color)
const Bool_t kTRUE
Definition: RtypesCore.h:91