Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rs801_HypoTestInverterOriginal.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roostats
3/// \notebook -js
4/// An example of using the HypoTestInverterOriginal class
5///
6/// \macro_image
7/// \macro_output
8/// \macro_code
9///
10/// \author Gregory Schott
11
12#include "RooRealVar.h"
13#include "RooConstVar.h"
14#include "RooProdPdf.h"
15#include "RooWorkspace.h"
16#include "RooDataSet.h"
17#include "RooPolynomial.h"
18#include "RooAddPdf.h"
19#include "RooExtendPdf.h"
20
25
26#include "TGraphErrors.h"
27
28using namespace RooFit;
29using namespace RooStats;
30
31void rs801_HypoTestInverterOriginal()
32{
33 // prepare the model
34 RooRealVar lumi("lumi", "luminosity", 1);
35 RooRealVar r("r", "cross-section ratio", 3.74, 0, 50);
36 RooFormulaVar ns("ns", "1*r*lumi", RooArgList(lumi, r));
37 RooRealVar nb("nb", "background yield", 1);
38 RooRealVar x("x", "dummy observable", 0, 1);
40 RooPolynomial flatPdf("flatPdf", "flat PDF", x, p0);
41 RooAddPdf totPdf("totPdf", "S+B model", RooArgList(flatPdf, flatPdf), RooArgList(ns, nb));
42 RooExtendPdf bkgPdf("bkgPdf", "B-only model", flatPdf, nb);
43 RooDataSet *data = totPdf.generate(x, 1);
44
45 // prepare the calculator
46 HybridCalculatorOriginal myhc(*data, totPdf, bkgPdf, 0, 0);
47 myhc.SetTestStatistic(2);
48 myhc.SetNumberOfToys(1000);
49 myhc.UseNuisance(false);
50
51 // run the hypothesis-test inversion
52 HypoTestInverterOriginal myInverter(myhc, r);
53 myInverter.SetTestSize(0.10);
54 myInverter.UseCLs(true);
55 // myInverter.RunFixedScan(5,1,6);
56 // scan for a 95% UL
57 myInverter.RunAutoScan(3., 5, myInverter.Size() / 2, 0.005);
58 // run an alternative autoscan algorithm
59 // myInverter.RunAutoScan(1,6,myInverter.Size()/2,0.005,1);
60 // myInverter.RunOnePoint(3.9);
61
62 HypoTestInverterResult *results = myInverter.GetInterval();
63
64 HypoTestInverterPlot myInverterPlot("myInverterPlot", "", results);
65 TGraphErrors *gr1 = myInverterPlot.MakePlot();
66 gr1->Draw("ALP");
67
68 double ulError = results->UpperLimitEstimatedError();
69
70 double upperLimit = results->UpperLimit();
71 std::cout << "The computed upper limit is: " << upperLimit << std::endl;
72 std::cout << "an estimated error on this upper limit is: " << ulError << std::endl;
73 // expected result: 4.10
74}
75int main()
76{
77 rs801_HypoTestInverterOriginal();
78}
ROOT::R::TRInterface & r
Definition Object.C:4
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:32
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:21
RooConstVar represent a constant real-valued object.
Definition RooConstVar.h:26
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:33
RooExtendPdf is a wrapper around an existing PDF that adds a parameteric extended likelihood term to ...
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
RooPolynomial implements a polynomial p.d.f of the form.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
HybridCalculatorOriginal class.
This class is now deprecated and to be replaced by the HypoTestInverter.
Class to plot a HypoTestInverterResult, the output of the HypoTestInverter calculator.
HypoTestInverterResult class holds the array of hypothesis test results and compute a confidence inte...
Double_t UpperLimitEstimatedError()
Estimate of lower limit error function evaluates only a rough error on the lower limit.
A TGraphErrors is a TGraph with error bars.
virtual void Draw(Option_t *chopt="")
Draw this graph with its current attributes.
Definition TGraph.cxx:769
int main()
RooConstVar & RooConst(Double_t val)
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Namespace for the RooStats classes.
Definition Asimov.h:19