ROOT logo
// @(#)root/roostats:$Id: SimpleLikelihoodRatioTestStat.h 39956 2011-06-24 21:06:20Z wouter $
// Author: Kyle Cranmer and Sven Kreiss    June 2010
/*************************************************************************
 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers.               *
 * All rights reserved.                                                  *
 *                                                                       *
 * For the licensing terms see $ROOTSYS/LICENSE.                         *
 * For the list of contributors see $ROOTSYS/README/CREDITS.             *
 *************************************************************************/

#ifndef ROOSTATS_SimpleLikelihoodRatioTestStat
#define ROOSTATS_SimpleLikelihoodRatioTestStat

//_________________________________________________
/*
 BEGIN_HTML
 <p>
 SimpleLikelihoodRatioTestStat: TestStatistic that returns -log(L[null] / L[alt]) where
 L is the likelihood.
 </p>
 END_HTML
 */
//

#ifndef ROOT_Rtypes
#include "Rtypes.h"
#endif

#ifndef ROO_NLL_VAR
#include "RooNLLVar.h"
#endif

#include "RooStats/TestStatistic.h"
#include "RooWorkspace.h"

namespace RooStats {

class SimpleLikelihoodRatioTestStat : public TestStatistic {

   public:

      //__________________________________
      SimpleLikelihoodRatioTestStat() :
         fNullPdf(NULL), fAltPdf(NULL)
      {
         // Constructor for proof. Do not use.
         fFirstEval = true;
         fNullParameters = NULL;
         fAltParameters = NULL;
	 fReuseNll=kFALSE ;
	 fNllNull=NULL ;
	 fNllAlt=NULL ;
      }

      //__________________________________
      SimpleLikelihoodRatioTestStat(
         RooAbsPdf& nullPdf,
         RooAbsPdf& altPdf
      ) :
         fFirstEval(true)
      {
         // Takes null and alternate parameters from PDF. Can be overridden.

         RooFit::MsgLevel msglevel = RooMsgService::instance().globalKillBelow();
         RooMsgService::instance().setGlobalKillBelow(RooFit::ERROR);
         w.import(nullPdf, RooFit::RecycleConflictNodes());
         w.import(altPdf, RooFit::RecycleConflictNodes());
         RooMsgService::instance().setGlobalKillBelow(msglevel);

         fNullPdf = w.pdf(nullPdf.GetName());
         fAltPdf = w.pdf(altPdf.GetName());

         fNullParameters = (RooArgSet*) fNullPdf->getVariables()->snapshot();
         fAltParameters = (RooArgSet*) fAltPdf->getVariables()->snapshot();
	 fReuseNll=kFALSE ;
	 fNllNull=NULL ;
	 fNllAlt=NULL ;
      }
      //__________________________________
      SimpleLikelihoodRatioTestStat(
         RooAbsPdf& nullPdf,
         RooAbsPdf& altPdf,
         const RooArgSet& nullParameters,
         const RooArgSet& altParameters
      ) :
         fFirstEval(true)
      {
         // Takes null and alternate parameters from values in nullParameters
         // and altParameters. Can be overridden.

         RooFit::MsgLevel msglevel = RooMsgService::instance().globalKillBelow();
         RooMsgService::instance().setGlobalKillBelow(RooFit::ERROR);
         w.import(nullPdf, RooFit::RecycleConflictNodes());
         w.import(altPdf, RooFit::RecycleConflictNodes());
         RooMsgService::instance().setGlobalKillBelow(msglevel);

         fNullPdf = w.pdf(nullPdf.GetName());
         fAltPdf = w.pdf(altPdf.GetName());

         fNullParameters = (RooArgSet*) nullParameters.snapshot();
         fAltParameters = (RooArgSet*) altParameters.snapshot();
	 fReuseNll=kFALSE ;
	 fNllNull=NULL ;
	 fNllAlt=NULL ;
      }

      //______________________________
      virtual ~SimpleLikelihoodRatioTestStat() {
         if (fNullParameters) delete fNullParameters;
         if (fAltParameters) delete fAltParameters;
	 if (fNllNull) delete fNllNull ;
	 if (fNllAlt) delete fNllAlt ;
      }

     static void setAlwaysReuseNLL(Bool_t flag) { fAlwaysReuseNll = flag ; }
     void setReuseNLL(Bool_t flag) { fReuseNll = flag ; }

