Logo ROOT   6.16/01
Reference Guide
StandardFeldmanCousinsDemo.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roostats
3/// \notebook -js
4/// Standard demo of the Feldman-Cousins calculator
5/// StandardFeldmanCousinsDemo
6///
7/// This is a standard demo that can be used with any ROOT file
8/// prepared in the standard way. You specify:
9/// - name for input ROOT file
10/// - name of workspace inside ROOT file that holds model and data
11/// - name of ModelConfig that specifies details for calculator tools
12/// - name of dataset
13///
14/// With default parameters the macro will attempt to run the
15/// standard hist2workspace example and read the ROOT file
16/// that it produces.
17///
18/// The actual heart of the demo is only about 10 lines long.
19///
20/// The FeldmanCousins tools is a classical frequentist calculation
21/// based on the Neyman Construction. The test statistic can be
22/// generalized for nuisance parameters by using the profile likelihood ratio.
23/// But unlike the ProfileLikelihoodCalculator, this tool explicitly
24/// builds the sampling distribution of the test statistic via toy Monte Carlo.
25///
26/// \macro_image
27/// \macro_output
28/// \macro_code
29///
30/// \author Kyle Cranmer
31
32#include "TFile.h"
33#include "TROOT.h"
34#include "TH1F.h"
35#include "TSystem.h"
36
37#include "RooWorkspace.h"
38#include "RooAbsData.h"
39
45
46
47using namespace RooFit;
48using namespace RooStats;
49
50void StandardFeldmanCousinsDemo(const char* infile = "",
51 const char* workspaceName = "combined",
52 const char* modelConfigName = "ModelConfig",
53 const char* dataName = "obsData"){
54
55 // -------------------------------------------------------
56 // First part is just to access a user-defined file
57 // or create the standard example file if it doesn't exist
58 const char* filename = "";
59 if (!strcmp(infile,"")) {
60 filename = "results/example_combined_GaussExample_model.root";
61 bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code
62 // if file does not exists generate with histfactory
63 if (!fileExist) {
64#ifdef _WIN32
65 cout << "HistFactory file cannot be generated on Windows - exit" << endl;
66 return;
67#endif
68 // Normally this would be run on the command line
69 cout <<"will run standard hist2workspace example"<<endl;
70 gROOT->ProcessLine(".! prepareHistFactory .");
71 gROOT->ProcessLine(".! hist2workspace config/example.xml");
72 cout <<"\n\n---------------------"<<endl;
73 cout <<"Done creating example input"<<endl;
74 cout <<"---------------------\n\n"<<endl;
75 }
76
77 }
78 else
79 filename = infile;
80
81 // Try to open the file
82 TFile *file = TFile::Open(filename);
83
84 // if input file was specified byt not found, quit
85 if(!file ){
86 cout <<"StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
87 return;
88 }
89
90
91 // -------------------------------------------------------
92 // Tutorial starts here
93 // -------------------------------------------------------
94
95 // get the workspace out of the file
96 RooWorkspace* w = (RooWorkspace*) file->Get(workspaceName);
97 if(!w){
98 cout <<"workspace not found" << endl;
99 return;
100 }
101
102 // get the modelConfig out of the file
103 ModelConfig* mc = (ModelConfig*) w->obj(modelConfigName);
104
105 // get the modelConfig out of the file
106 RooAbsData* data = w->data(dataName);
107
108 // make sure ingredients are found
109 if(!data || !mc){
110 w->Print();
111 cout << "data or ModelConfig was not found" <<endl;
112 return;
113 }
114
115 // -------------------------------------------------------
116 // create and use the FeldmanCousins tool
117 // to find and plot the 95% confidence interval
118 // on the parameter of interest as specified
119 // in the model config
120 FeldmanCousins fc(*data,*mc);
121 fc.SetConfidenceLevel(0.95); // 95% interval
122 //fc.AdditionalNToysFactor(0.1); // to speed up the result
123 fc.UseAdaptiveSampling(true); // speed it up a bit
124 fc.SetNBins(10); // set how many points per parameter of interest to scan
125 fc.CreateConfBelt(true); // save the information in the belt for plotting
126
127 // Since this tool needs to throw toy MC the PDF needs to be
128 // extended or the tool needs to know how many entries in a dataset
129 // per pseudo experiment.
130 // In the 'number counting form' where the entries in the dataset
131 // are counts, and not values of discriminating variables, the
132 // datasets typically only have one entry and the PDF is not
133 // extended.
134 if(!mc->GetPdf()->canBeExtended()){
135 if(data->numEntries()==1)
136 fc.FluctuateNumDataEntries(false);
137 else
138 cout <<"Not sure what to do about this model" <<endl;
139 }
140
141 // We can use PROOF to speed things along in parallel
142 // ProofConfig pc(*w, 1, "workers=4", kFALSE);
143 // ToyMCSampler* toymcsampler = (ToyMCSampler*) fc.GetTestStatSampler();
144 // toymcsampler->SetProofConfig(&pc); // enable proof
145
146
147 // Now get the interval
148 PointSetInterval* interval = fc.GetInterval();
149 ConfidenceBelt* belt = fc.GetConfidenceBelt();
150
151 // print out the interval on the first Parameter of Interest
152 RooRealVar* firstPOI = (RooRealVar*) mc->GetParametersOfInterest()->first();
153 cout << "\n95% interval on " <<firstPOI->GetName()<<" is : ["<<
154 interval->LowerLimit(*firstPOI) << ", "<<
155 interval->UpperLimit(*firstPOI) <<"] "<<endl;
156
157 // ---------------------------------------------
158 // No nice plots yet, so plot the belt by hand
159
160 // Ask the calculator which points were scanned
161 RooDataSet* parameterScan = (RooDataSet*) fc.GetPointsToScan();
162 RooArgSet* tmpPoint;
163
164 // make a histogram of parameter vs. threshold
165 TH1F* histOfThresholds = new TH1F("histOfThresholds","",
166 parameterScan->numEntries(),
167 firstPOI->getMin(),
168 firstPOI->getMax());
169
170 // loop through the points that were tested and ask confidence belt
171 // what the upper/lower thresholds were.
172 // For FeldmanCousins, the lower cut off is always 0
173 for(Int_t i=0; i<parameterScan->numEntries(); ++i){
174 tmpPoint = (RooArgSet*) parameterScan->get(i)->clone("temp");
175 double arMax = belt->GetAcceptanceRegionMax(*tmpPoint);
176 double arMin = belt->GetAcceptanceRegionMax(*tmpPoint);
177 double poiVal = tmpPoint->getRealValue(firstPOI->GetName()) ;
178 histOfThresholds->Fill(poiVal,arMax);
179 }
180 histOfThresholds->SetMinimum(0);
181 histOfThresholds->Draw();
182
183}
int Int_t
Definition: RtypesCore.h:41
#define gROOT
Definition: TROOT.h:410
R__EXTERN TSystem * gSystem
Definition: TSystem.h:540
static struct mg_connection * fc(struct mg_context *ctx)
Definition: civetweb.c:3728
RooAbsArg * first() const
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
virtual Int_t numEntries() const
Definition: RooAbsData.cxx:285
Bool_t canBeExtended() const
Definition: RooAbsPdf.h:230
virtual Double_t getMax(const char *name=0) const
virtual Double_t getMin(const char *name=0) const
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
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:472
virtual TObject * clone(const char *newname) const
Definition: RooArgSet.h:84
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
virtual const RooArgSet * get(Int_t index) const
Return RooArgSet with coordinates of event 'index'.
Definition: RooDataSet.cxx:995
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
ConfidenceBelt is a concrete implementation of the ConfInterval interface.
Double_t GetAcceptanceRegionMax(RooArgSet &, Double_t cl=-1., Double_t leftside=-1.)
The FeldmanCousins class (like the Feldman-Cousins technique) is essentially a specific configuration...
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:30
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return NULL if not existing)
Definition: ModelConfig.h:231
RooAbsPdf * GetPdf() const
get model PDF (return NULL if pdf has not been specified or does not exist)
Definition: ModelConfig.h:228
PointSetInterval is a concrete implementation of the ConfInterval interface.
Double_t UpperLimit(RooRealVar &param)
return upper limit on a given parameter
Double_t LowerLimit(RooRealVar &param)
return lower limit on a given parameter
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:43
RooAbsData * data(const char *name) const
Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found.
void Print(Option_t *opts=0) const
Print contents of the workspace.
TObject * obj(const char *name) const
Return any type of object (RooAbsArg, RooAbsData or generic object) with given name)
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseGeneralPurpose, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3975
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:571
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3251
virtual void SetMinimum(Double_t minimum=-1111)
Definition: TH1.h:395
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2974
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition: TSystem.cxx:1286
@(#)root/roostats:$Id$ Author: George Lewis, Kyle Cranmer
Definition: Asimov.h:20
Definition: file.py:1