Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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#include "RooArgSet.h"
15#include "RooAbsData.h"
17
19
20#include "TRef.h"
21
22#include <vector>
23#include <map>
24
25
26namespace RooStats {
27
28////////////////////////////////////////////////////////////////////////////////
29
31
32 typedef std::pair<double, double> AcceptanceCriteria; // defined by Confidence level, leftside tail probability
33 typedef std::map<Int_t, AcceptanceCriteria> LookupTable; // map ( Index, ( CL, leftside tail prob) )
34
35 public:
37
38 void Add(double cl, double leftside){
39 // add cl,leftside pair to lookup table
40 AcceptanceCriteria tmp(cl, leftside);
41
42 // should check to see if this is already in the map
43 if(GetLookupIndex(cl,leftside) >=0 ){
44 std::cout<< "SamplingSummaryLookup::Add, already in lookup table" << std::endl;
45 } else
46 fLookupTable[fLookupTable.size()]= tmp;
47 }
48
49 Int_t GetLookupIndex(double cl, double leftside){
50 // get index for cl,leftside pair
51 AcceptanceCriteria tmp(cl, leftside);
52
53 double tolerance = 1E-6; // some small number to protect floating point comparison. What is better way?
54 LookupTable::iterator it = fLookupTable.begin();
55 Int_t index = -1;
56 for(; it!= fLookupTable.end(); ++it) {
57 index++;
58 if( TMath::Abs( (*it).second.first - cl ) < tolerance &&
59 TMath::Abs( (*it).second.second - leftside ) < tolerance )
60 break; // exit loop, found
61 }
62
63 // check that it was found
64 if(index == (Int_t)fLookupTable.size())
65 index = -1;
66
67 return index;
68 }
69
71 if(index<0 || index>(Int_t)fLookupTable.size()) {
72 std::cout << "SamplingSummaryLookup::GetConfidenceLevel, index not in lookup table" << std::endl;
73 return -1;
74 }
75 return fLookupTable[index].first;
76 }
77
79 if(index<0 || index>(Int_t)fLookupTable.size()) {
80 std::cout << "SamplingSummaryLookup::GetLeftSideTailFraction, index not in lookup table" << std::endl;
81 return -1;
82 }
83 return fLookupTable[index].second;
84 }
85
86 private:
87 LookupTable fLookupTable; ///< map ( Index, ( CL, leftside tail prob) )
88
89 protected:
90 ClassDefOverride(SamplingSummaryLookup,1) // A simple class used by ConfidenceBelt
91 };
92
93
94 ///////////////////////////
95 class AcceptanceRegion : public TObject{
96 public:
98
99 AcceptanceRegion(Int_t lu, double ll, double ul){
100 fLookupIndex = lu;
101 fLowerLimit = ll;
102 fUpperLimit = ul;
103 }
105 double GetLowerLimit(){return fLowerLimit;}
106 double GetUpperLimit(){return fUpperLimit;}
107
108 private:
109 Int_t fLookupIndex; // want a small footprint reference to the RooArgSet for particular parameter point
110 double fLowerLimit; // lower limit on test statistic
111 double fUpperLimit; // upper limit on test statistic
112
113 protected:
114 ClassDefOverride(AcceptanceRegion,1) // A simple class for acceptance regions used for ConfidenceBelt
115
116 };
117
118////////////////////////////////////////////////////////////////////////////////
119
120 class SamplingSummary : public TObject {
121 public:
125 }
128 return (SamplingDistribution*) fSamplingDistribution.GetObject(); // dereference TRef
129 }
131
134 if( fAcceptanceRegions.count(index) !=0) {
135 std::cout << "SamplingSummary::AddAcceptanceRegion, need to implement merging protocol" << std::endl;
136 } else {
138 }
139 }
140
141 private:
142 Int_t fParameterPointIndex; // want a small footprint reference to the RooArgSet for particular parameter point
143 TRef fSamplingDistribution; // persistent pointer to a SamplingDistribution
144 std::map<Int_t, AcceptanceRegion> fAcceptanceRegions;
145
146 protected:
147 ClassDefOverride(SamplingSummary,1) // A summary of acceptance regions for confidence belt
148
149 };
150
151////////////////////////////////////////////////////////////////////////////////
152
153 class ConfidenceBelt : public TNamed {
154
155 private:
157 std::vector<SamplingSummary> fSamplingSummaries; // composite of several AcceptanceRegions
158 RooAbsData* fParameterPoints; // either a histogram (RooDataHist) or a tree (RooDataSet)
159
160
161 public:
162 // constructors,destructors
164 ConfidenceBelt(const char* name);
165 ConfidenceBelt(const char* name, const char* title);
166 ConfidenceBelt(const char* name, RooAbsData&);
167 ConfidenceBelt(const char* name, const char* title, RooAbsData&);
168 ~ConfidenceBelt() override;
169
170 /// add after creating a region
171 void AddAcceptanceRegion(RooArgSet&, AcceptanceRegion region, double cl=-1., double leftside=-1.);
172
173 /// add without creating a region, more useful for clients
174 void AddAcceptanceRegion(RooArgSet& point, Int_t dataSetIndex, double lower, double upper, double cl=-1., double leftside=-1.);
175
176 AcceptanceRegion* GetAcceptanceRegion(RooArgSet&, double cl=-1., double leftside=-1.);
177 double GetAcceptanceRegionMin(RooArgSet&, double cl=-1., double leftside=-1.);
178 double GetAcceptanceRegionMax(RooArgSet&, double cl=-1., double leftside=-1.);
179 std::vector<double> ConfidenceLevels() const ;
180
181 /// do we want it to return list of parameters
182 virtual RooArgSet* GetParameters() const;
183
184 /// check if parameters are correct. (dummy implementation to start)
185 bool CheckParameters(RooArgSet&) const ;
186
187 protected:
188 ClassDefOverride(ConfidenceBelt,1) // A confidence belt for the Neyman Construction
189
190 };
191}
192
193#endif
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
char name[80]
Definition TGX11.cxx:110
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
AcceptanceRegion(Int_t lu, double ll, double ul)
ConfidenceBelt is a concrete implementation of the ConfInterval interface.
std::vector< SamplingSummary > fSamplingSummaries
ConfidenceBelt()
Default constructor.
virtual RooArgSet * GetParameters() const
do we want it to return list of parameters
AcceptanceRegion * GetAcceptanceRegion(RooArgSet &, double cl=-1., double leftside=-1.)
Method to determine if a parameter point is in the interval.
double GetAcceptanceRegionMax(RooArgSet &, double cl=-1., double leftside=-1.)
void AddAcceptanceRegion(RooArgSet &, AcceptanceRegion region, double cl=-1., double leftside=-1.)
add after creating a region
SamplingSummaryLookup fSamplingSummaryLookup
std::vector< double > ConfidenceLevels() const
double GetAcceptanceRegionMin(RooArgSet &, double cl=-1., double leftside=-1.)
~ConfidenceBelt() override
Destructor.
bool CheckParameters(RooArgSet &) const
check if parameters are correct. (dummy implementation to start)
This class simply holds a sampling distribution of some test statistic.
void Add(double cl, double leftside)
Int_t GetLookupIndex(double cl, double leftside)
double GetLeftSideTailFraction(Int_t index)
std::map< Int_t, AcceptanceCriteria > LookupTable
std::pair< double, double > AcceptanceCriteria
double GetConfidenceLevel(Int_t index)
LookupTable fLookupTable
map ( Index, ( CL, leftside tail prob) )
std::map< Int_t, AcceptanceRegion > fAcceptanceRegions
SamplingDistribution * GetSamplingDistribution()
AcceptanceRegion & GetAcceptanceRegion(Int_t index=0)
void AddAcceptanceRegion(AcceptanceRegion &ar)
SamplingSummary(AcceptanceRegion &ar)
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
Mother of all ROOT objects.
Definition TObject.h:41
Persistent Reference link to a TObject A TRef is a lightweight object pointing to any TObject.
Definition TRef.h:32
TObject * GetObject() const
Return a pointer to the referenced object.
Definition TRef.cxx:377
Namespace for the RooStats classes.
Definition Asimov.h:19
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123