Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
SimpleLikelihoodRatioTestStat.h
Go to the documentation of this file.
1 // @(#)root/roostats:$Id$
2 // Author: Kyle Cranmer and Sven Kreiss June 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_SimpleLikelihoodRatioTestStat
12 #define ROOSTATS_SimpleLikelihoodRatioTestStat
13 
14 #include "Rtypes.h"
15 
16 #include "RooNLLVar.h"
17 
18 #include "RooRealVar.h"
19 
20 #include "RooStats/TestStatistic.h"
21 
22 namespace RooStats {
23 
24  class SimpleLikelihoodRatioTestStat : public TestStatistic {
25 
26  public:
27 
28  //__________________________________
29  SimpleLikelihoodRatioTestStat() :
30  fNullPdf(NULL), fAltPdf(NULL)
31  {
32  // Constructor for proof. Do not use.
33  fFirstEval = true;
34  fDetailedOutputEnabled = false;
35  fDetailedOutput = NULL;
36  fNullParameters = NULL;
37  fAltParameters = NULL;
38  fReuseNll=kFALSE ;
39  fNllNull=NULL ;
40  fNllAlt=NULL ;
41  }
42 
43  //__________________________________
44  SimpleLikelihoodRatioTestStat(
45  RooAbsPdf& nullPdf,
46  RooAbsPdf& altPdf
47  ) :
48  fFirstEval(true)
49  {
50  // Takes null and alternate parameters from PDF. Can be overridden.
51 
52  fNullPdf = &nullPdf;
53  fAltPdf = &altPdf;
54 
55  RooArgSet * allNullVars = fNullPdf->getVariables();
56  fNullParameters = (RooArgSet*) allNullVars->snapshot();
57  delete allNullVars;
58 
59  RooArgSet * allAltVars = fAltPdf->getVariables();
60  fAltParameters = (RooArgSet*) allAltVars->snapshot();
61  delete allAltVars;
62 
63  fDetailedOutputEnabled = false;
64  fDetailedOutput = NULL;
65 
66  fReuseNll=kFALSE ;
67  fNllNull=NULL ;
68  fNllAlt=NULL ;
69  }
70  //__________________________________
71  SimpleLikelihoodRatioTestStat(
72  RooAbsPdf& nullPdf,
73  RooAbsPdf& altPdf,
74  const RooArgSet& nullParameters,
75  const RooArgSet& altParameters
76  ) :
77  fFirstEval(true)
78  {
79  // Takes null and alternate parameters from values in nullParameters
80  // and altParameters. Can be overridden.
81  fNullPdf = &nullPdf;
82  fAltPdf = &altPdf;
83 
84  fNullParameters = (RooArgSet*) nullParameters.snapshot();
85  fAltParameters = (RooArgSet*) altParameters.snapshot();
86 
87  fDetailedOutputEnabled = false;
88  fDetailedOutput = NULL;
89 
90  fReuseNll=kFALSE ;
91  fNllNull=NULL ;
92  fNllAlt=NULL ;
93  }
94 
95  //______________________________
96  virtual ~SimpleLikelihoodRatioTestStat() {
97  if (fNullParameters) delete fNullParameters;
98  if (fAltParameters) delete fAltParameters;
99  if (fNllNull) delete fNllNull ;
100  if (fNllAlt) delete fNllAlt ;
101  if (fDetailedOutput) delete fDetailedOutput;
102  }
103 
104  static void SetAlwaysReuseNLL(Bool_t flag);
105 
106  void SetReuseNLL(Bool_t flag) { fReuseNll = flag ; }
107 
108  //_________________________________________
109  void SetNullParameters(const RooArgSet& nullParameters) {
110  if (fNullParameters) delete fNullParameters;
111  fFirstEval = true;
112  // if(fNullParameters) delete fNullParameters;
113  fNullParameters = (RooArgSet*) nullParameters.snapshot();
114  }
115 
116  //_________________________________________
117  void SetAltParameters(const RooArgSet& altParameters) {
118  if (fAltParameters) delete fAltParameters;
119  fFirstEval = true;
120  // if(fAltParameters) delete fAltParameters;
121  fAltParameters = (RooArgSet*) altParameters.snapshot();
122  }
123 
124  //______________________________
125  bool ParamsAreEqual() {
126  // this should be possible with RooAbsCollection
127  if (!fNullParameters->equals(*fAltParameters)) return false;
128 
129  RooAbsReal* null;
130  RooAbsReal* alt;
131 
132  TIterator* nullIt = fNullParameters->createIterator();
133  TIterator* altIt = fAltParameters->createIterator();
134  bool ret = true;
135  while ((null = (RooAbsReal*) nullIt->Next()) && (alt = (RooAbsReal*) altIt->Next())) {
136  if (null->getVal() != alt->getVal()) ret = false;
137  }
138  delete nullIt;
139  delete altIt;
140  return ret;
141  }
142 
143 
144  // set the conditional observables which will be used when creating the NLL
145  // so the pdf's will not be normalized on the conditional observables when computing the NLL
146  virtual void SetConditionalObservables(const RooArgSet& set) {fConditionalObs.removeAll(); fConditionalObs.add(set);}
147 
148  // set the global observables which will be used when creating the NLL
149  // so the constraint pdf's will be normalized correctly on the global observables when computing the NLL
150  virtual void SetGlobalObservables(const RooArgSet& set) {fGlobalObs.removeAll(); fGlobalObs.add(set);}
151 
152  //______________________________
153  virtual Double_t Evaluate(RooAbsData& data, RooArgSet& nullPOI);
154 
155  virtual void EnableDetailedOutput( bool e=true ) { fDetailedOutputEnabled = e; fDetailedOutput = NULL; }
156  virtual const RooArgSet* GetDetailedOutput(void) const { return fDetailedOutput; }
157 
158  virtual const TString GetVarName() const {
159  return "log(L(#mu_{1}) / L(#mu_{0}))";
160  }
161 
162  private:
163 
164  RooAbsPdf* fNullPdf;
165  RooAbsPdf* fAltPdf;
166  RooArgSet* fNullParameters;
167  RooArgSet* fAltParameters;
168  RooArgSet fConditionalObs;
169  RooArgSet fGlobalObs;
170  bool fFirstEval;
171 
172  bool fDetailedOutputEnabled;
173  RooArgSet* fDetailedOutput; //!
174 
175  RooAbsReal* fNllNull ; //! transient copy of the null NLL
176  RooAbsReal* fNllAlt ; //! transient copy of the alt NLL
177  static Bool_t fgAlwaysReuseNll ;
178  Bool_t fReuseNll ;
179 
180 
181  protected:
182  ClassDef(SimpleLikelihoodRatioTestStat,4)
183 };
184 
185 }
186 
187 #endif