Logo ROOT   6.07/09
Reference Guide
rs401c_FeldmanCousins.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roostats
3 /// \notebook
4 /// Produces an interval on the mean signal in a number counting
5 /// experiment with known background using the Feldman-Cousins technique.
6 ///
7 /// Using the RooStats FeldmanCousins tool with 200 bins
8 /// it takes 1 min and the interval is [0.2625, 10.6125]
9 /// with a step size of 0.075.
10 /// The interval in Feldman & Cousins's original paper is [.29, 10.81]
11 /// Phys.Rev.D57:3873-3889,1998.
12 ////
13 /// \macro_image
14 /// \macro_output
15 /// \macro_code
16 ///
17 /// \author Kyle Cranmer
18 
19 #include "RooGlobalFunc.h"
20 #include "RooStats/ConfInterval.h"
24 #include "RooStats/ModelConfig.h"
25 
26 #include "RooWorkspace.h"
27 #include "RooDataSet.h"
28 #include "RooRealVar.h"
29 #include "RooConstVar.h"
30 #include "RooAddition.h"
31 
32 #include "RooDataHist.h"
33 
34 #include "RooPoisson.h"
35 #include "RooPlot.h"
36 
37 #include "TCanvas.h"
38 #include "TTree.h"
39 #include "TH1F.h"
40 #include "TMarker.h"
41 #include "TStopwatch.h"
42 
43 #include <iostream>
44 
45 // use this order for safety on library loading
46 using namespace RooFit;
47 using namespace RooStats;
48 
49 
50 void rs401c_FeldmanCousins()
51 {
52  // to time the macro... about 30 s
53  TStopwatch t;
54  t.Start();
55 
56  // make a simple model
57  RooRealVar x("x","", 1,0,50);
58  RooRealVar mu("mu","", 2.5,0, 15); // with a limit on mu>=0
59  RooConstVar b("b","", 3.);
60  RooAddition mean("mean","",RooArgList(mu,b));
61  RooPoisson pois("pois", "", x, mean);
62  RooArgSet parameters(mu);
63 
64  // create a toy dataset
65  RooDataSet* data = pois.generate(RooArgSet(x), 1);
66  data->Print("v");
67 
68  TCanvas* dataCanvas = new TCanvas("dataCanvas");
69  RooPlot* frame = x.frame();
70  data->plotOn(frame);
71  frame->Draw();
72  dataCanvas->Update();
73 
74  RooWorkspace* w = new RooWorkspace();
75  ModelConfig modelConfig("poissonProblem",w);
76  modelConfig.SetPdf(pois);
77  modelConfig.SetParametersOfInterest(parameters);
78  modelConfig.SetObservables(RooArgSet(x));
79  w->Print();
80 
81  //////// show use of Feldman-Cousins
82  RooStats::FeldmanCousins fc(*data,modelConfig);
83  fc.SetTestSize(.05); // set size of test
84  fc.UseAdaptiveSampling(true);
85  fc.FluctuateNumDataEntries(false); // number counting analysis: dataset always has 1 entry with N events observed
86  fc.SetNBins(100); // number of points to test per parameter
87 
88  // use the Feldman-Cousins tool
89  PointSetInterval* interval = (PointSetInterval*)fc.GetInterval();
90 
91  // make a canvas for plots
92  TCanvas* intervalCanvas = new TCanvas("intervalCanvas");
93 
94  std::cout << "is this point in the interval? " <<
95  interval->IsInInterval(parameters) << std::endl;
96 
97  std::cout << "interval is ["<<
98  interval->LowerLimit(mu) << ", " <<
99  interval->UpperLimit(mu) << "]" << endl;
100 
101  // using 200 bins it takes 1 min and the answer is
102  // interval is [0.2625, 10.6125] with a step size of .075
103  // The interval in Feldman & Cousins's original paper is [.29, 10.81]
104  // Phys.Rev.D57:3873-3889,1998.
105 
106  // No dedicated plotting class yet, so do it by hand:
107 
108  RooDataHist* parameterScan = (RooDataHist*) fc.GetPointsToScan();
109  TH1F* hist = (TH1F*) parameterScan->createHistogram("mu",30);
110  hist->Draw();
111 
112 
113  RooArgSet* tmpPoint;
114  // loop over points to test
115  for(Int_t i=0; i<parameterScan->numEntries(); ++i){
116  // cout << "on parameter point " << i << " out of " << parameterScan->numEntries() << endl;
117  // get a parameter point from the list of points to test.
118  tmpPoint = (RooArgSet*) parameterScan->get(i)->clone("temp");
119 
120  TMarker* mark = new TMarker(tmpPoint->getRealValue("mu"), 1, 25);
121  if (interval->IsInInterval( *tmpPoint ) )
122  mark->SetMarkerColor(kBlue);
123  else
124  mark->SetMarkerColor(kRed);
125 
126  mark->Draw("s");
127  //delete tmpPoint;
128  // delete mark;
129  }
130  t.Stop();
131  t.Print();
132 
133 
134 }
Poisson pdf.
Definition: RooPoisson.h:19
void Print(Option_t *option="") const
Print the real and cpu time passed between the start and stop events.
Definition: TStopwatch.cxx:219
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:52
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Definition: TStopwatch.cxx:58
Definition: Rtypes.h:61
#define mark(osub)
Definition: triangle.c:1206
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooAbsData.h:157
virtual void Draw(Option_t *option="")
Draw this marker with its current attributes.
Definition: TMarker.cxx:139
tomato 1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:575
Manages Markers.
Definition: TMarker.h:31
int Int_t
Definition: RtypesCore.h:41
virtual const RooArgSet * get() const
Definition: RooDataHist.h:77
static struct mg_connection * fc(struct mg_context *ctx)
Definition: civetweb.c:1956
RooDataSet is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
void Stop()
Stop the stopwatch.
Definition: TStopwatch.cxx:77
Double_t x[n]
Definition: legend1.C:17
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition: TAttMarker.h:43
RooConstVar represent a constant real-valued object.
Definition: RooConstVar.h:25
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
virtual Bool_t IsInInterval(const RooArgSet &) const
check if parameter is in the interval
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
Double_t LowerLimit(RooRealVar &param)
return lower limit on a given parameter
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2853
TH1 * createHistogram(const char *name, const RooAbsRealLValue &xvar, 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
Create and fill a ROOT histogram TH1,TH2 or TH3 with the values of this dataset.
Definition: RooAbsData.cxx:656
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
Double_t UpperLimit(RooRealVar &param)
return upper limit on a given parameter
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
The FeldmanCousins class (like the Feldman-Cousins technique) is essentially a specific configuration...
The Canvas class.
Definition: TCanvas.h:41
Namespace for the RooStats classes.
Definition: Asimov.h:20
virtual TObject * clone(const char *newname) const
Definition: RooArgSet.h:82
PointSetInterval is a concrete implementation of the ConfInterval interface.
virtual Int_t numEntries() const
Return the number of bins.
void Print(Option_t *opts=0) const
Print contents of the workspace.
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 Double_t b
Definition: TRolke.cxx:630
Definition: Rtypes.h:61
RooAddition calculates the sum of a set of RooAbsReal terms, or when constructed with two sets...
Definition: RooAddition.h:26
virtual void Update()
Update canvas pad buffers.
Definition: TCanvas.cxx:2183
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
Double_t getRealValue(const char *name, Double_t defVal=0, Bool_t verbose=kFALSE) const
Get value of a RooAbsReal stored in set with given name.
Definition: RooArgSet.cxx:527
Stopwatch class.
Definition: TStopwatch.h:30