Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
OneSidedFrequentistUpperLimitWithBands.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roostats
3/// \notebook
4/// OneSidedFrequentistUpperLimitWithBands
5///
6/// This is a standard demo that can be used with any ROOT file
7/// prepared in the standard way. You specify:
8/// - name for input ROOT file
9/// - name of workspace inside ROOT file that holds model and data
10/// - name of ModelConfig that specifies details for calculator tools
11/// - name of dataset
12///
13/// With default parameters the macro will attempt to run the
14/// standard hist2workspace example and read the ROOT file
15/// that it produces.
16///
17/// The first ~100 lines define a new test statistic, then the main macro starts.
18/// You may want to control:
19/// ~~~{.cpp}
20/// double confidenceLevel=0.95;
21/// int nPointsToScan = 12;
22/// int nToyMC = 150;
23/// ~~~
24/// This uses a modified version of the profile likelihood ratio as
25/// a test statistic for upper limits (eg. test stat = 0 if muhat>mu).
26///
27/// Based on the observed data, one defines a set of parameter points
28/// to be tested based on the value of the parameter of interest
29/// and the conditional MLE (eg. profiled) values of the nuisance parameters.
30///
31/// At each parameter point, pseudo-experiments are generated using this
32/// fixed reference model and then the test statistic is evaluated.
33/// Note, the nuisance parameters are floating in the fits. For each point,
34/// the threshold that defines the 95% acceptance region is found. This
35/// forms a "Confidence Belt".
36///
37/// After constructing the confidence belt, one can find the confidence
38/// interval for any particular dataset by finding the intersection
39/// of the observed test statistic and the confidence belt. First
40/// this is done on the observed data to get an observed 1-sided upper limt.
41///
42/// Finally, there expected limit and bands (from background-only) are
43/// formed by generating background-only data and finding the upper limit.
44/// This is done by hand for now, will later be part of the RooStats tools.
45///
46/// On a technical note, this technique is NOT the Feldman-Cousins technique,
47/// because that is a 2-sided interval BY DEFINITION. However, like the
48/// Feldman-Cousins technique this is a Neyman-Construction. For technical
49/// reasons the easiest way to implement this right now is to use the
50/// FeldmanCousins tool and then change the test statistic that it is using.
51///
52/// Building the confidence belt can be computationally expensive. Once it is built,
53/// one could save it to a file and use it in a separate step.
54///
55/// We can use PROOF to speed things along in parallel, however,
56/// the test statistic has to be installed on the workers
57/// so either turn off PROOF or include the modified test statistic
58/// in your `$ROOTSYS/roofit/roostats/inc` directory,
59/// add the additional line to the LinkDef.h file,
60/// and recompile root.
61///
62/// Note, if you have a boundary on the parameter of interest (eg. cross-section)
63/// the threshold on the one-sided test statistic starts off very small because we
64/// are only including downward fluctuations. You can see the threshold in these printouts:
65/// ~~~{.cpp}
66/// [#0] PROGRESS:Generation -- generated toys: 500 / 999
67/// NeymanConstruction: Prog: 12/50 total MC = 39 this test stat = 0
68/// SigXsecOverSM=0.69 alpha_syst1=0.136515 alpha_syst3=0.425415 beta_syst2=1.08496 [-1e+30, 0.011215] in interval = 1
69/// ~~~
70/// this tells you the values of the parameters being used to generate the pseudo-experiments
71/// and the threshold in this case is 0.011215. One would expect for 95% that the threshold
72/// would be ~1.35 once the cross-section is far enough away from 0 that it is essentially
73/// unaffected by the boundary. As one reaches the last points in the scan, the
74/// threshold starts to get artificially high. This is because the range of the parameter in
75/// the fit is the same as the range in the scan. In the future, these should be independently
76/// controlled, but they are not now. As a result the ~50% of pseudo-experiments that have an
77/// upward fluctuation end up with muhat = muMax. Because of this, the upper range of the
78/// parameter should be well above the expected upper limit... but not too high or one will
79/// need a very large value of nPointsToScan to resolve the relevant region. This can be
80/// improved, but this is the first version of this script.
81///
82/// Important note: when the model includes external constraint terms, like a Gaussian
83/// constraint to a nuisance parameter centered around some nominal value there is
84/// a subtlety. The asymptotic results are all based on the assumption that all the
85/// measurements fluctuate... including the nominal values from auxiliary measurements.
86/// If these do not fluctuate, this corresponds to an "conditional ensemble". The
87/// result is that the distribution of the test statistic can become very non-chi^2.
88/// This results in thresholds that become very large. This can be seen in the following
89/// thought experiment. Say the model is
90/// \f$ Pois(N | s + b)G(b0|b,sigma) \f$
91/// where \f$ G(b0|b,sigma) \f$ is the external constraint and b0 is 100. If N is also 100
92/// then the profiled value of b given s is going to be some trade off between 100-s and b0.
93/// If sigma is \f$ \sqrt(N) \f$, then the profiled value of b is probably 100 - s/2 So for
94/// s=60 we are going to have a profiled value of b~70. Now when we generate pseudo-experiments
95/// for s=60, b=70 we will have N~130 and the average shat will be 30, not 60. In practice,
96/// this is only an issue for values of s that are very excluded. For values of s near the 95%
97/// limit this should not be a big effect. This can be avoided if the nominal values of the constraints also fluctuate,
98/// but that requires that those parameters are RooRealVars in the model.
99/// This version does not deal with this issue, but it will be addressed in a future version.
100///
101/// \macro_image
102/// \macro_output
103/// \macro_code
104///
105/// \authors Kyle Cranmer, Haichen Wang, Daniel Whiteson
106
107#include "TFile.h"
108#include "TROOT.h"
109#include "TH1F.h"
110#include "TCanvas.h"
111#include "TSystem.h"
112
113#include "RooWorkspace.h"
114#include "RooSimultaneous.h"
115#include "RooAbsData.h"
116
117#include "RooStats/ModelConfig.h"
122
125
126using namespace RooFit;
127using namespace RooStats;
128
129// -------------------------------------------------------
130// The actual macro
131
132void OneSidedFrequentistUpperLimitWithBands(const char *infile = "", const char *workspaceName = "combined",
133 const char *modelConfigName = "ModelConfig",
134 const char *dataName = "obsData")
135{
136
137 double confidenceLevel = 0.95;
138 int nPointsToScan = 12;
139 int nToyMC = 150;
140
141 // -------------------------------------------------------
142 // First part is just to access a user-defined file
143 // or create the standard example file if it doesn't exist
144 const char *filename = "";
145 if (!strcmp(infile, "")) {
146 filename = "results/example_combined_GaussExample_model.root";
147 bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code
148 // if file does not exists generate with histfactory
149 if (!fileExist) {
150 // Normally this would be run on the command line
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;
157 }
158
159 } else
160 filename = infile;
161
162 // Try to open the file
163 TFile *file = TFile::Open(filename);
164
165 // if input file was specified but not found, quit
166 if (!file) {
167 cout << "StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
168 return;
169 }
170
171 // -------------------------------------------------------
172 // Now get the data and workspace
173
174 // get the workspace out of the file
175 RooWorkspace *w = (RooWorkspace *)file->Get(workspaceName);
176 if (!w) {
177 cout << "workspace not found" << endl;
178 return;
179 }
180
181 // get the modelConfig out of the file
182 ModelConfig *mc = (ModelConfig *)w->obj(modelConfigName);
183
184 // get the modelConfig out of the file
185 RooAbsData *data = w->data(dataName);
186
187 // make sure ingredients are found
188 if (!data || !mc) {
189 w->Print();
190 cout << "data or ModelConfig was not found" << endl;
191 return;
192 }
193
194 // -------------------------------------------------------
195 // Now get the POI for convenience
196 // you may want to adjust the range of your POI
197
198 RooRealVar *firstPOI = (RooRealVar *)mc->GetParametersOfInterest()->first();
199 /* firstPOI->setMin(0);*/
200 /* firstPOI->setMax(10);*/
201
202 // --------------------------------------------
203 // Create and use the FeldmanCousins tool
204 // to find and plot the 95% confidence interval
205 // on the parameter of interest as specified
206 // in the model config
207 // REMEMBER, we will change the test statistic
208 // so this is NOT a Feldman-Cousins interval
209 FeldmanCousins fc(*data, *mc);
210 fc.SetConfidenceLevel(confidenceLevel);
211 fc.AdditionalNToysFactor(
212 0.5); // degrade/improve sampling that defines confidence belt: in this case makes the example faster
213 /* fc.UseAdaptiveSampling(true); // speed it up a bit, don't use for expected limits*/
214 fc.SetNBins(nPointsToScan); // set how many points per parameter of interest to scan
215 fc.CreateConfBelt(true); // save the information in the belt for plotting
216
217 // -------------------------------------------------------
218 // Feldman-Cousins is a unified limit by definition
219 // but the tool takes care of a few things for us like which values
220 // of the nuisance parameters should be used to generate toys.
221 // so let's just change the test statistic and realize this is
222 // no longer "Feldman-Cousins" but is a fully frequentist Neyman-Construction.
223 /* ProfileLikelihoodTestStatModified onesided(*mc->GetPdf());*/
224 /* fc.GetTestStatSampler()->SetTestStatistic(&onesided);*/
225 /* ((ToyMCSampler*) fc.GetTestStatSampler())->SetGenerateBinned(true); */
226 ToyMCSampler *toymcsampler = (ToyMCSampler *)fc.GetTestStatSampler();
227 ProfileLikelihoodTestStat *testStat = dynamic_cast<ProfileLikelihoodTestStat *>(toymcsampler->GetTestStatistic());
228 testStat->SetOneSided(true);
229
230 // Since this tool needs to throw toy MC the PDF needs to be
231 // extended or the tool needs to know how many entries in a dataset
232 // per pseudo experiment.
233 // In the 'number counting form' where the entries in the dataset
234 // are counts, and not values of discriminating variables, the
235 // datasets typically only have one entry and the PDF is not
236 // extended.
237 if (!mc->GetPdf()->canBeExtended()) {
238 if (data->numEntries() == 1)
239 fc.FluctuateNumDataEntries(false);
240 else
241 cout << "Not sure what to do about this model" << endl;
242 }
243
244 if (mc->GetGlobalObservables()) {
245 cout << "will use global observables for unconditional ensemble" << endl;
247 toymcsampler->SetGlobalObservables(*mc->GetGlobalObservables());
248 }
249
250 // Now get the interval
251 PointSetInterval *interval = fc.GetInterval();
252 ConfidenceBelt *belt = fc.GetConfidenceBelt();
253
254 // print out the interval on the first Parameter of Interest
255 cout << "\n95% interval on " << firstPOI->GetName() << " is : [" << interval->LowerLimit(*firstPOI) << ", "
256 << interval->UpperLimit(*firstPOI) << "] " << endl;
257
258 // get observed UL and value of test statistic evaluated there
259 RooArgSet tmpPOI(*firstPOI);
260 double observedUL = interval->UpperLimit(*firstPOI);
261 firstPOI->setVal(observedUL);
262 double obsTSatObsUL = fc.GetTestStatSampler()->EvaluateTestStatistic(*data, tmpPOI);
263
264 // Ask the calculator which points were scanned
265 RooDataSet *parameterScan = (RooDataSet *)fc.GetPointsToScan();
266 RooArgSet *tmpPoint;
267
268 // make a histogram of parameter vs. threshold
269 TH1F *histOfThresholds =
270 new TH1F("histOfThresholds", "", parameterScan->numEntries(), firstPOI->getMin(), firstPOI->getMax());
271 histOfThresholds->GetXaxis()->SetTitle(firstPOI->GetName());
272 histOfThresholds->GetYaxis()->SetTitle("Threshold");
273
274 // loop through the points that were tested and ask confidence belt
275 // what the upper/lower thresholds were.
276 // For FeldmanCousins, the lower cut off is always 0
277 for (Int_t i = 0; i < parameterScan->numEntries(); ++i) {
278 tmpPoint = (RooArgSet *)parameterScan->get(i)->clone("temp");
279 // cout <<"get threshold"<<endl;
280 double arMax = belt->GetAcceptanceRegionMax(*tmpPoint);
281 double poiVal = tmpPoint->getRealValue(firstPOI->GetName());
282 histOfThresholds->Fill(poiVal, arMax);
283 }
284 TCanvas *c1 = new TCanvas();
285 c1->Divide(2);
286 c1->cd(1);
287 histOfThresholds->SetMinimum(0);
288 histOfThresholds->Draw();
289 c1->cd(2);
290
291 // -------------------------------------------------------
292 // Now we generate the expected bands and power-constraint
293
294 // First: find parameter point for mu=0, with conditional MLEs for nuisance parameters
295 std::unique_ptr<RooAbsReal> nll{mc->GetPdf()->createNLL(*data)};
296 std::unique_ptr<RooAbsReal> profile{nll->createProfile(*mc->GetParametersOfInterest())};
297 firstPOI->setVal(0.);
298 profile->getVal(); // this will do fit and set nuisance parameters to profiled values
299 RooArgSet *poiAndNuisance = new RooArgSet();
300 if (mc->GetNuisanceParameters())
301 poiAndNuisance->add(*mc->GetNuisanceParameters());
302 poiAndNuisance->add(*mc->GetParametersOfInterest());
303 w->saveSnapshot("paramsToGenerateData", *poiAndNuisance);
304 RooArgSet *paramsToGenerateData = (RooArgSet *)poiAndNuisance->snapshot();
305 cout << "\nWill use these parameter points to generate pseudo data for bkg only" << endl;
306 paramsToGenerateData->Print("v");
307
308 RooArgSet unconditionalObs;
309 unconditionalObs.add(*mc->GetObservables());
310 unconditionalObs.add(*mc->GetGlobalObservables()); // comment this out for the original conditional ensemble
311
312 double CLb = 0;
313 double CLbinclusive = 0;
314
315 // Now we generate background only and find distribution of upper limits
316 TH1F *histOfUL = new TH1F("histOfUL", "", 100, 0, firstPOI->getMax());
317 histOfUL->GetXaxis()->SetTitle("Upper Limit (background only)");
318 histOfUL->GetYaxis()->SetTitle("Entries");
319 for (int imc = 0; imc < nToyMC; ++imc) {
320
321 // set parameters back to values for generating pseudo data
322 // cout << "\n get current nuis, set vals, print again" << endl;
323 w->loadSnapshot("paramsToGenerateData");
324 // poiAndNuisance->Print("v");
325
326 std::unique_ptr<RooDataSet> toyData;
327 // now generate a toy dataset
328 if (!mc->GetPdf()->canBeExtended()) {
329 if (data->numEntries() == 1)
330 toyData = std::unique_ptr<RooDataSet>{mc->GetPdf()->generate(*mc->GetObservables(), 1)};
331 else
332 cout << "Not sure what to do about this model" << endl;
333 } else {
334 // cout << "generating extended dataset"<<endl;
335 toyData = std::unique_ptr<RooDataSet>{mc->GetPdf()->generate(*mc->GetObservables(), Extended())};
336 }
337
338 // generate global observables
339 // need to be careful for simpdf
340 // RooDataSet* globalData = mc->GetPdf()->generate(*mc->GetGlobalObservables(),1);
341
342 RooSimultaneous *simPdf = dynamic_cast<RooSimultaneous *>(mc->GetPdf());
343 if (!simPdf) {
344 std::unique_ptr<RooDataSet> one{mc->GetPdf()->generate(*mc->GetGlobalObservables(), 1)};
345 const RooArgSet *values = one->get();
346 std::unique_ptr<RooArgSet> allVars{mc->GetPdf()->getVariables()};
347 allVars->assign(*values);
348 } else {
349
350 // try fix for sim pdf
351 for (auto const& tt : simPdf->indexCat()) {
352 auto const& catName = tt.first;
353
354 // Get pdf associated with state from simpdf
355 RooAbsPdf *pdftmp = simPdf->getPdf(catName.c_str());
356
357 // Generate only global variables defined by the pdf associated with this state
358 std::unique_ptr<RooArgSet> globtmp{pdftmp->getObservables(*mc->GetGlobalObservables())};
359 std::unique_ptr<RooDataSet> tmp{pdftmp->generate(*globtmp, 1)};
360
361 // Transfer values to output placeholder
362 globtmp->assign(*tmp->get(0));
363 }
364 }
365
366 // globalData->Print("v");
367 // unconditionalObs = *globalData->get();
368 // mc->GetGlobalObservables()->Print("v");
369 // delete globalData;
370 // cout << "toy data = " << endl;
371 // toyData->get()->Print("v");
372
373 // get test stat at observed UL in observed data
374 firstPOI->setVal(observedUL);
375 double toyTSatObsUL = fc.GetTestStatSampler()->EvaluateTestStatistic(*toyData, tmpPOI);
376 // toyData->get()->Print("v");
377 // cout <<"obsTSatObsUL " <<obsTSatObsUL << "toyTS " << toyTSatObsUL << endl;
378 if (obsTSatObsUL < toyTSatObsUL) // not sure about <= part yet
379 CLb += (1.) / nToyMC;
380 if (obsTSatObsUL <= toyTSatObsUL) // not sure about <= part yet
381 CLbinclusive += (1.) / nToyMC;
382
383 // loop over points in belt to find upper limit for this toy data
384 double thisUL = 0;
385 for (Int_t i = 0; i < parameterScan->numEntries(); ++i) {
386 tmpPoint = (RooArgSet *)parameterScan->get(i)->clone("temp");
387 double arMax = belt->GetAcceptanceRegionMax(*tmpPoint);
388 firstPOI->setVal(tmpPoint->getRealValue(firstPOI->GetName()));
389 // double thisTS = profile->getVal();
390 double thisTS = fc.GetTestStatSampler()->EvaluateTestStatistic(*toyData, tmpPOI);
391
392 // cout << "poi = " << firstPOI->getVal()
393 // << " max is " << arMax << " this profile = " << thisTS << endl;
394 // cout << "thisTS = " << thisTS<<endl;
395 if (thisTS <= arMax) {
396 thisUL = firstPOI->getVal();
397 } else {
398 break;
399 }
400 }
401
402 /*
403 // loop over points in belt to find upper limit for this toy data
404 double thisUL = 0;
405 for(Int_t i=0; i<histOfThresholds->GetNbinsX(); ++i){
406 tmpPoint = (RooArgSet*) parameterScan->get(i)->clone("temp");
407 cout <<"---------------- "<<i<<endl;
408 tmpPoint->Print("v");
409 cout << "from hist " << histOfThresholds->GetBinCenter(i+1) <<endl;
410 double arMax = histOfThresholds->GetBinContent(i+1);
411 // cout << " threshold from Hist = aMax " << arMax<<endl;
412 // double arMax2 = belt->GetAcceptanceRegionMax(*tmpPoint);
413 // cout << "from scan arMax2 = "<< arMax2 << endl; // not the same due to TH1F not TH1D
414 // cout << "scan - hist" << arMax2-arMax << endl;
415 firstPOI->setVal( histOfThresholds->GetBinCenter(i+1));
416 // double thisTS = profile->getVal();
417 double thisTS = fc.GetTestStatSampler()->EvaluateTestStatistic(*toyData,tmpPOI);
418
419 // cout << "poi = " << firstPOI->getVal()
420 // << " max is " << arMax << " this profile = " << thisTS << endl;
421 // cout << "thisTS = " << thisTS<<endl;
422
423 // NOTE: need to add a small epsilon term for single precision vs. double precision
424 if(thisTS<=arMax + 1e-7){
425 thisUL = firstPOI->getVal();
426 } else{
427 break;
428 }
429 }
430 */
431
432 histOfUL->Fill(thisUL);
433
434 // for few events, data is often the same, and UL is often the same
435 // cout << "thisUL = " << thisUL<<endl;
436 }
437 histOfUL->Draw();
438 c1->SaveAs("one-sided_upper_limit_output.pdf");
439
440 // if you want to see a plot of the sampling distribution for a particular scan point:
441 /*
442 SamplingDistPlot sampPlot;
443 int indexInScan = 0;
444 tmpPoint = (RooArgSet*) parameterScan->get(indexInScan)->clone("temp");
445 firstPOI->setVal( tmpPoint->getRealValue(firstPOI->GetName()) );
446 toymcsampler->SetParametersForTestStat(tmpPOI);
447 SamplingDistribution* samp = toymcsampler->GetSamplingDistribution(*tmpPoint);
448 sampPlot.AddSamplingDistribution(samp);
449 sampPlot.Draw();
450 */
451
452 // Now find bands and power constraint
453 Double_t *bins = histOfUL->GetIntegral();
454 TH1F *cumulative = (TH1F *)histOfUL->Clone("cumulative");
455 cumulative->SetContent(bins);
456 double band2sigDown, band1sigDown, bandMedian, band1sigUp, band2sigUp;
457 for (int i = 1; i <= cumulative->GetNbinsX(); ++i) {
458 if (bins[i] < RooStats::SignificanceToPValue(2))
459 band2sigDown = cumulative->GetBinCenter(i);
460 if (bins[i] < RooStats::SignificanceToPValue(1))
461 band1sigDown = cumulative->GetBinCenter(i);
462 if (bins[i] < 0.5)
463 bandMedian = cumulative->GetBinCenter(i);
464 if (bins[i] < RooStats::SignificanceToPValue(-1))
465 band1sigUp = cumulative->GetBinCenter(i);
466 if (bins[i] < RooStats::SignificanceToPValue(-2))
467 band2sigUp = cumulative->GetBinCenter(i);
468 }
469 cout << "-2 sigma band " << band2sigDown << endl;
470 cout << "-1 sigma band " << band1sigDown << " [Power Constraint)]" << endl;
471 cout << "median of band " << bandMedian << endl;
472 cout << "+1 sigma band " << band1sigUp << endl;
473 cout << "+2 sigma band " << band2sigUp << endl;
474
475 // print out the interval on the first Parameter of Interest
476 cout << "\nobserved 95% upper-limit " << interval->UpperLimit(*firstPOI) << endl;
477 cout << "CLb strict [P(toy>obs|0)] for observed 95% upper-limit " << CLb << endl;
478 cout << "CLb inclusive [P(toy>=obs|0)] for observed 95% upper-limit " << CLbinclusive << endl;
479}
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
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
#define gROOT
Definition TROOT.h:406
R__EXTERN TSystem * gSystem
Definition TSystem.h:561
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.
Definition RooAbsData.h:57
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
RooFit::OwningPtr< RooAbsReal > createNLL(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Construct representation of -log(L) of PDF with given dataset.
Definition RooAbsPdf.h:163
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition RooAbsPdf.h:218
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&,...
Definition RooAbsPdf.h:57
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.
Definition RooAbsReal.h:103
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
TObject * clone(const char *newname) const override
Definition RooArgSet.h:111
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:154
Container class to hold unbinned data.
Definition RooDataSet.h:34
const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
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...
Definition ModelConfig.h:35
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 &param)
return upper limit on a given parameter
double LowerLimit(RooRealVar &param)
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.
The Canvas class.
Definition TCanvas.h:23
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...
Definition TFile.h:53
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.
Definition TFile.cxx:4086
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:623
virtual Double_t GetBinCenter(Int_t bin) const
Return bin center for 1D histogram.
Definition TH1.cxx:9174
TAxis * GetXaxis()
Definition TH1.h:325
virtual Int_t GetNbinsX() const
Definition TH1.h:298
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3346
TAxis * GetYaxis()
Definition TH1.h:326
virtual void SetContent(const Double_t *content)
Replace bin contents by the contents of array content.
Definition TH1.cxx:8431
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3068
virtual void SetMinimum(Double_t minimum=-1111)
Definition TH1.h:406
virtual Double_t * GetIntegral()
Return a pointer to the array of bins integral.
Definition TH1.cxx:2588
TObject * Clone(const char *newname="") const override
Make a complete copy of the underlying object.
Definition TH1.cxx:2754
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition TSystem.cxx:1296
RooCmdArg Extended(bool flag=true)
return c1
Definition legend1.C:41
double nll(double pdf, double weight, int binnedL, int doBinOffset)
Definition MathFuncs.h:358
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64
Namespace for the RooStats classes.
Definition CodegenImpl.h:58
double SignificanceToPValue(double Z)
returns p-value corresponding to a 1-sided significance
auto * tt
Definition textangle.C:16