BATCalculator.h
1 /*
2  * Copyright (C) 2007-2018, the BAT core developer team
3  * All rights reserved.
4  *
5  * Original authors: Gregory Schott and Stefan A. Schmitz with inspiration from
6  * Kyle Cranmer, Lorenzo Moneta, Gregory Schott, and Wouter Verkerke.
7  *
8  * For the licensing terms see doc/COPYING.
9  * For documentation see http://mpp.mpg.de/bat
10  */
11 
12 #ifndef ROOSTATS_BATCalculator
13 #define ROOSTATS_BATCalculator
14 
15 #include "BCRooInterface.h"
16 
17 #include <TH1D.h>
18 #include <TNamed.h>
19 
20 #include <RooStats/IntervalCalculator.h>
21 #include <RooStats/MCMCInterval.h>
22 #include <RooStats/ModelConfig.h>
23 #include <RooStats/SimpleInterval.h>
24 
25 #include <RooAbsData.h>
26 #include <RooAbsPdf.h>
27 #include <RooAbsReal.h>
28 #include <RooArgSet.h>
29 #include <RooPlot.h>
30 
31 #include <map>
32 #include <vector>
33 
34 namespace RooStats
35 {
36 
39 class BATCalculator : public IntervalCalculator, public TNamed
40 {
41 
42 public:
43 
45  BATCalculator( );
46 
48  BATCalculator( RooAbsData& data,
49  RooAbsPdf& pdf,
50  RooArgSet& POI,
51  RooAbsPdf& prior,
52  RooArgSet* params = 0,
53  bool fillChain = false );
54 
55  BATCalculator( RooAbsData& data,
56  ModelConfig& model,
57  bool fillChain = false);
58 
60  virtual ~BATCalculator();
61 
62  RooPlot* GetPosteriorPlot1D() const;
63 
65  RooAbsPdf* GetPosteriorPdf1D() const;
66  RooAbsPdf* GetPosteriorPdf1D(const char* POIname) const;
67 
69  virtual SimpleInterval* GetInterval1D() const;
70  virtual SimpleInterval* GetInterval1D(const char* POIname) const;
71 
73  SimpleInterval* GetShortestInterval1D() const;
74  SimpleInterval* GetShortestInterval1D(const char* POIname, bool& checkConnected) const;
75 
77  Double_t GetOneSidedUperLim();
78 
79  virtual void SetData( RooAbsData& data )
80  { fData = &data; ClearAll(); }
81 
82  virtual void SetModel(const ModelConfig&) {}
83 
85  virtual void SetTestSize( Double_t size )
86  { fSize = size; fValidInterval = false; }
87 
89  void SetLeftSideTailFraction(Double_t leftSideFraction );
90 
92  virtual void SetConfidenceLevel( Double_t cl )
93  { SetTestSize(1. - cl); }
94 
96  virtual Double_t Size() const
97  { return fSize; }
100  { return fLeftSideFraction; }
101 
103  virtual Double_t ConfidenceLevel() const
104  { return 1. - fSize; }
105 
106  void SetBrfPrecision( double precision )
107  { fBrfPrecision = precision; }
108 
109  double GetBrfPrecision()
110  { return fBrfPrecision; }
111 
113  void SetnMCMC(int nMCMC)
114  { _nMCMC = nMCMC; }
115 
117  int GetnMCMC()
118  { return _nMCMC; }
119 
120  BCRooInterface* GetBCRooInterface() const
121  { return _myRooInterface; }
122 
123  RooStats::MarkovChain* GetRooStatsMarkovChain() const
124  { return _myRooInterface->GetRooStatsMarkovChain();}
125 
127  virtual MCMCInterval* GetInterval() const;
128 
131  { return fConnectedInterval; }
132 
134  std::vector<double> GetIntervalBorders1D()
135  { return _intervalBorders1D; }
136 
138  void SetNumBins(const char* parname, int nbins);
139 
141  void SetNumBins(int nbins);
142 
143  void CleanCalculatorForNewData()
144  { ClearAll(); }
145 
146 protected:
147 
148  void ClearAll() const;
149 
151 private:
152  RooAbsData* fData;
153  RooAbsPdf* fPdf;
154  const RooArgSet fPOI;
155  RooAbsPdf* fPrior;
156  const RooArgSet* fparams;
157  BCRooInterface* _myRooInterface;
158  mutable TH1D* _posteriorTH1D;
159 
160 
161  mutable RooAbsPdf* fProductPdf;
162  mutable RooAbsReal* fLogLike;
163  mutable RooAbsReal* fLikelihood;
164  mutable RooAbsReal* fIntegratedLikelihood;
165  mutable RooAbsPdf* fPosteriorPdf;
166  mutable Double_t fLower;
167  mutable Double_t fUpper;
168  double fBrfPrecision;
169  mutable Bool_t fValidInterval;
170  mutable Bool_t fValidMCMCInterval;
171  mutable bool fConnectedInterval;
172 
173  int _nMCMC;
174  double fSize;
175  double fLeftSideFraction; //
176  mutable std::vector<double> _intervalBorders1D;
177 
178 protected:
179 
180  ClassDef(BATCalculator, 1) // BATCalculator class
181 
182 };
183 }
184 
185 #endif
double GetLeftSideTailFraction()
Get left side tail fraction (only for 1D interval, not meaningful for shortest interval) ...
Definition: BATCalculator.h:99
BATCalculator()
constructor
Interface allowing to run BAT on a problem/data defined in a standard RooFit workspace format...
virtual ~BATCalculator()
destructor
bool GetConnected()
returns if last calculated shortest interval is connected (1 poi only)
virtual void SetConfidenceLevel(Double_t cl)
set the confidence level for the interval (eg. 0.95 for a 95% Confidence Interval) ...
Definition: BATCalculator.h:92
void SetnMCMC(int nMCMC)
set number of iterations per chain
Double_t GetOneSidedUperLim()
temporary solution ?
std::vector< double > GetIntervalBorders1D()
returns interval borders of shortest interval (1 poi only)
virtual void SetTestSize(Double_t size)
set the size of the test (rate of Type I error) ( Eg. 0.05 for a 95% Confidence Interval) ...
Definition: BATCalculator.h:85
RooAbsPdf * GetPosteriorPdf1D() const
return posterior pdf
int GetnMCMC()
return number of iterations per chain
SimpleInterval * GetShortestInterval1D() const
return SimpleInterval object: shortest interval (not necessarily connected, one poi only) ...
Adapter to call BAT from within RooStats to compute credibility intervals.
Definition: BATCalculator.h:39
virtual Double_t Size() const
Get the size of the test.
Definition: BATCalculator.h:96
virtual SimpleInterval * GetInterval1D() const
return SimpleInterval object: central interval (one poi only)
virtual MCMCInterval * GetInterval() const
returns MCMCInterval object
void SetLeftSideTailFraction(Double_t leftSideFraction)
set left side tail fraction (only for 1D interval, not meaningful for shortest interval) ...
RooStats::MarkovChain * GetRooStatsMarkovChain()
return the RooStats Markov Chain (empty if corresponding constructor option not set) ...
virtual Double_t ConfidenceLevel() const
Get the Confidence level for the test.
void SetNumBins(const char *parname, int nbins)
set the number of histogram bins for a specific parameter