66struct HypoTestOptions {
68 bool noSystematics =
false;
69 double nToysRatio = 4;
72 bool generateBinned =
false;
73 bool useProof =
false;
74 bool enableDetailedOutput =
false;
79void StandardHypoTestDemo(
const char *infile =
"",
const char *workspaceName =
"combined",
80 const char *modelSBName =
"ModelConfig",
const char *modelBName =
"",
81 const char *dataName =
"obsData",
int calcType = 0,
83 int ntoys = 5000,
bool useNC =
false,
const char *nuisPriorName = 0)
86 bool noSystematics = optHT.noSystematics;
87 double nToysRatio = optHT.nToysRatio;
88 double poiValue = optHT.poiValue;
89 int printLevel = optHT.printLevel;
90 bool generateBinned = optHT.generateBinned;
91 bool useProof = optHT.useProof;
92 bool enableDetOutput = optHT.enableDetailedOutput;
124 SimpleLikelihoodRatioTestStat::SetAlwaysReuseNLL(
true);
125 ProfileLikelihoodTestStat::SetAlwaysReuseNLL(
true);
126 RatioOfProfiledLikelihoodsTestStat::SetAlwaysReuseNLL(
true);
140 const char *filename =
"";
141 if (!strcmp(infile,
"")) {
142 filename =
"results/example_combined_GaussExample_model.root";
147 cout <<
"HistFactory file cannot be generated on Windows - exit" << endl;
151 cout <<
"will run standard hist2workspace example" << endl;
152 gROOT->ProcessLine(
".! prepareHistFactory .");
153 gROOT->ProcessLine(
".! hist2workspace config/example.xml");
154 cout <<
"\n\n---------------------" << endl;
155 cout <<
"Done creating example input" << endl;
156 cout <<
"---------------------\n\n" << endl;
167 cout <<
"StandardRooStatsDemoMacro: Input file " << filename <<
" is not found" << endl;
178 cout <<
"workspace not found" << endl;
190 if (!data || !sbModel) {
192 cout <<
"data or ModelConfig was not found" << endl;
202 if (nuisPar && nuisPar->
getSize() > 0) {
203 std::cout <<
"StandardHypoTestInvDemo"
204 <<
" - Switch off all systematics by setting them constant to their initial values" << std::endl;
215 Info(
"StandardHypoTestInvDemo",
"The background model %s does not exist", modelBName);
216 Info(
"StandardHypoTestInvDemo",
"Copy it from ModelConfig %s and set POI to zero", modelSBName);
222 double oldval = var->
getVal();
230 Info(
"StandardHypoTestDemo",
"Model %s has no snapshot - make one using model poi", modelSBName);
234 double oldval = var->
getVal();
263 if (testStatType == 3)
267 if (enableDetOutput) {
277 AsymptoticCalculator::SetPrintLevel(printLevel);
283 else if (calcType == 1)
285 else if (calcType == 2)
298 if (testStatType == 3)
300 if (testStatType != 2 && testStatType != 3)
301 Warning(
"StandardHypoTestDemo",
302 "Only the PL test statistic can be used with AsymptoticCalculator - use by default a two-sided PL");
309 nuisPdf = w->
pdf(nuisPriorName);
312 Info(
"StandardHypoTestDemo",
313 "No nuisance pdf given for the HybridCalculator - try to deduce pdf from the model");
322 Info(
"StandardHypoTestDemo",
323 "No nuisance pdf given - try to use %s that is defined as a prior pdf in the B model",
326 Error(
"StandardHypoTestDemo",
"Cannot run Hybrid calculator because no prior on the nuisance parameter is "
327 "specified or can be derived");
332 Info(
"StandardHypoTestDemo",
"Using as nuisance Pdf ... ");
339 Warning(
"StandardHypoTestDemo",
340 "Prior nuisance does not depend on nuisance parameters. They will be smeared in their full range");
353 if (sampler && (calcType == 0 || calcType == 1)) {
358 Warning(
"StandardHypoTestDemo",
"Pdf is extended: but number counting flag is set: ignore it ");
363 Info(
"StandardHypoTestDemo",
364 "Pdf is not extended: number of events to generate taken from observed data set is %d", nEvents);
367 Info(
"StandardHypoTestDemo",
"using a number counting pdf");
373 Info(
"StandardHypoTestDemo",
"Data set is weighted, nentries = %d and sum of weights = %8.1f but toy "
374 "generation is unbinned - it would be faster to set generateBinned to true\n",
387 if (testStatType == 0)
389 if (testStatType == 1)
391 if (testStatType == 2 || testStatType == 3)
410 std::cout <<
"Asymptotic results " << std::endl;
423 for (
int i = 0; i < 5; ++i) {
430 for (
int i = 0; i < 5; ++i) {
431 htExp.SetTestStatisticData(
q[i]);
433 std::cout <<
" Expected p -value and significance at " << sig <<
" sigma = " << htExp.NullPValue()
434 <<
" significance " << htExp.Significance() <<
" sigma " << std::endl;
438 for (
int i = 0; i < 5; ++i) {
442 std::cout <<
" Expected p -value and significance at " << sig <<
" sigma = " << pval <<
" significance "
448 bool writeResult = (calcType != 2);
450 if (enableDetOutput) {
452 Info(
"StandardHypoTestDemo",
"Detailed output will be written in output result file");
455 if (htr != NULL && writeResult) {
458 const char *calcTypeName = (calcType == 0) ?
"Freq" : (calcType == 1) ?
"Hybr" :
"Asym";
464 resultFileName +=
name;
466 TFile *fileOut =
new TFile(resultFileName,
"RECREATE");
468 Info(
"StandardHypoTestDemo",
"HypoTestResult has been written in the file %s", resultFileName.
Data());
void Info(const char *location, const char *msgfmt,...)
Use this function for informational messages.
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
void Warning(const char *location, const char *msgfmt,...)
Use this function in warning situations.
R__EXTERN TSystem * gSystem
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Given a set of possible observables, return the observables that this PDF depends on.
virtual void Print(Option_t *options=0) const
Print the object to the defaultPrintStream().
RooAbsArg * first() const
RooAbsData is the common abstract base class for binned and unbinned datasets.
virtual Double_t sumEntries() const =0
Return effective number of entries in dataset, i.e., sum all weights.
virtual Bool_t isWeighted() const
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
Bool_t canBeExtended() const
If true, PDF can provide extended likelihood term.
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
RooRealVar represents a variable that can be changed from the outside.
virtual void setVal(Double_t value)
Set value of variable to 'value'.
Hypothesis Test Calculator based on the asymptotic formulae for the profile likelihood ratio.
Does a frequentist hypothesis test.
Same purpose as HybridCalculatorOriginal, but different implementation.
Common base class for the Hypothesis Test Calculators.
TestStatSampler * GetTestStatSampler(void) const
Returns instance of TestStatSampler.
virtual HypoTestResult * GetHypoTest() const
inherited methods from HypoTestCalculator interface
This class provides the plots for the result of a study performed with any of the HypoTestCalculatorG...
HypoTestResult is a base class for results from hypothesis tests.
void SetBackgroundAsAlt(Bool_t l=kTRUE)
void SetPValueIsRightTail(Bool_t pr)
virtual Double_t AlternatePValue() const
Return p-value for alternate hypothesis.
virtual Double_t NullPValue() const
Return p-value for null hypothesis.
void Print(const Option_t *="") const
Print out some information about the results Note: use Alt/Null labels for the hypotheses here as the...
SamplingDistribution * GetAltDistribution(void) const
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
virtual void SetSnapshot(const RooArgSet &set)
Set parameter values for a particular hypothesis if using a common PDF by saving a snapshot in the wo...
virtual ModelConfig * Clone(const char *name="") const override
clone
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return NULL if not existing)
const RooArgSet * GetNuisanceParameters() const
get RooArgSet containing the nuisance parameters (return NULL if not existing)
const RooArgSet * GetObservables() const
get RooArgSet for observables (return NULL if not existing)
const RooArgSet * GetSnapshot() const
get RooArgSet for parameters for a particular hypothesis (return NULL if not existing)
RooAbsPdf * GetPdf() const
get model PDF (return NULL if pdf has not been specified or does not exist)
RooAbsPdf * GetPriorPdf() const
get parameters prior pdf (return NULL if not existing)
ProfileLikelihoodTestStat is an implementation of the TestStatistic interface that calculates the pro...
void SetPrintLevel(Int_t printlevel)
virtual void EnableDetailedOutput(bool e=true, bool withErrorsAndPulls=false)
void SetOneSidedDiscovery(Bool_t flag=true)
Holds configuration options for proof and proof-lite.
TestStatistic that returns the ratio of profiled likelihoods.
void SetSubtractMLE(bool subtract)
virtual void EnableDetailedOutput(bool e=true)
void SetLogYaxis(Bool_t ly)
changes plot to log scale on y axis
void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
This class simply holds a sampling distribution of some test statistic.
const std::vector< Double_t > & GetSamplingDistribution() const
Get test statistics values.
TestStatistic class that returns -log(L[null] / L[alt]) where L is the likelihood.
virtual void EnableDetailedOutput(bool e=true)
void SetNullParameters(const RooArgSet &nullParameters)
void SetAltParameters(const RooArgSet &altParameters)
ToyMCSampler is an implementation of the TestStatSampler interface.
void SetProofConfig(ProofConfig *pc=NULL)
virtual void SetTestStatistic(TestStatistic *testStatistic, unsigned int i)
void SetGenerateBinned(bool binned=true)
virtual void SetNEventsPerToy(const Int_t nevents)
Forces the generation of exactly n events even for extended PDFs.
The RooWorkspace is a persistable container for RooFit projects.
RooAbsData * data(const char *name) const
Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found.
void Print(Option_t *opts=0) const
Print contents of the workspace.
TObject * obj(const char *name) const
Return any type of object (RooAbsArg, RooAbsData or generic object) with given name)
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
void Close(Option_t *option="") override
Close a file.
virtual void SetName(const char *name)
Set the name of the TNamed.
virtual const char * GetName() const
Returns name of object.
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
TString & Replace(Ssiz_t pos, Ssiz_t n, const char *s)
const char * Data() const
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
double normal_cdf(double x, double sigma=1, double x0=0)
Cumulative distribution function of the normal (Gaussian) distribution (lower tail).
double normal_quantile_c(double z, double sigma)
Inverse ( ) of the cumulative distribution function of the upper tail of the normal (Gaussian) distri...
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Namespace for the RooStats classes.
bool SetAllConstant(const RooAbsCollection &coll, bool constant=true)
RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name)
void Quantiles(Int_t n, Int_t nprob, Double_t *x, Double_t *quantiles, Double_t *prob, Bool_t isSorted=kTRUE, Int_t *index=0, Int_t type=7)
Computes sample quantiles, corresponding to the given probabilities.