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
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 }
197};
198
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 }
207};
208
211 fChain(chain), fParam(param) {}
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 }
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 fHist = NULL;
389 }
390
391 if (fNumBurnInSteps >= fChain->Size()) {
393 "MCMCInterval::CreateHist: creation of histogram failed: " <<
394 "Number of burn-in steps (num steps to ignore) >= number of steps " <<
395 "in Markov chain." << endl;
396 return;
397 }
398
399 if (fDimension == 1)
400 fHist = new TH1F("posterior", "MCMC Posterior Histogram",
401 fAxes[0]->numBins(), fAxes[0]->getMin(), fAxes[0]->getMax());
402
403 else if (fDimension == 2)
404 fHist = new TH2F("posterior", "MCMC Posterior Histogram",
405 fAxes[0]->numBins(), fAxes[0]->getMin(), fAxes[0]->getMax(),
406 fAxes[1]->numBins(), fAxes[1]->getMin(), fAxes[1]->getMax());
407
408 else if (fDimension == 3)
409 fHist = new TH3F("posterior", "MCMC Posterior Histogram",
410 fAxes[0]->numBins(), fAxes[0]->getMin(), fAxes[0]->getMax(),
411 fAxes[1]->numBins(), fAxes[1]->getMin(), fAxes[1]->getMax(),
412 fAxes[2]->numBins(), fAxes[2]->getMin(), fAxes[2]->getMax());
413
414 else {
415 coutE(Eval) << "* Error in MCMCInterval::CreateHist() : " <<
416 "TH1* couldn't handle dimension: " << fDimension << endl;
417 return;
418 }
419
420 // Fill histogram
421 Int_t size = fChain->Size();
422 const RooArgSet* entry;
423 for (Int_t i = fNumBurnInSteps; i < size; i++) {
424 entry = fChain->Get(i);
425 if (fDimension == 1)
426 ((TH1F*)fHist)->Fill(entry->getRealValue(fAxes[0]->GetName()),
427 fChain->Weight());
428 else if (fDimension == 2)
429 ((TH2F*)fHist)->Fill(entry->getRealValue(fAxes[0]->GetName()),
430 entry->getRealValue(fAxes[1]->GetName()),
431 fChain->Weight());
432 else
433 ((TH3F*)fHist)->Fill(entry->getRealValue(fAxes[0]->GetName()),
434 entry->getRealValue(fAxes[1]->GetName()),
435 entry->getRealValue(fAxes[2]->GetName()),
436 fChain->Weight());
437 }
438
439 if (fDimension >= 1)
441 if (fDimension >= 2)
443 if (fDimension >= 3)
445}
446
447////////////////////////////////////////////////////////////////////////////////
448
450{
451 if (fAxes == NULL || fChain == NULL) {
452 coutE(InputArguments) << "* Error in MCMCInterval::CreateSparseHist(): "
453 << "Crucial data member was NULL." << endl;
454 coutE(InputArguments) << "Make sure to fully construct/initialize."
455 << endl;
456 return;
457 }
458 if (fSparseHist != NULL)
459 delete fSparseHist;
460
461 Double_t* min = new Double_t[fDimension];
462 Double_t* max = new Double_t[fDimension];
463 Int_t* bins = new Int_t[fDimension];
464 for (Int_t i = 0; i < fDimension; i++) {
465 min[i] = fAxes[i]->getMin();
466 max[i] = fAxes[i]->getMax();
467 bins[i] = fAxes[i]->numBins();
468 }
469 fSparseHist = new THnSparseF("posterior", "MCMC Posterior Histogram",
470 fDimension, bins, min, max);
471
472 delete[] min;
473 delete[] max;
474 delete[] bins;
475
476 // kbelasco: it appears we need to call Sumw2() just to get the
477 // histogram to keep a running total of the weight so that Getsumw doesn't
478 // just return 0
480
481 if (fNumBurnInSteps >= fChain->Size()) {
483 "MCMCInterval::CreateSparseHist: creation of histogram failed: " <<
484 "Number of burn-in steps (num steps to ignore) >= number of steps " <<
485 "in Markov chain." << endl;
486 }
487
488 // Fill histogram
489 Int_t size = fChain->Size();
490 const RooArgSet* entry;
492 for (Int_t i = fNumBurnInSteps; i < size; i++) {
493 entry = fChain->Get(i);
494 for (Int_t ii = 0; ii < fDimension; ii++)
495 x[ii] = entry->getRealValue(fAxes[ii]->GetName());
497 }
498 delete[] x;
499}
500
501////////////////////////////////////////////////////////////////////////////////
502
504{
505 if (fParameters.getSize() == 0 || fChain == NULL) {
506 coutE(Eval) << "* Error in MCMCInterval::CreateDataHist(): " <<
507 "Crucial data member was NULL or empty." << endl;
508 coutE(Eval) << "Make sure to fully construct/initialize." << endl;
509 return;
510 }
511
512 if (fNumBurnInSteps >= fChain->Size()) {
514 "MCMCInterval::CreateDataHist: creation of histogram failed: " <<
515 "Number of burn-in steps (num steps to ignore) >= number of steps " <<
516 "in Markov chain." << endl;
517 fDataHist = NULL;
518 return;
519 }
520
523}
524
525////////////////////////////////////////////////////////////////////////////////
526
528{
529 fVector.clear();
530 fVecWeight = 0;
531
532 if (fChain == NULL) {
533 coutE(InputArguments) << "* Error in MCMCInterval::CreateVector(): " <<
534 "Crucial data member (Markov chain) was NULL." << endl;
535 coutE(InputArguments) << "Make sure to fully construct/initialize."
536 << endl;
537 return;
538 }
539
540 if (fNumBurnInSteps >= fChain->Size()) {
542 "MCMCInterval::CreateVector: creation of vector failed: " <<
543 "Number of burn-in steps (num steps to ignore) >= number of steps " <<
544 "in Markov chain." << endl;
545 }
546
547 // Fill vector
549 fVector.resize(size);
550 Int_t i;
551 Int_t chainIndex;
552 for (i = 0; i < size; i++) {
553 chainIndex = i + fNumBurnInSteps;
554 fVector[i] = chainIndex;
555 fVecWeight += fChain->Weight(chainIndex);
556 }
557
558 stable_sort(fVector.begin(), fVector.end(),
560}
561
562////////////////////////////////////////////////////////////////////////////////
563
565{
567 fParameters.add(parameters);
569 if (fAxes != NULL)
570 delete[] fAxes;
573 Int_t n = 0;
574 TObject* obj;
575 while ((obj = it->Next()) != NULL) {
576 if (dynamic_cast<RooRealVar*>(obj) != NULL)
577 fAxes[n] = (RooRealVar*)obj;
578 else
579 coutE(Eval) << "* Error in MCMCInterval::SetParameters: " <<
580 obj->GetName() << " not a RooRealVar*" << std::endl;
581 n++;
582 }
583 delete it;
584}
585
586////////////////////////////////////////////////////////////////////////////////
587
589{
590 switch (fIntervalType) {
591 case kShortest:
593 break;
594 case kTailFraction:
596 break;
597 default:
598 coutE(InputArguments) << "MCMCInterval::DetermineInterval(): " <<
599 "Error: Interval type not set" << endl;
600 break;
601 }
602}
603
604////////////////////////////////////////////////////////////////////////////////
605
607{
608 if (fUseKeys)
610 else
612}
613
614////////////////////////////////////////////////////////////////////////////////
615
617{
618 if (fLeftSideTF < 0 || fLeftSideTF > 1) {
619 coutE(InputArguments) << "MCMCInterval::DetermineTailFractionInterval: "
620 << "Fraction must be in the range [0, 1]. "
621 << fLeftSideTF << "is not allowed." << endl;
622 return;
623 }
624
625 if (fDimension != 1) {
626 coutE(InputArguments) << "MCMCInterval::DetermineTailFractionInterval(): "
627 << "Error: Can only find a tail-fraction interval for 1-D intervals"
628 << endl;
629 return;
630 }
631
632 if (fAxes == NULL) {
633 coutE(InputArguments) << "MCMCInterval::DetermineTailFractionInterval(): "
634 << "Crucial data member was NULL." << endl;
635 coutE(InputArguments) << "Make sure to fully construct/initialize."
636 << endl;
637 return;
638 }
639
640 // kbelasco: fill in code here to find interval
641 //
642 // also make changes so that calling GetPosterior...() returns NULL
643 // when fIntervalType == kTailFraction, since there really
644 // is no posterior for this type of interval determination
645 if (fVector.size() == 0)
647
648 if (fVector.size() == 0 || fVecWeight == 0) {
649 // if size is still 0, then creation failed.
650 // if fVecWeight == 0, then there are no entries (indicates the same
651 // error as fVector.size() == 0 because that only happens when
652 // fNumBurnInSteps >= fChain->Size())
653 // either way, reset and return
654 fVector.clear();
655 fTFLower = -1.0 * RooNumber::infinity();
657 fTFConfLevel = 0.0;
658 fVecWeight = 0;
659 return;
660 }
661
662 RooRealVar* param = fAxes[0];
663
665 Double_t leftTailCutoff = fVecWeight * (1 - c) * fLeftSideTF;
666 Double_t rightTailCutoff = fVecWeight * (1 - c) * (1 - fLeftSideTF);
667 Double_t leftTailSum = 0;
668 Double_t rightTailSum = 0;
669
670 // kbelasco: consider changing these values to +infinity and -infinity
671 Double_t ll = param->getMin();
672 Double_t ul = param->getMax();
673
674 Double_t x;
675 Double_t w;
676
677 // save a lot of GetName() calls if compiler does not already optimize this
678 const char* name = param->GetName();
679
680 // find lower limit
681 Int_t i;
682 for (i = 0; i < (Int_t)fVector.size(); i++) {
684 w = fChain->Weight();
685
686 if (TMath::Abs(leftTailSum + w - leftTailCutoff) <
687 TMath::Abs(leftTailSum - leftTailCutoff)) {
688 // moving the lower limit to x would bring us closer to the desired
689 // left tail size
690 ll = x;
691 leftTailSum += w;
692 } else
693 break;
694 }
695
696 // find upper limit
697 for (i = (Int_t)fVector.size() - 1; i >= 0; i--) {
699 w = fChain->Weight();
700
701 if (TMath::Abs(rightTailSum + w - rightTailCutoff) <
702 TMath::Abs(rightTailSum - rightTailCutoff)) {
703 // moving the lower limit to x would bring us closer to the desired
704 // left tail size
705 ul = x;
706 rightTailSum += w;
707 } else
708 break;
709 }
710
711 fTFLower = ll;
712 fTFUpper = ul;
713 fTFConfLevel = 1 - (leftTailSum + rightTailSum) / fVecWeight;
714}
715
716////////////////////////////////////////////////////////////////////////////////
717
719{
720 if (fKeysPdf == NULL)
722
723 if (fKeysPdf == NULL) {
724 // if fKeysPdf is still NULL, then it means CreateKeysPdf failed
725 // so clear all the data members this function would normally determine
726 // and return
727 fFull = 0.0;
728 fKeysCutoff = -1;
729 fKeysConfLevel = 0.0;
730 return;
731 }
732
733 // now we have a keys pdf of the posterior
734
735 Double_t cutoff = 0.0;
736 fCutoffVar->setVal(cutoff);
738 Double_t full = integral->getVal(fParameters);
739 fFull = full;
740 delete integral;
741 if (full < 0.98) {
742 coutW(Eval) << "Warning: Integral of Keys PDF came out to " << full
743 << " instead of expected value 1. Will continue using this "
744 << "factor to normalize further integrals of this PDF." << endl;
745 }
746
747 // kbelasco: Is there a better way to set the search range?
748 // from 0 to max value of Keys
749 // kbelasco: how to get max value?
750 //Double_t max = product.maxVal(product.getMaxVal(fParameters));
751
752 Double_t volume = 1.0;
754 RooRealVar* var;
755 while ((var = (RooRealVar*)it->Next()) != NULL)
756 volume *= (var->getMax() - var->getMin());
757 delete it;
758
759 Double_t topCutoff = full / volume;
760 Double_t bottomCutoff = topCutoff;
761 Double_t confLevel = CalcConfLevel(topCutoff, full);
762 if (AcceptableConfLevel(confLevel)) {
763 fKeysConfLevel = confLevel;
764 fKeysCutoff = topCutoff;
765 return;
766 }
767 Bool_t changed = kFALSE;
768 // find high end of range
769 while (confLevel > fConfidenceLevel) {
770 topCutoff *= 2.0;
771 confLevel = CalcConfLevel(topCutoff, full);
772 if (AcceptableConfLevel(confLevel)) {
773 fKeysConfLevel = confLevel;
774 fKeysCutoff = topCutoff;
775 return;
776 }
777 changed = kTRUE;
778 }
779 if (changed) {
780 bottomCutoff = topCutoff / 2.0;
781 } else {
782 changed = kFALSE;
783 bottomCutoff /= 2.0;
784 confLevel = CalcConfLevel(bottomCutoff, full);
785 if (AcceptableConfLevel(confLevel)) {
786 fKeysConfLevel = confLevel;
787 fKeysCutoff = bottomCutoff;
788 return;
789 }
790 while (confLevel < fConfidenceLevel) {
791 bottomCutoff /= 2.0;
792 confLevel = CalcConfLevel(bottomCutoff, full);
793 if (AcceptableConfLevel(confLevel)) {
794 fKeysConfLevel = confLevel;
795 fKeysCutoff = bottomCutoff;
796 return;
797 }
798 changed = kTRUE;
799 }
800 if (changed) {
801 topCutoff = bottomCutoff * 2.0;
802 }
803 }
804
805 coutI(Eval) << "range set: [" << bottomCutoff << ", " << topCutoff << "]"
806 << endl;
807
808 cutoff = (topCutoff + bottomCutoff) / 2.0;
809 confLevel = CalcConfLevel(cutoff, full);
810
811 // need to use WithinDeltaFraction() because sometimes the integrating the
812 // posterior in this binary search seems to not have enough granularity to
813 // find an acceptable conf level (small no. of strange cases).
814 // WithinDeltaFraction causes the search to terminate when
815 // topCutoff is essentially equal to bottomCutoff (compared to the magnitude
816 // of their mean).
817 while (!AcceptableConfLevel(confLevel) &&
818 !WithinDeltaFraction(topCutoff, bottomCutoff)) {
819 if (confLevel > fConfidenceLevel)
820 bottomCutoff = cutoff;
821 else if (confLevel < fConfidenceLevel)
822 topCutoff = cutoff;
823 cutoff = (topCutoff + bottomCutoff) / 2.0;
824 coutI(Eval) << "cutoff range: [" << bottomCutoff << ", "
825 << topCutoff << "]" << endl;
826 confLevel = CalcConfLevel(cutoff, full);
827 }
828
829 fKeysCutoff = cutoff;
830 fKeysConfLevel = confLevel;
831}
832
833////////////////////////////////////////////////////////////////////////////////
834
836{
837 if (fUseSparseHist)
839 else
841}
842
843////////////////////////////////////////////////////////////////////////////////
844
846{
847 Long_t numBins;
848 if (fSparseHist == NULL)
850
851 if (fSparseHist == NULL) {
852 // if fSparseHist is still NULL, then CreateSparseHist failed
853 fHistCutoff = -1;
854 fHistConfLevel = 0.0;
855 return;
856 }
857
858 numBins = (Long_t)fSparseHist->GetNbins();
859
860 std::vector<Long_t> bins(numBins);
861 for (Int_t ibin = 0; ibin < numBins; ibin++)
862 bins[ibin] = (Long_t)ibin;
863 std::stable_sort(bins.begin(), bins.end(), CompareSparseHistBins(fSparseHist));
864
865 Double_t nEntries = fSparseHist->GetSumw();
866 Double_t sum = 0;
867 Double_t content;
868 Int_t i;
869 // see above note on indexing to understand numBins - 3
870 for (i = numBins - 1; i >= 0; i--) {
871 content = fSparseHist->GetBinContent(bins[i]);
872 if ((sum + content) / nEntries >= fConfidenceLevel) {
873 fHistCutoff = content;
874 if (fIsHistStrict) {
875 sum += content;
876 i--;
877 break;
878 } else {
879 i++;
880 break;
881 }
882 }
883 sum += content;
884 }
885
886 if (fIsHistStrict) {
887 // keep going to find the sum
888 for ( ; i >= 0; i--) {
889 content = fSparseHist->GetBinContent(bins[i]);
890 if (content == fHistCutoff)
891 sum += content;
892 else
893 break; // content must be < fHistCutoff
894 }
895 } else {
896 // backtrack to find the cutoff and sum
897 for ( ; i < numBins; i++) {
898 content = fSparseHist->GetBinContent(bins[i]);
899 if (content > fHistCutoff) {
900 fHistCutoff = content;
901 break;
902 } else // content == fHistCutoff
903 sum -= content;
904 if (i == numBins - 1)
905 // still haven't set fHistCutoff correctly yet, and we have no bins
906 // left, so set fHistCutoff to something higher than the tallest bin
907 fHistCutoff = content + 1.0;
908 }
909 }
910
911 fHistConfLevel = sum / nEntries;
912}
913
914////////////////////////////////////////////////////////////////////////////////
915
917{
918 Int_t numBins;
919 if (fDataHist == NULL)
921 if (fDataHist == NULL) {
922 // if fDataHist is still NULL, then CreateDataHist failed
923 fHistCutoff = -1;
924 fHistConfLevel = 0.0;
925 return;
926 }
927
928 numBins = fDataHist->numEntries();
929
930 std::vector<Int_t> bins(numBins);
931 for (Int_t ibin = 0; ibin < numBins; ibin++)
932 bins[ibin] = ibin;
933 std::stable_sort(bins.begin(), bins.end(), CompareDataHistBins(fDataHist));
934
935 Double_t nEntries = fDataHist->sum(kFALSE);
936 Double_t sum = 0;
937 Double_t content;
938 Int_t i;
939 for (i = numBins - 1; i >= 0; i--) {
940 fDataHist->get(bins[i]);
941 content = fDataHist->weight();
942 if ((sum + content) / nEntries >= fConfidenceLevel) {
943 fHistCutoff = content;
944 if (fIsHistStrict) {
945 sum += content;
946 i--;
947 break;
948 } else {
949 i++;
950 break;
951 }
952 }
953 sum += content;
954 }
955
956 if (fIsHistStrict) {
957 // keep going to find the sum
958 for ( ; i >= 0; i--) {
959 fDataHist->get(bins[i]);
960 content = fDataHist->weight();
961 if (content == fHistCutoff)
962 sum += content;
963 else
964 break; // content must be < fHistCutoff
965 }
966 } else {
967 // backtrack to find the cutoff and sum
968 for ( ; i < numBins; i++) {
969 fDataHist->get(bins[i]);
970 content = fDataHist->weight();
971 if (content > fHistCutoff) {
972 fHistCutoff = content;
973 break;
974 } else // content == fHistCutoff
975 sum -= content;
976 if (i == numBins - 1)
977 // still haven't set fHistCutoff correctly yet, and we have no bins
978 // left, so set fHistCutoff to something higher than the tallest bin
979 fHistCutoff = content + 1.0;
980 }
981 }
982
983 fHistConfLevel = sum / nEntries;
984}
985
986////////////////////////////////////////////////////////////////////////////////
987
989{
990 if (fIntervalType == kShortest) {
991 if (fUseKeys)
992 return fKeysConfLevel;
993 else
994 return fHistConfLevel;
995 } else if (fIntervalType == kTailFraction) {
996 return fTFConfLevel;
997 } else {
998 coutE(InputArguments) << "MCMCInterval::GetActualConfidenceLevel: "
999 << "not implemented for this type of interval. Returning 0." << endl;
1000 return 0;
1001 }
1002}
1003
1004////////////////////////////////////////////////////////////////////////////////
1005
1007{
1008 switch (fIntervalType) {
1009 case kShortest:
1010 return LowerLimitShortest(param);
1011 case kTailFraction:
1012 return LowerLimitTailFraction(param);
1013 default:
1014 coutE(InputArguments) << "MCMCInterval::LowerLimit(): " <<
1015 "Error: Interval type not set" << endl;
1016 return RooNumber::infinity();
1017 }
1018}
1019
1020////////////////////////////////////////////////////////////////////////////////
1021
1023{
1024 switch (fIntervalType) {
1025 case kShortest:
1026 return UpperLimitShortest(param);
1027 case kTailFraction:
1028 return UpperLimitTailFraction(param);
1029 default:
1030 coutE(InputArguments) << "MCMCInterval::UpperLimit(): " <<
1031 "Error: Interval type not set" << endl;
1032 return RooNumber::infinity();
1033 }
1034}
1035
1036////////////////////////////////////////////////////////////////////////////////
1037
1039{
1040 if (fTFLower == -1.0 * RooNumber::infinity())
1042
1043 return fTFLower;
1044}
1045
1046////////////////////////////////////////////////////////////////////////////////
1047
1049{
1052
1053 return fTFUpper;
1054}
1055
1056////////////////////////////////////////////////////////////////////////////////
1057
1059{
1060 if (fUseKeys)
1061 return LowerLimitByKeys(param);
1062 else
1063 return LowerLimitByHist(param);
1064}
1065
1066////////////////////////////////////////////////////////////////////////////////
1067
1069{
1070 if (fUseKeys)
1071 return UpperLimitByKeys(param);
1072 else
1073 return UpperLimitByHist(param);
1074}
1075
1076////////////////////////////////////////////////////////////////////////////////
1077/// Determine the lower limit for param on this interval
1078/// using the binned data set
1079
1081{
1082 if (fUseSparseHist)
1083 return LowerLimitBySparseHist(param);
1084 else
1085 return LowerLimitByDataHist(param);
1086}
1087
1088////////////////////////////////////////////////////////////////////////////////
1089/// Determine the upper limit for param on this interval
1090/// using the binned data set
1091
1093{
1094 if (fUseSparseHist)
1095 return UpperLimitBySparseHist(param);
1096 else
1097 return UpperLimitByDataHist(param);
1098}
1099
1100////////////////////////////////////////////////////////////////////////////////
1101/// Determine the lower limit for param on this interval
1102/// using the binned data set
1103
1105{
1106 if (fDimension != 1) {
1107 coutE(InputArguments) << "In MCMCInterval::LowerLimitBySparseHist: "
1108 << "Sorry, will not compute lower limit unless dimension == 1" << endl;
1109 return param.getMin();
1110 }
1111 if (fHistCutoff < 0)
1112 DetermineBySparseHist(); // this initializes fSparseHist
1113
1114 if (fHistCutoff < 0) {
1115 // if fHistCutoff < 0 still, then determination of interval failed
1116 coutE(Eval) << "In MCMCInterval::LowerLimitBySparseHist: "
1117 << "couldn't determine cutoff. Check that num burn in steps < num "
1118 << "steps in the Markov chain. Returning param.getMin()." << endl;
1119 return param.getMin();
1120 }
1121
1122 std::vector<Int_t> coord(fDimension);
1123 for (Int_t d = 0; d < fDimension; d++) {
1124 if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
1125 Long_t numBins = (Long_t)fSparseHist->GetNbins();
1126 Double_t lowerLimit = param.getMax();
1127 Double_t val;
1128 for (Long_t i = 0; i < numBins; i++) {
1129 if (fSparseHist->GetBinContent(i, &coord[0]) >= fHistCutoff) {
1130 val = fSparseHist->GetAxis(d)->GetBinCenter(coord[d]);
1131 if (val < lowerLimit)
1132 lowerLimit = val;
1133 }
1134 }
1135 return lowerLimit;
1136 }
1137 }
1138 return param.getMin();
1139}
1140
1141////////////////////////////////////////////////////////////////////////////////
1142/// Determine the lower limit for param on this interval
1143/// using the binned data set
1144
1146{
1147 if (fHistCutoff < 0)
1148 DetermineByDataHist(); // this initializes fDataHist
1149
1150 if (fHistCutoff < 0) {
1151 // if fHistCutoff < 0 still, then determination of interval failed
1152 coutE(Eval) << "In MCMCInterval::LowerLimitByDataHist: "
1153 << "couldn't determine cutoff. Check that num burn in steps < num "
1154 << "steps in the Markov chain. Returning param.getMin()." << endl;
1155 return param.getMin();
1156 }
1157
1158 for (Int_t d = 0; d < fDimension; d++) {
1159 if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
1160 Int_t numBins = fDataHist->numEntries();
1161 Double_t lowerLimit = param.getMax();
1162 Double_t val;
1163 for (Int_t i = 0; i < numBins; i++) {
1164 fDataHist->get(i);
1165 if (fDataHist->weight() >= fHistCutoff) {
1166 val = fDataHist->get()->getRealValue(param.GetName());
1167 if (val < lowerLimit)
1168 lowerLimit = val;
1169 }
1170 }
1171 return lowerLimit;
1172 }
1173 }
1174 return param.getMin();
1175}
1176
1177////////////////////////////////////////////////////////////////////////////////
1178/// Determine the upper limit for param on this interval
1179/// using the binned data set
1180
1182{
1183 if (fDimension != 1) {
1184 coutE(InputArguments) << "In MCMCInterval::UpperLimitBySparseHist: "
1185 << "Sorry, will not compute upper limit unless dimension == 1" << endl;
1186 return param.getMax();
1187 }
1188 if (fHistCutoff < 0)
1189 DetermineBySparseHist(); // this initializes fSparseHist
1190
1191 if (fHistCutoff < 0) {
1192 // if fHistCutoff < 0 still, then determination of interval failed
1193 coutE(Eval) << "In MCMCInterval::UpperLimitBySparseHist: "
1194 << "couldn't determine cutoff. Check that num burn in steps < num "
1195 << "steps in the Markov chain. Returning param.getMax()." << endl;
1196 return param.getMax();
1197 }
1198
1199 std::vector<Int_t> coord(fDimension);
1200 for (Int_t d = 0; d < fDimension; d++) {
1201 if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
1202 Long_t numBins = (Long_t)fSparseHist->GetNbins();
1203 Double_t upperLimit = param.getMin();
1204 Double_t val;
1205 for (Long_t i = 0; i < numBins; i++) {
1206 if (fSparseHist->GetBinContent(i, &coord[0]) >= fHistCutoff) {
1207 val = fSparseHist->GetAxis(d)->GetBinCenter(coord[d]);
1208 if (val > upperLimit)
1209 upperLimit = val;
1210 }
1211 }
1212 return upperLimit;
1213 }
1214 }
1215 return param.getMax();
1216}
1217
1218////////////////////////////////////////////////////////////////////////////////
1219/// Determine the upper limit for param on this interval
1220/// using the binned data set
1221
1223{
1224 if (fHistCutoff < 0)
1225 DetermineByDataHist(); // this initializes fDataHist
1226
1227 if (fHistCutoff < 0) {
1228 // if fHistCutoff < 0 still, then determination of interval failed
1229 coutE(Eval) << "In MCMCInterval::UpperLimitByDataHist: "
1230 << "couldn't determine cutoff. Check that num burn in steps < num "
1231 << "steps in the Markov chain. Returning param.getMax()." << endl;
1232 return param.getMax();
1233 }
1234
1235 for (Int_t d = 0; d < fDimension; d++) {
1236 if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
1237 Int_t numBins = fDataHist->numEntries();
1238 Double_t upperLimit = param.getMin();
1239 Double_t val;
1240 for (Int_t i = 0; i < numBins; i++) {
1241 fDataHist->get(i);
1242 if (fDataHist->weight() >= fHistCutoff) {
1243 val = fDataHist->get()->getRealValue(param.GetName());
1244 if (val > upperLimit)
1245 upperLimit = val;
1246 }
1247 }
1248 return upperLimit;
1249 }
1250 }
1251 return param.getMax();
1252}
1253
1254////////////////////////////////////////////////////////////////////////////////
1255/// Determine the lower limit for param on this interval
1256/// using the keys pdf
1257
1259{
1260 if (fKeysCutoff < 0)
1262
1263 if (fKeysDataHist == NULL)
1265
1266 if (fKeysCutoff < 0 || fKeysDataHist == NULL) {
1267 // failure in determination of cutoff and/or creation of histogram
1268 coutE(Eval) << "in MCMCInterval::LowerLimitByKeys(): "
1269 << "couldn't find lower limit, check that the number of burn in "
1270 << "steps < number of total steps in the Markov chain. Returning "
1271 << "param.getMin()" << endl;
1272 return param.getMin();
1273 }
1274
1275 for (Int_t d = 0; d < fDimension; d++) {
1276 if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
1277 Int_t numBins = fKeysDataHist->numEntries();
1278 Double_t lowerLimit = param.getMax();
1279 Double_t val;
1280 for (Int_t i = 0; i < numBins; i++) {
1281 fKeysDataHist->get(i);
1282 if (fKeysDataHist->weight() >= fKeysCutoff) {
1283 val = fKeysDataHist->get()->getRealValue(param.GetName());
1284 if (val < lowerLimit)
1285 lowerLimit = val;
1286 }
1287 }
1288 return lowerLimit;
1289 }
1290 }
1291 return param.getMin();
1292}
1293
1294////////////////////////////////////////////////////////////////////////////////
1295/// Determine the upper limit for param on this interval
1296/// using the keys pdf
1297
1299{
1300 if (fKeysCutoff < 0)
1302
1303 if (fKeysDataHist == NULL)
1305
1306 if (fKeysCutoff < 0 || fKeysDataHist == NULL) {
1307 // failure in determination of cutoff and/or creation of histogram
1308 coutE(Eval) << "in MCMCInterval::UpperLimitByKeys(): "
1309 << "couldn't find upper limit, check that the number of burn in "
1310 << "steps < number of total steps in the Markov chain. Returning "
1311 << "param.getMax()" << endl;
1312 return param.getMax();
1313 }
1314
1315 for (Int_t d = 0; d < fDimension; d++) {
1316 if (strcmp(fAxes[d]->GetName(), param.GetName()) == 0) {
1317 Int_t numBins = fKeysDataHist->numEntries();
1318 Double_t upperLimit = param.getMin();
1319 Double_t val;
1320 for (Int_t i = 0; i < numBins; i++) {
1321 fKeysDataHist->get(i);
1322 if (fKeysDataHist->weight() >= fKeysCutoff) {
1323 val = fKeysDataHist->get()->getRealValue(param.GetName());
1324 if (val > upperLimit)
1325 upperLimit = val;
1326 }
1327 }
1328 return upperLimit;
1329 }
1330 }
1331 return param.getMax();
1332}
1333
1334////////////////////////////////////////////////////////////////////////////////
1335/// Determine the approximate maximum value of the Keys PDF
1336
1338{
1339 if (fKeysCutoff < 0)
1341
1342 if (fKeysDataHist == NULL)
1344
1345 if (fKeysDataHist == NULL) {
1346 // failure in determination of cutoff and/or creation of histogram
1347 coutE(Eval) << "in MCMCInterval::KeysMax(): "
1348 << "couldn't find Keys max value, check that the number of burn in "
1349 << "steps < number of total steps in the Markov chain. Returning 0"
1350 << endl;
1351 return 0;
1352 }
1353
1354 Int_t numBins = fKeysDataHist->numEntries();
1355 Double_t max = 0;
1356 Double_t w;
1357 for (Int_t i = 0; i < numBins; i++) {
1358 fKeysDataHist->get(i);
1359 w = fKeysDataHist->weight();
1360 if (w > max)
1361 max = w;
1362 }
1363
1364 return max;
1365}
1366
1367////////////////////////////////////////////////////////////////////////////////
1368
1370{
1371 if (fHistCutoff < 0)
1373
1374 return fHistCutoff;
1375}
1376
1377////////////////////////////////////////////////////////////////////////////////
1378
1380{
1381 if (fKeysCutoff < 0)
1383
1384 // kbelasco: if fFull hasn't been set (because Keys creation failed because
1385 // fNumBurnInSteps >= fChain->Size()) then this will return infinity, which
1386 // seems ok to me since it will indicate error
1387
1388 return fKeysCutoff / fFull;
1389}
1390
1391////////////////////////////////////////////////////////////////////////////////
1392
1394{
1395 RooAbsReal* integral;
1396 Double_t confLevel;
1397 fCutoffVar->setVal(cutoff);
1399 confLevel = integral->getVal(fParameters) / full;
1400 coutI(Eval) << "cutoff = " << cutoff << ", conf = " << confLevel << endl;
1401 //cout << "tmp: cutoff = " << cutoff << ", conf = " << confLevel << endl;
1402 delete integral;
1403 return confLevel;
1404}
1405
1406////////////////////////////////////////////////////////////////////////////////
1407
1409{
1410 if(fConfidenceLevel == 0)
1411 coutE(InputArguments) << "Error in MCMCInterval::GetPosteriorHist: "
1412 << "confidence level not set " << endl;
1413 if (fHist == NULL)
1414 CreateHist();
1415
1416 if (fHist == NULL)
1417 // if fHist is still NULL, then CreateHist failed
1418 return NULL;
1419
1420 return (TH1*) fHist->Clone("MCMCposterior_hist");
1421}
1422
1423////////////////////////////////////////////////////////////////////////////////
1424
1426{
1427 if (fConfidenceLevel == 0)
1428 coutE(InputArguments) << "Error in MCMCInterval::GetPosteriorKeysPdf: "
1429 << "confidence level not set " << endl;
1430 if (fKeysPdf == NULL)
1431 CreateKeysPdf();
1432
1433 if (fKeysPdf == NULL)
1434 // if fKeysPdf is still NULL, then it means CreateKeysPdf failed
1435 return NULL;
1436
1437 return (RooNDKeysPdf*) fKeysPdf->Clone("MCMCPosterior_keys");
1438}
1439
1440////////////////////////////////////////////////////////////////////////////////
1441
1443{
1444 if (fConfidenceLevel == 0)
1445 coutE(InputArguments) << "MCMCInterval::GetPosteriorKeysProduct: "
1446 << "confidence level not set " << endl;
1447 if (fProduct == NULL) {
1448 CreateKeysPdf();
1450 }
1451
1452 if (fProduct == NULL)
1453 // if fProduct is still NULL, then it means CreateKeysPdf failed
1454 return NULL;
1455
1456 return (RooProduct*) fProduct->Clone("MCMCPosterior_keysproduct");
1457}
1458
1459////////////////////////////////////////////////////////////////////////////////
1460
1462{
1463 // returns list of parameters
1464 return new RooArgSet(fParameters);
1465}
1466
1467////////////////////////////////////////////////////////////////////////////////
1468
1470{
1471 return (TMath::Abs(confLevel - fConfidenceLevel) < fEpsilon);
1472}
1473
1474////////////////////////////////////////////////////////////////////////////////
1475
1477{
1478 return (TMath::Abs(a - b) < TMath::Abs(fDelta * (a + b)/2));
1479}
1480
1481////////////////////////////////////////////////////////////////////////////////
1482
1484{
1485 if (fAxes == NULL)
1486 return;
1487 if (fProduct == NULL)
1489 if (fProduct == NULL)
1490 // if fProduct still NULL, then creation failed
1491 return;
1492
1493 //RooAbsBinning** savedBinning = new RooAbsBinning*[fDimension];
1494 Int_t* savedBins = new Int_t[fDimension];
1495 Int_t i;
1496 Double_t numBins;
1497 RooRealVar* var;
1498
1499 // kbelasco: Note - the accuracy is only increased here if the binning for
1500 // each RooRealVar is uniform
1501
1502 // kbelasco: look into why saving the binnings and replacing them doesn't
1503 // work (replaces with 1 bin always).
1504 // Note: this code modifies the binning for the parameters (if they are
1505 // uniform) and sets them back to what they were. If the binnings are not
1506 // uniform, this code does nothing.
1507
1508 // first scan through fAxes to make sure all binnings are uniform, or else
1509 // we can't change the number of bins because there seems to be an error
1510 // when setting the binning itself rather than just the number of bins
1511 Bool_t tempChangeBinning = true;
1512 for (i = 0; i < fDimension; i++) {
1513 if (!fAxes[i]->getBinning(NULL, false, false).isUniform()) {
1514 tempChangeBinning = false;
1515 break;
1516 }
1517 }
1518
1519 // kbelasco: for 1 dimension this should be fine, but for more dimensions
1520 // the total number of bins in the histogram increases exponentially with
1521 // the dimension, so don't do this above 1-D for now.
1522 if (fDimension >= 2)
1523 tempChangeBinning = false;
1524
1525 if (tempChangeBinning) {
1526 // set high number of bins for high accuracy on lower/upper limit by keys
1527 for (i = 0; i < fDimension; i++) {
1528 var = fAxes[i];
1529 //savedBinning[i] = &var->getBinning("__binning_clone", false, true);
1530 savedBins[i] = var->getBinning(NULL, false, false).numBins();
1531 numBins = (var->getMax() - var->getMin()) / fEpsilon;
1532 var->setBins((Int_t)numBins);
1533 }
1534 }
1535
1536 fKeysDataHist = new RooDataHist("_productDataHist",
1537 "Keys PDF & Heaviside Product Data Hist", fParameters);
1539
1540 if (tempChangeBinning) {
1541 // set the binning back to normal
1542 for (i = 0; i < fDimension; i++)
1543 //fAxes[i]->setBinning(*savedBinning[i], NULL);
1544 //fAxes[i]->setBins(savedBinning[i]->numBins(), NULL);
1545 fAxes[i]->setBins(savedBins[i], NULL);
1546 }
1547
1548 //delete[] savedBinning;
1549 delete[] savedBins;
1550}
1551
1552////////////////////////////////////////////////////////////////////////////////
1553
1555{
1556 // check that the parameters are correct
1557
1558 if (parameterPoint.getSize() != fParameters.getSize() ) {
1559 coutE(Eval) << "MCMCInterval: size is wrong, parameters don't match" << std::endl;
1560 return kFALSE;
1561 }
1562 if ( ! parameterPoint.equals( fParameters ) ) {
1563 coutE(Eval) << "MCMCInterval: size is ok, but parameters don't match" << std::endl;
1564 return kFALSE;
1565 }
1566 return kTRUE;
1567}
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
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define coutI(a)
#define coutW(a)
#define coutE(a)
int Int_t
Definition RtypesCore.h:45
const Bool_t kFALSE
Definition RtypesCore.h:101
long Long_t
Definition RtypesCore.h:54
double Double_t
Definition RtypesCore.h:59
long long Long64_t
Definition RtypesCore.h:80
const Bool_t kTRUE
Definition RtypesCore.h:100
#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:81
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:64
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 ...
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:94
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:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:35
The RooDataHist is a container class to hold N-dimensional binned data.
Definition RooDataHist.h:45
Int_t getIndex(const RooAbsCollection &coord, Bool_t fast=false) const
Calculate bin number of the given coordinates.
double weight(std::size_t i) const
Return weight of i-th bin.
Int_t numEntries() const override
Return the number of bins.
const RooArgSet * get() const override
Get bin centre of current bin.
Definition RooDataHist.h:84
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:36
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:2741
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:149
Double_t GetSumw() const
Definition THnBase.h:211
TAxis * GetAxis(Int_t dim) const
Definition THnBase.h:130
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:41
virtual const char * GetName() const
Returns name of object.
Definition TObject.cxx:429
RooCmdArg SelectVars(const RooArgSet &vars)
RooCmdArg EventRange(Int_t nStart, Int_t nStop)
RooCmdArg NormSet(Args_t &&... argsOrArgSet)
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...
Definition Common.h:18
@ 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
RooDataHist * fDataHist
CompareDataHistBins(RooDataHist *hist)
CompareSparseHistBins(THnSparse *hist)
CompareVectorIndices(MarkovChain *chain, RooRealVar *param)
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345