ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rs603_HLFactoryElaborateExample.C
Go to the documentation of this file.
1 /////////////////////////////////////////////////////////////////////////
2 //
3 // 'High Level Factory Example' RooStats tutorial macro #602
4 // author: Danilo Piparo
5 // date August. 2009
6 //
7 // This tutorial shows an example of creating a combined
8 // model using the High Level model Factory.
9 //
10 //
11 /////////////////////////////////////////////////////////////////////////
12 
13 
14 
15 #include <fstream>
16 #include "TString.h"
17 #include "TFile.h"
18 #include "TROOT.h"
19 #include "RooGlobalFunc.h"
20 #include "RooWorkspace.h"
21 #include "RooRealVar.h"
22 #include "RooAbsPdf.h"
23 #include "RooDataSet.h"
24 #include "RooPlot.h"
25 #include "RooStats/HLFactory.h"
26 
27 
28 // use this order for safety on library loading
29 using namespace RooFit ;
30 using namespace RooStats ;
31 using namespace std;
32 
34 
35  // --- Prepare the 2 needed datacards for this example ---
36 
37  TString card_name("rs603_card_WsMaker.rs");
38  ofstream ofile(card_name);
39  ofile << "// The simplest card for combination\n\n";
40  ofile << "gauss1 = Gaussian(x[0,100],mean1[50,0,100],4);\n";
41  ofile << "flat1 = Polynomial(x,0);\n";
42  ofile << "sb_model1 = SUM(nsig1[120,0,300]*gauss1 , nbkg1[100,0,1000]*flat1);\n\n";
43  ofile << "echo In the middle!;\n\n";
44  ofile << "gauss2 = Gaussian(x,mean2[80,0,100],5);\n";
45  ofile << "flat2 = Polynomial(x,0);\n";
46  ofile << "sb_model2 = SUM(nsig2[90,0,400]*gauss2 , nbkg2[80,0,1000]*flat2);\n\n";
47  ofile << "echo At the end!;\n";
48  ofile.close();
49 
50  TString card_name2("rs603_card.rs");
51  ofstream ofile2(card_name2);
52  ofile2 << "// The simplest card for combination\n\n";
53  ofile2 << "gauss1 = Gaussian(x[0,100],mean1[50,0,100],4);\n";
54  ofile2 << "flat1 = Polynomial(x,0);\n";
55  ofile2 << "sb_model1 = SUM(nsig1[120,0,300]*gauss1 , nbkg1[100,0,1000]*flat1);\n\n";
56  ofile2 << "echo In the middle!;\n\n";
57  ofile2 << "gauss2 = Gaussian(x,mean2[80,0,100],5);\n";
58  ofile2 << "flat2 = Polynomial(x,0);\n";
59  ofile2 << "sb_model2 = SUM(nsig2[90,0,400]*gauss2 , nbkg2[80,0,1000]*flat2);\n\n";
60  ofile2 << "#include rs603_included_card.rs;\n\n";
61  ofile2 << "echo At the end!;\n";
62  ofile2.close();
63 
64  TString card_name3("rs603_included_card.rs");
65  ofstream ofile3(card_name3);
66  ofile3 << "echo Now reading the included file!;\n\n";
67  ofile3 << "echo Including datasets in a Workspace in a Root file...;\n";
68  ofile3 << "data1 = import(rs603_infile.root,\n";
69  ofile3 << " rs603_ws,\n";
70  ofile3 << " data1);\n\n";
71  ofile3 << "data2 = import(rs603_infile.root,\n";
72  ofile3 << " rs603_ws,\n";
73  ofile3 << " data2);\n";
74  ofile3.close();
75 
76 // --- Produce the two separate datasets into a WorkSpace ---
77 
78 HLFactory hlf("HLFactoryComplexExample",
79  "rs603_card_WsMaker.rs",
80  false);
81 
82 RooRealVar* x = static_cast<RooRealVar*>(hlf.GetWs()->arg("x"));
83 RooAbsPdf* pdf1 = hlf.GetWs()->pdf("sb_model1");
84 RooAbsPdf* pdf2 = hlf.GetWs()->pdf("sb_model2");
85 
86 RooWorkspace w("rs603_ws");
87 
88 RooDataSet* data1 = pdf1->generate(RooArgSet(*x),Extended());
89 data1->SetName("data1");
90 w.import(*data1);
91 
92 RooDataSet* data2 = pdf2->generate(RooArgSet(*x),Extended());
93 data2->SetName("data2");
94 w.import(*data2);
95 
96 // --- Write the WorkSpace into a rootfile ---
97 
98 TFile outfile("rs603_infile.root","RECREATE");
99 w.Write();
100 outfile.Close();
101 
102 cout << "-------------------------------------------------------------------\n"
103  << " Rootfile and Workspace prepared \n"
104  << "-------------------------------------------------------------------\n";
105 
106 
107 HLFactory hlf_2("HLFactoryElaborateExample",
108  "rs603_card.rs",
109  false);
110 
111 x = hlf_2.GetWs()->var("x");
112 pdf1 = hlf_2.GetWs()->pdf("sb_model1");
113 pdf2 = hlf_2.GetWs()->pdf("sb_model2");
114 
115 hlf_2.AddChannel("model1","sb_model1","flat1","data1");
116 hlf_2.AddChannel("model2","sb_model2","flat2","data2");
117 
118 RooDataSet* data = hlf_2.GetTotDataSet();
119 RooAbsPdf* pdf = hlf_2.GetTotSigBkgPdf();
120 RooCategory* thecat = hlf_2.GetTotCategory();
121 
122 // --- Perform extended ML fit of composite PDF to toy data ---
123 pdf->fitTo(*data) ;
124 
125 // --- Plot toy data and composite PDF overlaid ---
126 RooPlot* xframe = x->frame() ;
127 
128 data->plotOn(xframe);
129 thecat->setIndex(0);
130 pdf->plotOn(xframe,Slice(*thecat),ProjWData(*thecat,*data)) ;
131 
132 thecat->setIndex(1);
133 pdf->plotOn(xframe,Slice(*thecat),ProjWData(*thecat,*data)) ;
134 
135 gROOT->SetStyle("Plain");
136 
137 xframe->Draw();
138 
139 }
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition: TObject.cxx:823
void SetName(const char *name)
Change the name of this dataset into the given name.
int AddChannel(const char *label, const char *SigBkgPdfName, const char *BkgPdfName=0, const char *datasetName=0)
Add channel for the combination.
Definition: HLFactory.cxx:131
virtual Bool_t setIndex(Int_t index, Bool_t printError=kTRUE)
Set value by specifying the index code of the desired state.
#define gROOT
Definition: TROOT.h:344
RooDataSet * GetTotDataSet()
Get the combined dataset.
Definition: HLFactory.cxx:279
Basic string class.
Definition: TString.h:137
RooCmdArg ProjWData(const RooAbsData &projData, Bool_t binData=kFALSE)
RooCmdArg Extended(Bool_t flag=kTRUE)
Double_t x[n]
Definition: legend1.C:17
HLFactory is an High Level model Factory allows you to describe your models in a configuration file (...
Definition: HLFactory.h:41
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
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:626
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
void rs603_HLFactoryElaborateExample()
tuple outfile
Definition: mrt.py:21
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 RooCmdArg &arg9=RooCmdArg::none(), const RooCmdArg &arg10=RooCmdArg::none()) const
Plot (project) PDF on specified frame.
Definition: RooAbsPdf.h:105
tuple w
Definition: qtexample.py:51
RooPlot * frame(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
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
RooCategory represents a fundamental (non-derived) discrete value object.
Definition: RooCategory.h:25
RooRealVar * var(const char *name) const
Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found...
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
RooCategory * GetTotCategory()
Get the combined dataset.
Definition: HLFactory.cxx:332
RooWorkspace * GetWs()
Get the RooWorkspace containing the models and variables.
Definition: HLFactory.h:84
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
Bool_t import(const RooAbsArg &arg, 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())
Import a RooAbsArg object, e.g.
RooAbsArg * arg(const char *name) const
Return RooAbsArg with given name. A null pointer is returned if none is found.
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
virtual RooFitResult * fitTo(RooAbsData &data, 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())
Fit PDF to given dataset.
Definition: RooAbsPdf.cxx:1056
RooAbsPdf * GetTotSigBkgPdf()
Get the combined signal plus background pdf.
Definition: HLFactory.cxx:180
RooCmdArg Slice(const RooArgSet &sliceSet)
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:42
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:559