File indexing completed on 2025-12-16 10:29:44
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011 #ifndef ROOSTATS_AsymptoticCalculator
0012 #define ROOSTATS_AsymptoticCalculator
0013
0014 #include "RooStats/HypoTestCalculatorGeneric.h"
0015 #include "RooArgSet.h"
0016 #include "Rtypes.h"
0017
0018 class RooArgList;
0019 class RooCategory;
0020 class RooRealVar;
0021 class RooPoisson;
0022 class RooProdPdf;
0023
0024
0025 namespace RooStats {
0026
0027 class AsymptoticCalculator : public HypoTestCalculatorGeneric {
0028
0029 public:
0030 AsymptoticCalculator(
0031 RooAbsData &data,
0032 const ModelConfig &altModel,
0033 const ModelConfig &nullModel,
0034 bool nominalAsimov = false
0035 );
0036
0037
0038 bool Initialize() const;
0039
0040
0041 HypoTestResult *GetHypoTest() const override;
0042
0043
0044 static RooAbsData * MakeAsimovData( RooAbsData & data, const ModelConfig & model, const RooArgSet & poiValues, RooArgSet & globObs, const RooArgSet * genPoiValues = nullptr);
0045
0046
0047
0048 static RooAbsData * MakeAsimovData( const ModelConfig & model, const RooArgSet & allParamValues, RooArgSet & globObs);
0049
0050
0051
0052 static RooAbsData * GenerateAsimovData(const RooAbsPdf & pdf, const RooArgSet & observables );
0053
0054
0055 static double GetExpectedPValues(double pnull, double palt, double nsigma, bool usecls, bool oneSided = true );
0056
0057
0058 void SetOneSided(bool on) { fOneSided = on; }
0059
0060
0061
0062 void SetTwoSided() { fOneSided = false; fOneSidedDiscovery = false;}
0063
0064
0065 void SetOneSidedDiscovery(bool on) { fOneSidedDiscovery = on; }
0066
0067
0068 void SetNullModel(const ModelConfig &nullModel) override {
0069 HypoTestCalculatorGeneric::SetNullModel(nullModel);
0070 fIsInitialized = false;
0071 }
0072 void SetAlternateModel(const ModelConfig &altModel) override {
0073 HypoTestCalculatorGeneric::SetAlternateModel(altModel);
0074 fIsInitialized = false;
0075 }
0076 void SetData(RooAbsData &data) override {
0077 HypoTestCalculatorGeneric::SetData(data);
0078 fIsInitialized = false;
0079 }
0080
0081
0082 bool IsTwoSided() const { return (!fOneSided && !fOneSidedDiscovery); }
0083 bool IsOneSidedDiscovery() const { return fOneSidedDiscovery; }
0084
0085
0086
0087 void SetQTilde(bool on) { fUseQTilde = on; }
0088
0089
0090 const RooArgSet & GetBestFitPoi() const { return fBestFitPoi; }
0091
0092 const RooRealVar * GetMuHat() const { return dynamic_cast<RooRealVar*>(fBestFitPoi.first()); }
0093
0094 const RooArgSet & GetBestFitParams() const { return fBestFitPoi; }
0095
0096 static void SetPrintLevel(int level);
0097
0098 private:
0099 bool fOneSided;
0100 mutable bool fOneSidedDiscovery;
0101 bool fNominalAsimov;
0102 mutable bool fIsInitialized;
0103 mutable int fUseQTilde;
0104 mutable double fNLLObs;
0105 mutable double fNLLAsimov;
0106
0107 mutable RooAbsData * fAsimovData;
0108 mutable RooArgSet fAsimovGlobObs;
0109 mutable RooArgSet fBestFitPoi;
0110 mutable RooArgSet fBestFitParams;
0111
0112 ClassDefOverride(AsymptoticCalculator,0)
0113 };
0114 }
0115
0116 #endif