Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rs602_HLFactoryCombinationexample.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roostats
3 /// \notebook
4 /// High Level Factory: creation of a combined model
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 ///
10 /// \author Danilo Piparo
11 
12 #include <fstream>
13 #include "TString.h"
14 #include "TROOT.h"
15 #include "RooGlobalFunc.h"
16 #include "RooWorkspace.h"
17 #include "RooRealVar.h"
18 #include "RooAbsPdf.h"
19 #include "RooDataSet.h"
20 #include "RooPlot.h"
21 #include "RooStats/HLFactory.h"
22 
23 // use this order for safety on library loading
24 using namespace RooFit;
25 using namespace RooStats;
26 using namespace std;
27 
28 void rs602_HLFactoryCombinationexample()
29 {
30 
31  using namespace RooStats;
32  using namespace RooFit;
33 
34  // create a card
35  TString card_name("HLFavtoryCombinationexample.rs");
36  ofstream ofile(card_name);
37  ofile << "// The simplest card for combination\n\n"
38  << "gauss1 = Gaussian(x[0,100],mean1[50,0,100],4);\n"
39  << "flat1 = Polynomial(x,0);\n"
40  << "sb_model1 = SUM(nsig1[120,0,300]*gauss1 , nbkg1[100,0,1000]*flat1);\n"
41  << "gauss2 = Gaussian(x,mean2[80,0,100],5);\n"
42  << "flat2 = Polynomial(x,0);\n"
43  << "sb_model2 = SUM(nsig2[90,0,400]*gauss2 , nbkg2[80,0,1000]*flat2);\n";
44 
45  ofile.close();
46 
47  HLFactory hlf("HLFavtoryCombinationexample", card_name, false);
48 
49  hlf.AddChannel("model1", "sb_model1", "flat1");
50  hlf.AddChannel("model2", "sb_model2", "flat2");
51  auto pdf = hlf.GetTotSigBkgPdf();
52  auto thecat = hlf.GetTotCategory();
53  auto x = static_cast<RooRealVar *>(hlf.GetWs()->arg("x"));
54 
55  auto data = pdf->generate(RooArgSet(*x, *thecat), Extended());
56 
57  // --- Perform extended ML fit of composite PDF to toy data ---
58  pdf->fitTo(*data);
59 
60  // --- Plot toy data and composite PDF overlaid ---
61  auto xframe = x->frame();
62 
63  data->plotOn(xframe);
64  thecat->setIndex(0);
65  pdf->plotOn(xframe, Slice(*thecat), ProjWData(*thecat, *data));
66 
67  thecat->setIndex(1);
68  pdf->plotOn(xframe, Slice(*thecat), ProjWData(*thecat, *data));
69 
70  gROOT->SetStyle("Plain");
71  xframe->Draw();
72 }