117using std::cout, std::endl;
121void TwoSidedFrequentistUpperLimitWithBands(
const char *infile =
"",
const char *workspaceName =
"combined",
122 const char *modelConfigName =
"ModelConfig",
123 const char *dataName =
"obsData")
126 double confidenceLevel = 0.95;
129 double additionalToysFac = 0.5;
130 int nPointsToScan = 20;
136 const char *filename =
"";
137 if (!strcmp(infile,
"")) {
138 filename =
"results/example_combined_GaussExample_model.root";
139 bool fileExist = !
gSystem->AccessPathName(filename);
143 cout <<
"will run standard hist2workspace example" << endl;
144 gROOT->ProcessLine(
".! prepareHistFactory .");
145 gROOT->ProcessLine(
".! hist2workspace config/example.xml");
146 cout <<
"\n\n---------------------" << endl;
147 cout <<
"Done creating example input" << endl;
148 cout <<
"---------------------\n\n" << endl;
169 cout <<
"Found data and ModelConfig:" << endl;
187 fc.SetConfidenceLevel(confidenceLevel);
188 fc.AdditionalNToysFactor(additionalToysFac);
190 fc.SetNBins(nPointsToScan);
191 fc.CreateConfBelt(
true);
212 if (data->numEntries() == 1)
213 fc.FluctuateNumDataEntries(
false);
215 cout <<
"Not sure what to do about this model" << endl;
219 cout <<
"will use global observables for unconditional ensemble" << endl;
229 cout <<
"\n95% interval on " << firstPOI->
GetName() <<
" is : [" << interval->
LowerLimit(*firstPOI) <<
", "
230 << interval->
UpperLimit(*firstPOI) <<
"] " << endl;
234 double observedUL = interval->
UpperLimit(*firstPOI);
235 firstPOI->
setVal(observedUL);
236 double obsTSatObsUL = fc.GetTestStatSampler()->EvaluateTestStatistic(*data, tmpPOI);
243 TH1F *histOfThresholds =
256 histOfThresholds->
Fill(poiVal, arMax);
262 histOfThresholds->
Draw();
277 w->saveSnapshot(
"paramsToGenerateData", *poiAndNuisance);
279 cout <<
"\nWill use these parameter points to generate pseudo data for bkg only" << endl;
280 paramsToGenerateData->
Print(
"v");
287 double CLbinclusive = 0;
290 TH1F *histOfUL =
new TH1F(
"histOfUL",
"", 100, 0, firstPOI->
getMax());
293 for (
int imc = 0; imc < nToyMC; ++imc) {
297 w->loadSnapshot(
"paramsToGenerateData");
300 std::unique_ptr<RooDataSet> toyData;
303 if (data->numEntries() == 1)
306 cout <<
"Not sure what to do about this model" << endl;
316 allVars->assign(*values);
319 firstPOI->
setVal(observedUL);
320 double toyTSatObsUL = fc.GetTestStatSampler()->EvaluateTestStatistic(*toyData, tmpPOI);
323 if (obsTSatObsUL < toyTSatObsUL)
324 CLb += (1.) / nToyMC;
325 if (obsTSatObsUL <= toyTSatObsUL)
326 CLbinclusive += (1.) / nToyMC;
335 double thisTS = fc.GetTestStatSampler()->EvaluateTestStatistic(*toyData, tmpPOI);
340 if (thisTS <= arMax) {
341 thisUL = firstPOI->
getVal();
347 histOfUL->
Fill(thisUL);
353 c1->SaveAs(
"two-sided_upper_limit_output.pdf");
371 double band2sigDown = 0, band1sigDown = 0, bandMedian = 0, band1sigUp = 0, band2sigUp = 0;
372 for (
int i = 1; i <= cumulative->
GetNbinsX(); ++i) {
384 cout <<
"-2 sigma band " << band2sigDown << endl;
385 cout <<
"-1 sigma band " << band1sigDown <<
" [Power Constraint)]" << endl;
386 cout <<
"median of band " << bandMedian << endl;
387 cout <<
"+1 sigma band " << band1sigUp << endl;
388 cout <<
"+2 sigma band " << band2sigUp << endl;
391 cout <<
"\nobserved 95% upper-limit " << interval->
UpperLimit(*firstPOI) << endl;
392 cout <<
"CLb strict [P(toy>obs|0)] for observed 95% upper-limit " << CLb << endl;
393 cout <<
"CLb inclusive [P(toy>=obs|0)] for observed 95% upper-limit " << CLbinclusive << endl;
int Int_t
Signed integer 4 bytes (int).
double Double_t
Double 8 bytes.
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.
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.
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
RooFit::OwningPtr< RooAbsReal > createNLL(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Construct representation of -log(L) of PDF with given dataset.
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 RooFit::OwningPtr< RooDataSet > generateSimGlobal(const RooArgSet &whatVars, Int_t nEvents)
Special generator interface for generation of 'global observables' – for RooStats tools.
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.
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
TObject * clone(const char *newname=nullptr) const override
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'.
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...
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)
void Print(Option_t *option="") const override
overload the print method
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...
ToyMCSampler is an implementation of the TestStatSampler interface.
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 file, usually with extension .root, that stores data and code in the form of serialized objects in ...
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.
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...
RooStats::ModelConfig ModelConfig
Namespace for the RooStats classes.
double SignificanceToPValue(double Z)
returns p-value corresponding to a 1-sided significance