ROOT logo

From $ROOTSYS/tutorials/roostats/StandardFeldmanCousinsDemo.C

// Standard demo of the Feldman-Cousins calculator
/*
StandardFeldmanCousinsDemo 

Author: Kyle Cranmer
date: Dec. 2010

This is a standard demo that can be used with any ROOT file 
prepared in the standard way.  You specify:
 - name for input ROOT file
 - name of workspace inside ROOT file that holds model and data
 - name of ModelConfig that specifies details for calculator tools
 - name of dataset 

With default parameters the macro will attempt to run the 
standard hist2workspace example and read the ROOT file
that it produces.

The actual heart of the demo is only about 10 lines long.

The FeldmanCousins tools is a classical frequentist calculation
based on the Neyman Construction.  The test statistic can be
generalized for nuisance parameters by using the profile likeihood ratio.
But unlike the ProfileLikelihoodCalculator, this tool explicitly
builds the sampling distribution of the test statistic via toy Monte Carlo.
*/

#include "TFile.h"
#include "TROOT.h"
#include "TH1F.h"
#include "TSystem.h"

#include "RooWorkspace.h"
#include "RooAbsData.h"

#include "RooStats/ModelConfig.h"
#include "RooStats/FeldmanCousins.h"
#include "RooStats/ToyMCSampler.h"
#include "RooStats/PointSetInterval.h"
#include "RooStats/ConfidenceBelt.h"


using namespace RooFit;
using namespace RooStats;