      //_________________________________________
      void SetNullParameters(const RooArgSet& nullParameters) {
         if (fNullParameters) delete fNullParameters;
         fFirstEval = true;
         //      if(fNullParameters) delete fNullParameters;
         fNullParameters = (RooArgSet*) nullParameters.snapshot();
      }

      //_________________________________________
      void SetAltParameters(const RooArgSet& altParameters) {
         if (fAltParameters) delete fAltParameters;
         fFirstEval = true;
         //      if(fAltParameters) delete fAltParameters;
         fAltParameters = (RooArgSet*) altParameters.snapshot();
      }

      //______________________________
      bool ParamsAreEqual() {
         // this should be possible with RooAbsCollection
         if (!fNullParameters->equals(*fAltParameters)) return false;

         RooAbsReal* null;
         RooAbsReal* alt;

         TIterator* nullIt = fNullParameters->createIterator();
         TIterator* altIt = fAltParameters->createIterator();
         bool ret = true;
         while ((null = (RooAbsReal*) nullIt->Next()) && (alt = (RooAbsReal*) altIt->Next())) {
            if (null->getVal() != alt->getVal()) ret = false;
         }
         delete nullIt;
         delete altIt;
         return ret;
      }

      //______________________________
      virtual Double_t Evaluate(RooAbsData& data, RooArgSet& nullPOI) {

         if (fFirstEval && ParamsAreEqual()) {
            oocoutW(fNullParameters,InputArguments)
               << "Same RooArgSet used for null and alternate, so you must explicitly SetNullParameters and SetAlternateParameters or the likelihood ratio will always be 1."
               << endl;
         }
         fFirstEval = false;

         RooFit::MsgLevel msglevel = RooMsgService::instance().globalKillBelow();
         RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL);

	 Bool_t reuse = (fReuseNll || fAlwaysReuseNll) ;

	 Bool_t created = kFALSE ;
	 if (!fNllNull) {
	   fNllNull = (RooNLLVar*) fNullPdf->createNLL(data, RooFit::CloneData(kFALSE));
	   created = kTRUE ;
	 }
	 if (reuse && !created) {
	   fNllNull->setData(data, kFALSE) ;
	 }

         // make sure we set the variables attached to this nll
         RooArgSet* attachedSet = fNllNull->getVariables();
         *attachedSet = *fNullParameters;
         *attachedSet = nullPOI;
         double nullNLL = fNllNull->getVal();

         delete fNllNull ; fNllNull = NULL ;
         delete attachedSet;

	 created = kFALSE ;
	 if (!fNllAlt) {
	   fNllAlt = (RooNLLVar*) fAltPdf->createNLL(data, RooFit::CloneData(kFALSE));
	   created = kTRUE ;
	 }
	 if (reuse && !created) {
	   fNllAlt->setData(data, kFALSE) ;
	 }
         // make sure we set the variables attached to this nll
         attachedSet = fNllAlt->getVariables();
         *attachedSet = *fAltParameters;
         double altNLL = fNllAlt->getVal();

         delete fNllAlt ; fNllAlt = NULL ;
         delete attachedSet;

         RooMsgService::instance().setGlobalKillBelow(msglevel);
         return nullNLL - altNLL;
      }

      virtual const TString GetVarName() const {
         return "log(L(#mu_{1}) / L(#mu_{0}))";
      }

   private:

      RooWorkspace w;

      RooAbsPdf* fNullPdf;
      RooAbsPdf* fAltPdf;
      RooArgSet* fNullParameters;
      RooArgSet* fAltParameters;
      bool fFirstEval;

      RooNLLVar* fNllNull ;
      RooNLLVar* fNllAlt ;
      static Bool_t fAlwaysReuseNll ;
      Bool_t fReuseNll ;


   protected:
   ClassDef(SimpleLikelihoodRatioTestStat,1)
};

}

