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:
38
39 void Add(double cl, double leftside){
40 // add cl,leftside pair to lookup table
41 AcceptanceCriteria tmp(cl, leftside);
42
43 // should check to see if this is already in the map
44 if(GetLookupIndex(cl,leftside) >=0 ){
45 std::cout<< "SamplingSummaryLookup::Add, already in lookup table" << std::endl;
46 } else
47 fLookupTable[fLookupTable.size()]= tmp;
48 }
49
50 Int_t GetLookupIndex(double cl, double leftside){
51 // get index for cl,leftside pair
52 AcceptanceCriteria tmp(cl, leftside);
53
54 double tolerance = 1E-6; // some small number to protect floating point comparison. What is better way?
55 LookupTable::iterator it = fLookupTable.begin();
56 Int_t index = -1;
57 for(; it!= fLookupTable.end(); ++it) {
58 index++;
59 if( TMath::Abs( (*it).second.first - cl ) < tolerance &&
60 TMath::Abs( (*it).second.second - leftside ) < tolerance )
61 break; // exit loop, found
62 }
63
64 // check that it was found
65 if(index == (Int_t)fLookupTable.size())
66 index = -1;
67
68 return index;
69 }
70
72 if(index<0 || index>(Int_t)fLookupTable.size()) {
73 std::cout << "SamplingSummaryLookup::GetConfidenceLevel, index not in lookup table" << std::endl;
74 return -1;
75 }
76 return fLookupTable[index].first;
77 }
78
80 if(index<0 || index>(Int_t)fLookupTable.size()) {
81 std::cout << "SamplingSummaryLookup::GetLeftSideTailFraction, index not in lookup table" << std::endl;
82 return -1;
83 }
84 return fLookupTable[index].second;
85 }
86
87 private:
88 LookupTable fLookupTable; ///< map ( Index, ( CL, leftside tail prob) )
89
90 protected:
91 ClassDefOverride(SamplingSummaryLookup,1) // A simple class used by ConfidenceBelt
92 };
93
94
95 ///////////////////////////
96 class AcceptanceRegion : public TObject{
97 public:
99 ~AcceptanceRegion() override {}
100
101 AcceptanceRegion(Int_t lu, double ll, double ul){
102 fLookupIndex = lu;
103 fLowerLimit = ll;
104 fUpperLimit = ul;
105 }
107 double GetLowerLimit(){return fLowerLimit;}
108 double GetUpperLimit(){return fUpperLimit;}
109
110 private:
111 Int_t fLookupIndex; // want a small footprint reference to the RooArgSet for particular parameter point
112 double fLowerLimit; // lower limit on test statistic
113 double fUpperLimit; // upper limit on test statistic
114
115 protected:
116 ClassDefOverride(AcceptanceRegion,1) // A simple class for acceptance regions used for ConfidenceBelt
117
118 };
119
120////////////////////////////////////////////////////////////////////////////////
121
122 class SamplingSummary : public TObject {
123 public:
125 ~SamplingSummary() override {}
128 }
131 return (SamplingDistribution*) fSamplingDistribution.GetObject(); // dereference TRef
132 }
134
137 if( fAcceptanceRegions.count(index) !=0) {
138 std::cout << "SamplingSummary::AddAcceptanceRegion, need to implement merging protocol" << std::endl;
139 } else {
141 }
142 }
143
144 private:
145 Int_t fParameterPointIndex; // want a small footprint reference to the RooArgSet for particular parameter point
146 TRef fSamplingDistribution; // persistent pointer to a SamplingDistribution
147 std::map<Int_t, AcceptanceRegion> fAcceptanceRegions;
148
149 protected:
150 ClassDefOverride(SamplingSummary,1) // A summary of acceptance regions for confidence belt
151
152 };
153
154////////////////////////////////////////////////////////////////////////////////
155
156 class ConfidenceBelt : public TNamed {
157
158 private:
160 std::vector<SamplingSummary> fSamplingSummaries; // composite of several AcceptanceRegions
161 RooAbsData* fParameterPoints; // either a histogram (RooDataHist) or a tree (RooDataSet)
162
163
164 public:
165 // constructors,destructors
167 ConfidenceBelt(const char* name);
168 ConfidenceBelt(const char* name, const char* title);
169 ConfidenceBelt(const char* name, RooAbsData&);
170 ConfidenceBelt(const char* name, const char* title, RooAbsData&);
171 ~ConfidenceBelt() override;
172
173 /// add after creating a region
174 void AddAcceptanceRegion(RooArgSet&, AcceptanceRegion region, double cl=-1., double leftside=-1.);
175
176 /// add without creating a region, more useful for clients
177 void AddAcceptanceRegion(RooArgSet& point, Int_t dataSetIndex, double lower, double upper, double cl=-1., double leftside=-1.);
178
179 AcceptanceRegion* GetAcceptanceRegion(RooArgSet&, double cl=-1., double leftside=-1.);
180 double GetAcceptanceRegionMin(RooArgSet&, double cl=-1., double leftside=-1.);
181 double GetAcceptanceRegionMax(RooArgSet&, double cl=-1., double leftside=-1.);
182 std::vector<double> ConfidenceLevels() const ;
183
184 /// do we want it to return list of parameters
185 virtual RooArgSet* GetParameters() const;
186
187 /// check if parameters are correct. (dummy implementation to start)
188 bool CheckParameters(RooArgSet&) const ;
189
190 protected:
191 ClassDefOverride(ConfidenceBelt,1) // A confidence belt for the Neyman Construction
192
193 };
194}
195
196#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
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:59
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