129bool useProof =
false;
135void OneSidedFrequentistUpperLimitWithBands(
const char *infile =
"",
const char *workspaceName =
"combined",
136 const char *modelConfigName =
"ModelConfig",
137 const char *dataName =
"obsData")
140 double confidenceLevel = 0.95;
141 int nPointsToScan = 12;
148 if (!strcmp(infile,
"")) {
149 filename =
"results/example_combined_GaussExample_model.root";
154 cout <<
"will run standard hist2workspace example" << endl;
155 gROOT->ProcessLine(
".! prepareHistFactory .");
156 gROOT->ProcessLine(
".! hist2workspace config/example.xml");
157 cout <<
"\n\n---------------------" << endl;
158 cout <<
"Done creating example input" << endl;
159 cout <<
"---------------------\n\n" << endl;
170 cout <<
"StandardRooStatsDemoMacro: Input file " <<
filename <<
" is not found" << endl;
180 cout <<
"workspace not found" << endl;
193 cout <<
"data or ModelConfig was not found" << endl;
213 fc.SetConfidenceLevel(confidenceLevel);
214 fc.AdditionalNToysFactor(
217 fc.SetNBins(nPointsToScan);
218 fc.CreateConfBelt(
true);
241 if (
data->numEntries() == 1)
242 fc.FluctuateNumDataEntries(
false);
244 cout <<
"Not sure what to do about this model" << endl;
259 cout <<
"will use global observables for unconditional ensemble" << endl;
269 cout <<
"\n95% interval on " << firstPOI->
GetName() <<
" is : [" << interval->
LowerLimit(*firstPOI) <<
", "
270 << interval->
UpperLimit(*firstPOI) <<
"] " << endl;
274 double observedUL = interval->
UpperLimit(*firstPOI);
275 firstPOI->
setVal(observedUL);
276 double obsTSatObsUL = fc.GetTestStatSampler()->EvaluateTestStatistic(*
data, tmpPOI);
283 TH1F *histOfThresholds =
296 histOfThresholds->
Fill(poiVal, arMax);
302 histOfThresholds->
Draw();
317 w->saveSnapshot(
"paramsToGenerateData", *poiAndNuisance);
319 cout <<
"\nWill use these parameter points to generate pseudo data for bkg only" << endl;
320 paramsToGenerateData->
Print(
"v");
327 double CLbinclusive = 0;
330 TH1F *histOfUL =
new TH1F(
"histOfUL",
"", 100, 0, firstPOI->
getMax());
333 for (
int imc = 0; imc < nToyMC; ++imc) {
337 w->loadSnapshot(
"paramsToGenerateData");
340 std::unique_ptr<RooDataSet> toyData;
343 if (
data->numEntries() == 1)
346 cout <<
"Not sure what to do about this model" << endl;
361 allVars->assign(*values);
365 for (
auto const&
tt : simPdf->indexCat()) {
366 auto const& catName =
tt.first;
373 std::unique_ptr<RooDataSet>
tmp{pdftmp->
generate(*globtmp, 1)};
376 globtmp->assign(*
tmp->get(0));
388 firstPOI->
setVal(observedUL);
389 double toyTSatObsUL = fc.GetTestStatSampler()->EvaluateTestStatistic(*toyData, tmpPOI);
392 if (obsTSatObsUL < toyTSatObsUL)
393 CLb += (1.) / nToyMC;
394 if (obsTSatObsUL <= toyTSatObsUL)
395 CLbinclusive += (1.) / nToyMC;
404 double thisTS = fc.GetTestStatSampler()->EvaluateTestStatistic(*toyData, tmpPOI);
409 if (thisTS <= arMax) {
410 thisUL = firstPOI->
getVal();
446 histOfUL->
Fill(thisUL);
452 c1->SaveAs(
"one-sided_upper_limit_output.pdf");
470 double band2sigDown, band1sigDown, bandMedian, band1sigUp, band2sigUp;
471 for (
int i = 1; i <= cumulative->
GetNbinsX(); ++i) {
483 cout <<
"-2 sigma band " << band2sigDown << endl;
484 cout <<
"-1 sigma band " << band1sigDown <<
" [Power Constraint)]" << endl;
485 cout <<
"median of band " << bandMedian << endl;
486 cout <<
"+1 sigma band " << band1sigUp << endl;
487 cout <<
"+2 sigma band " << band2sigUp << endl;
490 cout <<
"\nobserved 95% upper-limit " << interval->
UpperLimit(*firstPOI) << endl;
491 cout <<
"CLb strict [P(toy>obs|0)] for observed 95% upper-limit " << CLb << endl;
492 cout <<
"CLb inclusive [P(toy>=obs|0)] for observed 95% upper-limit " << CLbinclusive << endl;
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
R__EXTERN TSystem * gSystem
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
RooFit::OwningPtr< RooArgSet > getVariables(bool stripDisconnected=true) const
Return RooArgSet with all variables (tree leaf nodes of expression tree)
double getRealValue(const char *name, double defVal=0.0, bool verbose=false) const
Get value of a RooAbsReal stored in set with given name.
Storage_t const & get() const
Const access to the underlying stl container.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
RooAbsArg * first() const
void Print(Option_t *options=nullptr) const override
This method must be overridden when a class wants to print itself.
Abstract base class for binned and unbinned datasets.
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
Abstract interface for all probability density functions.
RooFit::OwningPtr< RooAbsReal > createNLL(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Construct representation of -log(L) of PDF with given dataset.
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
RooFit::OwningPtr< RooDataSet > generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={})
See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,...
virtual double getMax(const char *name=nullptr) const
Get maximum of currently defined range.
virtual double getMin(const char *name=nullptr) const
Get minimum of currently defined range.
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
TObject * clone(const char *newname) const override
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Container class to hold unbinned data.
const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
Variable that can be changed from the outside.
void setVal(double value) override
Set value of variable to 'value'.
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
RooAbsPdf * getPdf(RooStringView catName) const
Return the p.d.f associated with the given index category name.
ConfidenceBelt is a concrete implementation of the ConfInterval interface.
double GetAcceptanceRegionMax(RooArgSet &, double cl=-1., double leftside=-1.)
The FeldmanCousins class (like the Feldman-Cousins technique) is essentially a specific configuration...
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
const RooArgSet * GetGlobalObservables() const
get RooArgSet for global observables (return nullptr if not existing)
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return nullptr if not existing)
const RooArgSet * GetNuisanceParameters() const
get RooArgSet containing the nuisance parameters (return nullptr if not existing)
const RooArgSet * GetObservables() const
get RooArgSet for observables (return nullptr if not existing)
RooAbsPdf * GetPdf() const
get model PDF (return nullptr if pdf has not been specified or does not exist)
PointSetInterval is a concrete implementation of the ConfInterval interface.
double UpperLimit(RooRealVar ¶m)
return upper limit on a given parameter
double LowerLimit(RooRealVar ¶m)
return lower limit on a given parameter
ProfileLikelihoodTestStat is an implementation of the TestStatistic interface that calculates the pro...
void SetOneSided(bool flag=true)
Holds configuration options for proof and proof-lite.
ToyMCSampler is an implementation of the TestStatSampler interface.
void SetProofConfig(ProofConfig *pc=nullptr)
calling with argument or nullptr deactivates proof
virtual TestStatistic * GetTestStatistic(unsigned int i) const
void SetGlobalObservables(const RooArgSet &o) override
specify the conditional observables
Persistable container for RooFit projects.
TObject * Get(const char *namecycle) override
Return pointer to object identified by namecycle.
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
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.
1-D histogram with a float per channel (see TH1 documentation)
virtual Double_t GetBinCenter(Int_t bin) const
Return bin center for 1D histogram.
virtual Int_t GetNbinsX() const
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
virtual void SetContent(const Double_t *content)
Replace bin contents by the contents of array content.
void Draw(Option_t *option="") override
Draw this histogram with options.
virtual void SetMinimum(Double_t minimum=-1111)
virtual Double_t * GetIntegral()
Return a pointer to the array of bins integral.
TObject * Clone(const char *newname="") const override
Make a complete copy of the underlying object.
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
const char * GetName() const override
Returns name of object.
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
RooCmdArg Extended(bool flag=true)
double nll(double pdf, double weight, int binnedL, int doBinOffset)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Namespace for the RooStats classes.
double SignificanceToPValue(double Z)
returns p-value corresponding to a 1-sided significance