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
129bool useProof = false; // flag to control whether to use Proof
130int nworkers = 0; // number of workers (default use all available cores)
131
132// -------------------------------------------------------
133// The actual macro
134
135void OneSidedFrequentistUpperLimitWithBands(const char *infile = "", const char *workspaceName = "combined",
136 const char *modelConfigName = "ModelConfig",
137 const char *dataName = "obsData")
138{
139
140 double confidenceLevel = 0.95;
141 int nPointsToScan = 12;
142 int nToyMC = 150;
143
144 // -------------------------------------------------------
145 // First part is just to access a user-defined file
146 // or create the standard example file if it doesn't exist
147 const char *filename = "";
148 if (!strcmp(infile, "")) {
149 filename = "results/example_combined_GaussExample_model.root";
150 bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code
151 // if file does not exists generate with histfactory
152 if (!fileExist) {
153 // Normally this would be run on the command line
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;
160 }
161
162 } else
163 filename = infile;
164
165 // Try to open the file
166 TFile *file = TFile::Open(filename);
167
168 // if input file was specified byt not found, quit
169 if (!file) {
170 cout << "StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
171 return;
172 }
173
174 // -------------------------------------------------------
175 // Now get the data and workspace
176
177 // get the workspace out of the file
178 RooWorkspace *w = (RooWorkspace *)file->Get(workspaceName);
179 if (!w) {
180 cout << "workspace not found" << endl;
181 return;
182 }
183
184 // get the modelConfig out of the file
185 ModelConfig *mc = (ModelConfig *)w->obj(modelConfigName);
186
187 // get the modelConfig out of the file
188 RooAbsData *data = w->data(dataName);
189
190 // make sure ingredients are found
191 if (!data || !mc) {
192 w->Print();
193 cout << "data or ModelConfig was not found" << endl;
194 return;
195 }
196
197 // -------------------------------------------------------
198 // Now get the POI for convenience
199 // you may want to adjust the range of your POI
200
201 RooRealVar *firstPOI = (RooRealVar *)mc->GetParametersOfInterest()->first();
202 /* firstPOI->setMin(0);*/
203 /* firstPOI->setMax(10);*/
204
205 // --------------------------------------------
206 // Create and use the FeldmanCousins tool
207 // to find and plot the 95% confidence interval
208 // on the parameter of interest as specified
209 // in the model config
210 // REMEMBER, we will change the test statistic
211 // so this is NOT a Feldman-Cousins interval
212 FeldmanCousins fc(*data, *mc);
213 fc.SetConfidenceLevel(confidenceLevel);
214 fc.AdditionalNToysFactor(
215 0.5); // degrade/improve sampling that defines confidence belt: in this case makes the example faster
216 /* fc.UseAdaptiveSampling(true); // speed it up a bit, don't use for expected limits*/
217 fc.SetNBins(nPointsToScan); // set how many points per parameter of interest to scan
218 fc.CreateConfBelt(true); // save the information in the belt for plotting
219
220 // -------------------------------------------------------
221 // Feldman-Cousins is a unified limit by definition
222 // but the tool takes care of a few things for us like which values
223 // of the nuisance parameters should be used to generate toys.
224 // so let's just change the test statistic and realize this is
225 // no longer "Feldman-Cousins" but is a fully frequentist Neyman-Construction.
226 /* ProfileLikelihoodTestStatModified onesided(*mc->GetPdf());*/
227 /* fc.GetTestStatSampler()->SetTestStatistic(&onesided);*/
228 /* ((ToyMCSampler*) fc.GetTestStatSampler())->SetGenerateBinned(true); */
229 ToyMCSampler *toymcsampler = (ToyMCSampler *)fc.GetTestStatSampler();
230 ProfileLikelihoodTestStat *testStat = dynamic_cast<ProfileLikelihoodTestStat *>(toymcsampler->GetTestStatistic());
231 testStat->SetOneSided(true);
232
233 // Since this tool needs to throw toy MC the PDF needs to be
234 // extended or the tool needs to know how many entries in a dataset
235 // per pseudo experiment.
236 // In the 'number counting form' where the entries in the dataset
237 // are counts, and not values of discriminating variables, the
238 // datasets typically only have one entry and the PDF is not
239 // extended.
240 if (!mc->GetPdf()->canBeExtended()) {
241 if (data->numEntries() == 1)
242 fc.FluctuateNumDataEntries(false);
243 else
244 cout << "Not sure what to do about this model" << endl;
245 }
246
247 // We can use PROOF to speed things along in parallel
248 // However, the test statistic has to be installed on the workers
249 // so either turn off PROOF or include the modified test statistic
250 // in your `$ROOTSYS/roofit/roostats/inc` directory,
251 // add the additional line to the LinkDef.h file,
252 // and recompile root.
253 if (useProof) {
254 ProofConfig pc(*w, nworkers, "", false);
255 toymcsampler->SetProofConfig(&pc); // enable proof
256 }
257
258 if (mc->GetGlobalObservables()) {
259 cout << "will use global observables for unconditional ensemble" << endl;
261 toymcsampler->SetGlobalObservables(*mc->GetGlobalObservables());
262 }
263
264 // Now get the interval
265 PointSetInterval *interval = fc.GetInterval();
266 ConfidenceBelt *belt = fc.GetConfidenceBelt();
267
268 // print out the interval on the first Parameter of Interest
269 cout << "\n95% interval on " << firstPOI->GetName() << " is : [" << interval->LowerLimit(*firstPOI) << ", "
270 << interval->UpperLimit(*firstPOI) << "] " << endl;
271
272 // get observed UL and value of test statistic evaluated there
273 RooArgSet tmpPOI(*firstPOI);
274 double observedUL = interval->UpperLimit(*firstPOI);
275 firstPOI->setVal(observedUL);
276 double obsTSatObsUL = fc.GetTestStatSampler()->EvaluateTestStatistic(*data, tmpPOI);
277
278 // Ask the calculator which points were scanned
279 RooDataSet *parameterScan = (RooDataSet *)fc.GetPointsToScan();
280 RooArgSet *tmpPoint;
281
282 // make a histogram of parameter vs. threshold
283 TH1F *histOfThresholds =
284 new TH1F("histOfThresholds", "", parameterScan->numEntries(), firstPOI->getMin(), firstPOI->getMax());
285 histOfThresholds->GetXaxis()->SetTitle(firstPOI->GetName());
286 histOfThresholds->GetYaxis()->SetTitle("Threshold");
287
288 // loop through the points that were tested and ask confidence belt
289 // what the upper/lower thresholds were.
290 // For FeldmanCousins, the lower cut off is always 0
291 for (Int_t i = 0; i < parameterScan->numEntries(); ++i) {
292 tmpPoint = (RooArgSet *)parameterScan->get(i)->clone("temp");
293 // cout <<"get threshold"<<endl;
294 double arMax = belt->GetAcceptanceRegionMax(*tmpPoint);
295 double poiVal = tmpPoint->getRealValue(firstPOI->GetName());
296 histOfThresholds->Fill(poiVal, arMax);
297 }
298 TCanvas *c1 = new TCanvas();
299 c1->Divide(2);
300 c1->cd(1);
301 histOfThresholds->SetMinimum(0);
302 histOfThresholds->Draw();
303 c1->cd(2);
304
305 // -------------------------------------------------------
306 // Now we generate the expected bands and power-constraint
307
308 // First: find parameter point for mu=0, with conditional MLEs for nuisance parameters
309 std::unique_ptr<RooAbsReal> nll{mc->GetPdf()->createNLL(*data)};
310 std::unique_ptr<RooAbsReal> profile{nll->createProfile(*mc->GetParametersOfInterest())};
311 firstPOI->setVal(0.);
312 profile->getVal(); // this will do fit and set nuisance parameters to profiled values
313 RooArgSet *poiAndNuisance = new RooArgSet();
314 if (mc->GetNuisanceParameters())
315 poiAndNuisance->add(*mc->GetNuisanceParameters());
316 poiAndNuisance->add(*mc->GetParametersOfInterest());
317 w->saveSnapshot("paramsToGenerateData", *poiAndNuisance);
318 RooArgSet *paramsToGenerateData = (RooArgSet *)poiAndNuisance->snapshot();
319 cout << "\nWill use these parameter points to generate pseudo data for bkg only" << endl;
320 paramsToGenerateData->Print("v");
321
322 RooArgSet unconditionalObs;
323 unconditionalObs.add(*mc->GetObservables());
324 unconditionalObs.add(*mc->GetGlobalObservables()); // comment this out for the original conditional ensemble
325
326 double CLb = 0;
327 double CLbinclusive = 0;
328
329 // Now we generate background only and find distribution of upper limits
330 TH1F *histOfUL = new TH1F("histOfUL", "", 100, 0, firstPOI->getMax());
331 histOfUL->GetXaxis()->SetTitle("Upper Limit (background only)");
332 histOfUL->GetYaxis()->SetTitle("Entries");
333 for (int imc = 0; imc < nToyMC; ++imc) {
334
335 // set parameters back to values for generating pseudo data
336 // cout << "\n get current nuis, set vals, print again" << endl;
337 w->loadSnapshot("paramsToGenerateData");
338 // poiAndNuisance->Print("v");
339
340 std::unique_ptr<RooDataSet> toyData;
341 // now generate a toy dataset
342 if (!mc->GetPdf()->canBeExtended()) {
343 if (data->numEntries() == 1)
344 toyData = std::unique_ptr<RooDataSet>{mc->GetPdf()->generate(*mc->GetObservables(), 1)};
345 else
346 cout << "Not sure what to do about this model" << endl;
347 } else {
348 // cout << "generating extended dataset"<<endl;
349 toyData = std::unique_ptr<RooDataSet>{mc->GetPdf()->generate(*mc->GetObservables(), Extended())};
350 }
351
352 // generate global observables
353 // need to be careful for simpdf
354 // RooDataSet* globalData = mc->GetPdf()->generate(*mc->GetGlobalObservables(),1);
355
356 RooSimultaneous *simPdf = dynamic_cast<RooSimultaneous *>(mc->GetPdf());
357 if (!simPdf) {
358 std::unique_ptr<RooDataSet> one{mc->GetPdf()->generate(*mc->GetGlobalObservables(), 1)};
359 const RooArgSet *values = one->get();
360 std::unique_ptr<RooArgSet> allVars{mc->GetPdf()->getVariables()};
361 allVars->assign(*values);
362 } else {
363
364 // try fix for sim pdf
365 for (auto const& tt : simPdf->indexCat()) {
366 auto const& catName = tt.first;
367
368 // Get pdf associated with state from simpdf
369 RooAbsPdf *pdftmp = simPdf->getPdf(catName.c_str());
370
371 // Generate only global variables defined by the pdf associated with this state
372 std::unique_ptr<RooArgSet> globtmp{pdftmp->getObservables(*mc->GetGlobalObservables())};
373 std::unique_ptr<RooDataSet> tmp{pdftmp->generate(*globtmp, 1)};
374
375 // Transfer values to output placeholder
376 globtmp->assign(*tmp->get(0));
377 }
378 }
379
380 // globalData->Print("v");
381 // unconditionalObs = *globalData->get();
382 // mc->GetGlobalObservables()->Print("v");
383 // delete globalData;
384 // cout << "toy data = " << endl;
385 // toyData->get()->Print("v");
386
387 // get test stat at observed UL in observed data
388 firstPOI->setVal(observedUL);
389 double toyTSatObsUL = fc.GetTestStatSampler()->EvaluateTestStatistic(*toyData, tmpPOI);
390 // toyData->get()->Print("v");
391 // cout <<"obsTSatObsUL " <<obsTSatObsUL << "toyTS " << toyTSatObsUL << endl;
392 if (obsTSatObsUL < toyTSatObsUL) // not sure about <= part yet
393 CLb += (1.) / nToyMC;
394 if (obsTSatObsUL <= toyTSatObsUL) // not sure about <= part yet
395 CLbinclusive += (1.) / nToyMC;
396
397 // loop over points in belt to find upper limit for this toy data
398 double thisUL = 0;
399 for (Int_t i = 0; i < parameterScan->numEntries(); ++i) {
400 tmpPoint = (RooArgSet *)parameterScan->get(i)->clone("temp");
401 double arMax = belt->GetAcceptanceRegionMax(*tmpPoint);
402 firstPOI->setVal(tmpPoint->getRealValue(firstPOI->GetName()));
403 // double thisTS = profile->getVal();
404 double thisTS = fc.GetTestStatSampler()->EvaluateTestStatistic(*toyData, tmpPOI);
405
406 // cout << "poi = " << firstPOI->getVal()
407 // << " max is " << arMax << " this profile = " << thisTS << endl;
408 // cout << "thisTS = " << thisTS<<endl;
409 if (thisTS <= arMax) {
410 thisUL = firstPOI->getVal();
411 } else {
412 break;
413 }
414 }
415
416 /*
417 // loop over points in belt to find upper limit for this toy data
418 double thisUL = 0;
419 for(Int_t i=0; i<histOfThresholds->GetNbinsX(); ++i){
420 tmpPoint = (RooArgSet*) parameterScan->get(i)->clone("temp");
421 cout <<"---------------- "<<i<<endl;
422 tmpPoint->Print("v");
423 cout << "from hist " << histOfThresholds->GetBinCenter(i+1) <<endl;
424 double arMax = histOfThresholds->GetBinContent(i+1);
425 // cout << " threshold from Hist = aMax " << arMax<<endl;
426 // double arMax2 = belt->GetAcceptanceRegionMax(*tmpPoint);
427 // cout << "from scan arMax2 = "<< arMax2 << endl; // not the same due to TH1F not TH1D
428 // cout << "scan - hist" << arMax2-arMax << endl;
429 firstPOI->setVal( histOfThresholds->GetBinCenter(i+1));
430 // double thisTS = profile->getVal();
431 double thisTS = fc.GetTestStatSampler()->EvaluateTestStatistic(*toyData,tmpPOI);
432
433 // cout << "poi = " << firstPOI->getVal()
434 // << " max is " << arMax << " this profile = " << thisTS << endl;
435 // cout << "thisTS = " << thisTS<<endl;
436
437 // NOTE: need to add a small epsilon term for single precision vs. double precision
438 if(thisTS<=arMax + 1e-7){
439 thisUL = firstPOI->getVal();
440 } else{
441 break;
442 }
443 }
444 */
445
446 histOfUL->Fill(thisUL);
447
448 // for few events, data is often the same, and UL is often the same
449 // cout << "thisUL = " << thisUL<<endl;
450 }
451 histOfUL->Draw();
452 c1->SaveAs("one-sided_upper_limit_output.pdf");
453
454 // if you want to see a plot of the sampling distribution for a particular scan point:
455 /*
456 SamplingDistPlot sampPlot;
457 int indexInScan = 0;
458 tmpPoint = (RooArgSet*) parameterScan->get(indexInScan)->clone("temp");
459 firstPOI->setVal( tmpPoint->getRealValue(firstPOI->GetName()) );
460 toymcsampler->SetParametersForTestStat(tmpPOI);
461 SamplingDistribution* samp = toymcsampler->GetSamplingDistribution(*tmpPoint);
462 sampPlot.AddSamplingDistribution(samp);
463 sampPlot.Draw();
464 */
465
466 // Now find bands and power constraint
467 Double_t *bins = histOfUL->GetIntegral();
468 TH1F *cumulative = (TH1F *)histOfUL->Clone("cumulative");
469 cumulative->SetContent(bins);
470 double band2sigDown, band1sigDown, bandMedian, band1sigUp, band2sigUp;
471 for (int i = 1; i <= cumulative->GetNbinsX(); ++i) {
472 if (bins[i] < RooStats::SignificanceToPValue(2))
473 band2sigDown = cumulative->GetBinCenter(i);
474 if (bins[i] < RooStats::SignificanceToPValue(1))
475 band1sigDown = cumulative->GetBinCenter(i);
476 if (bins[i] < 0.5)
477 bandMedian = cumulative->GetBinCenter(i);
478 if (bins[i] < RooStats::SignificanceToPValue(-1))
479 band1sigUp = cumulative->GetBinCenter(i);
480 if (bins[i] < RooStats::SignificanceToPValue(-2))
481 band2sigUp = cumulative->GetBinCenter(i);
482 }
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;
488
489 // print out the interval on the first Parameter of Interest
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;
493}
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 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
#define gROOT
Definition TROOT.h:406
R__EXTERN TSystem * gSystem
Definition TSystem.h:555
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:219
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:55
TObject * clone(const char *newname) const override
Definition RooArgSet.h:148
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:191
Container class to hold unbinned data.
Definition RooDataSet.h:57
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...
Holds configuration options for proof and proof-lite.
Definition ProofConfig.h:45
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.
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:4089
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:621
virtual Double_t GetBinCenter(Int_t bin) const
Return bin center for 1D histogram.
Definition TH1.cxx:9109
TAxis * GetXaxis()
Definition TH1.h:324
virtual Int_t GetNbinsX() const
Definition TH1.h:297
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3344
TAxis * GetYaxis()
Definition TH1.h:325
virtual void SetContent(const Double_t *content)
Replace bin contents by the contents of array content.
Definition TH1.cxx:8366
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3066
virtual void SetMinimum(Double_t minimum=-1111)
Definition TH1.h:404
virtual Double_t * GetIntegral()
Return a pointer to the array of bins integral.
Definition TH1.cxx:2586
TObject * Clone(const char *newname="") const override
Make a complete copy of the underlying object.
Definition TH1.cxx:2752
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:353
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
Namespace for the RooStats classes.
Definition Asimov.h:19
double SignificanceToPValue(double Z)
returns p-value corresponding to a 1-sided significance
auto * tt
Definition textangle.C:16