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