Logo ROOT   6.12/07
Reference Guide
AsymptoticCalculator.h
Go to the documentation of this file.
1 // @(#)root/roostats:$Id$
2 // Author: Sven Kreiss, Kyle Cranmer Nov 2010
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_AsymptoticCalculator
12 #define ROOSTATS_AsymptoticCalculator
13 
15 #include "RooArgSet.h"
16 #include "Rtypes.h"
17 
18 class RooArgList;
19 class RooCategory;
20 class RooRealVar;
21 class RooPoisson;
22 class RooProdPdf;
23 
24 
25 namespace RooStats {
26 
28 
29  public:
31  RooAbsData &data, // need to pass non-const since RooAbsPdf::fitTo takes a non-const data set
32  const ModelConfig &altModel,
33  const ModelConfig &nullModel,
34  bool nominalAsimov = false
35  );
36  // HypoTestCalculatorGeneric(data, altModel, nullModel, 0)
37  // {
38  // }
39 
41  }
42 
43  /// initialize the calculator by performing a global fit and make the Asimov data set
44  bool Initialize() const;
45 
46  /// re-implement HypoTest computation using the asymptotic
47  virtual HypoTestResult *GetHypoTest() const;
48 
49  /// make the asimov data from the ModelConfig and list of poi - return data set and snapshot of global obs
50  /// poiValues is the snapshot of POI used for finding the best nuisance parameter values (conditioned at these values)
51  /// genPoiValues is optionally a different set of POI values used for generating. By default the same POI are used for generating and for finding the nuisance parameters
52  static RooAbsData * MakeAsimovData( RooAbsData & data, const ModelConfig & model, const RooArgSet & poiValues, RooArgSet & globObs, const RooArgSet * genPoiValues = 0);
53 
54 
55  /// make a nominal asimov data from the ModelConfig and parameter values
56  /// The parameter values (including the nuisance) could be given from a fit to data or be at the nominal values
57  static RooAbsData * MakeAsimovData( const ModelConfig & model, const RooArgSet & allParamValues, RooArgSet & globObs);
58 
59 
60 
61  static RooAbsData * GenerateAsimovData(const RooAbsPdf & pdf, const RooArgSet & observables );
62 
63  /// function given the null and the alt p value - return the expected one given the N - sigma value
64  static double GetExpectedPValues(double pnull, double palt, double nsigma, bool usecls, bool oneSided = true );
65 
66  /// set test statistic for one sided (upper limits)
67  void SetOneSided(bool on) { fOneSided = on; }
68 
69  /// set the test statistics for two sided (in case of upper limits
70  /// for discovery does not make really sense)
71  void SetTwoSided() { fOneSided = false; fOneSidedDiscovery = false;}
72 
73  /// set the test statistics for one-sided discovery
74  void SetOneSidedDiscovery(bool on) { fOneSidedDiscovery = on; }
75 
76  /// re-implementation of setters since they are needed to re-initialize the calculator
77  virtual void SetNullModel(const ModelConfig &nullModel) {
79  fIsInitialized = false;
80  }
81  virtual void SetAlternateModel(const ModelConfig &altModel) {
83  fIsInitialized = false;
84  }
85  virtual void SetData(RooAbsData &data) {
87  fIsInitialized = false;
88  }
89 
90 
91  bool IsTwoSided() const { return (!fOneSided && !fOneSidedDiscovery); }
92  bool IsOneSidedDiscovery() const { return fOneSidedDiscovery; }
93 
94 
95  /// set using of qtilde, by default is controlled if RoORealVar is limited or not
96  void SetQTilde(bool on) { fUseQTilde = on; }
97 
98  /// return snapshot of the best fit parameter
99  const RooArgSet & GetBestFitPoi() const { return fBestFitPoi; }
100  /// return best fit parameter (firs of poi)
101  const RooRealVar * GetMuHat() const { return dynamic_cast<RooRealVar*>(fBestFitPoi.first()); }
102  /// return best fit value for all parameters
103  const RooArgSet & GetBestFitParams() const { return fBestFitPoi; }
104 
105  static void SetPrintLevel(int level);
106 
107  protected:
108  // // configure TestStatSampler for the Null run
109  // int PreNullHook(RooArgSet *parameterPoint, double obsTestStat) const;
110 
111  // // configure TestStatSampler for the Alt run
112  // int PreAltHook(RooArgSet *parameterPoint, double obsTestStat) const;
113 
114 
115  static RooAbsData * GenerateAsimovDataSinglePdf(const RooAbsPdf & pdf, const RooArgSet & obs, const RooRealVar & weightVar,
116  RooCategory * channelCat = 0);
117 
118  static RooAbsData * GenerateCountingAsimovData(RooAbsPdf & pdf, const RooArgSet & obs, const RooRealVar & weightVar,
119  RooCategory * channelCat = 0);
120 
121 
122  static void FillBins(const RooAbsPdf & pdf, const RooArgList &obs, RooAbsData & data, int &index, double
123  &binVolume, int &ibin);
124 
125  static double EvaluateNLL(RooAbsPdf & pdf, RooAbsData& data, const RooArgSet * condObs, const RooArgSet * globObs, const RooArgSet *poiSet = 0 );
126 
127  static bool SetObsToExpected(RooAbsPdf &pdf, const RooArgSet &obs);
128  static bool SetObsToExpected(RooProdPdf &prod, const RooArgSet &obs);
129 
130  protected:
132 
133  private:
134 
135  bool fOneSided; // for one sided PL test statistic (upper limits)
136  mutable bool fOneSidedDiscovery; // for one sided PL test statistic (for discovery)
137  bool fNominalAsimov; // make Asimov at nominal parameter values
138  mutable bool fIsInitialized; //! flag to check if calculator is initialized
139  mutable int fUseQTilde; // flag to indicate if using qtilde or not (-1 (default based on RooRealVar)), 0 false, 1 (true)
140  static int fgPrintLevel; // control print level (0 minimal, 1 normal, 2 debug)
141  mutable double fNLLObs;
142  mutable double fNLLAsimov;
143 
144  mutable RooAbsData * fAsimovData; // asimov data set
145  mutable RooArgSet fAsimovGlobObs; // snapshot of Asimov global observables
146  mutable RooArgSet fBestFitPoi; // snapshot of best fitted POI values
147  mutable RooArgSet fBestFitParams; // snapshot of all best fitted Parameter values
148 
149 
150  };
151 }
152 
153 #endif
static RooAbsData * GenerateAsimovData(const RooAbsPdf &pdf, const RooArgSet &observables)
generate the asimov data for the observables (not the global ones) need to deal with the case of a si...
Poisson pdf.
Definition: RooPoisson.h:19
virtual void SetNullModel(const ModelConfig &nullModel)
re-implementation of setters since they are needed to re-initialize the calculator ...
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:30
static double EvaluateNLL(RooAbsPdf &pdf, RooAbsData &data, const RooArgSet *condObs, const RooArgSet *globObs, const RooArgSet *poiSet=0)
static double GetExpectedPValues(double pnull, double palt, double nsigma, bool usecls, bool oneSided=true)
function given the null and the alt p value - return the expected one given the N - sigma value ...
RooProdPdf is an efficient implementation of a product of PDFs of the form.
Definition: RooProdPdf.h:31
HypoTestResult is a base class for results from hypothesis tests.
const RooArgSet & GetBestFitParams() const
return best fit value for all parameters
AsymptoticCalculator(RooAbsData &data, const ModelConfig &altModel, const ModelConfig &nullModel, bool nominalAsimov=false)
constructor for asymptotic calculator from Data set and ModelConfig
static void SetPrintLevel(int level)
set print level (static function)
static RooAbsData * GenerateCountingAsimovData(RooAbsPdf &pdf, const RooArgSet &obs, const RooRealVar &weightVar, RooCategory *channelCat=0)
generate counting Asimov data for the case when the pdf cannot be extended assume pdf is a RooPoisson...
Common base class for the Hypothesis Test Calculators.
static RooAbsData * MakeAsimovData(RooAbsData &data, const ModelConfig &model, const RooArgSet &poiValues, RooArgSet &globObs, const RooArgSet *genPoiValues=0)
make the asimov data from the ModelConfig and list of poi - return data set and snapshot of global ob...
static void FillBins(const RooAbsPdf &pdf, const RooArgList &obs, RooAbsData &data, int &index, double &binVolume, int &ibin)
fill bins by looping recursively on observables
#define ClassDef(name, id)
Definition: Rtypes.h:320
const RooArgSet & GetBestFitPoi() const
return snapshot of the best fit parameter
void SetQTilde(bool on)
set using of qtilde, by default is controlled if RoORealVar is limited or not
virtual HypoTestResult * GetHypoTest() const
re-implement HypoTest computation using the asymptotic
virtual void SetData(RooAbsData &data)
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
bool Initialize() const
initialize the calculator by performing a global fit and make the Asimov data set ...
void SetTwoSided()
set the test statistics for two sided (in case of upper limits for discovery does not make really sen...
RooAbsArg * first() const
static bool SetObsToExpected(RooAbsPdf &pdf, const RooArgSet &obs)
set observed value to the expected one works for Gaussian, Poisson or LogNormal assumes mean paramete...
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
RooCategory represents a fundamental (non-derived) discrete value object.
Definition: RooCategory.h:24
virtual void SetAlternateModel(const ModelConfig &altModel)
Namespace for the RooStats classes.
Definition: Asimov.h:20
int fUseQTilde
flag to check if calculator is initialized
void SetOneSidedDiscovery(bool on)
set the test statistics for one-sided discovery
void SetOneSided(bool on)
set test statistic for one sided (upper limits)
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
virtual void SetData(RooAbsData &data)
virtual void SetAlternateModel(const ModelConfig &altModel)
virtual void SetNullModel(const ModelConfig &nullModel)
Hypothesis Test Calculator based on the asymptotic formulae for the profile likelihood ratio...
const RooRealVar * GetMuHat() const
return best fit parameter (firs of poi)
static RooAbsData * GenerateAsimovDataSinglePdf(const RooAbsPdf &pdf, const RooArgSet &obs, const RooRealVar &weightVar, RooCategory *channelCat=0)
compute the asimov data set for an observable of a pdf use the number of bins sets in the observables...