Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-12-16 10:29:46

0001 // @(#)root/roostats:$Id$
0002 // Author: Kyle Cranmer and Sven Kreiss    June 2010
0003 /*************************************************************************
0004  * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers.               *
0005  * All rights reserved.                                                  *
0006  *                                                                       *
0007  * For the licensing terms see $ROOTSYS/LICENSE.                         *
0008  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
0009  *************************************************************************/
0010 
0011 #ifndef ROOSTATS_SimpleLikelihoodRatioTestStat
0012 #define ROOSTATS_SimpleLikelihoodRatioTestStat
0013 
0014 #include "Rtypes.h"
0015 
0016 #include "RooAbsPdf.h"
0017 #include "RooRealVar.h"
0018 
0019 #include "RooStats/TestStatistic.h"
0020 
0021 namespace RooStats {
0022 
0023    class SimpleLikelihoodRatioTestStat : public TestStatistic {
0024 
0025    public:
0026 
0027       /// Constructor for proof, now removed. Do not use.
0028       SimpleLikelihoodRatioTestStat() = default;
0029 
0030       /// Takes null and alternate parameters from PDF. Can be overridden.
0031       SimpleLikelihoodRatioTestStat(RooAbsPdf &nullPdf, RooAbsPdf &altPdf)
0032          : fNullPdf(&nullPdf),
0033            fAltPdf(&altPdf)
0034       {
0035          std::unique_ptr<RooArgSet> allNullVars{fNullPdf->getVariables()};
0036          fNullParameters = allNullVars->snapshot();
0037 
0038          std::unique_ptr<RooArgSet> allAltVars{fAltPdf->getVariables()};
0039          fAltParameters = allAltVars->snapshot();
0040       }
0041 
0042       /// Takes null and alternate parameters from values in nullParameters
0043       /// and altParameters. Can be overridden.
0044       SimpleLikelihoodRatioTestStat(RooAbsPdf &nullPdf, RooAbsPdf &altPdf, const RooArgSet &nullParameters,
0045                                     const RooArgSet &altParameters)
0046          : fNullPdf(&nullPdf),
0047            fAltPdf(&altPdf),
0048            fNullParameters(nullParameters.snapshot()),
0049            fAltParameters(altParameters.snapshot())
0050       {
0051       }
0052 
0053       ~SimpleLikelihoodRatioTestStat() override {
0054          if (fNullParameters) delete fNullParameters;
0055          if (fAltParameters) delete fAltParameters;
0056       }
0057 
0058       static void SetAlwaysReuseNLL(bool flag);
0059 
0060       void SetReuseNLL(bool flag) { fReuseNll = flag ; }
0061 
0062       void SetNullParameters(const RooArgSet& nullParameters) {
0063          if (fNullParameters) delete fNullParameters;
0064          fFirstEval = true;
0065          fNullParameters = (RooArgSet*) nullParameters.snapshot();
0066       }
0067 
0068       void SetAltParameters(const RooArgSet& altParameters) {
0069          if (fAltParameters) delete fAltParameters;
0070          fFirstEval = true;
0071          fAltParameters = (RooArgSet*) altParameters.snapshot();
0072       }
0073 
0074       /// this should be possible with RooAbsCollection
0075       bool ParamsAreEqual() {
0076          if (!fNullParameters->equals(*fAltParameters)) return false;
0077 
0078          bool ret = true;
0079 
0080          for (auto nullIt = fNullParameters->begin(), altIt = fAltParameters->begin();
0081               nullIt != fNullParameters->end() && altIt != fAltParameters->end(); ++nullIt, ++altIt) {
0082             RooAbsReal *null = static_cast<RooAbsReal *>(*nullIt);
0083             RooAbsReal *alt = static_cast<RooAbsReal *>(*altIt);
0084             if (null->getVal() != alt->getVal())
0085                ret = false;
0086          }
0087 
0088          return ret;
0089       }
0090 
0091 
0092       /// set the conditional observables which will be used when creating the NLL
0093       /// so the pdf's will not be normalized on the conditional observables when computing the NLL
0094       void SetConditionalObservables(const RooArgSet& set) override {fConditionalObs.removeAll(); fConditionalObs.add(set);}
0095 
0096       /// set the global observables which will be used when creating the NLL
0097       /// so the constraint pdf's will be normalized correctly on the global observables when computing the NLL
0098       void SetGlobalObservables(const RooArgSet& set) override {fGlobalObs.removeAll(); fGlobalObs.add(set);}
0099 
0100       double Evaluate(RooAbsData& data, RooArgSet& nullPOI) override;
0101 
0102       virtual void EnableDetailedOutput( bool e=true ) { fDetailedOutputEnabled = e; fDetailedOutput = nullptr; }
0103       const RooArgSet* GetDetailedOutput(void) const override { return fDetailedOutput.get(); }
0104 
0105       const TString GetVarName() const override {
0106          return "log(L(#mu_{1}) / L(#mu_{0}))";
0107       }
0108 
0109    private:
0110 
0111       RooAbsPdf* fNullPdf = nullptr;
0112       RooAbsPdf* fAltPdf = nullptr;
0113       RooArgSet* fNullParameters = nullptr;
0114       RooArgSet* fAltParameters = nullptr;
0115       RooArgSet fConditionalObs;
0116       RooArgSet fGlobalObs;
0117       bool fFirstEval = true;
0118 
0119       bool fDetailedOutputEnabled = false;
0120       std::unique_ptr<RooArgSet> fDetailedOutput; ///<!
0121 
0122       std::unique_ptr<RooAbsReal> fNllNull; ///<! transient copy of the null NLL
0123       std::unique_ptr<RooAbsReal> fNllAlt;  ///<!  transient copy of the alt NLL
0124       static bool fgAlwaysReuseNll ;
0125       bool fReuseNll = false;
0126 
0127 
0128    ClassDefOverride(SimpleLikelihoodRatioTestStat,4)
0129 };
0130 
0131 }
0132 
0133 #endif