void StandardFeldmanCousinsDemo(const char* infile = "",
		      const char* workspaceName = "combined",
		      const char* modelConfigName = "ModelConfig",
		      const char* dataName = "obsData"){

  /////////////////////////////////////////////////////////////
  // First part is just to access a user-defined file 
  // or create the standard example file if it doesn't exist
  ////////////////////////////////////////////////////////////
  const char* filename = "";
  if (!strcmp(infile,"")) {
    filename = "results/example_combined_GaussExample_model.root";
     bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code
     // if file does not exists generate with histfactory
     if (!fileExist) {
#ifdef _WIN32
        cout << "HistFactory file cannot be generated on Windows - exit" << endl;
        return;
#endif
        // Normally this would be run on the command line
        cout <<"will run standard hist2workspace example"<<endl;
        gROOT->ProcessLine(".! prepareHistFactory .");
        gROOT->ProcessLine(".! hist2workspace config/example.xml");
        cout <<"\n\n---------------------"<<endl;
        cout <<"Done creating example input"<<endl;
        cout <<"---------------------\n\n"<<endl;
     }

  }
  else
    filename = infile;
   
  // Try to open the file
  TFile *file = TFile::Open(filename);

  // if input file was specified byt not found, quit
  if(!file ){
    cout <<"StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
    return;
  } 

  
  /////////////////////////////////////////////////////////////
  // Tutorial starts here
  ////////////////////////////////////////////////////////////

  // get the workspace out of the file
  RooWorkspace* w = (RooWorkspace*) file->Get(workspaceName);
  if(!w){
    cout <<"workspace not found" << endl;
    return;
  }

  // get the modelConfig out of the file
  ModelConfig* mc = (ModelConfig*) w->obj(modelConfigName);

  // get the modelConfig out of the file
  RooAbsData* data = w->data(dataName);

  // make sure ingredients are found
  if(!data || !mc){
    w->Print();
    cout << "data or ModelConfig was not found" <<endl;
    return;
  }

  /////////////////////////////////////////////
  // create and use the FeldmanCousins tool
  // to find and plot the 95% confidence interval
  // on the parameter of interest as specified
  // in the model config
  FeldmanCousins fc(*data,*mc);
  fc.SetConfidenceLevel(0.95); // 95% interval
  //fc.AdditionalNToysFactor(0.1); // to speed up the result 
  fc.UseAdaptiveSampling(true); // speed it up a bit
  fc.SetNBins(10); // set how many points per parameter of interest to scan
  fc.CreateConfBelt(true); // save the information in the belt for plotting

  // Since this tool needs to throw toy MC the PDF needs to be
  // extended or the tool needs to know how many entries in a dataset
  // per pseudo experiment.  
  // In the 'number counting form' where the entries in the dataset
  // are counts, and not values of discriminating variables, the
  // datasets typically only have one entry and the PDF is not
  // extended.  
  if(!mc->GetPdf()->canBeExtended()){
    if(data->numEntries()==1)     
      fc.FluctuateNumDataEntries(false);
    else
      cout <<"Not sure what to do about this model" <<endl;
  }

  // We can use PROOF to speed things along in parallel
  //  ProofConfig pc(*w, 1, "workers=4", kFALSE);
  //  ToyMCSampler*  toymcsampler = (ToyMCSampler*) fc.GetTestStatSampler();
  //  toymcsampler->SetProofConfig(&pc);	// enable proof


  // Now get the interval
  PointSetInterval* interval = fc.GetInterval();
  ConfidenceBelt* belt = fc.GetConfidenceBelt();
 
  // print out the iterval on the first Parameter of Interest
  RooRealVar* firstPOI = (RooRealVar*) mc->GetParametersOfInterest()->first();
  cout << "\n95% interval on " <<firstPOI->GetName()<<" is : ["<<
    interval->LowerLimit(*firstPOI) << ", "<<
    interval->UpperLimit(*firstPOI) <<"] "<<endl;

  //////////////////////////////////////////////
  // No nice plots yet, so plot the belt by hand
  
  // Ask the calculator which points were scanned
  RooDataSet* parameterScan = (RooDataSet*) fc.GetPointsToScan();
  RooArgSet* tmpPoint;

  // make a histogram of parameter vs. threshold
  TH1F* histOfThresholds = new TH1F("histOfThresholds","",
				    parameterScan->numEntries(),
				    firstPOI->getMin(),
				    firstPOI->getMax());

  // loop through the points that were tested and ask confidence belt
  // what the upper/lower thresholds were.
  // For FeldmanCousins, the lower cut off is always 0
  for(Int_t i=0; i<parameterScan->numEntries(); ++i){
    tmpPoint = (RooArgSet*) parameterScan->get(i)->clone("temp");
    double arMax = belt->GetAcceptanceRegionMax(*tmpPoint);
    double arMin = belt->GetAcceptanceRegionMax(*tmpPoint);
    double poiVal = tmpPoint->getRealValue(firstPOI->GetName()) ;
    histOfThresholds->Fill(poiVal,arMax);
  }
  histOfThresholds->SetMinimum(0);
  histOfThresholds->Draw();

}
 StandardFeldmanCousinsDemo.C:1
 StandardFeldmanCousinsDemo.C:2
 StandardFeldmanCousinsDemo.C:3
 StandardFeldmanCousinsDemo.C:4
 StandardFeldmanCousinsDemo.C:5
 StandardFeldmanCousinsDemo.C:6
 StandardFeldmanCousinsDemo.C:7
 StandardFeldmanCousinsDemo.C:8
 StandardFeldmanCousinsDemo.C:9
 StandardFeldmanCousinsDemo.C:10
 StandardFeldmanCousinsDemo.C:11
 StandardFeldmanCousinsDemo.C:12
 StandardFeldmanCousinsDemo.C:13
 StandardFeldmanCousinsDemo.C:14
 StandardFeldmanCousinsDemo.C:15
 StandardFeldmanCousinsDemo.C:16
 StandardFeldmanCousinsDemo.C:17
 StandardFeldmanCousinsDemo.C:18
 StandardFeldmanCousinsDemo.C:19
 StandardFeldmanCousinsDemo.C:20
 StandardFeldmanCousinsDemo.C:21
 StandardFeldmanCousinsDemo.C:22
 StandardFeldmanCousinsDemo.C:23
 StandardFeldmanCousinsDemo.C:24
 StandardFeldmanCousinsDemo.C:25
 StandardFeldmanCousinsDemo.C:26
 StandardFeldmanCousinsDemo.C:27
 StandardFeldmanCousinsDemo.C:28
 StandardFeldmanCousinsDemo.C:29
 StandardFeldmanCousinsDemo.C:30
 StandardFeldmanCousinsDemo.C:31
 StandardFeldmanCousinsDemo.C:32
 StandardFeldmanCousinsDemo.C:33
 StandardFeldmanCousinsDemo.C:34
 StandardFeldmanCousinsDemo.C:35
 StandardFeldmanCousinsDemo.C:36
 StandardFeldmanCousinsDemo.C:37
 StandardFeldmanCousinsDemo.C:38
 StandardFeldmanCousinsDemo.C:39
 StandardFeldmanCousinsDemo.C:40
 StandardFeldmanCousinsDemo.C:41
 StandardFeldmanCousinsDemo.C:42
 StandardFeldmanCousinsDemo.C:43
 StandardFeldmanCousinsDemo.C:44
 StandardFeldmanCousinsDemo.C:45
 StandardFeldmanCousinsDemo.C:46
 StandardFeldmanCousinsDemo.C:47
 StandardFeldmanCousinsDemo.C:48
 StandardFeldmanCousinsDemo.C:49
 StandardFeldmanCousinsDemo.C:50
 StandardFeldmanCousinsDemo.C:51
 StandardFeldmanCousinsDemo.C:52
 StandardFeldmanCousinsDemo.C:53
 StandardFeldmanCousinsDemo.C:54
 StandardFeldmanCousinsDemo.C:55
 StandardFeldmanCousinsDemo.C:56
 StandardFeldmanCousinsDemo.C:57
 StandardFeldmanCousinsDemo.C:58
 StandardFeldmanCousinsDemo.C:59
 StandardFeldmanCousinsDemo.C:60
 StandardFeldmanCousinsDemo.C:61
 StandardFeldmanCousinsDemo.C:62
 StandardFeldmanCousinsDemo.C:63
 StandardFeldmanCousinsDemo.C:64
 StandardFeldmanCousinsDemo.C:65
 StandardFeldmanCousinsDemo.C:66
 StandardFeldmanCousinsDemo.C:67
 StandardFeldmanCousinsDemo.C:68
 StandardFeldmanCousinsDemo.C:69
 StandardFeldmanCousinsDemo.C:70
 StandardFeldmanCousinsDemo.C:71
 StandardFeldmanCousinsDemo.C:72
 StandardFeldmanCousinsDemo.C:73
 StandardFeldmanCousinsDemo.C:74
 StandardFeldmanCousinsDemo.C:75
 StandardFeldmanCousinsDemo.C:76
 StandardFeldmanCousinsDemo.C:77
 StandardFeldmanCousinsDemo.C:78
 StandardFeldmanCousinsDemo.C:79
 StandardFeldmanCousinsDemo.C:80
 StandardFeldmanCousinsDemo.C:81
 StandardFeldmanCousinsDemo.C:82
 StandardFeldmanCousinsDemo.C:83
 StandardFeldmanCousinsDemo.C:84
 StandardFeldmanCousinsDemo.C:85
 StandardFeldmanCousinsDemo.C:86
 StandardFeldmanCousinsDemo.C:87
 StandardFeldmanCousinsDemo.C:88
 StandardFeldmanCousinsDemo.C:89
 StandardFeldmanCousinsDemo.C:90
 StandardFeldmanCousinsDemo.C:91
 StandardFeldmanCousinsDemo.C:92
 StandardFeldmanCousinsDemo.C:93
 StandardFeldmanCousinsDemo.C:94
 StandardFeldmanCousinsDemo.C:95
 StandardFeldmanCousinsDemo.C:96
 StandardFeldmanCousinsDemo.C:97
 StandardFeldmanCousinsDemo.C:98
 StandardFeldmanCousinsDemo.C:99
 StandardFeldmanCousinsDemo.C:100
 StandardFeldmanCousinsDemo.C:101
 StandardFeldmanCousinsDemo.C:102
 StandardFeldmanCousinsDemo.C:103
 StandardFeldmanCousinsDemo.C:104
 StandardFeldmanCousinsDemo.C:105
 StandardFeldmanCousinsDemo.C:106
 StandardFeldmanCousinsDemo.C:107
 StandardFeldmanCousinsDemo.C:108
 StandardFeldmanCousinsDemo.C:109
 StandardFeldmanCousinsDemo.C:110
 StandardFeldmanCousinsDemo.C:111
 StandardFeldmanCousinsDemo.C:112
 StandardFeldmanCousinsDemo.C:113
 StandardFeldmanCousinsDemo.C:114
 StandardFeldmanCousinsDemo.C:115
 StandardFeldmanCousinsDemo.C:116
 StandardFeldmanCousinsDemo.C:117
 StandardFeldmanCousinsDemo.C:118
 StandardFeldmanCousinsDemo.C:119
 StandardFeldmanCousinsDemo.C:120
 StandardFeldmanCousinsDemo.C:121
 StandardFeldmanCousinsDemo.C:122
 StandardFeldmanCousinsDemo.C:123
 StandardFeldmanCousinsDemo.C:124
 StandardFeldmanCousinsDemo.C:125
 StandardFeldmanCousinsDemo.C:126
 StandardFeldmanCousinsDemo.C:127
 StandardFeldmanCousinsDemo.C:128
 StandardFeldmanCousinsDemo.C:129
 StandardFeldmanCousinsDemo.C:130
 StandardFeldmanCousinsDemo.C:131
 StandardFeldmanCousinsDemo.C:132
 StandardFeldmanCousinsDemo.C:133
 StandardFeldmanCousinsDemo.C:134
 StandardFeldmanCousinsDemo.C:135
 StandardFeldmanCousinsDemo.C:136
 StandardFeldmanCousinsDemo.C:137
 StandardFeldmanCousinsDemo.C:138
 StandardFeldmanCousinsDemo.C:139
 StandardFeldmanCousinsDemo.C:140
 StandardFeldmanCousinsDemo.C:141
 StandardFeldmanCousinsDemo.C:142
 StandardFeldmanCousinsDemo.C:143
 StandardFeldmanCousinsDemo.C:144
 StandardFeldmanCousinsDemo.C:145
 StandardFeldmanCousinsDemo.C:146
 StandardFeldmanCousinsDemo.C:147
 StandardFeldmanCousinsDemo.C:148
 StandardFeldmanCousinsDemo.C:149
 StandardFeldmanCousinsDemo.C:150
 StandardFeldmanCousinsDemo.C:151
 StandardFeldmanCousinsDemo.C:152
 StandardFeldmanCousinsDemo.C:153
 StandardFeldmanCousinsDemo.C:154
 StandardFeldmanCousinsDemo.C:155
 StandardFeldmanCousinsDemo.C:156
 StandardFeldmanCousinsDemo.C:157
 StandardFeldmanCousinsDemo.C:158
 StandardFeldmanCousinsDemo.C:159
 StandardFeldmanCousinsDemo.C:160
 StandardFeldmanCousinsDemo.C:161
 StandardFeldmanCousinsDemo.C:162
 StandardFeldmanCousinsDemo.C:163
 StandardFeldmanCousinsDemo.C:164
 StandardFeldmanCousinsDemo.C:165
 StandardFeldmanCousinsDemo.C:166
 StandardFeldmanCousinsDemo.C:167
 StandardFeldmanCousinsDemo.C:168
 StandardFeldmanCousinsDemo.C:169
 StandardFeldmanCousinsDemo.C:170
 StandardFeldmanCousinsDemo.C:171
 StandardFeldmanCousinsDemo.C:172
 StandardFeldmanCousinsDemo.C:173
 StandardFeldmanCousinsDemo.C:174
 StandardFeldmanCousinsDemo.C:175
 StandardFeldmanCousinsDemo.C:176
 StandardFeldmanCousinsDemo.C:177
 StandardFeldmanCousinsDemo.C:178
 StandardFeldmanCousinsDemo.C:179
 StandardFeldmanCousinsDemo.C:180
 StandardFeldmanCousinsDemo.C:181