ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rs401c_FeldmanCousins.C
Go to the documentation of this file.
1 /////////////////////////////////////////////////////////////////////////
2 //
3 // Produces an interval on the mean signal in a number counting
4 // experiment with known background using the Feldman-Cousins technique.
5 // date Jan. 2009
6 // updated June 2010
7 //
8 // Using the RooStats FeldmanCousins tool with 200 bins
9 // it takes 1 min and the interval is [0.2625, 10.6125]
10 // with a step size of 0.075.
11 // The interval in Feldman & Cousins's original paper is [.29, 10.81]
12 // Phys.Rev.D57:3873-3889,1998.
13 /////////////////////////////////////////////////////////////////////////
14 
15 #include "RooGlobalFunc.h"
16 #include "RooStats/ConfInterval.h"
20 #include "RooStats/ModelConfig.h"
21 
22 #include "RooWorkspace.h"
23 #include "RooDataSet.h"
24 #include "RooRealVar.h"
25 #include "RooConstVar.h"
26 #include "RooAddition.h"
27 
28 #include "RooDataHist.h"
29 
30 #include "RooPoisson.h"
31 #include "RooPlot.h"
32 
33 #include "TCanvas.h"
34 #include "TTree.h"
35 #include "TH1F.h"
36 #include "TMarker.h"
37 #include "TStopwatch.h"
38 
39 #include <iostream>
40 
41 // use this order for safety on library loading
42 using namespace RooFit ;
43 using namespace RooStats ;
44 
45 
47 {
48  // to time the macro... about 30 s
49  TStopwatch t;
50  t.Start();
51 
52  // make a simple model
53  RooRealVar x("x","", 1,0,50);
54  RooRealVar mu("mu","", 2.5,0, 15); // with a limit on mu>=0
55  RooConstVar b("b","", 3.);
56  RooAddition mean("mean","",RooArgList(mu,b));
57  RooPoisson pois("pois", "", x, mean);
58  RooArgSet parameters(mu);
59 
60  // create a toy dataset
61  RooDataSet* data = pois.generate(RooArgSet(x), 1);
62  data->Print("v");
63 
64  TCanvas* dataCanvas = new TCanvas("dataCanvas");
65  RooPlot* frame = x.frame();
66  data->plotOn(frame);
67  frame->Draw();
68  dataCanvas->Update();
69 
70  RooWorkspace* w = new RooWorkspace();
71  ModelConfig modelConfig("poissonProblem",w);
72  modelConfig.SetPdf(pois);
73  modelConfig.SetParametersOfInterest(parameters);
74  modelConfig.SetObservables(RooArgSet(x));
75  w->Print();
76 
77  //////// show use of Feldman-Cousins
78  RooStats::FeldmanCousins fc(*data,modelConfig);
79  fc.SetTestSize(.05); // set size of test
80  fc.UseAdaptiveSampling(true);
81  fc.FluctuateNumDataEntries(false); // number counting analysis: dataset always has 1 entry with N events observed
82  fc.SetNBins(100); // number of points to test per parameter
83 
84  // use the Feldman-Cousins tool
86 
87  // make a canvas for plots
88  TCanvas* intervalCanvas = new TCanvas("intervalCanvas");
89 
90  std::cout << "is this point in the interval? " <<
91  interval->IsInInterval(parameters) << std::endl;
92 
93  std::cout << "interval is ["<<
94  interval->LowerLimit(mu) << ", " <<
95  interval->UpperLimit(mu) << "]" << endl;
96 
97  // using 200 bins it takes 1 min and the answer is
98  // interval is [0.2625, 10.6125] with a step size of .075
99  // The interval in Feldman & Cousins's original paper is [.29, 10.81]
100  // Phys.Rev.D57:3873-3889,1998.
101 
102  // No dedicated plotting class yet, so do it by hand:
103 
104  RooDataHist* parameterScan = (RooDataHist*) fc.GetPointsToScan();
105  TH1F* hist = (TH1F*) parameterScan->createHistogram("mu",30);
106  hist->Draw();
107 
108 
109  RooArgSet* tmpPoint;
110  // loop over points to test
111  for(Int_t i=0; i<parameterScan->numEntries(); ++i){
112  // cout << "on parameter point " << i << " out of " << parameterScan->numEntries() << endl;
113  // get a parameter point from the list of points to test.
114  tmpPoint = (RooArgSet*) parameterScan->get(i)->clone("temp");
115 
116  TMarker* mark = new TMarker(tmpPoint->getRealValue("mu"), 1, 25);
117  if (interval->IsInInterval( *tmpPoint ) )
118  mark->SetMarkerColor(kBlue);
119  else
120  mark->SetMarkerColor(kRed);
121 
122  mark->Draw("s");
123  //delete tmpPoint;
124  // delete mark;
125  }
126  t.Stop();
127  t.Print();
128 
129 
130 }
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:217
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:56
Definition: Rtypes.h:61
#define mark(osub)
Definition: triangle.c:1206
virtual void SetTestSize(Double_t size)
set the size of the test (rate of Type I error) ( Eg. 0.05 for a 95% Confidence Interval) ...
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
Manages Markers.
Definition: TMarker.h:40
int Int_t
Definition: RtypesCore.h:41
virtual void SetObservables(const RooArgSet &set)
specify the observables
Definition: ModelConfig.h:156
virtual const RooArgSet * get() const
Definition: RooDataHist.h:77
RooAbsData * GetPointsToScan()
static struct mg_connection * fc(struct mg_context *ctx)
Definition: civetweb.c:839
RooDataSet is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
void Stop()
Stop the stopwatch.
Definition: TStopwatch.cxx:75
Double_t x[n]
Definition: legend1.C:17
virtual void SetMarkerColor(Color_t mcolor=1)
Definition: TAttMarker.h:51
virtual void SetPdf(const RooAbsPdf &pdf)
Set the Pdf, add to the the workspace if not already there.
Definition: ModelConfig.h:97
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:626
virtual Bool_t IsInInterval(const RooArgSet &) const
check if parameter is in the interval
void rs401c_FeldmanCousins()
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
TThread * t[5]
Definition: threadsh1.C:13
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
void SetNBins(Int_t bins)
void FluctuateNumDataEntries(bool flag=true)
tuple w
Definition: qtexample.py:51
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:735
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
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...
virtual TObject * clone(const char *newname) const
Definition: RooArgSet.h:82
PointSetInterval is a concrete implementation of the ConfInterval interface.
virtual PointSetInterval * GetInterval() const
Main interface to get a ConfInterval (will be a PointSetInterval)
virtual Int_t numEntries() const
Return the number of bins.
void Print(Option_t *opts=0) const
Print contents of the workspace.
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 void SetParametersOfInterest(const RooArgSet &set)
Definition: ModelConfig.h:115
Definition: Rtypes.h:61
RooAddition calculates the sum of a set of RooAbsReal terms, or when constructed with two sets...
Definition: RooAddition.h:26
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:42
void UseAdaptiveSampling(bool flag=true)
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