Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
FeldmanCousins.cxx
Go to the documentation of this file.
1 // @(#)root/roostats:$Id$
2 // Author: Kyle Cranmer January 2009
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 
14 
15 #include "RooStats/RooStatsUtils.h"
16 
18 
19 #include "RooStats/ModelConfig.h"
20 
22 
25 #include "RooStats/RooStatsUtils.h"
26 
27 #include "RooDataSet.h"
28 #include "RooDataHist.h"
29 #include "RooGlobalFunc.h"
30 #include "RooMsgService.h"
31 #include "TFile.h"
32 #include "TTree.h"
33 
34 /** \class RooStats::FeldmanCousins
35  \ingroup Roostats
36 
37 The FeldmanCousins class (like the Feldman-Cousins technique) is essentially a
38 specific configuration of the more general NeymanConstruction. It is a concrete
39 implementation of the IntervalCalculator interface that, which uses the
40 NeymanConstruction in a particular way. As the name suggests, it returns a
41 ConfidenceInterval. In particular, it produces a RooStats::PointSetInterval,
42 which is a concrete implementation of the ConfInterval interface.
43 
44 The Neyman Construction is not a uniquely defined statistical technique, it
45 requires that one specify an ordering rule or ordering principle, which is
46 usually encoded by choosing a specific test statistic and limits of integration
47 (corresponding to upper/lower/central limits). As a result, this class must be
48 configured with the corresponding information before it can produce an interval.
49 
50 In the case of the Feldman-Cousins approach, the ordering principle is the
51 likelihood ratio -- motivated by the Neyman-Pearson lemma. When nuisance
52 parameters are involved, the profile likelihood ratio is the natural
53 generalization. One may either choose to perform the construction over the full
54 space of the nuisance parameters, or restrict the nuisance parameters to their
55 conditional MLE (eg. profiled values).
56 
57 */
58 
59 ClassImp(RooStats::FeldmanCousins); ;
60 
61 using namespace RooFit;
62 using namespace RooStats;
63 using namespace std;
64 
65 ////////////////////////////////////////////////////////////////////////////////
66 /// standard constructor
67 
68 FeldmanCousins::FeldmanCousins(RooAbsData& data, ModelConfig& model) :
69  fSize(0.05),
70  fModel(model),
71  fData(data),
72  fTestStatSampler(0),
73  fPointsToTest(0),
74  fPOIToTest(0),
75  fConfBelt(0),
76  fAdaptiveSampling(false),
77  fAdditionalNToysFactor(1.),
78  fNbins(10),
79  fFluctuateData(true),
80  fDoProfileConstruction(true),
81  fSaveBeltToFile(false),
82  fCreateBelt(false)
83 {
84 }
85 
86 ////////////////////////////////////////////////////////////////////////////////
87 /// destructor
88 
89 FeldmanCousins::~FeldmanCousins() {
90  if(fPointsToTest) delete fPointsToTest;
91  if(fPOIToTest) delete fPOIToTest;
92  if(fTestStatSampler) delete fTestStatSampler;
93 }
94 
95 ////////////////////////////////////////////////////////////////////////////////
96 /// set the model
97 
98 void FeldmanCousins::SetModel(const ModelConfig & model) {
99  fModel = model;
100 }
101 
102 ////////////////////////////////////////////////////////////////////////////////
103 
104 TestStatSampler* FeldmanCousins::GetTestStatSampler() const{
105  if(!fTestStatSampler)
106  this->CreateTestStatSampler();
107  return fTestStatSampler;
108 }
109 
110 ////////////////////////////////////////////////////////////////////////////////
111 /// specify the Test Statistic and create a ToyMC test statistic sampler
112 
113 void FeldmanCousins::CreateTestStatSampler() const{
114  // use the profile likelihood ratio as the test statistic
115  ProfileLikelihoodTestStat* testStatistic = new ProfileLikelihoodTestStat(*fModel.GetPdf());
116 
117  // create the ToyMC test statistic sampler
118  fTestStatSampler = new ToyMCSampler(*testStatistic,int(fAdditionalNToysFactor*50./fSize)) ;
119  fTestStatSampler->SetParametersForTestStat(*fModel.GetParametersOfInterest() );
120  if(fModel.GetObservables())
121  fTestStatSampler->SetObservables(*fModel.GetObservables());
122  fTestStatSampler->SetPdf(*fModel.GetPdf());
123 
124  if(!fAdaptiveSampling){
125  ooccoutP(&fModel,Generation) << "FeldmanCousins: ntoys per point = " << (int) (fAdditionalNToysFactor*50./fSize) << endl;
126  } else{
127  ooccoutP(&fModel,Generation) << "FeldmanCousins: ntoys per point: adaptive" << endl;
128  }
129  if(fFluctuateData){
130  ooccoutP(&fModel,Generation) << "FeldmanCousins: nEvents per toy will fluctuate about expectation" << endl;
131  } else{
132  ooccoutP(&fModel,Generation) << "FeldmanCousins: nEvents per toy will not fluctuate, will always be " << fData.numEntries() << endl;
133  fTestStatSampler->SetNEventsPerToy(fData.numEntries());
134  }
135 }
136 
137 ////////////////////////////////////////////////////////////////////////////////
138 /// specify the parameter points to perform the construction.
139 /// allow ability to profile on some nuisance parameters
140 
141 void FeldmanCousins::CreateParameterPoints() const{
142  // get ingredients
143  RooAbsPdf* pdf = fModel.GetPdf();
144  if (!pdf ){
145  ooccoutE(&fModel,Generation) << "FeldmanCousins: ModelConfig has no PDF" << endl;
146  return;
147  }
148 
149  // get list of all parameters
150  RooArgSet* parameters = new RooArgSet(*fModel.GetParametersOfInterest());
151  if(fModel.GetNuisanceParameters())
152  parameters->add(*fModel.GetNuisanceParameters());
153 
154 
155  if( fModel.GetNuisanceParameters() && ! fModel.GetParametersOfInterest()->equals(*parameters) && fDoProfileConstruction) {
156  // if parameters include nuisance parameters, do profile construction
157  ooccoutP(&fModel,Generation) << "FeldmanCousins: Model has nuisance parameters, will do profile construction" << endl;
158 
159  // set nbins for the POI
160  TIter it2 = fModel.GetParametersOfInterest()->createIterator();
161  RooRealVar *myarg2;
162  while ((myarg2 = dynamic_cast<RooRealVar*>(it2.Next()))) {
163  myarg2->setBins(fNbins);
164  }
165 
166  // get dataset for POI scan
167  // RooDataHist* parameterScan = NULL;
168  RooAbsData* parameterScan = NULL;
169  if(fPOIToTest)
170  parameterScan = fPOIToTest;
171  else
172  parameterScan = new RooDataHist("parameterScan", "", *fModel.GetParametersOfInterest());
173 
174 
175  ooccoutP(&fModel,Generation) << "FeldmanCousins: # points to test = " << parameterScan->numEntries() << endl;
176  // make profile construction
177  RooFit::MsgLevel previous = RooMsgService::instance().globalKillBelow();
178  RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL) ;
179  RooAbsReal* nll = pdf->createNLL(fData,RooFit::CloneData(false));
180  RooAbsReal* profile = nll->createProfile(*fModel.GetParametersOfInterest());
181 
182  RooDataSet* profileConstructionPoints = new RooDataSet("profileConstruction",
183  "profileConstruction",
184  *parameters);
185 
186 
187  for(Int_t i=0; i<parameterScan->numEntries(); ++i){
188  // here's where we figure out the profiled value of nuisance parameters
189  *parameters = *parameterScan->get(i);
190  profile->getVal();
191  profileConstructionPoints->add(*parameters);
192  }
193  RooMsgService::instance().setGlobalKillBelow(previous) ;
194  delete profile;
195  delete nll;
196  if(!fPOIToTest) delete parameterScan;
197 
198  // done
199  fPointsToTest = profileConstructionPoints;
200 
201  } else{
202  // Do full construction
203  ooccoutP(&fModel,Generation) << "FeldmanCousins: Model has no nuisance parameters" << endl;
204 
205  TIter it = parameters->createIterator();
206  RooRealVar *myarg;
207  while ((myarg = dynamic_cast<RooRealVar*>(it.Next()))) {
208  myarg->setBins(fNbins);
209  }
210 
211  RooDataHist* parameterScan = new RooDataHist("parameterScan", "", *parameters);
212  ooccoutP(&fModel,Generation) << "FeldmanCousins: # points to test = " << parameterScan->numEntries() << endl;
213 
214  fPointsToTest = parameterScan;
215  }
216 
217  delete parameters;
218 
219 }
220 
221 ////////////////////////////////////////////////////////////////////////////////
222 /// Main interface to get a RooStats::ConfInterval.
223 /// It constructs a RooStats::PointSetInterval.
224 
225 PointSetInterval* FeldmanCousins::GetInterval() const {
226  // local variables
227  // RooAbsData* data = fData; //fWS->data(fDataName);
228 
229  // fill in implied variables given data
230  fModel.GuessObsAndNuisance(fData);
231 
232  // create the test statistic sampler (private data member fTestStatSampler)
233  if(!fTestStatSampler)
234  this->CreateTestStatSampler();
235 
236  fTestStatSampler->SetObservables(*fModel.GetObservables());
237 
238  if(!fFluctuateData)
239  fTestStatSampler->SetNEventsPerToy(fData.numEntries());
240 
241  // create parameter points to perform construction (private data member fPointsToTest)
242  this->CreateParameterPoints();
243 
244  // Create a Neyman Construction
245  NeymanConstruction nc(fData,fModel);
246  // configure it
247  nc.SetTestStatSampler(*fTestStatSampler);
248  nc.SetTestSize(fSize); // set size of test
249  // nc.SetParameters( fModel.GetParametersOfInterest);
250  nc.SetParameterPointsToTest( *fPointsToTest );
251  nc.SetLeftSideTailFraction(0.); // part of definition of Feldman-Cousins
252  nc.SetData(fData);
253  nc.UseAdaptiveSampling(fAdaptiveSampling);
254  nc.AdditionalNToysFactor(fAdditionalNToysFactor);
255  nc.SaveBeltToFile(fSaveBeltToFile);
256  nc.CreateConfBelt(fCreateBelt);
257  fConfBelt = nc.GetConfidenceBelt();
258  // use it
259  return nc.GetInterval();
260 }