Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
MCMCInterval.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Authors: Kevin Belasco 17/06/2009
3// Authors: Kyle Cranmer 17/06/2009
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/** \class RooStats::MCMCInterval
13 \ingroup Roostats
14
15 MCMCInterval is a concrete implementation of the RooStats::ConfInterval
16 interface. It takes as input Markov Chain of data points in the parameter
17 space generated by Monte Carlo using the Metropolis algorithm. From the Markov
18 Chain, the confidence interval can be determined in two ways:
19
20#### Using a Kernel-Estimated PDF: (not the default method)
21
22 A RooNDKeysPdf is constructed from the data set using adaptive kernel width.
23 With this RooNDKeysPdf F, we then integrate over the most likely domain in the
24 parameter space (tallest points in the posterior RooNDKeysPdf) until the target
25 confidence level is reached within an acceptable neighborhood as defined by
26 SetEpsilon(). More specifically: we calculate the following for different
27 cutoff values C until we reach the target confidence level: \f$\int_{ F >= C } F
28 d{normset} \f$.
29 Important note: this is not the default method because of a bug in constructing
30 the RooNDKeysPdf from a weighted data set. Configure to use this method by
31 calling SetUseKeys(true), and the data set will be interpreted without weights.
32
33
34#### Using a binned data set: (the default method)
35
36 This is the binned analog of the continuous integrative method that uses the
37 kernel-estimated PDF. The points in the Markov Chain are put into a binned
38 data set and the interval is then calculated by adding the heights of the bins
39 in decreasing order until the desired level of confidence has been reached.
40 Note that this means the actual confidence level is >= the confidence level
41 prescribed by the client (unless the user calls SetHistStrict(kFALSE)). This
42 method is the default but may not remain as such in future releases, so you may
43 wish to explicitly configure to use this method by calling SetUseKeys(false)
44
45
46 These are not the only ways for the confidence interval to be determined, and
47 other possibilities are being considered being added, especially for the
48 1-dimensional case.
49
50
51 One can ask an MCMCInterval for the lower and upper limits on a specific
52 parameter of interest in the interval. Note that this works better for some
53 distributions (ones with exactly one local maximum) than others, and sometimes
54 has little value.
55*/
56
57
58#include "Rtypes.h"
59
60#include "TMath.h"
61
64#include "RooStats/Heaviside.h"
65#include "RooDataHist.h"
66#include "RooNDKeysPdf.h"
67#include "RooProduct.h"
69#include "RooRealVar.h"
70#include "RooArgList.h"
71#include "TIterator.h"
72#include "TH1.h"
73#include "TH1F.h"
74#include "TH2F.h"
75#include "TH3F.h"
76#include "RooMsgService.h"
77#include "RooGlobalFunc.h"
78#include "TObject.h"
79#include "THnSparse.h"
80#include "RooNumber.h"
81
82#include <cstdlib>
83#include <string>
84#include <algorithm>
85
87
88using namespace RooFit;
89using namespace RooStats;
90using namespace std;
91
92static const Double_t DEFAULT_EPSILON = 0.01;
93static const Double_t DEFAULT_DELTA = 10e-6;
94
95////////////////////////////////////////////////////////////////////////////////
96
99{
100 fConfidenceLevel = 0.0;
101 fHistConfLevel = 0.0;
102 fKeysConfLevel = 0.0;
103 fTFConfLevel = 0.0;
104 fFull = 0.0;
105 fChain = NULL;
106 fAxes = NULL;
107 fDataHist = NULL;
108 fSparseHist = NULL;
109 fVector.clear();
110 fKeysPdf = NULL;
111 fProduct = NULL;
112 fHeaviside = NULL;
113 fKeysDataHist = NULL;
114 fCutoffVar = NULL;
115 fHist = NULL;
116 fNumBurnInSteps = 0;
117 fHistCutoff = -1;
118 fKeysCutoff = -1;
119 fDimension = 1;
126 fTFLower = -1.0 * RooNumber::infinity();
128 fVecWeight = 0;
129 fLeftSideTF = -1;
130}
131
132////////////////////////////////////////////////////////////////////////////////
133
135 const RooArgSet& parameters, MarkovChain& chain) : ConfInterval(name)
136{
137 fNumBurnInSteps = 0;
138 fConfidenceLevel = 0.0;
139 fHistConfLevel = 0.0;
140 fKeysConfLevel = 0.0;
141 fTFConfLevel = 0.0;
142 fFull = 0.0;
143 fAxes = NULL;
144 fChain = &chain;
145 fDataHist = NULL;
146 fSparseHist = NULL;
147 fVector.clear();
148 fKeysPdf = NULL;
149 fProduct = NULL;
150 fHeaviside = NULL;
151 fKeysDataHist = NULL;
152 fCutoffVar = NULL;
153 fHist = NULL;
154 fHistCutoff = -1;
155 fKeysCutoff = -1;
160 SetParameters(parameters);
163 fTFLower = -1.0 * RooNumber::infinity();
165 fVecWeight = 0;
166 fLeftSideTF = -1;
167}
168
169////////////////////////////////////////////////////////////////////////////////
170
172{
173 // destructor
174 delete[] fAxes;
175 delete fHist;
176 delete fChain;
177 // kbelasco: check here for memory management errors
178 delete fDataHist;
179 delete fSparseHist;
180 delete fKeysPdf;
181 delete fProduct;
182 delete fHeaviside;
183 delete fKeysDataHist;
184 delete fCutoffVar;
185}
186
187struct CompareDataHistBins {
188 CompareDataHistBins(RooDataHist* hist) : fDataHist(hist) {}
189 bool operator() (Int_t bin1 , Int_t bin2) {
190 fDataHist->get(bin1);
191 Double_t n1 = fDataHist->weight();
192 fDataHist->get(bin2);
193 Double_t n2 = fDataHist->weight();
194 return (n1 < n2);
195 }
196 RooDataHist* fDataHist;
197};
198
199struct CompareSparseHistBins {
200 CompareSparseHistBins(THnSparse* hist) : fSparseHist(hist) {}
201 bool operator() (Long_t bin1, Long_t bin2) {
202 Double_t n1 = fSparseHist->GetBinContent(bin1);
203 Double_t n2 = fSparseHist->GetBinContent(bin2);
204 return (n1 < n2);
205 }
206 THnSparse* fSparseHist;
207};
208
209struct CompareVectorIndices {
210 CompareVectorIndices(MarkovChain* chain, RooRealVar* param) :
211 fChain(chain), fParam(param) {}
212 bool operator() (Int_t i, Int_t j) {
213 Double_t xi = fChain->Get(i)->getRealValue(fParam->GetName());
214 Double_t xj = fChain->Get(j)->getRealValue(fParam->GetName());
215 return (xi < xj);
216 }
217 MarkovChain* fChain;
218 RooRealVar* fParam;
219};
220
221////////////////////////////////////////////////////////////////////////////////
222/// kbelasco: for this method, consider running DetermineInterval() if
223/// fKeysPdf==NULL, fSparseHist==NULL, fDataHist==NULL, or fVector.size()==0
224/// rather than just returning false. Though this should not be an issue
225/// because nobody should be able to get an MCMCInterval that has their interval
226/// or posterior representation NULL/empty since they should only get this
227/// through the MCMCCalculator
228
230{
231 if (fIntervalType == kShortest) {
232 if (fUseKeys) {
233 if (fKeysPdf == NULL)
234 return false;
235
236 // evaluate keyspdf at point and return whether >= cutoff
237 RooStats::SetParameters(&point, const_cast<RooArgSet *>(&fParameters));
239 } else {
240 if (fUseSparseHist) {
241 if (fSparseHist == NULL)
242 return false;
243
244 // evaluate sparse hist at bin where point lies and return
245 // whether >= cutoff
247 const_cast<RooArgSet*>(&fParameters));
248 Long_t bin;
249 // kbelasco: consider making x static
251 for (Int_t i = 0; i < fDimension; i++)
252 x[i] = fAxes[i]->getVal();
253 bin = fSparseHist->GetBin(x, kFALSE);
255 delete[] x;
256 return (weight >= (Double_t)fHistCutoff);
257 } else {
258 if (fDataHist == NULL)
259 return false;
260
261 // evaluate data hist at bin where point lies and return whether
262 // >= cutoff
263 Int_t bin;
264 bin = fDataHist->getIndex(point);
265 fDataHist->get(bin);
266 return (fDataHist->weight() >= (Double_t)fHistCutoff);
267 }
268 }
269 } else if (fIntervalType == kTailFraction) {
270 if (fVector.size() == 0)
271 return false;
272
273 // return whether value of point is within the range
274 Double_t x = point.getRealValue(fAxes[0]->GetName());
275 if (fTFLower <= x && x <= fTFUpper)
276 return true;
277
278 return false;
279 }
280
281 coutE(InputArguments) << "Error in MCMCInterval::IsInInterval: "
282 << "Interval type not set. Returning false." << endl;
283 return false;
284}
285
286////////////////////////////////////////////////////////////////////////////////
287
289{
290 fConfidenceLevel = cl;
292}
293
294// kbelasco: update this or just take it out
295// kbelasco: consider keeping this around but changing the implementation
296// to set the number of bins for each RooRealVar and then recreating the
297// histograms
298//void MCMCInterval::SetNumBins(Int_t numBins)
299//{
300// if (numBins > 0) {
301// fPreferredNumBins = numBins;
302// for (Int_t d = 0; d < fDimension; d++)
303// fNumBins[d] = numBins;
304// }
305// else {
306// coutE(Eval) << "* Error in MCMCInterval::SetNumBins: " <<
307// "Negative number of bins given: " << numBins << endl;
308// return;
309// }
310//
311// // If the histogram already exists, recreate it with the new bin numbers
312// if (fHist != NULL)
313// CreateHist();
314//}
315
316////////////////////////////////////////////////////////////////////////////////
317
319{
320 Int_t size = axes.getSize();
321 if (size != fDimension) {
322 coutE(InputArguments) << "* Error in MCMCInterval::SetAxes: " <<
323 "number of variables in axes (" << size <<
324 ") doesn't match number of parameters ("
325 << fDimension << ")" << endl;
326 return;
327 }
328 for (Int_t i = 0; i < size; i++)
329 fAxes[i] = (RooRealVar*)axes.at(i);
330}
331
332////////////////////////////////////////////////////////////////////////////////
333
335{
336 // kbelasco: check here for memory leak. does RooNDKeysPdf use
337 // the RooArgList passed to it or does it make a clone?
338 // also check for memory leak from chain, does RooNDKeysPdf clone that?
339 if (fAxes == NULL || fParameters.getSize() == 0) {
340 coutE(InputArguments) << "Error in MCMCInterval::CreateKeysPdf: "
341 << "parameters have not been set." << endl;
342 return;
343 }
344
345 if (fNumBurnInSteps >= fChain->Size()) {
347 "MCMCInterval::CreateKeysPdf: creation of Keys PDF failed: " <<
348 "Number of burn-in steps (num steps to ignore) >= number of steps " <<
349 "in Markov chain." << endl;
350 delete fKeysPdf;
351 delete fCutoffVar;
352 delete fHeaviside;
353 delete fProduct;
354 fKeysPdf = NULL;
355 fCutoffVar = NULL;
356 fHeaviside = NULL;
357 fProduct = NULL;
358 return;
359 }
360
363 RooArgList* paramsList = new RooArgList();
364 for (Int_t i = 0; i < fDimension; i++)
365 paramsList->add(*fAxes[i]);
366
367 // kbelasco: check for memory leaks with chain. who owns it? does
368 // RooNDKeysPdf take ownership?
369 fKeysPdf = new RooNDKeysPdf("keysPDF", "Keys PDF", *paramsList, *chain, "a");
370 fCutoffVar = new RooRealVar("cutoff", "cutoff", 0);
371 fHeaviside = new Heaviside("heaviside", "Heaviside", *fKeysPdf, *fCutoffVar);
372 fProduct = new RooProduct("product", "Keys PDF & Heaviside Product",
374}
375
376////////////////////////////////////////////////////////////////////////////////
377
379{
380 if (fAxes == NULL || fChain == NULL) {
381 coutE(Eval) << "* Error in MCMCInterval::CreateHist(): " <<
382 "Crucial data member was NULL." << endl;
383 coutE(Eval) << "Make sure to fully construct/initialize." << endl;
384 return;
385 }
386 if (fHist != NULL)
387 delete fHist;
388
389 if (fNumBurnInSteps >= fChain->Size()) {
391 "MCMCInterval::CreateHist: creation of histogram failed: " <<
392 "Number of burn-in steps (num steps to ignore) >= number of steps " <<
393 "in Markov chain." << endl;
394 fHist = NULL;
395 return;
396 }
397
398 if (fDimension == 1)
399 fHist = new TH1F("posterior", "MCMC Posterior Histogram",
400 fAxes[0]->numBins(), fAxes[0]->getMin(), fAxes[0]->getMax());
401
402 else if (fDimension == 2)
403 fHist = new TH2F("posterior", "MCMC Posterior Histogram",
404 fAxes[0]->numBins(), fAxes[0]->getMin(), fAxes[0]->getMax(),
405 fAxes[1]->numBins(), fAxes[1]->getMin(), fAxes[1]->getMax());
406
407 else if (fDimension == 3)
408 fHist = new TH3F("posterior", "MCMC Posterior Histogram",
409 fAxes[0]->numBins(), fAxes[0]->getMin(), fAxes[0]->getMax(),
410 fAxes[1]->numBins(), fAxes[1]->getMin(), fAxes[1]->getMax(),
411 fAxes[2]->numBins(), fAxes[2]->getMin(), fAxes[2]->getMax());
412
413 else {
414 coutE(Eval) << "* Error in MCMCInterval::CreateHist() : " <<
415 "TH1* couldn't handle dimension: " << fDimension << endl;
416 return;
417 }
418
419 // Fill histogram
420 Int_t size = fChain->Size();
421 const RooArgSet* entry;
422 for (Int_t i = fNumBurnInSteps; i < size; i++) {
423 entry = fChain->Get(i);
424 if (fDimension == 1)
425 ((TH1F*)fHist)->Fill(entry->getRealValue(fAxes[0]->GetName()),
426 fChain->Weight());
427 else if (fDimension == 2)
428 ((TH2F*)fHist)->Fill(entry->getRealValue(fAxes[0]->GetName()),
429 entry->getRealValue(fAxes[1]->GetName()),
430 fChain->Weight());
431 else
432 ((TH3F*)fHist)->Fill(entry->getRealValue(fAxes[0]->GetName()),
433 entry->getRealValue(fAxes[1]->GetName()),
434 entry->getRealValue(fAxes[2]->GetName()),
435 fChain->Weight());
436 }
437
438 if (fDimension >= 1)
440 if (fDimension >= 2)
442 if (fDimension >= 3)
444}
445
446////////////////////////////////////////////////////////////////////////////////
447
449{
450 if (fAxes == NULL || fChain == NULL) {
451 coutE(InputArguments) << "* Error in MCMCInterval::CreateSparseHist(): "
452 << "Crucial data member was NULL." << endl;
453 coutE(InputArguments) << "Make sure to fully construct/initialize."
454 << endl;
455 return;
456 }
457 if (fSparseHist != NULL)
458 delete fSparseHist;
459
460 Double_t* min = new Double_t[fDimension];
461 Double_t* max = new Double_t[fDimension];
462 Int_t* bins = new Int_t[fDimension];
463 for (Int_t i = 0; i < fDimension; i++) {
464 min[i] = fAxes[i]->getMin();
465 max[i] = fAxes[i]->getMax();
466 bins[i] = fAxes[i]->numBins();
467 }
468 fSparseHist = new THnSparseF("posterior", "MCMC Posterior Histogram",
469 fDimension, bins, min, max);
470
471 delete[] min;
472 delete[] max;
473 delete[] bins;
474
475 // kbelasco: it appears we need to call Sumw2() just to get the
476 // histogram to keep a running total of the weight so that Getsumw doesn't
477 // just return 0
479
480 if (fNumBurnInSteps >= fChain->Size()) {
482 "MCMCInterval::CreateSparseHist: creation of histogram failed: " <<
483 "Number of burn-in steps (num steps to ignore) >= number of steps " <<
484 "in Markov chain." << endl;
485 }
486
487 // Fill histogram
488 Int_t size = fChain->Size();
489 const RooArgSet* entry;
491 for (Int_t i = fNumBurnInSteps; i < size; i++) {
492 entry = fChain->Get(i);
493 for (Int_t ii = 0; ii < fDimension; ii++)
494 x[ii] = entry->getRealValue(fAxes[ii]->GetName());
496 }
497 delete[] x;
498}
499
500////////////////////////////////////////////////////////////////////////////////
501
503{
504 if (fParameters.getSize() == 0 || fChain == NULL) {
505 coutE(Eval) << "* Error in MCMCInterval::CreateDataHist(): " <<
506 "Crucial data member was NULL or empty." << endl;
507 coutE(Eval) << "Make sure to fully construct/initialize." << endl;
508 return;
509 }
510
511 if (fNumBurnInSteps >= fChain->Size()) {
513 "MCMCInterval::CreateDataHist: creation of histogram failed: " <<
514 "Number of burn-in steps (num steps to ignore) >= number of steps " <<
515 "in Markov chain." << endl;
516 fDataHist = NULL;
517 return;
518 }
519
522}
523
524////////////////////////////////////////////////////////////////////////////////
525
527{
528 fVector.clear();
529 fVecWeight = 0;
530
531 if (fChain == NULL) {
532 coutE(InputArguments) << "* Error in MCMCInterval::CreateVector(): " <<
533 "Crucial data member (Markov chain) was NULL." << endl;
534 coutE(InputArguments) << "Make sure to fully construct/initialize."
535 << endl;
536 return;
537 }
538
539 if (fNumBurnInSteps >= fChain->Size()) {
541 "MCMCInterval::CreateVector: creation of vector failed: " <<
542 "Number of burn-in steps (num steps to ignore) >= number of steps " <<
543 "in Markov chain." << endl;
544 }
545
546 // Fill vector
547 Int_t size = fChain->Size() - fNumBurnInSteps;
548 fVector.resize(size);
549 Int_t i;
550 Int_t chainIndex;
551 for (i = 0; i < size; i++) {
552 chainIndex = i + fNumBurnInSteps;
553 fVector[i] = chainIndex;
554 fVecWeight += fChain->Weight(chainIndex);
555 }
556
557 stable_sort(fVector.begin(), fVector.end(),
558 CompareVectorIndices(fChain, param));
559}
560
561////////////////////////////////////////////////////////////////////////////////
562
564{
566 fParameters.add(parameters);
568 if (fAxes != NULL)
569 delete[] fAxes;
572 Int_t n = 0;
573 TObject* obj;
574 while ((obj = it->Next()) != NULL) {
575 if (dynamic_cast<RooRealVar*>(obj) != NULL)
576 fAxes[n] = (RooRealVar*)obj;
577 else
578 coutE(Eval) << "* Error in MCMCInterval::SetParameters: " <<
579 obj->GetName() << " not a RooRealVar*" << std::endl;
580 n++;
581 }
582 delete it;
583}
584
585////////////////////////////////////////////////////////////////////////////////
586
588{
589 switch (fIntervalType) {
590 case kShortest:
592 break;
593 case kTailFraction:
595 break;
596 default:
597 coutE(InputArguments) << "MCMCInterval::DetermineInterval(): " <<
598 "Error: Interval type not set" << endl;
599 break;
600 }
601}
602
603////////////////////////////////////////////////////////////////////////////////
604
606{
607 if (fUseKeys)
609 else
611}
612
613////////////////////////////////////////////////////////////////////////////////
614
616{
617 if (fLeftSideTF < 0 || fLeftSideTF > 1) {
618 coutE(InputArguments) << "MCMCInterval::DetermineTailFractionInterval: "
619 << "Fraction must be in the range [0, 1]. "
620 << fLeftSideTF << "is not allowed." << endl;
621 return;
622 }
623
624 if (fDimension != 1) {
625 coutE(InputArguments) << "MCMCInterval::DetermineTailFractionInterval(): "
626 << "Error: Can only find a tail-fraction interval for 1-D intervals"
627 << endl;
628 return;
629 }
630
631 if (fAxes == NULL) {
632 coutE(InputArguments) << "MCMCInterval::DetermineTailFractionInterval(): "
633 << "Crucial data member was NULL." << endl;
634 coutE(InputArguments) << "Make sure to fully construct/initialize."
635 << endl;
636 return;
637 }
638
639 // kbelasco: fill in code here to find interval
640 //
641 // also make changes so that calling GetPosterior...() returns NULL
642 // when fIntervalType == kTailFraction, since there really
643 // is no posterior for this type of interval determination
644 if (fVector.size() == 0)
646
647 if (fVector.size() == 0 || fVecWeight == 0) {
648 // if size is still 0, then creation failed.
649 // if fVecWeight == 0, then there are no entries (indicates the same
650 // error as fVector.size() == 0 because that only happens when
651 // fNumBurnInSteps >= fChain->Size())
652 // either way, reset and return
653 fVector.clear();
654 fTFLower = -1.0 * RooNumber::infinity();
656 fTFConfLevel = 0.0;
657 fVecWeight = 0;
658 return;
659 }
660
661 RooRealVar* param = fAxes[0];
662
664 Double_t leftTailCutoff = fVecWeight * (1 - c) * fLeftSideTF;
665 Double_t rightTailCutoff = fVecWeight * (1 - c) * (1 - fLeftSideTF);
666 Double_t leftTailSum = 0;
667 Double_t rightTailSum = 0;
668
669 // kbelasco: consider changing these values to +infinity and -infinity
670 Double_t ll = param->getMin();
671 Double_t ul = param->getMax();
672
673 Double_t x;
674 Double_t w;
675
676 // save a lot of GetName() calls if compiler does not already optimize this
677 const char* name = param->GetName();
678
679 // find lower limit
680 Int_t i;
681 for (i = 0; i < (Int_t)fVector.size(); i++) {
683 w = fChain->Weight();
684
685 if (TMath::Abs(leftTailSum + w - leftTailCutoff) <
686 TMath::Abs(leftTailSum - leftTailCutoff)) {
687 // moving the lower limit to x would bring us closer to the desired
688 // left tail size
689 ll = x;
690 leftTailSum += w;
691 } else
692 break;
693 }
694
695 // find upper limit
696 for (i = (Int_t)fVector.size() - 1; i >= 0; i--) {
698 w = fChain->Weight();
699
700 if (TMath::Abs(rightTailSum + w - rightTailCutoff) <
701 TMath::Abs(rightTailSum - rightTailCutoff)) {
702 // moving the lower limit to x would bring us closer to the desired
703 // left tail size
704 ul = x;
705 rightTailSum += w;
706 } else
707 break;
708 }
709
710 fTFLower = ll;
711 fTFUpper = ul;
712 fTFConfLevel = 1 - (leftTailSum + rightTailSum) / fVecWeight;
713}
714
715////////////////////////////////////////////////////////////////////////////////
716
718{
719 if (fKeysPdf == NULL)
721
722 if (fKeysPdf == NULL) {
723 // if fKeysPdf is still NULL, then it means CreateKeysPdf failed
724 // so clear all the data members this function would normally determine
725 // and return
726 fFull = 0.0;
727 fKeysCutoff = -1;
728 fKeysConfLevel = 0.0;
729 return;
730 }
731
732 // now we have a keys pdf of the posterior
733
734 Double_t cutoff = 0.0;
735 fCutoffVar->setVal(cutoff);
737 Double_t full = integral->getVal(fParameters);
738 fFull = full;
739 delete integral;
740 if (full < 0.98) {
741 coutW(Eval) << "Warning: Integral of Keys PDF came out to " << full
742 << " instead of expected value 1. Will continue using this "
743 << "factor to normalize further integrals of this PDF." << endl;
744 }
745
746 // kbelasco: Is there a better way to set the search range?
747 // from 0 to max value of Keys
748 // kbelasco: how to get max value?
749 //Double_t max = product.maxVal(product.getMaxVal(fParameters));
750
751 Double_t volume = 1.0;
753 RooRealVar* var;
754 while ((var = (RooRealVar*)it->Next()) != NULL)
755 volume *= (var->getMax() - var->getMin());
756 delete it;
757
758 Double_t topCutoff = full / volume;
759 Double_t bottomCutoff = topCutoff;
760 Double_t confLevel = CalcConfLevel(topCutoff, full);
761 if (AcceptableConfLevel(confLevel)) {
762 fKeysConfLevel = confLevel;
763 fKeysCutoff = topCutoff;
764 return;
765 }
766 Bool_t changed = kFALSE;
767 // find high end of range
768 while (confLevel > fConfidenceLevel) {
769 topCutoff *= 2.0;
770 confLevel = CalcConfLevel(topCutoff, full);
771 if (AcceptableConfLevel(confLevel)) {
772 fKeysConfLevel = confLevel;
773 fKeysCutoff = topCutoff;
774 return;
775 }
776 changed = kTRUE;
777 }
778 if (changed) {
779 bottomCutoff = topCutoff / 2.0;
780 } else {
781 changed = kFALSE;
782 bottomCutoff /= 2.0;
783 confLevel = CalcConfLevel(bottomCutoff, full);
784 if (AcceptableConfLevel(confLevel)) {
785 fKeysConfLevel = confLevel;
786 fKeysCutoff = bottomCutoff;
787 return;
788 }
789 while (confLevel < fConfidenceLevel) {
790 bottomCutoff /= 2.0;
791 confLevel = CalcConfLevel(bottomCutoff, full);
792 if (AcceptableConfLevel(confLevel)) {
793 fKeysConfLevel = confLevel;
794 fKeysCutoff = bottomCutoff;
795 return;
796 }
797 changed = kTRUE;
798 }
799 if (changed) {
800 topCutoff = bottomCutoff * 2.0;
801 }
802 }
803
804 coutI(Eval) << "range set: [" << bottomCutoff << ", " << topCutoff << "]"
805 << endl;
806
807 cutoff = (topCutoff + bottomCutoff) / 2.0;
808 confLevel = CalcConfLevel(cutoff, full);
809
810 // need to use WithinDeltaFraction() because sometimes the integrating the
811 // posterior in this binary search seems to not have enough granularity to
812 // find an acceptable conf level (small no. of strange cases).
813 // WithinDeltaFraction causes the search to terminate when
814 // topCutoff is essentially equal to bottomCutoff (compared to the magnitude
815 // of their mean).
816 while (!AcceptableConfLevel(confLevel) &&
817 !WithinDeltaFraction(topCutoff, bottomCutoff)) {
818 if (confLevel > fConfidenceLevel)
819 bottomCutoff = cutoff;
820 else if (confLevel < fConfidenceLevel)
821 topCutoff = cutoff;
822 cutoff = (topCutoff + bottomCutoff) / 2.0;
823 coutI(Eval) << "cutoff range: [" << bottomCutoff << ", "
824 << topCutoff << "]" << endl;
825 confLevel = CalcConfLevel(cutoff, full);
826 }
827
828 fKeysCutoff = cutoff;
829 fKeysConfLevel = confLevel;
830}
831
832////////////////////////////////////////////////////////////////////////////////
833
835{
836 if (fUseSparseHist)
838 else
840}
841
842////////////////////////////////////////////////////////////////////////////////
843
845{
846 Long_t numBins;
847 if (fSparseHist == NULL)
849
850 if (fSparseHist == NULL) {
851 // if fSparseHist is still NULL, then CreateSparseHist failed
852 fHistCutoff = -1;
853 fHistConfLevel = 0.0;
854 return;
855 }
856
857 numBins = (Long_t)fSparseHist->GetNbins();
858
859 std::vector<Long_t> bins(numBins);
860 for (Int_t ibin = 0; ibin < numBins; ibin++)
861 bins[ibin] = (Long_t)ibin;
862 std::stable_sort(bins.begin(), bins.end(), CompareSparseHistBins(fSparseHist));
863
864 Double_t nEntries = fSparseHist->GetSumw();
865 Double_t sum = 0;
866 Double_t content;
867 Int_t i;
868 // see above note on indexing to understand numBins - 3
869 for (i = numBins - 1; i >= 0; i--) {
870 content = fSparseHist->GetBinContent(bins[i]);
871 if ((sum + content) / nEntries >= fConfidenceLevel) {
872 fHistCutoff = content;
873 if (fIsHistStrict) {
874 sum += content;
875 i--;
876 break;
877 } else {
878 i++;
879 break;
880 }
881 }
882 sum += content;
883 }
884
885 if (fIsHistStrict) {
886 // keep going to find the sum
887 for ( ; i >= 0; i--) {
888 content = fSparseHist->GetBinContent(bins[i]);
889 if (content == fHistCutoff)
890 sum += content;
891 else
892 break; // content must be < fHistCutoff
893 }
894 } else {
895 // backtrack to find the cutoff and sum
896 for ( ; i < numBins; i++) {
897 content = fSparseHist->GetBinContent(bins[i]);
898 if (content > fHistCutoff) {
899 fHistCutoff = content;
900 break;
901 } else // content == fHistCutoff
902 sum -= content;
903 if (i == numBins - 1)
904 // still haven't set fHistCutoff correctly yet, and we have no bins
905 // left, so set fHistCutoff to something higher than the tallest bin
906 fHistCutoff = content + 1.0;
907 }
908 }
909
910 fHistConfLevel = sum / nEntries;
911}
912
913////////////////////////////////////////////////////////////////////////////////
914
916{
917 Int_t numBins;
918 if (fDataHist == NULL)
920 if (fDataHist == NULL) {
921 // if fDataHist is still NULL, then CreateDataHist failed
922 fHistCutoff = -1;
923 fHistConfLevel = 0.0;
924 return;
925 }
926
927 numBins = fDataHist->numEntries();
928
929 std::vector<Int_t> bins(numBins);
930 for (Int_t ibin = 0; ibin < numBins; ibin++)
931 bins[ibin] = ibin;
932 std::stable_sort(bins.begin(), bins.end(), CompareDataHistBins(fDataHist));
933
934 Double_t nEntries = fDataHist->sum(kFALSE);
935 Double_t sum = 0;
936 Double_t content;
937 Int_t i;
938 for (i = numBins - 1; i >= 0; i--) {
939 fDataHist->get(bins[i]);
940 content = fDataHist->weight();
941 if ((sum + content) / nEntries >= fConfidenceLevel) {
942 fHistCutoff = content;
943 if (fIsHistStrict) {
944 sum += content;
945 i--;
946 break;
947 } else {
948 i++;
949 break;
950 }
951 }
952 sum += content;
953 }
954
955 if (fIsHistStrict) {
956 // keep going to find the sum
957 for ( ; i >= 0; i--) {
958 fDataHist->get(bins[i]);
959 content = fDataHist->weight();
960 if (content == fHistCutoff)
961 sum += content;
962 else
963 break; // content must be < fHistCutoff
964 }
965 } else {
966 // backtrack to find the cutoff and sum
967 for ( ; i < numBins; i++) {
968 fDataHist->get(bins[i]);
969 content = fDataHist->weight();
970 if (content > fHistCutoff) {
971 fHistCutoff = content;
972 break;
973 } else // content == fHistCutoff
974 sum -= content;
975 if (i == numBins - 1)
976 // still haven't set fHistCutoff correctly yet, and we have no bins
977 // left, so set fHistCutoff to something higher than the tallest bin
978 fHistCutoff = content + 1.0;
979 }
980 }
981
982 fHistConfLevel = sum / nEntries;
983}
984
985////////////////////////////////////////////////////////////////////////////////
986
988{
989 if (fIntervalType == kShortest) {
990 if (fUseKeys)
991 return fKeysConfLevel;
992 else
993 return fHistConfLevel;
994 } else if (fIntervalType == kTailFraction) {
995 return fTFConfLevel;
996 } else {
997 coutE(InputArguments) << "MCMCInterval::GetActualConfidenceLevel: "
998 << "not implemented for this type of interval. Returning 0." << endl;
999 return 0;
1000 }
1001}
1002
1003////////////////////////////////////////////////////////////////////////////////
1004
1006{
1007 switch (fIntervalType) {
1008 case kShortest:
1009 return LowerLimitShortest(param);
1010 case kTailFraction:
1011 return LowerLimitTailFraction(param);
1012 default:
1013 coutE(InputArguments) << "MCMCInterval::LowerLimit(): " <<
1014 "Error: Interval type not set" << endl;
1015 return RooNumber::infinity();
1016 }
1017}
1018
1019////////////////////////////////////////////////////////////////////////////////
1020
1022{
1023 switch (fIntervalType) {
1024 case kShortest:
1025 return UpperLimitShortest(param);
1026 case kTailFraction:
1027 return UpperLimitTailFraction(param);
1028 default:
1029 coutE(InputArguments) << "MCMCInterval::UpperLimit(): " <<
1030 "Error: Interval type not set" << endl;
1031 return RooNumber::infinity();
1032 }
1033}
1034
1035////////////////////////////////////////////////////////////////////////////////
1036
1038{
1039 if (fTFLower == -1.0 * RooNumber::infinity())
1041
1042 return fTFLower;
1043}
1044
1045////////////////////////////////////////////////////////////////////////////////
1046
1048{
1051
1052 return fTFUpper;
1053}
1054
1055////////////////////////////////////////////////////////////////////////////////
1056
1058{
1059 if (fUseKeys)
1060 return LowerLimitByKeys(param);
1061 else
1062 return LowerLimitByHist(param);
1063}
1064
1065////////////////////////////////////////////////////////////////////////////////
1066
1068{
1069 if (fUseKeys)
1070 return UpperLimitByKeys(param);
1071 else
1072 return UpperLimitByHist(param);
1073}
1074
1075////////////////////////////////////////////////////////////////////////////////
1076/// Determine the lower limit for param on this interval
1077/// using the binned data set
1078
1080{
1081 if (fUseSparseHist)
1082 return LowerLimitBySparseHist(param);
1083 else
1084 return LowerLimitByDataHist(param);
1085}
1086
1087////////////////////////////////////////////////////////////////////////////////
1088/// Determine the upper limit for param on this interval
1089/// using the binned data set
1090
1092{
1093 if (fUseSparseHist)
1094 return UpperLimitBySparseHist(param);
1095 else
1096 return UpperLimitByDataHist(param);
1097}
1098
1099////////////////////////////////////////////////////////////////////////////////
1100/// Determine the lower limit for param on this interval
1101/// using the binned data set
1102
1104{
1105 if (fDimension != 1) {
1106 coutE(InputArguments) << "In MCMCInterval::LowerLimitBySparseHist: "
1107 << "Sorry, will not compute lower limit unless dimension == 1" << endl;
1108 return param.getMin();
1109 }
1110 if (fHistCutoff < 0)
1111 DetermineBySparseHist(); // this initializes fSparseHist
1112
1113 if (fHistCutoff < 0) {
1114 // if fHistCutoff < 0 still, then determination of interval failed
1115 coutE(Eval) << "In MCMCInterval::LowerLimitBySparseHist: "
1116 << "couldn't determine cutoff. Check that num burn in steps < num "
1117 << "steps in the Markov chain. Returning param.getMin()." << endl;
1118 return param.getMin();
1119 }
1120
1121 std::vector<Int_t> coord(fDimension);
1122 for (Int_t d = 0; d < fDimension; d++) {
1123 if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
1124 Long_t numBins = (Long_t)fSparseHist->GetNbins();
1125 Double_t lowerLimit = param.getMax();
1126 Double_t val;
1127 for (Long_t i = 0; i < numBins; i++) {
1128 if (fSparseHist->GetBinContent(i, &coord[0]) >= fHistCutoff) {
1129 val = fSparseHist->GetAxis(d)->GetBinCenter(coord[d]);
1130 if (val < lowerLimit)
1131 lowerLimit = val;
1132 }
1133 }
1134 return lowerLimit;
1135 }
1136 }
1137 return param.getMin();
1138}
1139
1140////////////////////////////////////////////////////////////////////////////////
1141/// Determine the lower limit for param on this interval
1142/// using the binned data set
1143
1145{
1146 if (fHistCutoff < 0)
1147 DetermineByDataHist(); // this initializes fDataHist
1148
1149 if (fHistCutoff < 0) {
1150 // if fHistCutoff < 0 still, then determination of interval failed
1151 coutE(Eval) << "In MCMCInterval::LowerLimitByDataHist: "
1152 << "couldn't determine cutoff. Check that num burn in steps < num "
1153 << "steps in the Markov chain. Returning param.getMin()." << endl;
1154 return param.getMin();
1155 }
1156
1157 for (Int_t d = 0; d < fDimension; d++) {
1158 if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
1159 Int_t numBins = fDataHist->numEntries();
1160 Double_t lowerLimit = param.getMax();
1161 Double_t val;
1162 for (Int_t i = 0; i < numBins; i++) {
1163 fDataHist->get(i);
1164 if (fDataHist->weight() >= fHistCutoff) {
1165 val = fDataHist->get()->getRealValue(param.GetName());
1166 if (val < lowerLimit)
1167 lowerLimit = val;
1168 }
1169 }
1170 return lowerLimit;
1171 }
1172 }
1173 return param.getMin();
1174}
1175
1176////////////////////////////////////////////////////////////////////////////////
1177/// Determine the upper limit for param on this interval
1178/// using the binned data set
1179
1181{
1182 if (fDimension != 1) {
1183 coutE(InputArguments) << "In MCMCInterval::UpperLimitBySparseHist: "
1184 << "Sorry, will not compute upper limit unless dimension == 1" << endl;
1185 return param.getMax();
1186 }
1187 if (fHistCutoff < 0)
1188 DetermineBySparseHist(); // this initializes fSparseHist
1189
1190 if (fHistCutoff < 0) {
1191 // if fHistCutoff < 0 still, then determination of interval failed
1192 coutE(Eval) << "In MCMCInterval::UpperLimitBySparseHist: "
1193 << "couldn't determine cutoff. Check that num burn in steps < num "
1194 << "steps in the Markov chain. Returning param.getMax()." << endl;
1195 return param.getMax();
1196 }
1197
1198 std::vector<Int_t> coord(fDimension);
1199 for (Int_t d = 0; d < fDimension; d++) {
1200 if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
1201 Long_t numBins = (Long_t)fSparseHist->GetNbins();
1202 Double_t upperLimit = param.getMin();
1203 Double_t val;
1204 for (Long_t i = 0; i < numBins; i++) {
1205 if (fSparseHist->GetBinContent(i, &coord[0]) >= fHistCutoff) {
1206 val = fSparseHist->GetAxis(d)->GetBinCenter(coord[d]);
1207 if (val > upperLimit)
1208 upperLimit = val;
1209 }
1210 }
1211 return upperLimit;
1212 }
1213 }
1214 return param.getMax();
1215}
1216
1217////////////////////////////////////////////////////////////////////////////////
1218/// Determine the upper limit for param on this interval
1219/// using the binned data set
1220
1222{
1223 if (fHistCutoff < 0)
1224 DetermineByDataHist(); // this initializes fDataHist
1225
1226 if (fHistCutoff < 0) {
1227 // if fHistCutoff < 0 still, then determination of interval failed
1228 coutE(Eval) << "In MCMCInterval::UpperLimitByDataHist: "
1229 << "couldn't determine cutoff. Check that num burn in steps < num "
1230 << "steps in the Markov chain. Returning param.getMax()." << endl;
1231 return param.getMax();
1232 }
1233
1234 for (Int_t d = 0; d < fDimension; d++) {
1235 if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
1236 Int_t numBins = fDataHist->numEntries();
1237 Double_t upperLimit = param.getMin();
1238 Double_t val;
1239 for (Int_t i = 0; i < numBins; i++) {
1240 fDataHist->get(i);
1241 if (fDataHist->weight() >= fHistCutoff) {
1242 val = fDataHist->get()->getRealValue(param.GetName());
1243 if (val > upperLimit)
1244 upperLimit = val;
1245 }
1246 }
1247 return upperLimit;
1248 }
1249 }
1250 return param.getMax();
1251}
1252
1253////////////////////////////////////////////////////////////////////////////////
1254/// Determine the lower limit for param on this interval
1255/// using the keys pdf
1256
1258{
1259 if (fKeysCutoff < 0)
1261
1262 if (fKeysDataHist == NULL)
1264
1265 if (fKeysCutoff < 0 || fKeysDataHist == NULL) {
1266 // failure in determination of cutoff and/or creation of histogram
1267 coutE(Eval) << "in MCMCInterval::LowerLimitByKeys(): "
1268 << "couldn't find lower limit, check that the number of burn in "
1269 << "steps < number of total steps in the Markov chain. Returning "
1270 << "param.getMin()" << endl;
1271 return param.getMin();
1272 }
1273
1274 for (Int_t d = 0; d < fDimension; d++) {
1275 if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
1276 Int_t numBins = fKeysDataHist->numEntries();
1277 Double_t lowerLimit = param.getMax();
1278 Double_t val;
1279 for (Int_t i = 0; i < numBins; i++) {
1280 fKeysDataHist->get(i);
1281 if (fKeysDataHist->weight() >= fKeysCutoff) {
1282 val = fKeysDataHist->get()->getRealValue(param.GetName());
1283 if (val < lowerLimit)
1284 lowerLimit = val;
1285 }
1286 }
1287 return lowerLimit;
1288 }
1289 }
1290 return param.getMin();
1291}
1292
1293////////////////////////////////////////////////////////////////////////////////
1294/// Determine the upper limit for param on this interval
1295/// using the keys pdf
1296
1298{
1299 if (fKeysCutoff < 0)
1301
1302 if (fKeysDataHist == NULL)
1304
1305 if (fKeysCutoff < 0 || fKeysDataHist == NULL) {
1306 // failure in determination of cutoff and/or creation of histogram
1307 coutE(Eval) << "in MCMCInterval::UpperLimitByKeys(): "
1308 << "couldn't find upper limit, check that the number of burn in "
1309 << "steps < number of total steps in the Markov chain. Returning "
1310 << "param.getMax()" << endl;
1311 return param.getMax();
1312 }
1313
1314 for (Int_t d = 0; d < fDimension; d++) {
1315 if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
1316 Int_t numBins = fKeysDataHist->numEntries();
1317 Double_t upperLimit = param.getMin();
1318 Double_t val;
1319 for (Int_t i = 0; i < numBins; i++) {
1320 fKeysDataHist->get(i);
1321 if (fKeysDataHist->weight() >= fKeysCutoff) {
1322 val = fKeysDataHist->get()->getRealValue(param.GetName());
1323 if (val > upperLimit)
1324 upperLimit = val;
1325 }
1326 }
1327 return upperLimit;
1328 }
1329 }
1330 return param.getMax();
1331}
1332
1333////////////////////////////////////////////////////////////////////////////////
1334/// Determine the approximate maximum value of the Keys PDF
1335
1337{
1338 if (fKeysCutoff < 0)
1340
1341 if (fKeysDataHist == NULL)
1343
1344 if (fKeysDataHist == NULL) {
1345 // failure in determination of cutoff and/or creation of histogram
1346 coutE(Eval) << "in MCMCInterval::KeysMax(): "
1347 << "couldn't find Keys max value, check that the number of burn in "
1348 << "steps < number of total steps in the Markov chain. Returning 0"
1349 << endl;
1350 return 0;
1351 }
1352
1353 Int_t numBins = fKeysDataHist->numEntries();
1354 Double_t max = 0;
1355 Double_t w;
1356 for (Int_t i = 0; i < numBins; i++) {
1357 fKeysDataHist->get(i);
1358 w = fKeysDataHist->weight();
1359 if (w > max)
1360 max = w;
1361 }
1362
1363 return max;
1364}
1365
1366////////////////////////////////////////////////////////////////////////////////
1367
1369{
1370 if (fHistCutoff < 0)
1372
1373 return fHistCutoff;
1374}
1375
1376////////////////////////////////////////////////////////////////////////////////
1377
1379{
1380 if (fKeysCutoff < 0)
1382
1383 // kbelasco: if fFull hasn't been set (because Keys creation failed because
1384 // fNumBurnInSteps >= fChain->Size()) then this will return infinity, which
1385 // seems ok to me since it will indicate error
1386
1387 return fKeysCutoff / fFull;
1388}
1389
1390////////////////////////////////////////////////////////////////////////////////
1391
1393{
1394 RooAbsReal* integral;
1395 Double_t confLevel;
1396 fCutoffVar->setVal(cutoff);
1398 confLevel = integral->getVal(fParameters) / full;
1399 coutI(Eval) << "cutoff = " << cutoff << ", conf = " << confLevel << endl;
1400 //cout << "tmp: cutoff = " << cutoff << ", conf = " << confLevel << endl;
1401 delete integral;
1402 return confLevel;
1403}
1404
1405////////////////////////////////////////////////////////////////////////////////
1406
1408{
1409 if(fConfidenceLevel == 0)
1410 coutE(InputArguments) << "Error in MCMCInterval::GetPosteriorHist: "
1411 << "confidence level not set " << endl;
1412 if (fHist == NULL)
1413 CreateHist();
1414
1415 if (fHist == NULL)
1416 // if fHist is still NULL, then CreateHist failed
1417 return NULL;
1418
1419 return (TH1*) fHist->Clone("MCMCposterior_hist");
1420}
1421
1422////////////////////////////////////////////////////////////////////////////////
1423
1425{
1426 if (fConfidenceLevel == 0)
1427 coutE(InputArguments) << "Error in MCMCInterval::GetPosteriorKeysPdf: "
1428 << "confidence level not set " << endl;
1429 if (fKeysPdf == NULL)
1430 CreateKeysPdf();
1431
1432 if (fKeysPdf == NULL)
1433 // if fKeysPdf is still NULL, then it means CreateKeysPdf failed
1434 return NULL;
1435
1436 return (RooNDKeysPdf*) fKeysPdf->Clone("MCMCPosterior_keys");
1437}
1438
1439////////////////////////////////////////////////////////////////////////////////
1440
1442{
1443 if (fConfidenceLevel == 0)
1444 coutE(InputArguments) << "MCMCInterval::GetPosteriorKeysProduct: "
1445 << "confidence level not set " << endl;
1446 if (fProduct == NULL) {
1447 CreateKeysPdf();
1449 }
1450
1451 if (fProduct == NULL)
1452 // if fProduct is still NULL, then it means CreateKeysPdf failed
1453 return NULL;
1454
1455 return (RooProduct*) fProduct->Clone("MCMCPosterior_keysproduct");
1456}
1457
1458////////////////////////////////////////////////////////////////////////////////
1459
1461{
1462 // returns list of parameters
1463 return new RooArgSet(fParameters);
1464}
1465
1466////////////////////////////////////////////////////////////////////////////////
1467
1469{
1470 return (TMath::Abs(confLevel - fConfidenceLevel) < fEpsilon);
1471}
1472
1473////////////////////////////////////////////////////////////////////////////////
1474
1476{
1477 return (TMath::Abs(a - b) < TMath::Abs(fDelta * (a + b)/2));
1478}
1479
1480////////////////////////////////////////////////////////////////////////////////
1481
1483{
1484 if (fAxes == NULL)
1485 return;
1486 if (fProduct == NULL)
1488 if (fProduct == NULL)
1489 // if fProduct still NULL, then creation failed
1490 return;
1491
1492 //RooAbsBinning** savedBinning = new RooAbsBinning*[fDimension];
1493 Int_t* savedBins = new Int_t[fDimension];
1494 Int_t i;
1495 Double_t numBins;
1496 RooRealVar* var;
1497
1498 // kbelasco: Note - the accuracy is only increased here if the binning for
1499 // each RooRealVar is uniform
1500
1501 // kbelasco: look into why saving the binnings and replacing them doesn't
1502 // work (replaces with 1 bin always).
1503 // Note: this code modifies the binning for the parameters (if they are
1504 // uniform) and sets them back to what they were. If the binnings are not
1505 // uniform, this code does nothing.
1506
1507 // first scan through fAxes to make sure all binnings are uniform, or else
1508 // we can't change the number of bins because there seems to be an error
1509 // when setting the binning itself rather than just the number of bins
1510 Bool_t tempChangeBinning = true;
1511 for (i = 0; i < fDimension; i++) {
1512 if (!fAxes[i]->getBinning(NULL, false, false).isUniform()) {
1513 tempChangeBinning = false;
1514 break;
1515 }
1516 }
1517
1518 // kbelasco: for 1 dimension this should be fine, but for more dimensions
1519 // the total number of bins in the histogram increases exponentially with
1520 // the dimension, so don't do this above 1-D for now.
1521 if (fDimension >= 2)
1522 tempChangeBinning = false;
1523
1524 if (tempChangeBinning) {
1525 // set high number of bins for high accuracy on lower/upper limit by keys
1526 for (i = 0; i < fDimension; i++) {
1527 var = fAxes[i];
1528 //savedBinning[i] = &var->getBinning("__binning_clone", false, true);
1529 savedBins[i] = var->getBinning(NULL, false, false).numBins();
1530 numBins = (var->getMax() - var->getMin()) / fEpsilon;
1531 var->setBins((Int_t)numBins);
1532 }
1533 }
1534
1535 fKeysDataHist = new RooDataHist("_productDataHist",
1536 "Keys PDF & Heaviside Product Data Hist", fParameters);
1538
1539 if (tempChangeBinning) {
1540 // set the binning back to normal
1541 for (i = 0; i < fDimension; i++)
1542 //fAxes[i]->setBinning(*savedBinning[i], NULL);
1543 //fAxes[i]->setBins(savedBinning[i]->numBins(), NULL);
1544 fAxes[i]->setBins(savedBins[i], NULL);
1545 }
1546
1547 //delete[] savedBinning;
1548 delete[] savedBins;
1549}
1550
1551////////////////////////////////////////////////////////////////////////////////
1552
1554{
1555 // check that the parameters are correct
1556
1557 if (parameterPoint.getSize() != fParameters.getSize() ) {
1558 coutE(Eval) << "MCMCInterval: size is wrong, parameters don't match" << std::endl;
1559 return kFALSE;
1560 }
1561 if ( ! parameterPoint.equals( fParameters ) ) {
1562 coutE(Eval) << "MCMCInterval: size is ok, but parameters don't match" << std::endl;
1563 return kFALSE;
1564 }
1565 return kTRUE;
1566}
static const Double_t DEFAULT_DELTA
static const Double_t DEFAULT_EPSILON
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
#define e(i)
Definition RSha256.hxx:103
#define coutI(a)
#define coutW(a)
#define coutE(a)
int Int_t
Definition RtypesCore.h:45
const Bool_t kFALSE
Definition RtypesCore.h:92
long Long_t
Definition RtypesCore.h:54
double Double_t
Definition RtypesCore.h:59
long long Long64_t
Definition RtypesCore.h:73
const Bool_t kTRUE
Definition RtypesCore.h:91
#define ClassImp(name)
Definition Rtypes.h:364
char name[80]
Definition TGX11.cxx:110
THnSparseT< TArrayF > THnSparseF
Definition THnSparse.h:220
TRObject operator()(const T1 &t1) const
virtual TObject * Clone(const char *newname=0) const
Make a clone of an object using the Streamer facility.
Definition RooAbsArg.h:84
Int_t numBins() const
Return number of bins.
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
Int_t getSize() const
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
Double_t getRealValue(const char *name, Double_t defVal=0, Bool_t verbose=kFALSE) const
Get value of a RooAbsReal stored in set with given name.
Bool_t equals(const RooAbsCollection &otherColl) const
Check if this and other collection have identically-named contents.
TIterator * createIterator(Bool_t dir=kIterForward) const
TIterator-style iteration over contained elements.
virtual Double_t getMax(const char *name=0) const
Get maximum of currently defined range.
virtual Int_t numBins(const char *rangeName=0) const
virtual Double_t getMin(const char *name=0) const
Get miniminum of currently defined range.
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:61
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
RooAbsReal * createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables listed in ...
RooDataHist * fillDataHist(RooDataHist *hist, const RooArgSet *nset, Double_t scaleFactor, Bool_t correctForBinVolume=kFALSE, Bool_t showProgress=kFALSE) const
Fill a RooDataHist with values sampled from this function at the bin centers.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:21
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:70
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Add element to non-owning set.
The RooDataHist is a container class to hold N-dimensional binned data.
Definition RooDataHist.h:37
double weight(std::size_t i) const
Return weight of i-th bin.
Definition RooDataHist.h:98
Int_t getIndex(const RooArgSet &coord, Bool_t fast=false) const
Calculate bin number of the given coordinates.
Int_t numEntries() const override
Return the number of bins.
const RooArgSet * get() const override
Get bin centre of current bin.
Definition RooDataHist.h:74
Double_t sum(bool correctForBinSize, bool inverseCorr=false) const
Return the sum of the weights of all bins in the histogram.
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:33
Generic N-dimensional implementation of a kernel estimation p.d.f.
static Double_t infinity()
Return internal infinity representation.
Definition RooNumber.cxx:49
A RooProduct represents the product of a given set of RooAbsReal objects.
Definition RooProduct.h:29
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
void setBins(Int_t nBins, const char *name=0)
Create a uniform binning under name 'name' for this variable.
const RooAbsBinning & getBinning(const char *name=0, Bool_t verbose=kTRUE, Bool_t createOnTheFly=kFALSE) const
Return binning definition with name.
virtual void setVal(Double_t value)
Set value of variable to 'value'.
ConfInterval is an interface class for a generic interval in the RooStats framework.
Represents the Heaviside function.
Definition Heaviside.h:18
MCMCInterval is a concrete implementation of the RooStats::ConfInterval interface.
virtual void CreateDataHist()
virtual Double_t LowerLimitBySparseHist(RooRealVar &param)
determine lower limit using histogram
virtual void DetermineByDataHist()
virtual Double_t UpperLimitByKeys(RooRealVar &param)
determine upper limit in the shortest interval by using keys pdf
virtual void DetermineShortestInterval()
virtual void CreateVector(RooRealVar *param)
virtual Double_t UpperLimit(RooRealVar &param)
get the highest value of param that is within the confidence interval
RooDataHist * fDataHist
virtual Bool_t IsInInterval(const RooArgSet &point) const
determine whether this point is in the confidence interval
Bool_t CheckParameters(const RooArgSet &point) const
check if parameters are correct. (dummy implementation to start)
virtual Double_t UpperLimitShortest(RooRealVar &param)
get the upper limit of param in the confidence interval Note that this works better for some distribu...
enum IntervalType fIntervalType
virtual Double_t UpperLimitBySparseHist(RooRealVar &param)
determine upper limit using histogram
virtual void SetAxes(RooArgList &axes)
Set which parameters go on which axis.
virtual void DetermineByKeys()
virtual void DetermineBySparseHist()
virtual TH1 * GetPosteriorHist()
set the number of bins to use (same for all axes, for now) virtual void SetNumBins(Int_t numBins);
virtual Double_t LowerLimitByKeys(RooRealVar &param)
determine lower limit in the shortest interval by using keys pdf
virtual Double_t GetHistCutoff()
get the cutoff bin height for being considered in the confidence interval
virtual Double_t UpperLimitTailFraction(RooRealVar &param)
determine upper limit of the lower confidence interval
MCMCInterval(const char *name=0)
default constructor
virtual void DetermineTailFractionInterval()
virtual void CreateSparseHist()
virtual void DetermineInterval()
Bool_t AcceptableConfLevel(Double_t confLevel)
virtual Double_t GetActualConfidenceLevel()
virtual Double_t GetKeysPdfCutoff() { return fKeysCutoff; }
virtual void SetConfidenceLevel(Double_t cl)
set the desired confidence level (see GetActualConfidenceLevel()) Note: calling this function trigger...
virtual Double_t LowerLimitShortest(RooRealVar &param)
get the lower limit of param in the shortest confidence interval Note that this works better for some...
virtual Double_t LowerLimit(RooRealVar &param)
get the lowest value of param that is within the confidence interval
virtual Double_t GetKeysPdfCutoff()
get the cutoff RooNDKeysPdf value for being considered in the confidence interval
std::vector< Int_t > fVector
virtual Double_t LowerLimitTailFraction(RooRealVar &param)
determine lower limit of the lower confidence interval
virtual RooProduct * GetPosteriorKeysProduct()
Get a clone of the (keyspdf * heaviside) product of the posterior.
virtual Double_t UpperLimitByHist(RooRealVar &param)
determine upper limit using histogram
virtual Double_t UpperLimitByDataHist(RooRealVar &param)
determine upper limit using histogram
virtual Double_t LowerLimitByHist(RooRealVar &param)
determine lower limit using histogram
virtual RooNDKeysPdf * GetPosteriorKeysPdf()
Get a clone of the keys pdf of the posterior.
virtual void CreateKeysPdf()
virtual Double_t LowerLimitByDataHist(RooRealVar &param)
determine lower limit using histogram
virtual void CreateHist()
RooDataHist * fKeysDataHist
virtual Double_t CalcConfLevel(Double_t cutoff, Double_t full)
virtual void DetermineByHist()
virtual void SetParameters(const RooArgSet &parameters)
Set the parameters of interest for this interval and change other internal data members accordingly.
virtual void CreateKeysDataHist()
RooNDKeysPdf * fKeysPdf
Double_t GetKeysMax()
Determine the approximate maximum value of the Keys PDF.
Bool_t WithinDeltaFraction(Double_t a, Double_t b)
virtual RooArgSet * GetParameters() const
return a set containing the parameters of this interval the caller owns the returned RooArgSet*
Stores the steps in a Markov Chain of points.
Definition MarkovChain.h:30
virtual Double_t Weight() const
get the weight of the current (last indexed) entry
virtual RooDataHist * GetAsDataHist(RooArgSet *whichVars=NULL) const
get this MarkovChain as a RooDataHist whose entries contain the values of whichVars.
virtual const RooArgSet * Get(Int_t i) const
get the entry at position i
Definition MarkovChain.h:51
virtual RooDataSet * GetAsDataSet(RooArgSet *whichVars=NULL) const
get this MarkovChain as a RooDataSet whose entries contain the values of whichVars.
virtual Int_t Size() const
get the number of steps in the chain
Definition MarkovChain.h:49
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition TAxis.cxx:478
1-D histogram with a float per channel (see TH1 documentation)}
Definition TH1.h:575
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
TAxis * GetZaxis()
Definition TH1.h:322
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition TH1.h:320
TObject * Clone(const char *newname=0) const
Make a complete copy of the underlying object.
Definition TH1.cxx:2740
TAxis * GetYaxis()
Definition TH1.h:321
2-D histogram with a float per channel (see TH1 documentation)}
Definition TH2.h:251
3-D histogram with a float per channel (see TH1 documentation)}
Definition TH3.h:268
Long64_t Fill(const Double_t *x, Double_t w=1.)
Definition THnBase.h:144
Double_t GetSumw() const
Definition THnBase.h:184
TAxis * GetAxis(Int_t dim) const
Definition THnBase.h:125
Efficient multidimensional histogram.
Definition THnSparse.h:36
Long64_t GetBin(const Int_t *idx) const
Definition THnSparse.h:95
Double_t GetBinContent(const Int_t *idx) const
Forwards to THnBase::GetBinContent() overload.
Definition THnSparse.h:120
void Sumw2()
Enable calculation of errors.
Long64_t GetNbins() const
Definition THnSparse.h:92
Iterator abstract base class.
Definition TIterator.h:30
virtual TObject * Next()=0
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
Mother of all ROOT objects.
Definition TObject.h:37
virtual const char * GetName() const
Returns name of object.
Definition TObject.cxx:359
RooCmdArg NormSet(const RooArgSet &nset)
RooCmdArg SelectVars(const RooArgSet &vars)
RooCmdArg EventRange(Int_t nStart, Int_t nStop)
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
@ InputArguments
Namespace for the RooStats classes.
Definition Asimov.h:19
void SetParameters(const RooArgSet *desiredVals, RooArgSet *paramsToChange)
Short_t Abs(Short_t d)
Definition TMathBase.h:120
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345