Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
BayesianCalculator.h
Go to the documentation of this file.
1 // @(#)root/roostats:$Id$
2 // Author: Kyle Cranmer, Lorenzo Moneta, Gregory Schott, Wouter Verkerke
3 /*************************************************************************
4  * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
5  * All rights reserved. *
6  * *
7  * For the licensing terms see $ROOTSYS/LICENSE. *
8  * For the list of contributors see $ROOTSYS/README/CREDITS. *
9  *************************************************************************/
10 
11 #ifndef ROOSTATS_BayesianCalculator
12 #define ROOSTATS_BayesianCalculator
13 
14 #include "TNamed.h"
15 
16 #include "Math/IFunctionfwd.h"
17 
18 #include "RooArgSet.h"
19 
21 
23 
24 class RooAbsData;
25 class RooAbsPdf;
26 class RooPlot;
27 class RooAbsReal;
28 class TF1;
29 class TH1;
30 
31 
32 namespace RooStats {
33 
34  class ModelConfig;
35  class SimpleInterval;
36 
37  class BayesianCalculator : public IntervalCalculator, public TNamed {
38 
39  public:
40 
41  // constructor
42  BayesianCalculator( );
43 
44  BayesianCalculator( RooAbsData& data,
45  RooAbsPdf& pdf,
46  const RooArgSet& POI,
47  RooAbsPdf& priorPdf,
48  const RooArgSet* nuisanceParameters = 0 );
49 
50  BayesianCalculator( RooAbsData& data,
51  ModelConfig& model );
52 
53  // destructor
54  virtual ~BayesianCalculator();
55 
56  // get the plot with option to get it normalized
57  RooPlot* GetPosteriorPlot(bool norm = false, double precision = 0.01) const;
58 
59  // return posterior pdf (object is managed by the user)
60  RooAbsPdf* GetPosteriorPdf() const;
61  // return posterior function (object is managed by the BayesianCalculator class)
62  RooAbsReal* GetPosteriorFunction() const;
63 
64  // return the approximate posterior as histogram (TH1 object). Note the object is managed by the BayesianCalculator class
65  TH1 * GetPosteriorHistogram() const;
66 
67  // compute the interval. By Default a central interval is computed
68  // By using SetLeftTileFraction can control if central/ upper/lower interval
69  // For shortest interval use SetShortestInterval(true)
70  virtual SimpleInterval* GetInterval() const ;
71 
72  virtual void SetData( RooAbsData & data ) {
73  fData = &data;
74  ClearAll();
75  }
76 
77 
78  // set the model via the ModelConfig
79  virtual void SetModel( const ModelConfig& model );
80 
81  // specify the parameters of interest in the interval
82  virtual void SetParameters(const RooArgSet& set) { fPOI.removeAll(); fPOI.add(set); }
83 
84  // specify the nuisance parameters (eg. the rest of the parameters)
85  virtual void SetNuisanceParameters(const RooArgSet& set) {fNuisanceParameters.removeAll(); fNuisanceParameters.add(set);}
86 
87  // Set only the Prior Pdf
88  virtual void SetPriorPdf(RooAbsPdf& pdf) { fPriorPdf = &pdf; }
89 
90  // set the conditional observables which will be used when creating the NLL
91  // so the pdf's will not be normalized on the conditional observables when computing the NLL
92  virtual void SetConditionalObservables(const RooArgSet& set) {fConditionalObs.removeAll(); fConditionalObs.add(set);}
93 
94  // set the global observables which will be used when creating the NLL
95  // so the constraint pdf's will be normalized correctly on the global observables when computing the NLL
96  virtual void SetGlobalObservables(const RooArgSet& set) {fGlobalObs.removeAll(); fGlobalObs.add(set);}
97 
98  // set the size of the test (rate of Type I error) ( Eg. 0.05 for a 95% Confidence Interval)
99  virtual void SetTestSize( Double_t size ) {
100  fSize = size;
101  fValidInterval = false;
102  }
103  // set the confidence level for the interval (eg. 0.95 for a 95% Confidence Interval)
104  virtual void SetConfidenceLevel( Double_t cl ) { SetTestSize(1.-cl); }
105  // Get the size of the test (eg. rate of Type I error)
106  virtual Double_t Size() const { return fSize; }
107  // Get the Confidence level for the test
108  virtual Double_t ConfidenceLevel() const { return 1.-fSize; }
109 
110  // set the fraction of probability content on the left tail
111  // Central limits use 0.5 (default case)
112  // for upper limits it is 0 and 1 for lower limit
113  // For shortest intervals a negative value (i.e. -1) must be given
114  void SetLeftSideTailFraction(Double_t leftSideFraction ) {fLeftSideFraction = leftSideFraction;}
115 
116  // set the Bayesian calculator to compute the shortest interval (default is central interval)
117  // to switch off SetLeftSideTailFraction to the right value
118  void SetShortestInterval() { fLeftSideFraction = -1; }
119 
120  // set the precision of the Root Finder
121  void SetBrfPrecision( double precision ) { fBrfPrecision = precision; }
122 
123  // use directly the approximate posterior function obtained by binning it in nbins
124  // by default the cdf is used by integrating the posterior
125  // if a value of nbin <= 0 the cdf function will be used
126  void SetScanOfPosterior(int nbin = 100) { fNScanBins = nbin; }
127 
128  // set the number of iterations when running a MC integration algorithm
129  // If not set use default algorithmic values
130  // In case of ToyMC sampling of the nuisance the value is 100
131  // In case of using the GSL MCintegrations types the default value is
132  // defined in ROOT::Math::IntegratorMultiDimOptions::DefaultNCalls()
133  virtual void SetNumIters(Int_t numIters) { fNumIterations = numIters; }
134 
135  // set the integration type (possible type are) :
136  void SetIntegrationType(const char * type);
137 
138  // return the mode (most probable value of the posterior function)
139  double GetMode() const;
140 
141  // force the nuisance pdf when using the toy mc sampling
142  void ForceNuisancePdf(RooAbsPdf & pdf) { fNuisancePdf = &pdf; }
143 
144  protected:
145 
146  void ClearAll() const;
147 
148  void ApproximatePosterior() const;
149 
150  void ComputeIntervalFromApproxPosterior(double c1, double c2) const;
151 
152  void ComputeIntervalFromCdf(double c1, double c2) const;
153 
154  void ComputeIntervalUsingRooFit(double c1, double c2) const;
155 
156  void ComputeShortestInterval() const;
157 
158  private:
159 
160  // plan to replace the above: return a SimpleInterval integrating
161  // over all other parameters except the one specified as argument
162  // virtual SimpleInterval* GetInterval( RooRealVar* parameter ) const { return 0; }
163 
164  RooAbsData* fData; // data set
165  RooAbsPdf* fPdf; // model pdf (could contain the nuisance pdf as constraint term)
166  RooArgSet fPOI; // POI
167  RooAbsPdf* fPriorPdf; // prior pdf (typically for the POI)
168  RooAbsPdf* fNuisancePdf; // nuisance pdf (needed when using nuisance sampling technique)
169  RooArgSet fNuisanceParameters; // nuisance parameters
170  RooArgSet fConditionalObs ; // conditional observables
171  RooArgSet fGlobalObs; // global observables
172 
173  mutable RooAbsPdf* fProductPdf; // internal pointer to model * prior
174  mutable RooAbsReal* fLogLike; // internal pointer to log likelihood function
175  mutable RooAbsReal* fLikelihood; // internal pointer to likelihood function
176  mutable RooAbsReal* fIntegratedLikelihood; // integrated likelihood function, i.e - unnormalized posterior function
177  mutable RooAbsPdf* fPosteriorPdf; // normalized (on the poi) posterior pdf
178  mutable ROOT::Math::IGenFunction * fPosteriorFunction; // function representing the posterior
179  mutable TF1 * fApproxPosterior; // TF1 representing the scanned posterior function
180  mutable Double_t fLower; // computer lower interval bound
181  mutable Double_t fUpper; // upper interval bound
182  mutable Double_t fNLLMin; // minimum value of Nll
183  double fSize; // size used for getting the interval
184  double fLeftSideFraction; // fraction of probability content on left side of interval
185  double fBrfPrecision; // root finder precision
186  mutable int fNScanBins; // number of bins to scan, if = -1 no scan is done (default)
187  int fNumIterations; // number of iterations (when using ToyMC)
188  mutable Bool_t fValidInterval;
189 
190  TString fIntegrationType;
191 
192  protected:
193 
194  ClassDef(BayesianCalculator,2) // BayesianCalculator class
195 
196  };
197 }
198 
199 #endif