18 #ifndef ROOSTATS_MCMCIntervalPlot
33 #ifndef ROOT_TObjArray
63 #ifndef ROO_GLOBAL_FUNC
81 using namespace RooStats;
83 MCMCIntervalPlot::MCMCIntervalPlot()
87 fPosteriorHist =
NULL;
88 fPosteriorKeysPdf =
NULL;
89 fPosteriorKeysProduct =
NULL;
103 fPosteriorHistHistCopy =
NULL;
104 fPosteriorHistTFCopy =
NULL;
107 MCMCIntervalPlot::MCMCIntervalPlot(MCMCInterval& interval)
109 SetMCMCInterval(interval);
110 fPosteriorHist =
NULL;
111 fPosteriorKeysPdf =
NULL;
112 fPosteriorKeysProduct =
NULL;
125 fPosteriorHistHistCopy =
NULL;
126 fPosteriorHistTFCopy =
NULL;
129 MCMCIntervalPlot::~MCMCIntervalPlot()
140 delete fPosteriorKeysPdf;
141 delete fPosteriorKeysProduct;
150 void MCMCIntervalPlot::SetMCMCInterval(MCMCInterval& interval)
152 fInterval = &interval;
153 fDimension = fInterval->GetDimension();
154 fParameters = fInterval->GetParameters();
159 DrawInterval(options);
162 void MCMCIntervalPlot::DrawPosterior(
const Option_t* options)
164 if (fInterval->GetUseKeys())
165 DrawPosteriorKeysPdf(options);
167 DrawPosteriorHist(options);
170 void* MCMCIntervalPlot::DrawPosteriorHist(
const Option_t* ,
173 if (fPosteriorHist ==
NULL)
174 fPosteriorHist = fInterval->GetPosteriorHist();
176 if (fPosteriorHist ==
NULL) {
178 <<
"Couldn't get posterior histogram." << endl;
192 fPosteriorHist->Scale(1/fPosteriorHist->GetBinContent(
193 fPosteriorHist->GetMaximumBin()));
198 fPosteriorHist->SetTitle(title);
200 fPosteriorHist->SetTitle(GetTitle());
204 return (
void*)fPosteriorHist;
207 void* MCMCIntervalPlot::DrawPosteriorKeysPdf(
const Option_t* options)
209 if (fPosteriorKeysPdf ==
NULL)
210 fPosteriorKeysPdf = fInterval->GetPosteriorKeysPdf();
212 if (fPosteriorKeysPdf ==
NULL) {
214 <<
"Couldn't get posterior Keys PDF." << endl;
221 if (fDimension == 1) {
226 <<
"Invalid parameter" << endl;
238 }
else if (fDimension == 2) {
242 TH2F* keysHist = (TH2F*)fPosteriorKeysPdf->createHistogram(
246 Form(
"MCMC histogram of posterior Keys PDF for %s, %s",
249 keysHist->SetTitle(GetTitle());
251 keysHist->Draw(options);
258 void MCMCIntervalPlot::DrawInterval(
const Option_t* options)
260 switch (fInterval->GetIntervalType()) {
261 case MCMCInterval::kShortest:
262 DrawShortestInterval(options);
264 case MCMCInterval::kTailFraction:
265 DrawTailFractionInterval(options);
269 "Interval type not supported" << endl;
274 void MCMCIntervalPlot::DrawShortestInterval(
const Option_t* options)
276 if (fInterval->GetUseKeys())
277 DrawKeysPdfInterval(options);
279 DrawHistInterval(options);
282 void MCMCIntervalPlot::DrawKeysPdfInterval(
const Option_t* options)
287 if (fDimension == 1) {
295 Double_t height = fInterval->GetKeysMax();
298 Double_t ul = fInterval->UpperLimitByKeys(*p);
299 Double_t ll = fInterval->LowerLimitByKeys(*p);
301 if (frame !=
NULL && fPosteriorKeysPdf !=
NULL) {
309 fPosteriorKeysPdf->plotOn(frame,
320 fPosteriorKeysPdf->plotOn(frame,
324 frame->
Draw(options);
327 TLine* llLine =
new TLine(ll, 0, ll, height);
328 TLine* ulLine =
new TLine(ul, 0, ul, height);
329 llLine->SetLineColor(fLineColor);
330 ulLine->SetLineColor(fLineColor);
331 llLine->SetLineWidth(fLineWidth);
332 ulLine->SetLineWidth(fLineWidth);
333 llLine->Draw(options);
334 ulLine->Draw(options);
335 }
else if (fDimension == 2) {
336 if (fPosteriorKeysPdf ==
NULL)
337 fPosteriorKeysPdf = fInterval->GetPosteriorKeysPdf();
339 if (fPosteriorKeysPdf ==
NULL) {
341 <<
"Couldn't get posterior Keys PDF." << endl;
348 TH2F* contHist = (TH2F*)fPosteriorKeysPdf->createHistogram(
358 contHist->SetTitle(
NULL);
360 contHist->SetStats(
kFALSE);
365 Double_t cutoff = fInterval->GetKeysPdfCutoff();
366 contHist->SetContour(1, &cutoff);
367 contHist->SetLineColor(fLineColor);
368 contHist->SetLineWidth(fLineWidth);
369 contHist->Draw(tmpOpt.
Data());
373 <<
" Sorry: " << fDimension <<
"-D plots not currently supported" << endl;
377 void MCMCIntervalPlot::DrawHistInterval(
const Option_t* options)
382 if (fDimension == 1) {
385 Double_t ul = fInterval->UpperLimitByHist(*p);
386 Double_t ll = fInterval->LowerLimitByHist(*p);
391 TH1F* hist = (TH1F*)DrawPosteriorHist(options,
NULL,
false);
392 if (hist ==
NULL)
return;
396 hist->SetTitle(GetTitle());
397 hist->GetYaxis()->SetTitle(
Form(
"Posterior for parameter %s",
400 TH1F* copy = (TH1F*)hist->Clone(
Form(
"%s_copy", hist->GetTitle()));
401 Double_t histCutoff = fInterval->GetHistCutoff();
404 Int_t nBins = copy->GetNbinsX();
406 for (i = 1; i <= nBins; i++) {
408 height = copy->GetBinContent(i);
409 if (height < histCutoff)
410 copy->SetBinContent(i, 0);
413 hist->Scale(1/hist->GetBinContent(hist->GetMaximumBin()));
414 copy->Scale(1/copy->GetBinContent(hist->GetMaximumBin()));
416 copy->SetFillStyle(1001);
417 copy->SetFillColor(fShadeColor);
421 fPosteriorHistHistCopy = copy;
423 TLine* llLine =
new TLine(ll, 0, ll, 1);
424 TLine* ulLine =
new TLine(ul, 0, ul, 1);
425 llLine->SetLineColor(fLineColor);
426 ulLine->SetLineColor(fLineColor);
427 llLine->SetLineWidth(fLineWidth);
428 ulLine->SetLineWidth(fLineWidth);
429 llLine->Draw(options);
430 ulLine->Draw(options);
432 }
else if (fDimension == 2) {
433 if (fPosteriorHist ==
NULL)
434 fPosteriorHist = fInterval->GetPosteriorHist();
436 if (fPosteriorHist ==
NULL) {
438 <<
"Couldn't get posterior histogram." << endl;
450 fPosteriorHist->SetTitle(GetTitle());
452 fPosteriorHist->SetTitle(
NULL);
455 fPosteriorHist->SetStats(
kFALSE);
460 Double_t cutoff = fInterval->GetHistCutoff();
461 fPosteriorHist->SetContour(1, &cutoff);
462 fPosteriorHist->SetLineColor(fLineColor);
463 fPosteriorHist->SetLineWidth(fLineWidth);
464 fPosteriorHist->Draw(tmpOpt.
Data());
467 <<
" Sorry: " << fDimension <<
"-D plots not currently supported" << endl;
471 void MCMCIntervalPlot::DrawTailFractionInterval(
const Option_t* options)
476 if (fDimension == 1) {
480 Double_t ul = fInterval->UpperLimitTailFraction(*p);
481 Double_t ll = fInterval->LowerLimitTailFraction(*p);
483 TH1F* hist = (TH1F*)DrawPosteriorHist(options,
NULL,
false);
484 if (hist ==
NULL)
return;
488 hist->SetTitle(GetTitle());
489 hist->GetYaxis()->SetTitle(
Form(
"Posterior for parameter %s",
492 TH1F* copy = (TH1F*)hist->Clone(
Form(
"%s_copy", hist->GetTitle()));
495 Int_t nBins = copy->GetNbinsX();
497 for (i = 1; i <= nBins; i++) {
499 center = copy->GetBinCenter(i);
500 if (center < ll || center > ul)
501 copy->SetBinContent(i, 0);
504 hist->Scale(1/hist->GetBinContent(hist->GetMaximumBin()));
505 copy->Scale(1/copy->GetBinContent(hist->GetMaximumBin()));
507 copy->SetFillStyle(1001);
508 copy->SetFillColor(fShadeColor);
513 TLine* llLine =
new TLine(ll, 0, ll, 1);
514 TLine* ulLine =
new TLine(ul, 0, ul, 1);
515 llLine->SetLineColor(fLineColor);
516 ulLine->SetLineColor(fLineColor);
517 llLine->SetLineWidth(fLineWidth);
518 ulLine->SetLineWidth(fLineWidth);
519 llLine->Draw(options);
520 ulLine->Draw(options);
523 <<
" Sorry: " << fDimension <<
"-D plots not currently supported"
528 void* MCMCIntervalPlot::DrawPosteriorKeysProduct(
const Option_t* options)
530 if (fPosteriorKeysProduct ==
NULL)
531 fPosteriorKeysProduct = fInterval->GetPosteriorKeysProduct();
533 if (fPosteriorKeysProduct ==
NULL) {
535 <<
"Couldn't get posterior Keys product." << endl;
544 if (fDimension == 1) {
546 if (!frame)
return NULL;
548 frame->
SetTitle(
Form(
"Posterior Keys PDF * Heaviside product for %s",
553 fPosteriorKeysProduct->plotOn(frame,
555 frame->
Draw(options);
557 }
else if (fDimension == 2) {
560 TH2F* productHist = (TH2F*)fPosteriorKeysProduct->createHistogram(
564 Form(
"MCMC Posterior Keys Product Hist. for %s, %s",
567 productHist->SetTitle(GetTitle());
568 productHist->Draw(options);
577 const MarkovChain* markovChain = fInterval->GetChain();
582 burnInSteps = fInterval->GetNumBurnInSteps();
590 if (burnInSteps > 0) {
591 burnInX =
new Double_t[burnInSteps];
592 burnInY =
new Double_t[burnInSteps];
597 for (
Int_t i = burnInSteps; i < size; i++) {
602 for (
Int_t i = 0; i < burnInSteps; i++) {
613 TGraph* walk =
new TGraph(size - burnInSteps, x, y);
615 walk->SetTitle(
Form(
"2-D Scatter Plot of Markov chain for %s, %s",
618 walk->SetTitle(GetTitle());
621 walk->GetXaxis()->SetTitle(xVar.
GetName());
623 walk->GetYaxis()->SetTitle(yVar.
GetName());
624 walk->SetLineColor(
kGray);
625 walk->SetMarkerStyle(6);
627 walk->Draw(
"A,L,P,same");
629 TGraph* burnIn =
NULL;
630 if (burnInX !=
NULL && burnInY !=
NULL) {
631 burnIn =
new TGraph(burnInSteps - 1, burnInX, burnInY);
632 burnIn->SetLineColor(
kPink);
633 burnIn->SetMarkerStyle(6);
634 burnIn->SetMarkerColor(
kPink);
635 burnIn->Draw(
"L,P,same");
638 TGraph*
first =
new TGraph(1, &firstX, &firstY);
639 first->SetLineColor(
kGreen);
640 first->SetMarkerStyle(3);
641 first->SetMarkerSize(2);
642 first->SetMarkerColor(
kGreen);
643 first->Draw(
"L,P,same");
648 if (burnInX !=
NULL)
delete [] burnInX;
649 if (burnInY !=
NULL)
delete [] burnInY;
655 void MCMCIntervalPlot::DrawParameterVsTime(
RooRealVar& param)
657 const MarkovChain* markovChain = fInterval->GetChain();
659 Int_t numEntries = 2 * size;
665 for (
Int_t i = 0; i < size; i++) {
669 value[2*i + 1] = val;
678 TGraph* paramGraph =
new TGraph(numEntries, time, value);
680 paramGraph->SetTitle(
Form(
"%s vs. time in Markov chain",param.
GetName()));
682 paramGraph->SetTitle(GetTitle());
683 paramGraph->GetXaxis()->SetTitle(
"Time (discrete steps)");
684 paramGraph->GetYaxis()->SetTitle(param.
GetName());
686 paramGraph->Draw(
"A,L,same");
692 void MCMCIntervalPlot::DrawNLLVsTime()
694 const MarkovChain* markovChain = fInterval->GetChain();
696 Int_t numEntries = 2 * size;
702 for (
Int_t i = 0; i < size; i++) {
703 nll = markovChain->
NLL(i);
706 nllValue[2*i + 1] = nll;
715 TGraph* nllGraph =
new TGraph(numEntries, time, nllValue);
717 nllGraph->SetTitle(
"NLL value vs. time in Markov chain");
719 nllGraph->SetTitle(GetTitle());
720 nllGraph->GetXaxis()->SetTitle(
"Time (discrete steps)");
721 nllGraph->GetYaxis()->SetTitle(
"NLL (-log(likelihood))");
723 nllGraph->Draw(
"A,L,same");
729 void MCMCIntervalPlot::DrawNLLHist(
const Option_t* options)
731 if (fNLLHist ==
NULL) {
732 const MarkovChain* markovChain = fInterval->GetChain();
736 for (
Int_t i = 0; i < size; i++)
737 if (markovChain->
NLL(i) > maxNLL)
738 maxNLL = markovChain->
NLL(i);
740 fNLLHist =
new TH1F(
"mcmc_nll_hist",
"MCMC NLL Histogram",
743 Bool_t isEmpty = (title.CompareTo(
"") == 0);
745 fNLLHist->SetTitle(GetTitle());
746 fNLLHist->GetXaxis()->SetTitle(
"-log(likelihood)");
747 for (
Int_t i = 0; i < size; i++)
748 fNLLHist->Fill(markovChain->
NLL(i), markovChain->
Weight());
750 fNLLHist->Draw(options);
753 void MCMCIntervalPlot::DrawWeightHist(
const Option_t* options)
755 if (fWeightHist ==
NULL) {
756 const MarkovChain* markovChain = fInterval->GetChain();
760 for (
Int_t i = 0; i < size; i++)
761 if (markovChain->
Weight(i) > maxWeight)
762 maxWeight = markovChain->
Weight(i);
763 fWeightHist =
new TH1F(
"mcmc_weight_hist",
"MCMC Weight Histogram",
764 (
Int_t)(maxWeight + 1), 0, maxWeight * 1.02);
765 for (
Int_t i = 0; i < size; i++)
766 fWeightHist->Fill(markovChain->
Weight(i));
768 fWeightHist->Draw(options);
773 // 3-d plot of the parameter points
775 // also plot the points in the markov chain
776 RooDataSet* markovChainData = ((MCMCInterval*)mcmcint)->GetChainAsDataSet();
778 TTree& chain = ((RooTreeDataStore*) markovChainData->store())->tree();
779 chain.SetMarkerStyle(6);
780 chain.SetMarkerColor(kRed);
781 chain.Draw("s:ratioSigEff:ratioBkgEff","","box"); // 3-d box proporional to posterior
783 // the points used in the profile construction
784 TTree& parameterScan = ((RooTreeDataStore*) fc.GetPointsToScan()->store())->tree();
785 parameterScan.SetMarkerStyle(24);
786 parameterScan.Draw("s:ratioSigEff:ratioBkgEff","","same");
788 chain.SetMarkerStyle(6);
789 chain.SetMarkerColor(kRed);
790 //chain.Draw("s:ratioSigEff:ratioBkgEff", "_MarkovChain_local_nll","box");
791 //chain.Draw("_MarkovChain_local_nll");
ClassImp(RooStats::MCMCIntervalPlot)
RooCmdArg DrawOption(const char *opt)
virtual Int_t numBins(const char *rangeName=0) const
This class provides simple and straightforward utilities to plot a MCMCInterval object.
virtual Double_t NLL(Int_t i) const
get the NLL value of entry at position i
virtual Double_t getMin(const char *name=0) const
void SetTitle(const char *name)
Set the title of the RooPlot to 'title'.
const char * Data() const
RooCmdArg Range(const char *rangeName, Bool_t adjustNorm=kTRUE)
virtual Int_t getBins(const char *name=0) const
TString & Append(const char *cs)
RooRealVar represents a fundamental (non-derived) real valued object.
TPaveLabel title(3, 27.1, 15, 28.7,"ROOT Environment and Tools")
char * Form(const char *fmt,...)
virtual const char * GetName() const
Returns name of object.
RooPlot * frame(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
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Stores the steps in a Markov Chain of points.
RooAbsArg * at(Int_t idx) const
RooCmdArg FillColor(Color_t color)
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
RooCmdArg Normalization(Double_t scaleFactor)
RooCmdArg Scaling(Bool_t flag)
virtual Double_t getMax(const char *name=0) const
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
virtual Int_t Size() const
get the number of steps in the chain
virtual void SetTitle(const char *title="")
Change (i.e. set) the title of the TNamed.
virtual const RooArgSet * Get(Int_t i) const
get the entry at position i
int CompareTo(const char *cs, ECaseCompare cmp=kExact) const
Compare a string to char *cs2.
virtual Double_t Weight() const
get the weight of the current (last indexed) entry
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
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.