Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ConfidenceBelt.cxx
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/*****************************************************************************
12 * Project: RooStats
13 * Package: RooFit/RooStats
14 * @(#)root/roofit/roostats:$Id$
15 * Original Author: Kyle Cranmer
16 * Kyle Cranmer, Lorenzo Moneta, Gregory Schott, Wouter Verkerke
17 *
18 *****************************************************************************/
19
20/** \class RooStats::ConfidenceBelt
21 \ingroup Roostats
22
23ConfidenceBelt is a concrete implementation of the ConfInterval interface.
24It implements simple general purpose interval of arbitrary dimensions and shape.
25It does not assume the interval is connected.
26It uses either a RooDataSet (eg. a list of parameter points in the interval) or
27a RooDataHist (eg. a Histogram-like object for small regions of the parameter space) to
28store the interval.
29
30*/
31
33
34#include "RooDataSet.h"
35#include "RooDataHist.h"
36
38
39#include <TMath.h>
40
41
42using namespace RooStats;
43using std::vector;
44
45////////////////////////////////////////////////////////////////////////////////
46/// Alternate constructor
47
52
53////////////////////////////////////////////////////////////////////////////////
54/// Alternate constructor
55
56ConfidenceBelt::ConfidenceBelt(const char* name, const char* title) :
57 TNamed(name,title)
58{
59}
60
61////////////////////////////////////////////////////////////////////////////////
62/// Alternate constructor
63
65 TNamed(name,name), fParameterPoints(static_cast<RooAbsData*>(data.Clone("PointsToTestForBelt")))
66{
67}
68
69////////////////////////////////////////////////////////////////////////////////
70/// Alternate constructor
71
72ConfidenceBelt::ConfidenceBelt(const char* name, const char* title, RooAbsData& data) :
73 TNamed(name,title), fParameterPoints(static_cast<RooAbsData*>(data.Clone("PointsToTestForBelt")))
74{
75}
76
77////////////////////////////////////////////////////////////////////////////////
78
80 if(cl>0 || leftside > 0) std::cout <<"using default cl, leftside for now" << std::endl;
82 return (region) ? region->GetLowerLimit() : TMath::QuietNaN();
83}
84
85////////////////////////////////////////////////////////////////////////////////
86
88 if(cl>0 || leftside > 0) std::cout <<"using default cl, leftside for now" << std::endl;
90 return (region) ? region->GetUpperLimit() : TMath::QuietNaN();
91}
92
93////////////////////////////////////////////////////////////////////////////////
94
99
100////////////////////////////////////////////////////////////////////////////////
101
103 double lower, double upper,
104 double cl, double leftside){
105 if(cl>0 || leftside > 0) std::cout <<"using default cl, leftside for now" << std::endl;
106
107 RooDataSet* tree = dynamic_cast<RooDataSet*>( fParameterPoints );
108 RooDataHist* hist = dynamic_cast<RooDataHist*>(fParameterPoints );
109
110 // std::cout << "add: " << tree << " " << hist << std::endl;
111
112 if( !this->CheckParameters(parameterPoint) )
113 std::cout << "problem with parameters" << std::endl;
114
116 // std::cout << "lookup index = " << luIndex << std::endl;
117 if(luIndex <0 ) {
120 std::cout << "lookup index = " << luIndex << std::endl;
121 }
123
124 if( hist ) {
125 // need a way to get index for given point
126 // Can do this by setting hist's internal parameters to desired values
127 // need a better way
128 // RooStats::SetParameters(&parameterPoint, const_cast<RooArgSet*>(hist->get()));
129 // int index = hist->calcTreeIndex(); // get index
130 int index = hist->getIndex(parameterPoint); // get index
131 // std::cout << "hist index = " << index << std::endl;
132
133 // allocate memory if necessary. numEntries is overkill?
134 if((Int_t)fSamplingSummaries.size() <= index) {
135 fSamplingSummaries.reserve( hist->numEntries() );
136 fSamplingSummaries.resize( hist->numEntries() );
137 }
138
139 // set the region for this point (check for duplicate?)
141 }
142 else if( tree ){
143 // int index = tree->getIndex(parameterPoint);
144 int index = dsIndex;
145 // std::cout << "tree index = " << index << std::endl;
146
147 // allocate memory if necessary. numEntries is overkill?
148 if((Int_t)fSamplingSummaries.size() <= index){
149 fSamplingSummaries.reserve( tree->numEntries() );
150 fSamplingSummaries.resize( tree->numEntries() );
151 }
152
153 // set the region for this point (check for duplicate?)
155 }
156}
157
158////////////////////////////////////////////////////////////////////////////////
159
161 double cl, double leftside){
162 if(cl>0 || leftside > 0) std::cout <<"using default cl, leftside for now" << std::endl;
163
164 RooDataSet* tree = dynamic_cast<RooDataSet*>( fParameterPoints );
165 RooDataHist* hist = dynamic_cast<RooDataHist*>(fParameterPoints );
166
167 if( !this->CheckParameters(parameterPoint) )
168 std::cout << "problem with parameters" << std::endl;
169
170
171 if( hist ) {
172 // need a way to get index for given point
173 // Can do this by setting hist's internal parameters to desired values
174 // need a better way
175 // RooStats::SetParameters(&parameterPoint, const_cast<RooArgSet*>(hist->get()));
176 // int index = hist->calcTreeIndex(); // get index
177 int index = hist->getIndex(parameterPoint); // get index
178
179 // allocate memory if necessary. numEntries is overkill?
180 if((Int_t)fSamplingSummaries.size() < index) fSamplingSummaries.reserve( hist->numEntries() );
181
182 // set the region for this point (check for duplicate?)
184 }
185 else if( tree ){
186 tree->add( parameterPoint ); // assume it's unique for now
187 int index = tree->numEntries() - 1; //check that last point added has index nEntries -1
188 // allocate memory if necessary. numEntries is overkill?
189 if((Int_t)fSamplingSummaries.size() < index) fSamplingSummaries.reserve( tree->numEntries() );
190
191 // set the region for this point (check for duplicate?)
193 }
194}
195
196////////////////////////////////////////////////////////////////////////////////
197/// Method to determine if a parameter point is in the interval
198
200{
201 if(cl>0 || leftside > 0) std::cout <<"using default cl, leftside for now" << std::endl;
202
203 RooDataSet* tree = dynamic_cast<RooDataSet*>( fParameterPoints );
204 RooDataHist* hist = dynamic_cast<RooDataHist*>(fParameterPoints );
205
206 if( !this->CheckParameters(parameterPoint) ){
207 std::cout << "problem with parameters" << std::endl;
208 return nullptr;
209 }
210
211 if( hist ) {
212 // need a way to get index for given point
213 // Can do this by setting hist's internal parameters to desired values
214 // need a better way
215 // RooStats::SetParameters(&parameterPoint, const_cast<RooArgSet*>(hist->get()));
216 // Int_t index = hist->calcTreeIndex(); // get index
217 int index = hist->getIndex(parameterPoint); // get index
218 if (index >= (int)fSamplingSummaries.size())
219 throw std::runtime_error("ConfidenceBelt::GetAcceptanceRegion: Sampling summaries are not filled yet. Switch on NeymanConstruction::CreateConfBelt() or FeldmanCousins::CreateConfBelt().");
220
221 return &(fSamplingSummaries[index].GetAcceptanceRegion());
222 }
223 else if( tree ){
224 // need a way to get index for given point
225 // RooStats::SetParameters(&parameterPoint, tree->get()); // set tree's parameters to desired values
226 Int_t index = 0; //need something like tree->calcTreeIndex();
227 const RooArgSet* thisPoint = nullptr;
228 for(index=0; index<tree->numEntries(); ++index){
229 thisPoint = tree->get(index);
230 bool samePoint = true;
232 if (samePoint == false)
233 break;
234 if(myarg->getVal() != thisPoint->getRealValue(myarg->GetName()))
235 samePoint = false;
236 }
237 if(samePoint)
238 break;
239 }
240
241 if (index >= static_cast<int>(fSamplingSummaries.size()))
242 throw std::runtime_error("ConfidenceBelt::GetAcceptanceRegion: Sampling summaries are not filled yet. Switch on NeymanConstruction::CreateConfBelt() or FeldmanCousins::CreateConfBelt().");
243
244 return &(fSamplingSummaries[index].GetAcceptanceRegion());
245 }
246 else {
247 std::cout << "dataset is not initialized properly" << std::endl;
248 }
249
250 return nullptr;
251
252}
253
254////////////////////////////////////////////////////////////////////////////////
255/// returns list of parameters
256
258{
259 return new RooArgSet(*(fParameterPoints->get()));
260}
261
262////////////////////////////////////////////////////////////////////////////////
263
265{
266 if (parameterPoint.size() != fParameterPoints->get()->size() ) {
267 std::cout << "size is wrong, parameters don't match" << std::endl;
268 return false;
269 }
270 if ( ! parameterPoint.equals( *(fParameterPoints->get() ) ) ) {
271 std::cout << "size is ok, but parameters don't match" << std::endl;
272 return false;
273 }
274 return true;
275}
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
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
virtual const RooArgSet * get() const
Definition RooAbsData.h:101
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
Int_t getIndex(const RooAbsCollection &coord, bool fast=false) const
Calculate bin number of the given coordinates.
Container class to hold unbinned data.
Definition RooDataSet.h:34
const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
void add(const RooArgSet &row, double weight, double weightError)
Add one ore more rows of data.
std::vector< SamplingSummary > fSamplingSummaries
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.)
bool CheckParameters(RooArgSet &) const
check if parameters are correct. (dummy implementation to start)
void Add(double cl, double leftside)
Int_t GetLookupIndex(double cl, double leftside)
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
Namespace for the RooStats classes.
Definition CodegenImpl.h:58
Double_t QuietNaN()
Returns a quiet NaN as defined by IEEE 754.
Definition TMath.h:906