ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ConfidenceBelt.h
Go to the documentation of this file.
1 // @(#)root/roostats:$Id$
2 // Author: Kyle Cranmer, Lorenzo Moneta, Gregory Schott, Wouter Verkerke
3 /*************************************************************************
4  * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
5  * All rights reserved. *
6  * *
7  * For the licensing terms see $ROOTSYS/LICENSE. *
8  * For the list of contributors see $ROOTSYS/README/CREDITS. *
9  *************************************************************************/
10 
11 #ifndef RooStats_ConfidenceBelt
12 #define RooStats_ConfidenceBelt
13 
14 #ifndef ROO_ARG_SET
15 #include "RooArgSet.h"
16 #endif
17 #ifndef ROO_TREE_DATA
18 #include "RooAbsData.h"
19 #endif
20 #ifndef RooStats_ConfInterval
21 #include "RooStats/ConfInterval.h"
22 #endif
23 
25 
26 #include "TRef.h"
27 
28 #include <vector>
29 #include <map>
30 
31 
32 namespace RooStats {
33 
34  ///////////////////////////
35  class SamplingSummaryLookup : public TObject {
36 
37  typedef std::pair<Double_t, Double_t> AcceptanceCriteria; // defined by Confidence level, leftside tail probability
38  typedef std::map<Int_t, AcceptanceCriteria> LookupTable; // map ( Index, ( CL, leftside tail prob) )
39 
40  public:
43 
44  void Add(Double_t cl, Double_t leftside){
45  // add cl,leftside pair to lookup table
46  AcceptanceCriteria tmp(cl, leftside);
47 
48  // should check to see if this is already in the map
49  if(GetLookupIndex(cl,leftside) >=0 ){
50  std::cout<< "SamplingSummaryLookup::Add, already in lookup table" << std::endl;
51  } else
52  fLookupTable[fLookupTable.size()]= tmp;
53  }
54 
56  // get index for cl,leftside pair
57  AcceptanceCriteria tmp(cl, leftside);
58 
59  Double_t tolerance = 1E-6; // some small number to protect floating point comparison. What is better way?
60  LookupTable::iterator it = fLookupTable.begin();
61  Int_t index = -1;
62  for(; it!= fLookupTable.end(); ++it) {
63  index++;
64  if( TMath::Abs( (*it).second.first - cl ) < tolerance &&
65  TMath::Abs( (*it).second.second - leftside ) < tolerance )
66  break; // exit loop, found
67  }
68 
69  // check that it was found
70  if(index == (Int_t)fLookupTable.size())
71  index = -1;
72 
73  return index;
74  }
75 
77  if(index<0 || index>(Int_t)fLookupTable.size()) {
78  std::cout << "SamplingSummaryLookup::GetConfidenceLevel, index not in lookup table" << std::endl;
79  return -1;
80  }
81  return fLookupTable[index].first;
82  }
83 
85  if(index<0 || index>(Int_t)fLookupTable.size()) {
86  std::cout << "SamplingSummaryLookup::GetLeftSideTailFraction, index not in lookup table" << std::endl;
87  return -1;
88  }
89  return fLookupTable[index].second;
90  }
91 
92  private:
93  LookupTable fLookupTable; // map ( Index, ( CL, leftside tail prob) )
94 
95  protected:
96  ClassDef(SamplingSummaryLookup,1) // A simple class used by ConfidenceBelt
97  };
98 
99 
100  ///////////////////////////
101  class AcceptanceRegion : public TObject{
102  public:
104  virtual ~AcceptanceRegion() {}
105 
107  fLookupIndex = lu;
108  fLowerLimit = ll;
109  fUpperLimit = ul;
110  }
114 
115  private:
116  Int_t fLookupIndex; // want a small footprint reference to the RooArgSet for particular parameter point
117  Double_t fLowerLimit; // lower limit on test statistic
118  Double_t fUpperLimit; // upper limit on test statistic
119 
120  protected:
121  ClassDef(AcceptanceRegion,1) // A simple class for acceptance regions used for ConfidenceBelt
122 
123  };
124 
125 
126  ///////////////////////////
127  class SamplingSummary : public TObject {
128  public:
130  virtual ~SamplingSummary() {}
133  }
136  return (SamplingDistribution*) fSamplingDistribution.GetObject(); // dereference TRef
137  }
139 
141  Int_t index = ar.GetLookupIndex();
142  if( fAcceptanceRegions.count(index) !=0) {
143  std::cout << "SamplingSummary::AddAcceptanceRegion, need to implement merging protocol" << std::endl;
144  } else {
145  fAcceptanceRegions[index]=ar;
146  }
147  }
148 
149  private:
150  Int_t fParameterPointIndex; // want a small footprint reference to the RooArgSet for particular parameter point
151  TRef fSamplingDistribution; // persistent pointer to a SamplingDistribution
152  std::map<Int_t, AcceptanceRegion> fAcceptanceRegions;
153 
154  protected:
155  ClassDef(SamplingSummary,1) // A summary of acceptance regions for confidence belt
156 
157  };
158 
159  /////////////////////////////////////////
160 
161 
162  /**
163 
164  \ingroup Roostats
165 
166  ConfidenceBelt is a concrete implementation of the ConfInterval interface.
167  It implements simple general purpose interval of arbitrary dimensions and shape.
168  It does not assume the interval is connected.
169  It uses either a RooDataSet (eg. a list of parameter points in the interval) or
170  a RooDataHist (eg. a Histogram-like object for small regions of the parameter space) to
171  store the interval.
172 
173  */
174 
175 
176  class ConfidenceBelt : public TNamed {
177 
178  private:
180  std::vector<SamplingSummary> fSamplingSummaries; // composite of several AcceptanceRegions
181  RooAbsData* fParameterPoints; // either a histogram (RooDataHist) or a tree (RooDataSet)
182 
183 
184  public:
185  // constructors,destructors
186  ConfidenceBelt();
187  ConfidenceBelt(const char* name);
188  ConfidenceBelt(const char* name, const char* title);
189  ConfidenceBelt(const char* name, RooAbsData&);
190  ConfidenceBelt(const char* name, const char* title, RooAbsData&);
191  virtual ~ConfidenceBelt();
192 
193  // add after creating a region
194  void AddAcceptanceRegion(RooArgSet&, AcceptanceRegion region, Double_t cl=-1., Double_t leftside=-1.);
195 
196  // add without creating a region, more useful for clients
197  void AddAcceptanceRegion(RooArgSet& point, Int_t dataSetIndex, Double_t lower, Double_t upper, Double_t cl=-1., Double_t leftside=-1.);
198 
202  std::vector<Double_t> ConfidenceLevels() const ;
203 
204  // Method to return lower limit on a given parameter
205  // Double_t LowerLimit(RooRealVar& param) ; // could provide, but misleading?
206  // Double_t UpperLimit(RooRealVar& param) ; // could provide, but misleading?
207 
208  // do we want it to return list of parameters
209  virtual RooArgSet* GetParameters() const;
210 
211  // check if parameters are correct. (dummy implementation to start)
213 
214  protected:
215  ClassDef(ConfidenceBelt,1) // A confidence belt for the Neyman Construction
216 
217  };
218 }
219 
220 #endif
SamplingSummary(AcceptanceRegion &ar)
ConfidenceBelt()
Default constructor.
AcceptanceRegion * GetAcceptanceRegion(RooArgSet &, Double_t cl=-1., Double_t leftside=-1.)
Method to determine if a parameter point is in the interval.
void Add(Double_t cl, Double_t leftside)
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
virtual ~ConfidenceBelt()
Destructor.
static const float upper
Definition: main.cpp:49
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
void AddAcceptanceRegion(RooArgSet &, AcceptanceRegion region, Double_t cl=-1., Double_t leftside=-1.)
std::pair< Double_t, Double_t > AcceptanceCriteria
Persistent Reference link to a TObject A TRef is a lightweight object pointing to any TObject...
Definition: TRef.h:34
Double_t GetConfidenceLevel(Int_t index)
#define ClassDef(name, id)
Definition: Rtypes.h:254
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:33
ConfidenceBelt is a concrete implementation of the ConfInterval interface.
virtual RooArgSet * GetParameters() const
returns list of parameters
std::map< Int_t, AcceptanceCriteria > LookupTable
SamplingSummaryLookup fSamplingSummaryLookup
AcceptanceRegion(Int_t lu, Double_t ll, Double_t ul)
Bool_t CheckParameters(RooArgSet &) const
TPaveLabel title(3, 27.1, 15, 28.7,"ROOT Environment and Tools")
Int_t GetLookupIndex(Double_t cl, Double_t leftside)
std::map< Int_t, AcceptanceRegion > fAcceptanceRegions
std::vector< Double_t > ConfidenceLevels() const
AcceptanceRegion & GetAcceptanceRegion(Int_t index=0)
Double_t E()
Definition: TMath.h:54
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
This class simply holds a sampling distribution of some test statistic.
double Double_t
Definition: RtypesCore.h:55
Double_t GetLeftSideTailFraction(Int_t index)
#define name(a, b)
Definition: linkTestLib0.cpp:5
Mother of all ROOT objects.
Definition: TObject.h:58
TObject * GetObject() const
Return a pointer to the referenced object.
Definition: TRef.cxx:377
TArrow ar(9, 23, 9, 21.6, 0.015,"|>")
Double_t GetAcceptanceRegionMax(RooArgSet &, Double_t cl=-1., Double_t leftside=-1.)
std::vector< SamplingSummary > fSamplingSummaries
void AddAcceptanceRegion(AcceptanceRegion &ar)
SamplingDistribution * GetSamplingDistribution()
static const float lower
Definition: main.cpp:48
Double_t GetAcceptanceRegionMin(RooArgSet &, Double_t cl=-1., Double_t leftside=-1.)