#endif
 SimpleLikelihoodRatioTestStat.h:1
 SimpleLikelihoodRatioTestStat.h:2
 SimpleLikelihoodRatioTestStat.h:3
 SimpleLikelihoodRatioTestStat.h:4
 SimpleLikelihoodRatioTestStat.h:5
 SimpleLikelihoodRatioTestStat.h:6
 SimpleLikelihoodRatioTestStat.h:7
 SimpleLikelihoodRatioTestStat.h:8
 SimpleLikelihoodRatioTestStat.h:9
 SimpleLikelihoodRatioTestStat.h:10
 SimpleLikelihoodRatioTestStat.h:11
 SimpleLikelihoodRatioTestStat.h:12
 SimpleLikelihoodRatioTestStat.h:13
 SimpleLikelihoodRatioTestStat.h:14
 SimpleLikelihoodRatioTestStat.h:15
 SimpleLikelihoodRatioTestStat.h:16
 SimpleLikelihoodRatioTestStat.h:17
 SimpleLikelihoodRatioTestStat.h:18
 SimpleLikelihoodRatioTestStat.h:19
 SimpleLikelihoodRatioTestStat.h:20
 SimpleLikelihoodRatioTestStat.h:21
 SimpleLikelihoodRatioTestStat.h:22
 SimpleLikelihoodRatioTestStat.h:23
 SimpleLikelihoodRatioTestStat.h:24
 SimpleLikelihoodRatioTestStat.h:25
 SimpleLikelihoodRatioTestStat.h:26
 SimpleLikelihoodRatioTestStat.h:27
 SimpleLikelihoodRatioTestStat.h:28
 SimpleLikelihoodRatioTestStat.h:29
 SimpleLikelihoodRatioTestStat.h:30
 SimpleLikelihoodRatioTestStat.h:31
 SimpleLikelihoodRatioTestStat.h:32
 SimpleLikelihoodRatioTestStat.h:33
 SimpleLikelihoodRatioTestStat.h:34
 SimpleLikelihoodRatioTestStat.h:35
 SimpleLikelihoodRatioTestStat.h:36
 SimpleLikelihoodRatioTestStat.h:37
 SimpleLikelihoodRatioTestStat.h:38
 SimpleLikelihoodRatioTestStat.h:39
 SimpleLikelihoodRatioTestStat.h:40
 SimpleLikelihoodRatioTestStat.h:41
 SimpleLikelihoodRatioTestStat.h:42
 SimpleLikelihoodRatioTestStat.h:43
 SimpleLikelihoodRatioTestStat.h:44
 SimpleLikelihoodRatioTestStat.h:45
 SimpleLikelihoodRatioTestStat.h:46
 SimpleLikelihoodRatioTestStat.h:47
 SimpleLikelihoodRatioTestStat.h:48
 SimpleLikelihoodRatioTestStat.h:49
 SimpleLikelihoodRatioTestStat.h:50
 SimpleLikelihoodRatioTestStat.h:51
 SimpleLikelihoodRatioTestStat.h:52
 SimpleLikelihoodRatioTestStat.h:53
 SimpleLikelihoodRatioTestStat.h:54
 SimpleLikelihoodRatioTestStat.h:55
 SimpleLikelihoodRatioTestStat.h:56
 SimpleLikelihoodRatioTestStat.h:57
 SimpleLikelihoodRatioTestStat.h:58
 SimpleLikelihoodRatioTestStat.h:59
 SimpleLikelihoodRatioTestStat.h:60
 SimpleLikelihoodRatioTestStat.h:61
 SimpleLikelihoodRatioTestStat.h:62
 SimpleLikelihoodRatioTestStat.h:63
 SimpleLikelihoodRatioTestStat.h:64
 SimpleLikelihoodRatioTestStat.h:65
 SimpleLikelihoodRatioTestStat.h:66
 SimpleLikelihoodRatioTestStat.h:67
 SimpleLikelihoodRatioTestStat.h:68
 SimpleLikelihoodRatioTestStat.h:69
 SimpleLikelihoodRatioTestStat.h:70
 SimpleLikelihoodRatioTestStat.h:71
 SimpleLikelihoodRatioTestStat.h:72
 SimpleLikelihoodRatioTestStat.h:73
 SimpleLikelihoodRatioTestStat.h:74
 SimpleLikelihoodRatioTestStat.h:75
 SimpleLikelihoodRatioTestStat.h:76
 SimpleLikelihoodRatioTestStat.h:77
 SimpleLikelihoodRatioTestStat.h:78
 SimpleLikelihoodRatioTestStat.h:79
 SimpleLikelihoodRatioTestStat.h:80
 SimpleLikelihoodRatioTestStat.h:81
 SimpleLikelihoodRatioTestStat.h:82
 SimpleLikelihoodRatioTestStat.h:83
 SimpleLikelihoodRatioTestStat.h:84
 SimpleLikelihoodRatioTestStat.h:85
 SimpleLikelihoodRatioTestStat.h:86
 SimpleLikelihoodRatioTestStat.h:87
 SimpleLikelihoodRatioTestStat.h:88
 SimpleLikelihoodRatioTestStat.h:89
 SimpleLikelihoodRatioTestStat.h:90
 SimpleLikelihoodRatioTestStat.h:91
 SimpleLikelihoodRatioTestStat.h:92
 SimpleLikelihoodRatioTestStat.h:93
 SimpleLikelihoodRatioTestStat.h:94
 SimpleLikelihoodRatioTestStat.h:95
 SimpleLikelihoodRatioTestStat.h:96
 SimpleLikelihoodRatioTestStat.h:97
 SimpleLikelihoodRatioTestStat.h:98
 SimpleLikelihoodRatioTestStat.h:99
 SimpleLikelihoodRatioTestStat.h:100
 SimpleLikelihoodRatioTestStat.h:101
 SimpleLikelihoodRatioTestStat.h:102
 SimpleLikelihoodRatioTestStat.h:103
 SimpleLikelihoodRatioTestStat.h:104
 SimpleLikelihoodRatioTestStat.h:105
 SimpleLikelihoodRatioTestStat.h:106
 SimpleLikelihoodRatioTestStat.h:107
 SimpleLikelihoodRatioTestStat.h:108
 SimpleLikelihoodRatioTestStat.h:109
 SimpleLikelihoodRatioTestStat.h:110
 SimpleLikelihoodRatioTestStat.h:111
 SimpleLikelihoodRatioTestStat.h:112
 SimpleLikelihoodRatioTestStat.h:113
 SimpleLikelihoodRatioTestStat.h:114
 SimpleLikelihoodRatioTestStat.h:115
 SimpleLikelihoodRatioTestStat.h:116
 SimpleLikelihoodRatioTestStat.h:117
 SimpleLikelihoodRatioTestStat.h:118
 SimpleLikelihoodRatioTestStat.h:119
 SimpleLikelihoodRatioTestStat.h:120
 SimpleLikelihoodRatioTestStat.h:121
 SimpleLikelihoodRatioTestStat.h:122
 SimpleLikelihoodRatioTestStat.h:123
 SimpleLikelihoodRatioTestStat.h:124
 SimpleLikelihoodRatioTestStat.h:125
 SimpleLikelihoodRatioTestStat.h:126
 SimpleLikelihoodRatioTestStat.h:127
 SimpleLikelihoodRatioTestStat.h:128
 SimpleLikelihoodRatioTestStat.h:129
 SimpleLikelihoodRatioTestStat.h:130
 SimpleLikelihoodRatioTestStat.h:131
 SimpleLikelihoodRatioTestStat.h:132
 SimpleLikelihoodRatioTestStat.h:133
 SimpleLikelihoodRatioTestStat.h:134
 SimpleLikelihoodRatioTestStat.h:135
 SimpleLikelihoodRatioTestStat.h:136
 SimpleLikelihoodRatioTestStat.h:137
 SimpleLikelihoodRatioTestStat.h:138
 SimpleLikelihoodRatioTestStat.h:139
 SimpleLikelihoodRatioTestStat.h:140
 SimpleLikelihoodRatioTestStat.h:141
 SimpleLikelihoodRatioTestStat.h:142
 SimpleLikelihoodRatioTestStat.h:143
 SimpleLikelihoodRatioTestStat.h:144
 SimpleLikelihoodRatioTestStat.h:145
 SimpleLikelihoodRatioTestStat.h:146
 SimpleLikelihoodRatioTestStat.h:147
 SimpleLikelihoodRatioTestStat.h:148
 SimpleLikelihoodRatioTestStat.h:149
 SimpleLikelihoodRatioTestStat.h:150
 SimpleLikelihoodRatioTestStat.h:151
 SimpleLikelihoodRatioTestStat.h:152
 SimpleLikelihoodRatioTestStat.h:153
 SimpleLikelihoodRatioTestStat.h:154
 SimpleLikelihoodRatioTestStat.h:155
 SimpleLikelihoodRatioTestStat.h:156
 SimpleLikelihoodRatioTestStat.h:157
 SimpleLikelihoodRatioTestStat.h:158
 SimpleLikelihoodRatioTestStat.h:159
 SimpleLikelihoodRatioTestStat.h:160
 SimpleLikelihoodRatioTestStat.h:161
 SimpleLikelihoodRatioTestStat.h:162
 SimpleLikelihoodRatioTestStat.h:163
 SimpleLikelihoodRatioTestStat.h:164
 SimpleLikelihoodRatioTestStat.h:165
 SimpleLikelihoodRatioTestStat.h:166
 SimpleLikelihoodRatioTestStat.h:167
 SimpleLikelihoodRatioTestStat.h:168
 SimpleLikelihoodRatioTestStat.h:169
 SimpleLikelihoodRatioTestStat.h:170
 SimpleLikelihoodRatioTestStat.h:171
 SimpleLikelihoodRatioTestStat.h:172
 SimpleLikelihoodRatioTestStat.h:173
 SimpleLikelihoodRatioTestStat.h:174
 SimpleLikelihoodRatioTestStat.h:175
 SimpleLikelihoodRatioTestStat.h:176
 SimpleLikelihoodRatioTestStat.h:177
 SimpleLikelihoodRatioTestStat.h:178
 SimpleLikelihoodRatioTestStat.h:179
 SimpleLikelihoodRatioTestStat.h:180
 SimpleLikelihoodRatioTestStat.h:181
 SimpleLikelihoodRatioTestStat.h:182
 SimpleLikelihoodRatioTestStat.h:183
 SimpleLikelihoodRatioTestStat.h:184
 SimpleLikelihoodRatioTestStat.h:185
 SimpleLikelihoodRatioTestStat.h:186
 SimpleLikelihoodRatioTestStat.h:187
 SimpleLikelihoodRatioTestStat.h:188
 SimpleLikelihoodRatioTestStat.h:189
 SimpleLikelihoodRatioTestStat.h:190
 SimpleLikelihoodRatioTestStat.h:191
 SimpleLikelihoodRatioTestStat.h:192
 SimpleLikelihoodRatioTestStat.h:193
 SimpleLikelihoodRatioTestStat.h:194
 SimpleLikelihoodRatioTestStat.h:195
 SimpleLikelihoodRatioTestStat.h:196
 SimpleLikelihoodRatioTestStat.h:197
 SimpleLikelihoodRatioTestStat.h:198
 SimpleLikelihoodRatioTestStat.h:199
 SimpleLikelihoodRatioTestStat.h:200
 SimpleLikelihoodRatioTestStat.h:201
 SimpleLikelihoodRatioTestStat.h:202
 SimpleLikelihoodRatioTestStat.h:203
 SimpleLikelihoodRatioTestStat.h:204
 SimpleLikelihoodRatioTestStat.h:205
 SimpleLikelihoodRatioTestStat.h:206
 SimpleLikelihoodRatioTestStat.h:207
 SimpleLikelihoodRatioTestStat.h:208
 SimpleLikelihoodRatioTestStat.h:209
 SimpleLikelihoodRatioTestStat.h:210
 SimpleLikelihoodRatioTestStat.h:211
 SimpleLikelihoodRatioTestStat.h:212
 SimpleLikelihoodRatioTestStat.h:213
 SimpleLikelihoodRatioTestStat.h:214
 SimpleLikelihoodRatioTestStat.h:215
 SimpleLikelihoodRatioTestStat.h:216
 SimpleLikelihoodRatioTestStat.h:217
 SimpleLikelihoodRatioTestStat.h:218
 SimpleLikelihoodRatioTestStat.h:219
 SimpleLikelihoodRatioTestStat.h:220
 SimpleLikelihoodRatioTestStat.h:221
 SimpleLikelihoodRatioTestStat.h:222
 SimpleLikelihoodRatioTestStat.h:223
 SimpleLikelihoodRatioTestStat.h:224
 SimpleLikelihoodRatioTestStat.h:225
 SimpleLikelihoodRatioTestStat.h:226
 SimpleLikelihoodRatioTestStat.h:227
 SimpleLikelihoodRatioTestStat.h:228
 SimpleLikelihoodRatioTestStat.h:229
 SimpleLikelihoodRatioTestStat.h:230
 SimpleLikelihoodRatioTestStat.h:231
 SimpleLikelihoodRatioTestStat.h:232