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/// theshold 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#ifdef _WIN32
154 cout << "HistFactory file cannot be generated on Windows - exit" << endl;
155 return;
156#endif
157 // Normally this would be run on the command line
158 cout << "will run standard hist2workspace example" << endl;
159 gROOT->ProcessLine(".! prepareHistFactory .");
160 gROOT->ProcessLine(".! hist2workspace config/example.xml");
161 cout << "\n\n---------------------" << endl;
162 cout << "Done creating example input" << endl;
163 cout << "---------------------\n\n" << endl;
164 }
165
166 } else
167 filename = infile;
168
169 // Try to open the file
170 TFile *file = TFile::Open(filename);
171
172 // if input file was specified byt not found, quit
173 if (!file) {
174 cout << "StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
175 return;
176 }
177
178 // -------------------------------------------------------
179 // Now get the data and workspace
180
181 // get the workspace out of the file
182 RooWorkspace *w = (RooWorkspace *)file->Get(workspaceName);
183 if (!w) {
184 cout << "workspace not found" << endl;
185 return;
186 }
187
188 // get the modelConfig out of the file
189 ModelConfig *mc = (ModelConfig *)w->obj(modelConfigName);
190
191 // get the modelConfig out of the file
192 RooAbsData *data = w->data(dataName);
193
194 // make sure ingredients are found
195 if (!data || !mc) {
196 w->Print();
197 cout << "data or ModelConfig was not found" << endl;
198 return;
199 }
200
201 // -------------------------------------------------------
202 // Now get the POI for convenience
203 // you may want to adjust the range of your POI
204
205 RooRealVar *firstPOI = (RooRealVar *)mc->GetParametersOfInterest()->first();
206 /* firstPOI->setMin(0);*/
207 /* firstPOI->setMax(10);*/
208
209 // --------------------------------------------
210 // Create and use the FeldmanCousins tool
211 // to find and plot the 95% confidence interval
212 // on the parameter of interest as specified
213 // in the model config
214 // REMEMBER, we will change the test statistic
215 // so this is NOT a Feldman-Cousins interval
216 FeldmanCousins fc(*data, *mc);
217 fc.SetConfidenceLevel(confidenceLevel);
218 fc.AdditionalNToysFactor(
219 0.5); // degrade/improve sampling that defines confidence belt: in this case makes the example faster
220 /* fc.UseAdaptiveSampling(true); // speed it up a bit, don't use for expected limits*/
221 fc.SetNBins(nPointsToScan); // set how many points per parameter of interest to scan
222 fc.CreateConfBelt(true); // save the information in the belt for plotting
223
224 // -------------------------------------------------------
225 // Feldman-Cousins is a unified limit by definition
226 // but the tool takes care of a few things for us like which values
227 // of the nuisance parameters should be used to generate toys.
228 // so let's just change the test statistic and realize this is
229 // no longer "Feldman-Cousins" but is a fully frequentist Neyman-Construction.
230 /* ProfileLikelihoodTestStatModified onesided(*mc->GetPdf());*/
231 /* fc.GetTestStatSampler()->SetTestStatistic(&onesided);*/
232 /* ((ToyMCSampler*) fc.GetTestStatSampler())->SetGenerateBinned(true); */
233 ToyMCSampler *toymcsampler = (ToyMCSampler *)fc.GetTestStatSampler();
234 ProfileLikelihoodTestStat *testStat = dynamic_cast<ProfileLikelihoodTestStat *>(toymcsampler->GetTestStatistic());
235 testStat->SetOneSided(true);
236
237 // Since this tool needs to throw toy MC the PDF needs to be
238 // extended or the tool needs to know how many entries in a dataset
239 // per pseudo experiment.
240 // In the 'number counting form' where the entries in the dataset
241 // are counts, and not values of discriminating variables, the
242 // datasets typically only have one entry and the PDF is not
243 // extended.
244 if (!mc->GetPdf()->canBeExtended()) {
245 if (data->numEntries() == 1)
246 fc.FluctuateNumDataEntries(false);
247 else
248 cout << "Not sure what to do about this model" << endl;
249 }
250
251 // We can use PROOF to speed things along in parallel
252 // However, the test statistic has to be installed on the workers
253 // so either turn off PROOF or include the modified test statistic
254 // in your `$ROOTSYS/roofit/roostats/inc` directory,
255 // add the additional line to the LinkDef.h file,
256 // and recompile root.
257 if (useProof) {
258 ProofConfig pc(*w, nworkers, "", false);
259 toymcsampler->SetProofConfig(&pc); // enable proof
260 }
261
262 if (mc->GetGlobalObservables()) {
263 cout << "will use global observables for unconditional ensemble" << endl;
265 toymcsampler->SetGlobalObservables(*mc->GetGlobalObservables());
266 }
267
268 // Now get the interval
269 PointSetInterval *interval = fc.GetInterval();
270 ConfidenceBelt *belt = fc.GetConfidenceBelt();
271
272 // print out the interval on the first Parameter of Interest
273 cout << "\n95% interval on " << firstPOI->GetName() << " is : [" << interval->LowerLimit(*firstPOI) << ", "
274 << interval->UpperLimit(*firstPOI) << "] " << endl;
275
276 // get observed UL and value of test statistic evaluated there
277 RooArgSet tmpPOI(*firstPOI);
278 double observedUL = interval->UpperLimit(*firstPOI);
279 firstPOI->setVal(observedUL);
280 double obsTSatObsUL = fc.GetTestStatSampler()->EvaluateTestStatistic(*data, tmpPOI);
281
282 // Ask the calculator which points were scanned
283 RooDataSet *parameterScan = (RooDataSet *)fc.GetPointsToScan();
284 RooArgSet *tmpPoint;
285
286 // make a histogram of parameter vs. threshold
287 TH1F *histOfThresholds =
288 new TH1F("histOfThresholds", "", parameterScan->numEntries(), firstPOI->getMin(), firstPOI->getMax());
289 histOfThresholds->GetXaxis()->SetTitle(firstPOI->GetName());
290 histOfThresholds->GetYaxis()->SetTitle("Threshold");
291
292 // loop through the points that were tested and ask confidence belt
293 // what the upper/lower thresholds were.
294 // For FeldmanCousins, the lower cut off is always 0
295 for (Int_t i = 0; i < parameterScan->numEntries(); ++i) {
296 tmpPoint = (RooArgSet *)parameterScan->get(i)->clone("temp");
297 // cout <<"get threshold"<<endl;
298 double arMax = belt->GetAcceptanceRegionMax(*tmpPoint);
299 double poiVal = tmpPoint->getRealValue(firstPOI->GetName());
300 histOfThresholds->Fill(poiVal, arMax);
301 }
302 TCanvas *c1 = new TCanvas();
303 c1->Divide(2);
304 c1->cd(1);
305 histOfThresholds->SetMinimum(0);
306 histOfThresholds->Draw();
307 c1->cd(2);
308
309 // -------------------------------------------------------
310 // Now we generate the expected bands and power-constraint
311
312 // First: find parameter point for mu=0, with conditional MLEs for nuisance parameters
313 RooAbsReal *nll = mc->GetPdf()->createNLL(*data);
314 RooAbsReal *profile = nll->createProfile(*mc->GetParametersOfInterest());
315 firstPOI->setVal(0.);
316 profile->getVal(); // this will do fit and set nuisance parameters to profiled values
317 RooArgSet *poiAndNuisance = new RooArgSet();
318 if (mc->GetNuisanceParameters())
319 poiAndNuisance->add(*mc->GetNuisanceParameters());
320 poiAndNuisance->add(*mc->GetParametersOfInterest());
321 w->saveSnapshot("paramsToGenerateData", *poiAndNuisance);
322 RooArgSet *paramsToGenerateData = (RooArgSet *)poiAndNuisance->snapshot();
323 cout << "\nWill use these parameter points to generate pseudo data for bkg only" << endl;
324 paramsToGenerateData->Print("v");
325
326 RooArgSet unconditionalObs;
327 unconditionalObs.add(*mc->GetObservables());
328 unconditionalObs.add(*mc->GetGlobalObservables()); // comment this out for the original conditional ensemble
329
330 double CLb = 0;
331 double CLbinclusive = 0;
332
333 // Now we generate background only and find distribution of upper limits
334 TH1F *histOfUL = new TH1F("histOfUL", "", 100, 0, firstPOI->getMax());
335 histOfUL->GetXaxis()->SetTitle("Upper Limit (background only)");
336 histOfUL->GetYaxis()->SetTitle("Entries");
337 for (int imc = 0; imc < nToyMC; ++imc) {
338
339 // set parameters back to values for generating pseudo data
340 // cout << "\n get current nuis, set vals, print again" << endl;
341 w->loadSnapshot("paramsToGenerateData");
342 // poiAndNuisance->Print("v");
343
344 RooDataSet *toyData = 0;
345 // now generate a toy dataset
346 if (!mc->GetPdf()->canBeExtended()) {
347 if (data->numEntries() == 1)
348 toyData = mc->GetPdf()->generate(*mc->GetObservables(), 1);
349 else
350 cout << "Not sure what to do about this model" << endl;
351 } else {
352 // cout << "generating extended dataset"<<endl;
353 toyData = mc->GetPdf()->generate(*mc->GetObservables(), Extended());
354 }
355
356 // generate global observables
357 // need to be careful for simpdf
358 // RooDataSet* globalData = mc->GetPdf()->generate(*mc->GetGlobalObservables(),1);
359
360 RooSimultaneous *simPdf = dynamic_cast<RooSimultaneous *>(mc->GetPdf());
361 if (!simPdf) {
362 RooDataSet *one = mc->GetPdf()->generate(*mc->GetGlobalObservables(), 1);
363 const RooArgSet *values = one->get();
364 RooArgSet *allVars = mc->GetPdf()->getVariables();
365 *allVars = *values;
366 delete allVars;
367 delete values;
368 delete one;
369 } else {
370
371 // try fix for sim pdf
372 TIterator *iter = simPdf->indexCat().typeIterator();
373 RooCatType *tt = NULL;
374 while ((tt = (RooCatType *)iter->Next())) {
375
376 // Get pdf associated with state from simpdf
377 RooAbsPdf *pdftmp = simPdf->getPdf(tt->GetName());
378
379 // Generate only global variables defined by the pdf associated with this state
380 RooArgSet *globtmp = pdftmp->getObservables(*mc->GetGlobalObservables());
381 RooDataSet *tmp = pdftmp->generate(*globtmp, 1);
382
383 // Transfer values to output placeholder
384 *globtmp = *tmp->get(0);
385
386 // Cleanup
387 delete globtmp;
388 delete tmp;
389 }
390 }
391
392 // globalData->Print("v");
393 // unconditionalObs = *globalData->get();
394 // mc->GetGlobalObservables()->Print("v");
395 // delete globalData;
396 // cout << "toy data = " << endl;
397 // toyData->get()->Print("v");
398
399 // get test stat at observed UL in observed data
400 firstPOI->setVal(observedUL);
401 double toyTSatObsUL = fc.GetTestStatSampler()->EvaluateTestStatistic(*toyData, tmpPOI);
402 // toyData->get()->Print("v");
403 // cout <<"obsTSatObsUL " <<obsTSatObsUL << "toyTS " << toyTSatObsUL << endl;
404 if (obsTSatObsUL < toyTSatObsUL) // not sure about <= part yet
405 CLb += (1.) / nToyMC;
406 if (obsTSatObsUL <= toyTSatObsUL) // not sure about <= part yet
407 CLbinclusive += (1.) / nToyMC;
408
409 // loop over points in belt to find upper limit for this toy data
410 double thisUL = 0;
411 for (Int_t i = 0; i < parameterScan->numEntries(); ++i) {
412 tmpPoint = (RooArgSet *)parameterScan->get(i)->clone("temp");
413 double arMax = belt->GetAcceptanceRegionMax(*tmpPoint);
414 firstPOI->setVal(tmpPoint->getRealValue(firstPOI->GetName()));
415 // double thisTS = profile->getVal();
416 double thisTS = fc.GetTestStatSampler()->EvaluateTestStatistic(*toyData, tmpPOI);
417
418 // cout << "poi = " << firstPOI->getVal()
419 // << " max is " << arMax << " this profile = " << thisTS << endl;
420 // cout << "thisTS = " << thisTS<<endl;
421 if (thisTS <= arMax) {
422 thisUL = firstPOI->getVal();
423 } else {
424 break;
425 }
426 }
427
428 /*
429 // loop over points in belt to find upper limit for this toy data
430 double thisUL = 0;
431 for(Int_t i=0; i<histOfThresholds->GetNbinsX(); ++i){
432 tmpPoint = (RooArgSet*) parameterScan->get(i)->clone("temp");
433 cout <<"---------------- "<<i<<endl;
434 tmpPoint->Print("v");
435 cout << "from hist " << histOfThresholds->GetBinCenter(i+1) <<endl;
436 double arMax = histOfThresholds->GetBinContent(i+1);
437 // cout << " threhold from Hist = aMax " << arMax<<endl;
438 // double arMax2 = belt->GetAcceptanceRegionMax(*tmpPoint);
439 // cout << "from scan arMax2 = "<< arMax2 << endl; // not the same due to TH1F not TH1D
440 // cout << "scan - hist" << arMax2-arMax << endl;
441 firstPOI->setVal( histOfThresholds->GetBinCenter(i+1));
442 // double thisTS = profile->getVal();
443 double thisTS = fc.GetTestStatSampler()->EvaluateTestStatistic(*toyData,tmpPOI);
444
445 // cout << "poi = " << firstPOI->getVal()
446 // << " max is " << arMax << " this profile = " << thisTS << endl;
447 // cout << "thisTS = " << thisTS<<endl;
448
449 // NOTE: need to add a small epsilon term for single precision vs. double precision
450 if(thisTS<=arMax + 1e-7){
451 thisUL = firstPOI->getVal();
452 } else{
453 break;
454 }
455 }
456 */
457
458 histOfUL->Fill(thisUL);
459
460 // for few events, data is often the same, and UL is often the same
461 // cout << "thisUL = " << thisUL<<endl;
462
463 delete toyData;
464 }
465 histOfUL->Draw();
466 c1->SaveAs("one-sided_upper_limit_output.pdf");
467
468 // if you want to see a plot of the sampling distribution for a particular scan point:
469 /*
470 SamplingDistPlot sampPlot;
471 int indexInScan = 0;
472 tmpPoint = (RooArgSet*) parameterScan->get(indexInScan)->clone("temp");
473 firstPOI->setVal( tmpPoint->getRealValue(firstPOI->GetName()) );
474 toymcsampler->SetParametersForTestStat(tmpPOI);
475 SamplingDistribution* samp = toymcsampler->GetSamplingDistribution(*tmpPoint);
476 sampPlot.AddSamplingDistribution(samp);
477 sampPlot.Draw();
478 */
479
480 // Now find bands and power constraint
481 Double_t *bins = histOfUL->GetIntegral();
482 TH1F *cumulative = (TH1F *)histOfUL->Clone("cumulative");
483 cumulative->SetContent(bins);
484 double band2sigDown, band1sigDown, bandMedian, band1sigUp, band2sigUp;
485 for (int i = 1; i <= cumulative->GetNbinsX(); ++i) {
486 if (bins[i] < RooStats::SignificanceToPValue(2))
487 band2sigDown = cumulative->GetBinCenter(i);
488 if (bins[i] < RooStats::SignificanceToPValue(1))
489 band1sigDown = cumulative->GetBinCenter(i);
490 if (bins[i] < 0.5)
491 bandMedian = cumulative->GetBinCenter(i);
492 if (bins[i] < RooStats::SignificanceToPValue(-1))
493 band1sigUp = cumulative->GetBinCenter(i);
494 if (bins[i] < RooStats::SignificanceToPValue(-2))
495 band2sigUp = cumulative->GetBinCenter(i);
496 }
497 cout << "-2 sigma band " << band2sigDown << endl;
498 cout << "-1 sigma band " << band1sigDown << " [Power Constraint)]" << endl;
499 cout << "median of band " << bandMedian << endl;
500 cout << "+1 sigma band " << band1sigUp << endl;
501 cout << "+2 sigma band " << band2sigUp << endl;
502
503 // print out the interval on the first Parameter of Interest
504 cout << "\nobserved 95% upper-limit " << interval->UpperLimit(*firstPOI) << endl;
505 cout << "CLb strict [P(toy>obs|0)] for observed 95% upper-limit " << CLb << endl;
506 cout << "CLb inclusive [P(toy>=obs|0)] for observed 95% upper-limit " << CLbinclusive << endl;
507
508 delete profile;
509 delete nll;
510}
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
#define gROOT
Definition TROOT.h:406
R__EXTERN TSystem * gSystem
Definition TSystem.h:559
static struct mg_connection * fc(struct mg_context *ctx)
Definition civetweb.c:3728
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Given a set of possible observables, return the observables that this PDF depends on.
Definition RooAbsArg.h:313
RooArgSet * getVariables(Bool_t stripDisconnected=kTRUE) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
TIterator * typeIterator() const
Double_t getRealValue(const char *name, Double_t defVal=0, Bool_t verbose=kFALSE) const
Get value of a RooAbsReal stored in set with given name.
RooAbsArg * first() const
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:49
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
virtual RooAbsReal * createNLL(RooAbsData &data, const RooLinkedList &cmdList)
Construct representation of -log(L) of PDFwith given dataset.
Bool_t canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition RooAbsPdf.h:238
RooDataSet * generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none())
See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,...
Definition RooAbsPdf.h:58
virtual Double_t getMax(const char *name=0) const
Get maximum of currently defined range.
virtual Double_t getMin(const char *name=0) const
Get miniminum of currently defined range.
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:61
virtual RooAbsReal * createProfile(const RooArgSet &paramsOfInterest)
Create a RooProfileLL object that eliminates all nuisance parameters in the present function.
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
TObject * clone(const char *newname) const override
Definition RooArgSet.h:83
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:118
Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Add element to non-owning set.
RooCatType is an auxilary class for RooAbsCategory and defines a a single category state.
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:33
virtual const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
virtual void setVal(Double_t value)
Set value of variable to 'value'.
RooSimultaneous facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
const RooAbsCategoryLValue & indexCat() const
RooAbsPdf * getPdf(const char *catName) const
Return the p.d.f associated with the given index category name.
ConfidenceBelt is a concrete implementation of the ConfInterval interface.
Double_t GetAcceptanceRegionMax(RooArgSet &, Double_t cl=-1., Double_t 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:30
const RooArgSet * GetGlobalObservables() const
get RooArgSet for global observables (return NULL if not existing)
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)
RooAbsPdf * GetPdf() const
get model PDF (return NULL if pdf has not been specified or does not exist)
PointSetInterval is a concrete implementation of the ConfInterval interface.
Double_t UpperLimit(RooRealVar &param)
return upper limit on a given parameter
Double_t 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:46
ToyMCSampler is an implementation of the TestStatSampler interface.
void SetProofConfig(ProofConfig *pc=NULL)
virtual TestStatistic * GetTestStatistic(unsigned int i) const
virtual void SetGlobalObservables(const RooArgSet &o)
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.
Bool_t saveSnapshot(const char *name, const char *paramNames)
Save snapshot of values and attributes (including "Constant") of given parameters.
Bool_t loadSnapshot(const char *name)
Load the values and attributes of the parameters in the snapshot saved with the given name.
TObject * obj(const char *name) const
Return any type of object (RooAbsArg, RooAbsData or generic object) with given name)
The Canvas class.
Definition TCanvas.h:23
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition TFile.h:54
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:3997
1-D histogram with a float per channel (see TH1 documentation)}
Definition TH1.h:575
virtual Double_t GetBinCenter(Int_t bin) const
Return bin center for 1D histogram.
Definition TH1.cxx:8981
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition TH1.h:320
TObject * Clone(const char *newname=0) const
Make a complete copy of the underlying object.
Definition TH1.cxx:2740
virtual Int_t GetNbinsX() const
Definition TH1.h:296
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3350
TAxis * GetYaxis()
Definition TH1.h:321
virtual void SetContent(const Double_t *content)
Replace bin contents by the contents of array content.
Definition TH1.cxx:8246
virtual void SetMinimum(Double_t minimum=-1111)
Definition TH1.h:399
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition TH1.cxx:3073
virtual Double_t * GetIntegral()
Return a pointer to the array of bins integral.
Definition TH1.cxx:2578
Iterator abstract base class.
Definition TIterator.h:30
virtual TObject * Next()=0
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
virtual const char * GetName() const
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:1294
RooCmdArg Extended(Bool_t flag=kTRUE)
return c1
Definition legend1.C:41
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Namespace for the RooStats classes.
Definition Asimov.h:19
Double_t SignificanceToPValue(Double_t Z)
returns p-value corresponding to a 1-sided significance
Definition file.py:1
auto * tt
Definition textangle.C:16