11 #ifndef ROOSTATS_SimpleLikelihoodRatioTestStat 
   12 #define ROOSTATS_SimpleLikelihoodRatioTestStat 
   24    class SimpleLikelihoodRatioTestStat : 
public TestStatistic {
 
   29       SimpleLikelihoodRatioTestStat() :
 
   30          fNullPdf(NULL), fAltPdf(NULL)
 
   34         fDetailedOutputEnabled = 
false;
 
   35         fDetailedOutput = NULL;
 
   36          fNullParameters = NULL;
 
   37          fAltParameters = NULL;
 
   44       SimpleLikelihoodRatioTestStat(
 
   55          RooArgSet * allNullVars = fNullPdf->getVariables();
 
   56          fNullParameters = (RooArgSet*) allNullVars->snapshot();
 
   59          RooArgSet * allAltVars = fAltPdf->getVariables();
 
   60          fAltParameters = (RooArgSet*) allAltVars->snapshot();
 
   63          fDetailedOutputEnabled = 
false;
 
   64          fDetailedOutput = NULL;
 
   71       SimpleLikelihoodRatioTestStat(
 
   74          const RooArgSet& nullParameters,
 
   75          const RooArgSet& altParameters
 
   84          fNullParameters = (RooArgSet*) nullParameters.snapshot();
 
   85          fAltParameters = (RooArgSet*) altParameters.snapshot();
 
   87          fDetailedOutputEnabled = 
false;
 
   88          fDetailedOutput = NULL;
 
   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;
 
  104       static void SetAlwaysReuseNLL(Bool_t flag);
 
  106       void SetReuseNLL(Bool_t flag) { fReuseNll = flag ; }
 
  109       void SetNullParameters(
const RooArgSet& nullParameters) {
 
  110          if (fNullParameters) 
delete fNullParameters;
 
  113          fNullParameters = (RooArgSet*) nullParameters.snapshot();
 
  117       void SetAltParameters(
const RooArgSet& altParameters) {
 
  118          if (fAltParameters) 
delete fAltParameters;
 
  121          fAltParameters = (RooArgSet*) altParameters.snapshot();
 
  125       bool ParamsAreEqual() {
 
  127          if (!fNullParameters->equals(*fAltParameters)) 
return false;
 
  132          TIterator* nullIt = fNullParameters->createIterator();
 
  133          TIterator* altIt = fAltParameters->createIterator();
 
  135          while ((null = (RooAbsReal*) nullIt->Next()) && (alt = (RooAbsReal*) altIt->Next())) {
 
  136             if (null->getVal() != alt->getVal()) ret = 
false;
 
  146       virtual void SetConditionalObservables(
const RooArgSet& set) {fConditionalObs.removeAll(); fConditionalObs.add(set);}
 
  150       virtual void SetGlobalObservables(
const RooArgSet& set) {fGlobalObs.removeAll(); fGlobalObs.add(set);}
 
  153       virtual Double_t Evaluate(RooAbsData& data, RooArgSet& nullPOI);
 
  155       virtual void EnableDetailedOutput( 
bool e=
true ) { fDetailedOutputEnabled = e; fDetailedOutput = NULL; }
 
  156       virtual const RooArgSet* GetDetailedOutput(
void)
 const { 
return fDetailedOutput; }
 
  158       virtual const TString GetVarName()
 const {
 
  159          return "log(L(#mu_{1}) / L(#mu_{0}))";
 
  166       RooArgSet* fNullParameters;
 
  167       RooArgSet* fAltParameters;
 
  168       RooArgSet fConditionalObs;
 
  169       RooArgSet fGlobalObs;
 
  172       bool fDetailedOutputEnabled;
 
  173       RooArgSet* fDetailedOutput; 
 
  175       RooAbsReal* fNllNull ;  
 
  176       RooAbsReal* fNllAlt ; 
 
  177       static Bool_t fgAlwaysReuseNll ;
 
  182    ClassDef(SimpleLikelihoodRatioTestStat,4)