Logo ROOT   6.14/05
Reference Guide
StandardHypoTestInvDemo.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roostats
3 /// \notebook
4 /// Standard tutorial macro for performing an inverted hypothesis test for computing an interval
5 ///
6 /// This macro will perform a scan of the p-values for computing the interval or limit
7 ///
8 /// Usage:
9 ///
10 /// ~~~{.cpp}
11 /// root>.L StandardHypoTestInvDemo.C
12 /// root> StandardHypoTestInvDemo("fileName","workspace name","S+B modelconfig name","B model name","data set name",calculator type, test statistic type, use CLS,
13 /// number of points, xmin, xmax, number of toys, use number counting)
14 ///
15 /// type = 0 Freq calculator
16 /// type = 1 Hybrid calculator
17 /// type = 2 Asymptotic calculator
18 /// type = 3 Asymptotic calculator using nominal Asimov data sets (not using fitted parameter values but nominal ones)
19 ///
20 /// testStatType = 0 LEP
21 /// = 1 Tevatron
22 /// = 2 Profile Likelihood two sided
23 /// = 3 Profile Likelihood one sided (i.e. = 0 if mu < mu_hat)
24 /// = 4 Profile Likelihood signed ( pll = -pll if mu < mu_hat)
25 /// = 5 Max Likelihood Estimate as test statistic
26 /// = 6 Number of observed event as test statistic
27 /// ~~~
28 ///
29 /// \macro_image
30 /// \macro_output
31 /// \macro_code
32 ///
33 /// \author Lorenzo Moneta
34 
35 #include "TFile.h"
36 #include "RooWorkspace.h"
37 #include "RooAbsPdf.h"
38 #include "RooRealVar.h"
39 #include "RooDataSet.h"
40 #include "RooStats/ModelConfig.h"
41 #include "RooRandom.h"
42 #include "TGraphErrors.h"
43 #include "TGraphAsymmErrors.h"
44 #include "TCanvas.h"
45 #include "TLine.h"
46 #include "TROOT.h"
47 #include "TSystem.h"
48 
52 #include "RooStats/ToyMCSampler.h"
53 #include "RooStats/HypoTestPlot.h"
54 
61 
65 
66 #include <cassert>
67 
68 using namespace RooFit;
69 using namespace RooStats;
70 using namespace std;
71 
72 // structure defining the options
73 struct HypoTestInvOptions {
74 
75 bool plotHypoTestResult = true; // plot test statistic result at each point
76 bool writeResult = true; // write HypoTestInverterResult in a file
77 TString resultFileName; // file with results (by default is built automatically using the workspace input file name)
78 bool optimize = true; // optimize evaluation of test statistic
79 bool useVectorStore = true; // convert data to use new roofit data store
80 bool generateBinned = false; // generate binned data sets
81 bool noSystematics = false; // force all systematics to be off (i.e. set all nuisance parameters as constat
82  // to their nominal values)
83 double nToysRatio = 2; // ratio Ntoys S+b/ntoysB
84 double maxPOI = -1; // max value used of POI (in case of auto scan)
85 bool useProof = false; // use Proof Lite when using toys (for freq or hybrid)
86 int nworkers = 0; // number of worker for ProofLite (default use all available cores)
87 bool enableDetailedOutput = false; // enable detailed output with all fit information for each toys (output will be written in result file)
88 bool rebuild = false; // re-do extra toys for computing expected limits and rebuild test stat
89  // distributions (N.B this requires much more CPU (factor is equivalent to nToyToRebuild)
90 int nToyToRebuild = 100; // number of toys used to rebuild
91 int rebuildParamValues=0; // = 0 do a profile of all the parameters on the B (alt snapshot) before performing a rebuild operation (default)
92  // = 1 use initial workspace parameters with B snapshot values
93  // = 2 use all initial workspace parameters with B
94  // Otherwise the rebuild will be performed using
95 int initialFit = -1; // do a first fit to the model (-1 : default, 0 skip fit, 1 do always fit)
96 int randomSeed = -1; // random seed (if = -1: use default value, if = 0 always random )
97  // NOTE: Proof uses automatically a random seed
98 
99 int nAsimovBins = 0; // number of bins in observables used for Asimov data sets (0 is the default and it is given by workspace, typically is 100)
100 
101 bool reuseAltToys = false; // reuse same toys for alternate hypothesis (if set one gets more stable bands)
102 double confLevel = 0.95; // confidence level value
103 
104 
105 
106 std::string minimizerType = ""; // minimizer type (default is what is in ROOT::Math::MinimizerOptions::DefaultMinimizerType()
107 std::string massValue = ""; // extra string to tag output file of result
108 int printLevel = 0; // print level for debugging PL test statistics and calculators
109 
110 bool useNLLOffset = false; // use NLL offset when fitting (this increase stability of fits)
111 
112 };
113 
114 
115 HypoTestInvOptions optHTInv;
116 
117 // internal class to run the inverter and more
118 
119 namespace RooStats {
120 
121  class HypoTestInvTool{
122 
123  public:
124  HypoTestInvTool();
125  ~HypoTestInvTool(){};
126 
128  RunInverter(RooWorkspace * w,
129  const char * modelSBName, const char * modelBName,
130  const char * dataName,
131  int type, int testStatType,
132  bool useCLs,
133  int npoints, double poimin, double poimax, int ntoys,
134  bool useNumberCounting = false,
135  const char * nuisPriorName = 0);
136 
137 
138 
139  void
140  AnalyzeResult( HypoTestInverterResult * r,
141  int calculatorType,
142  int testStatType,
143  bool useCLs,
144  int npoints,
145  const char * fileNameBase = 0 );
146 
147  void SetParameter(const char * name, const char * value);
148  void SetParameter(const char * name, bool value);
149  void SetParameter(const char * name, int value);
150  void SetParameter(const char * name, double value);
151 
152  private:
153 
154  bool mPlotHypoTestResult;
155  bool mWriteResult;
156  bool mOptimize;
157  bool mUseVectorStore;
158  bool mGenerateBinned;
159  bool mUseProof;
160  bool mRebuild;
161  bool mReuseAltToys;
162  bool mEnableDetOutput;
163  int mNWorkers;
164  int mNToyToRebuild;
165  int mRebuildParamValues;
166  int mPrintLevel;
167  int mInitialFit;
168  int mRandomSeed;
169  double mNToysRatio;
170  double mMaxPoi;
171  int mAsimovBins;
172  std::string mMassValue;
173  std::string mMinimizerType; // minimizer type (default is what is in ROOT::Math::MinimizerOptions::DefaultMinimizerType()
174  TString mResultFileName;
175  };
176 
177 } // end namespace RooStats
178 
179 RooStats::HypoTestInvTool::HypoTestInvTool() : mPlotHypoTestResult(true),
180  mWriteResult(false),
181  mOptimize(true),
182  mUseVectorStore(true),
183  mGenerateBinned(false),
184  mUseProof(false),
185  mEnableDetOutput(false),
186  mRebuild(false),
187  mReuseAltToys(false),
188  mNWorkers(4),
189  mNToyToRebuild(100),
190  mRebuildParamValues(0),
191  mPrintLevel(0),
192  mInitialFit(-1),
193  mRandomSeed(-1),
194  mNToysRatio(2),
195  mMaxPoi(-1),
196  mAsimovBins(0),
197  mMassValue(""),
198  mMinimizerType(""),
199  mResultFileName() {
200 }
201 
202 
203 
204 void
205 RooStats::HypoTestInvTool::SetParameter(const char * name, bool value){
206  //
207  // set boolean parameters
208  //
209 
210  std::string s_name(name);
211 
212  if (s_name.find("PlotHypoTestResult") != std::string::npos) mPlotHypoTestResult = value;
213  if (s_name.find("WriteResult") != std::string::npos) mWriteResult = value;
214  if (s_name.find("Optimize") != std::string::npos) mOptimize = value;
215  if (s_name.find("UseVectorStore") != std::string::npos) mUseVectorStore = value;
216  if (s_name.find("GenerateBinned") != std::string::npos) mGenerateBinned = value;
217  if (s_name.find("UseProof") != std::string::npos) mUseProof = value;
218  if (s_name.find("EnableDetailedOutput") != std::string::npos) mEnableDetOutput = value;
219  if (s_name.find("Rebuild") != std::string::npos) mRebuild = value;
220  if (s_name.find("ReuseAltToys") != std::string::npos) mReuseAltToys = value;
221 
222  return;
223 }
224 
225 
226 
227 void
228 RooStats::HypoTestInvTool::SetParameter(const char * name, int value){
229  //
230  // set integer parameters
231  //
232 
233  std::string s_name(name);
234 
235  if (s_name.find("NWorkers") != std::string::npos) mNWorkers = value;
236  if (s_name.find("NToyToRebuild") != std::string::npos) mNToyToRebuild = value;
237  if (s_name.find("RebuildParamValues") != std::string::npos) mRebuildParamValues = value;
238  if (s_name.find("PrintLevel") != std::string::npos) mPrintLevel = value;
239  if (s_name.find("InitialFit") != std::string::npos) mInitialFit = value;
240  if (s_name.find("RandomSeed") != std::string::npos) mRandomSeed = value;
241  if (s_name.find("AsimovBins") != std::string::npos) mAsimovBins = value;
242 
243  return;
244 }
245 
246 
247 
248 void
249 RooStats::HypoTestInvTool::SetParameter(const char * name, double value){
250  //
251  // set double precision parameters
252  //
253 
254  std::string s_name(name);
255 
256  if (s_name.find("NToysRatio") != std::string::npos) mNToysRatio = value;
257  if (s_name.find("MaxPOI") != std::string::npos) mMaxPoi = value;
258 
259  return;
260 }
261 
262 
263 
264 void
265 RooStats::HypoTestInvTool::SetParameter(const char * name, const char * value){
266  //
267  // set string parameters
268  //
269 
270  std::string s_name(name);
271 
272  if (s_name.find("MassValue") != std::string::npos) mMassValue.assign(value);
273  if (s_name.find("MinimizerType") != std::string::npos) mMinimizerType.assign(value);
274  if (s_name.find("ResultFileName") != std::string::npos) mResultFileName = value;
275 
276  return;
277 }
278 
279 
280 
281 void
282 StandardHypoTestInvDemo(const char * infile = 0,
283  const char * wsName = "combined",
284  const char * modelSBName = "ModelConfig",
285  const char * modelBName = "",
286  const char * dataName = "obsData",
287  int calculatorType = 0,
288  int testStatType = 0,
289  bool useCLs = true ,
290  int npoints = 6,
291  double poimin = 0,
292  double poimax = 5,
293  int ntoys=1000,
294  bool useNumberCounting = false,
295  const char * nuisPriorName = 0){
296 /*
297 
298  Other Parameter to pass in tutorial
299  apart from standard for filename, ws, modelconfig and data
300 
301  type = 0 Freq calculator
302  type = 1 Hybrid calculator
303  type = 2 Asymptotic calculator
304  type = 3 Asymptotic calculator using nominal Asimov data sets (not using fitted parameter values but nominal ones)
305 
306  testStatType = 0 LEP
307  = 1 Tevatron
308  = 2 Profile Likelihood
309  = 3 Profile Likelihood one sided (i.e. = 0 if mu < mu_hat)
310  = 4 Profiel Likelihood signed ( pll = -pll if mu < mu_hat)
311  = 5 Max Likelihood Estimate as test statistic
312  = 6 Number of observed event as test statistic
313 
314  useCLs scan for CLs (otherwise for CLs+b)
315 
316  npoints: number of points to scan , for autoscan set npoints = -1
317 
318  poimin,poimax: min/max value to scan in case of fixed scans
319  (if min > max, try to find automatically)
320 
321  ntoys: number of toys to use
322 
323  useNumberCounting: set to true when using number counting events
324 
325  nuisPriorName: name of prior for the nuisance. This is often expressed as constraint term in the global model
326  It is needed only when using the HybridCalculator (type=1)
327  If not given by default the prior pdf from ModelConfig is used.
328 
329  extra options are available as global parameters of the macro. They major ones are:
330 
331  plotHypoTestResult plot result of tests at each point (TS distributions) (default is true)
332  useProof use Proof (default is true)
333  writeResult write result of scan (default is true)
334  rebuild rebuild scan for expected limits (require extra toys) (default is false)
335  generateBinned generate binned data sets for toys (default is false) - be careful not to activate with
336  a too large (>=3) number of observables
337  nToyRatio ratio of S+B/B toys (default is 2)
338 
339 
340 */
341 
342 
343 
344  TString filename(infile);
345  if (filename.IsNull()) {
346  filename = "results/example_combined_GaussExample_model.root";
347  bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code
348  // if file does not exists generate with histfactory
349  if (!fileExist) {
350 #ifdef _WIN32
351  cout << "HistFactory file cannot be generated on Windows - exit" << endl;
352  return;
353 #endif
354  // Normally this would be run on the command line
355  cout <<"will run standard hist2workspace example"<<endl;
356  gROOT->ProcessLine(".! prepareHistFactory .");
357  gROOT->ProcessLine(".! hist2workspace config/example.xml");
358  cout <<"\n\n---------------------"<<endl;
359  cout <<"Done creating example input"<<endl;
360  cout <<"---------------------\n\n"<<endl;
361  }
362 
363  }
364  else
365  filename = infile;
366 
367  // Try to open the file
368  TFile *file = TFile::Open(filename);
369 
370  // if input file was specified byt not found, quit
371  if(!file ){
372  cout <<"StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
373  return;
374  }
375 
376 
377 
378  HypoTestInvTool calc;
379 
380  // set parameters
381  calc.SetParameter("PlotHypoTestResult", optHTInv.plotHypoTestResult);
382  calc.SetParameter("WriteResult", optHTInv.writeResult);
383  calc.SetParameter("Optimize", optHTInv.optimize);
384  calc.SetParameter("UseVectorStore", optHTInv.useVectorStore);
385  calc.SetParameter("GenerateBinned", optHTInv.generateBinned);
386  calc.SetParameter("NToysRatio", optHTInv.nToysRatio);
387  calc.SetParameter("MaxPOI", optHTInv.maxPOI);
388  calc.SetParameter("UseProof", optHTInv.useProof);
389  calc.SetParameter("EnableDetailedOutput", optHTInv.enableDetailedOutput);
390  calc.SetParameter("NWorkers", optHTInv.nworkers);
391  calc.SetParameter("Rebuild", optHTInv.rebuild);
392  calc.SetParameter("ReuseAltToys", optHTInv.reuseAltToys);
393  calc.SetParameter("NToyToRebuild", optHTInv.nToyToRebuild);
394  calc.SetParameter("RebuildParamValues", optHTInv.rebuildParamValues);
395  calc.SetParameter("MassValue", optHTInv.massValue.c_str());
396  calc.SetParameter("MinimizerType", optHTInv.minimizerType.c_str());
397  calc.SetParameter("PrintLevel", optHTInv.printLevel);
398  calc.SetParameter("InitialFit", optHTInv.initialFit);
399  calc.SetParameter("ResultFileName", optHTInv.resultFileName);
400  calc.SetParameter("RandomSeed", optHTInv.randomSeed);
401  calc.SetParameter("AsimovBins", optHTInv.nAsimovBins);
402 
403  // enable offset for all roostats
404  if (optHTInv.useNLLOffset) RooStats::UseNLLOffset(true);
405 
406  RooWorkspace * w = dynamic_cast<RooWorkspace*>( file->Get(wsName) );
408  std::cout << w << "\t" << filename << std::endl;
409  if (w != NULL) {
410  r = calc.RunInverter(w, modelSBName, modelBName,
411  dataName, calculatorType, testStatType, useCLs,
412  npoints, poimin, poimax,
413  ntoys, useNumberCounting, nuisPriorName );
414  if (!r) {
415  std::cerr << "Error running the HypoTestInverter - Exit " << std::endl;
416  return;
417  }
418  }
419  else {
420  // case workspace is not present look for the inverter result
421  std::cout << "Reading an HypoTestInverterResult with name " << wsName << " from file " << filename << std::endl;
422  r = dynamic_cast<HypoTestInverterResult*>( file->Get(wsName) ); //
423  if (!r) {
424  std::cerr << "File " << filename << " does not contain a workspace or an HypoTestInverterResult - Exit "
425  << std::endl;
426  file->ls();
427  return;
428  }
429  }
430 
431  calc.AnalyzeResult( r, calculatorType, testStatType, useCLs, npoints, infile );
432 
433  return;
434 }
435 
436 
437 
438 void
439 RooStats::HypoTestInvTool::AnalyzeResult( HypoTestInverterResult * r,
440  int calculatorType,
441  int testStatType,
442  bool useCLs,
443  int npoints,
444  const char * fileNameBase ){
445 
446  // analyze result produced by the inverter, optionally save it in a file
447 
448 
449  double lowerLimit = 0;
450  double llError = 0;
451 #if defined ROOT_SVN_VERSION && ROOT_SVN_VERSION >= 44126
452  if (r->IsTwoSided()) {
453  lowerLimit = r->LowerLimit();
454  llError = r->LowerLimitEstimatedError();
455  }
456 #else
457  lowerLimit = r->LowerLimit();
458  llError = r->LowerLimitEstimatedError();
459 #endif
460 
461  double upperLimit = r->UpperLimit();
462  double ulError = r->UpperLimitEstimatedError();
463 
464  //std::cout << "DEBUG : [ " << lowerLimit << " , " << upperLimit << " ] " << std::endl;
465 
466  if (lowerLimit < upperLimit*(1.- 1.E-4) && lowerLimit != 0)
467  std::cout << "The computed lower limit is: " << lowerLimit << " +/- " << llError << std::endl;
468  std::cout << "The computed upper limit is: " << upperLimit << " +/- " << ulError << std::endl;
469 
470 
471  // compute expected limit
472  std::cout << "Expected upper limits, using the B (alternate) model : " << std::endl;
473  std::cout << " expected limit (median) " << r->GetExpectedUpperLimit(0) << std::endl;
474  std::cout << " expected limit (-1 sig) " << r->GetExpectedUpperLimit(-1) << std::endl;
475  std::cout << " expected limit (+1 sig) " << r->GetExpectedUpperLimit(1) << std::endl;
476  std::cout << " expected limit (-2 sig) " << r->GetExpectedUpperLimit(-2) << std::endl;
477  std::cout << " expected limit (+2 sig) " << r->GetExpectedUpperLimit(2) << std::endl;
478 
479 
480  // detailed output
481  if (mEnableDetOutput) {
482  mWriteResult=true;
483  Info("StandardHypoTestInvDemo","detailed output will be written in output result file");
484  }
485 
486 
487  // write result in a file
488  if (r != NULL && mWriteResult) {
489 
490  // write to a file the results
491  const char * calcType = (calculatorType == 0) ? "Freq" : (calculatorType == 1) ? "Hybr" : "Asym";
492  const char * limitType = (useCLs) ? "CLs" : "Cls+b";
493  const char * scanType = (npoints < 0) ? "auto" : "grid";
494  if (mResultFileName.IsNull()) {
495  mResultFileName = TString::Format("%s_%s_%s_ts%d_",calcType,limitType,scanType,testStatType);
496  //strip the / from the filename
497  if (mMassValue.size()>0) {
498  mResultFileName += mMassValue.c_str();
499  mResultFileName += "_";
500  }
501 
502  TString name = fileNameBase;
503  name.Replace(0, name.Last('/')+1, "");
504  mResultFileName += name;
505  }
506 
507  // get (if existing) rebuilt UL distribution
508  TString uldistFile = "RULDist.root";
509  TObject * ulDist = 0;
510  bool existULDist = !gSystem->AccessPathName(uldistFile);
511  if (existULDist) {
512  TFile * fileULDist = TFile::Open(uldistFile);
513  if (fileULDist) ulDist= fileULDist->Get("RULDist");
514  }
515 
516 
517  TFile * fileOut = new TFile(mResultFileName,"RECREATE");
518  r->Write();
519  if (ulDist) ulDist->Write();
520  Info("StandardHypoTestInvDemo","HypoTestInverterResult has been written in the file %s",mResultFileName.Data());
521 
522  fileOut->Close();
523  }
524 
525 
526  // plot the result ( p values vs scan points)
527  std::string typeName = "";
528  if (calculatorType == 0 )
529  typeName = "Frequentist";
530  if (calculatorType == 1 )
531  typeName = "Hybrid";
532  else if (calculatorType == 2 || calculatorType == 3) {
533  typeName = "Asymptotic";
534  mPlotHypoTestResult = false;
535  }
536 
537  const char * resultName = r->GetName();
538  TString plotTitle = TString::Format("%s CL Scan for workspace %s",typeName.c_str(),resultName);
539  HypoTestInverterPlot *plot = new HypoTestInverterPlot("HTI_Result_Plot",plotTitle,r);
540 
541  // plot in a new canvas with style
542  TString c1Name = TString::Format("%s_Scan",typeName.c_str());
543  TCanvas * c1 = new TCanvas(c1Name);
544  c1->SetLogy(false);
545 
546  plot->Draw("CLb 2CL"); // plot all and Clb
547 
548  // if (useCLs)
549  // plot->Draw("CLb 2CL"); // plot all and Clb
550  // else
551  // plot->Draw(""); // plot all and Clb
552 
553  const int nEntries = r->ArraySize();
554 
555  // plot test statistics distributions for the two hypothesis
556  if (mPlotHypoTestResult) {
557  TCanvas * c2 = new TCanvas("c2");
558  if (nEntries > 1) {
559  int ny = TMath::CeilNint(TMath::Sqrt(nEntries));
560  int nx = TMath::CeilNint(double(nEntries)/ny);
561  c2->Divide( nx,ny);
562  }
563  for (int i=0; i<nEntries; i++) {
564  if (nEntries > 1) c2->cd(i+1);
565  SamplingDistPlot * pl = plot->MakeTestStatPlot(i);
566  pl->SetLogYaxis(true);
567  pl->Draw();
568  }
569  }
570  gPad = c1;
571 
572 }
573 
574 
575 
576 // internal routine to run the inverter
578 RooStats::HypoTestInvTool::RunInverter(RooWorkspace * w,
579  const char * modelSBName, const char * modelBName,
580  const char * dataName, int type, int testStatType,
581  bool useCLs, int npoints, double poimin, double poimax,
582  int ntoys,
583  bool useNumberCounting,
584  const char * nuisPriorName ){
585 
586  std::cout << "Running HypoTestInverter on the workspace " << w->GetName() << std::endl;
587 
588  w->Print();
589 
590 
591  RooAbsData * data = w->data(dataName);
592  if (!data) {
593  Error("StandardHypoTestDemo","Not existing data %s",dataName);
594  return 0;
595  }
596  else
597  std::cout << "Using data set " << dataName << std::endl;
598 
599  if (mUseVectorStore) {
601  data->convertToVectorStore() ;
602  }
603 
604 
605  // get models from WS
606  // get the modelConfig out of the file
607  ModelConfig* bModel = (ModelConfig*) w->obj(modelBName);
608  ModelConfig* sbModel = (ModelConfig*) w->obj(modelSBName);
609 
610  if (!sbModel) {
611  Error("StandardHypoTestDemo","Not existing ModelConfig %s",modelSBName);
612  return 0;
613  }
614  // check the model
615  if (!sbModel->GetPdf()) {
616  Error("StandardHypoTestDemo","Model %s has no pdf ",modelSBName);
617  return 0;
618  }
619  if (!sbModel->GetParametersOfInterest()) {
620  Error("StandardHypoTestDemo","Model %s has no poi ",modelSBName);
621  return 0;
622  }
623  if (!sbModel->GetObservables()) {
624  Error("StandardHypoTestInvDemo","Model %s has no observables ",modelSBName);
625  return 0;
626  }
627  if (!sbModel->GetSnapshot() ) {
628  Info("StandardHypoTestInvDemo","Model %s has no snapshot - make one using model poi",modelSBName);
629  sbModel->SetSnapshot( *sbModel->GetParametersOfInterest() );
630  }
631 
632  // case of no systematics
633  // remove nuisance parameters from model
634  if (optHTInv.noSystematics) {
635  const RooArgSet * nuisPar = sbModel->GetNuisanceParameters();
636  if (nuisPar && nuisPar->getSize() > 0) {
637  std::cout << "StandardHypoTestInvDemo" << " - Switch off all systematics by setting them constant to their initial values" << std::endl;
638  RooStats::SetAllConstant(*nuisPar);
639  }
640  if (bModel) {
641  const RooArgSet * bnuisPar = bModel->GetNuisanceParameters();
642  if (bnuisPar)
643  RooStats::SetAllConstant(*bnuisPar);
644  }
645  }
646 
647  if (!bModel || bModel == sbModel) {
648  Info("StandardHypoTestInvDemo","The background model %s does not exist",modelBName);
649  Info("StandardHypoTestInvDemo","Copy it from ModelConfig %s and set POI to zero",modelSBName);
650  bModel = (ModelConfig*) sbModel->Clone();
651  bModel->SetName(TString(modelSBName)+TString("_with_poi_0"));
652  RooRealVar * var = dynamic_cast<RooRealVar*>(bModel->GetParametersOfInterest()->first());
653  if (!var) return 0;
654  double oldval = var->getVal();
655  var->setVal(0);
656  bModel->SetSnapshot( RooArgSet(*var) );
657  var->setVal(oldval);
658  }
659  else {
660  if (!bModel->GetSnapshot() ) {
661  Info("StandardHypoTestInvDemo","Model %s has no snapshot - make one using model poi and 0 values ",modelBName);
662  RooRealVar * var = dynamic_cast<RooRealVar*>(bModel->GetParametersOfInterest()->first());
663  if (var) {
664  double oldval = var->getVal();
665  var->setVal(0);
666  bModel->SetSnapshot( RooArgSet(*var) );
667  var->setVal(oldval);
668  }
669  else {
670  Error("StandardHypoTestInvDemo","Model %s has no valid poi",modelBName);
671  return 0;
672  }
673  }
674  }
675 
676  // check model has global observables when there are nuisance pdf
677  // for the hybrid case the globals are not needed
678  if (type != 1 ) {
679  bool hasNuisParam = (sbModel->GetNuisanceParameters() && sbModel->GetNuisanceParameters()->getSize() > 0);
680  bool hasGlobalObs = (sbModel->GetGlobalObservables() && sbModel->GetGlobalObservables()->getSize() > 0);
681  if (hasNuisParam && !hasGlobalObs ) {
682  // try to see if model has nuisance parameters first
683  RooAbsPdf * constrPdf = RooStats::MakeNuisancePdf(*sbModel,"nuisanceConstraintPdf_sbmodel");
684  if (constrPdf) {
685  Warning("StandardHypoTestInvDemo","Model %s has nuisance parameters but no global observables associated",sbModel->GetName());
686  Warning("StandardHypoTestInvDemo","\tThe effect of the nuisance parameters will not be treated correctly ");
687  }
688  }
689  }
690 
691  // save all initial parameters of the model including the global observables
692  RooArgSet initialParameters;
693  RooArgSet * allParams = sbModel->GetPdf()->getParameters(*data);
694  allParams->snapshot(initialParameters);
695  delete allParams;
696 
697  // run first a data fit
698 
699  const RooArgSet * poiSet = sbModel->GetParametersOfInterest();
700  RooRealVar *poi = (RooRealVar*)poiSet->first();
701 
702  std::cout << "StandardHypoTestInvDemo : POI initial value: " << poi->GetName() << " = " << poi->getVal() << std::endl;
703 
704  // fit the data first (need to use constraint )
705  TStopwatch tw;
706 
707  bool doFit = mInitialFit;
708  if (testStatType == 0 && mInitialFit == -1) doFit = false; // case of LEP test statistic
709  if (type == 3 && mInitialFit == -1) doFit = false; // case of Asymptoticcalculator with nominal Asimov
710  double poihat = 0;
711 
712  if (mMinimizerType.size()==0) mMinimizerType = ROOT::Math::MinimizerOptions::DefaultMinimizerType();
713  else
715 
716  Info("StandardHypoTestInvDemo","Using %s as minimizer for computing the test statistic",
718 
719  if (doFit) {
720 
721  // do the fit : By doing a fit the POI snapshot (for S+B) is set to the fit value
722  // and the nuisance parameters nominal values will be set to the fit value.
723  // This is relevant when using LEP test statistics
724 
725  Info( "StandardHypoTestInvDemo"," Doing a first fit to the observed data ");
726  RooArgSet constrainParams;
727  if (sbModel->GetNuisanceParameters() ) constrainParams.add(*sbModel->GetNuisanceParameters());
728  RooStats::RemoveConstantParameters(&constrainParams);
729  tw.Start();
730  RooFitResult * fitres = sbModel->GetPdf()->fitTo(*data,InitialHesse(false), Hesse(false),
731  Minimizer(mMinimizerType.c_str(),"Migrad"), Strategy(0), PrintLevel(mPrintLevel), Constrain(constrainParams), Save(true), Offset(RooStats::IsNLLOffset()) );
732  if (fitres->status() != 0) {
733  Warning("StandardHypoTestInvDemo","Fit to the model failed - try with strategy 1 and perform first an Hesse computation");
734  fitres = sbModel->GetPdf()->fitTo(*data,InitialHesse(true), Hesse(false),Minimizer(mMinimizerType.c_str(),"Migrad"), Strategy(1), PrintLevel(mPrintLevel+1), Constrain(constrainParams),
735  Save(true), Offset(RooStats::IsNLLOffset()) );
736  }
737  if (fitres->status() != 0)
738  Warning("StandardHypoTestInvDemo"," Fit still failed - continue anyway.....");
739 
740 
741  poihat = poi->getVal();
742  std::cout << "StandardHypoTestInvDemo - Best Fit value : " << poi->GetName() << " = "
743  << poihat << " +/- " << poi->getError() << std::endl;
744  std::cout << "Time for fitting : "; tw.Print();
745 
746  //save best fit value in the poi snapshot
747  sbModel->SetSnapshot(*sbModel->GetParametersOfInterest());
748  std::cout << "StandardHypoTestInvo: snapshot of S+B Model " << sbModel->GetName()
749  << " is set to the best fit value" << std::endl;
750 
751  }
752 
753  // print a message in case of LEP test statistics because it affects result by doing or not doing a fit
754  if (testStatType == 0) {
755  if (!doFit)
756  Info("StandardHypoTestInvDemo","Using LEP test statistic - an initial fit is not done and the TS will use the nuisances at the model value");
757  else
758  Info("StandardHypoTestInvDemo","Using LEP test statistic - an initial fit has been done and the TS will use the nuisances at the best fit value");
759  }
760 
761 
762  // build test statistics and hypotest calculators for running the inverter
763 
764  SimpleLikelihoodRatioTestStat slrts(*sbModel->GetPdf(),*bModel->GetPdf());
765 
766  // null parameters must includes snapshot of poi plus the nuisance values
767  RooArgSet nullParams(*sbModel->GetSnapshot());
768  if (sbModel->GetNuisanceParameters()) nullParams.add(*sbModel->GetNuisanceParameters());
769  if (sbModel->GetSnapshot()) slrts.SetNullParameters(nullParams);
770  RooArgSet altParams(*bModel->GetSnapshot());
771  if (bModel->GetNuisanceParameters()) altParams.add(*bModel->GetNuisanceParameters());
772  if (bModel->GetSnapshot()) slrts.SetAltParameters(altParams);
773  if (mEnableDetOutput) slrts.EnableDetailedOutput();
774 
775  // ratio of profile likelihood - need to pass snapshot for the alt
777  ropl(*sbModel->GetPdf(), *bModel->GetPdf(), bModel->GetSnapshot());
778  ropl.SetSubtractMLE(false);
779  if (testStatType == 11) ropl.SetSubtractMLE(true);
780  ropl.SetPrintLevel(mPrintLevel);
781  ropl.SetMinimizer(mMinimizerType.c_str());
782  if (mEnableDetOutput) ropl.EnableDetailedOutput();
783 
784  ProfileLikelihoodTestStat profll(*sbModel->GetPdf());
785  if (testStatType == 3) profll.SetOneSided(true);
786  if (testStatType == 4) profll.SetSigned(true);
787  profll.SetMinimizer(mMinimizerType.c_str());
788  profll.SetPrintLevel(mPrintLevel);
789  if (mEnableDetOutput) profll.EnableDetailedOutput();
790 
791  profll.SetReuseNLL(mOptimize);
792  slrts.SetReuseNLL(mOptimize);
793  ropl.SetReuseNLL(mOptimize);
794 
795  if (mOptimize) {
796  profll.SetStrategy(0);
797  ropl.SetStrategy(0);
799  }
800 
801  if (mMaxPoi > 0) poi->setMax(mMaxPoi); // increase limit
802 
803  MaxLikelihoodEstimateTestStat maxll(*sbModel->GetPdf(),*poi);
804  NumEventsTestStat nevtts;
805 
806  AsymptoticCalculator::SetPrintLevel(mPrintLevel);
807 
808  // create the HypoTest calculator class
809  HypoTestCalculatorGeneric * hc = 0;
810  if (type == 0) hc = new FrequentistCalculator(*data, *bModel, *sbModel);
811  else if (type == 1) hc = new HybridCalculator(*data, *bModel, *sbModel);
812  // else if (type == 2 ) hc = new AsymptoticCalculator(*data, *bModel, *sbModel, false, mAsimovBins);
813  // else if (type == 3 ) hc = new AsymptoticCalculator(*data, *bModel, *sbModel, true, mAsimovBins); // for using Asimov data generated with nominal values
814  else if (type == 2 ) hc = new AsymptoticCalculator(*data, *bModel, *sbModel, false );
815  else if (type == 3 ) hc = new AsymptoticCalculator(*data, *bModel, *sbModel, true ); // for using Asimov data generated with nominal values
816  else {
817  Error("StandardHypoTestInvDemo","Invalid - calculator type = %d supported values are only :\n\t\t\t 0 (Frequentist) , 1 (Hybrid) , 2 (Asymptotic) ",type);
818  return 0;
819  }
820 
821  // set the test statistic
822  TestStatistic * testStat = 0;
823  if (testStatType == 0) testStat = &slrts;
824  if (testStatType == 1 || testStatType == 11) testStat = &ropl;
825  if (testStatType == 2 || testStatType == 3 || testStatType == 4) testStat = &profll;
826  if (testStatType == 5) testStat = &maxll;
827  if (testStatType == 6) testStat = &nevtts;
828 
829  if (testStat == 0) {
830  Error("StandardHypoTestInvDemo","Invalid - test statistic type = %d supported values are only :\n\t\t\t 0 (SLR) , 1 (Tevatron) , 2 (PLR), 3 (PLR1), 4(MLE)",testStatType);
831  return 0;
832  }
833 
834 
835  ToyMCSampler *toymcs = (ToyMCSampler*)hc->GetTestStatSampler();
836  if (toymcs && (type == 0 || type == 1) ) {
837  // look if pdf is number counting or extended
838  if (sbModel->GetPdf()->canBeExtended() ) {
839  if (useNumberCounting) Warning("StandardHypoTestInvDemo","Pdf is extended: but number counting flag is set: ignore it ");
840  }
841  else {
842  // for not extended pdf
843  if (!useNumberCounting ) {
844  int nEvents = data->numEntries();
845  Info("StandardHypoTestInvDemo","Pdf is not extended: number of events to generate taken from observed data set is %d",nEvents);
846  toymcs->SetNEventsPerToy(nEvents);
847  }
848  else {
849  Info("StandardHypoTestInvDemo","using a number counting pdf");
850  toymcs->SetNEventsPerToy(1);
851  }
852  }
853 
854  toymcs->SetTestStatistic(testStat);
855 
856  if (data->isWeighted() && !mGenerateBinned) {
857  Info("StandardHypoTestInvDemo","Data set is weighted, nentries = %d and sum of weights = %8.1f but toy generation is unbinned - it would be faster to set mGenerateBinned to true\n",data->numEntries(), data->sumEntries());
858  }
859  toymcs->SetGenerateBinned(mGenerateBinned);
860 
861  toymcs->SetUseMultiGen(mOptimize);
862 
863  if (mGenerateBinned && sbModel->GetObservables()->getSize() > 2) {
864  Warning("StandardHypoTestInvDemo","generate binned is activated but the number of observable is %d. Too much memory could be needed for allocating all the bins",sbModel->GetObservables()->getSize() );
865  }
866 
867  // set the random seed if needed
868  if (mRandomSeed >= 0) RooRandom::randomGenerator()->SetSeed(mRandomSeed);
869 
870  }
871 
872  // specify if need to re-use same toys
873  if (mReuseAltToys) {
874  hc->UseSameAltToys();
875  }
876 
877  if (type == 1) {
878  HybridCalculator *hhc = dynamic_cast<HybridCalculator*> (hc);
879  assert(hhc);
880 
881  hhc->SetToys(ntoys,ntoys/mNToysRatio); // can use less ntoys for b hypothesis
882 
883  // remove global observables from ModelConfig (this is probably not needed anymore in 5.32)
884  bModel->SetGlobalObservables(RooArgSet() );
885  sbModel->SetGlobalObservables(RooArgSet() );
886 
887 
888  // check for nuisance prior pdf in case of nuisance parameters
889  if (bModel->GetNuisanceParameters() || sbModel->GetNuisanceParameters() ) {
890 
891  // fix for using multigen (does not work in this case)
892  toymcs->SetUseMultiGen(false);
893  ToyMCSampler::SetAlwaysUseMultiGen(false);
894 
895  RooAbsPdf * nuisPdf = 0;
896  if (nuisPriorName) nuisPdf = w->pdf(nuisPriorName);
897  // use prior defined first in bModel (then in SbModel)
898  if (!nuisPdf) {
899  Info("StandardHypoTestInvDemo","No nuisance pdf given for the HybridCalculator - try to deduce pdf from the model");
900  if (bModel->GetPdf() && bModel->GetObservables() )
901  nuisPdf = RooStats::MakeNuisancePdf(*bModel,"nuisancePdf_bmodel");
902  else
903  nuisPdf = RooStats::MakeNuisancePdf(*sbModel,"nuisancePdf_sbmodel");
904  }
905  if (!nuisPdf ) {
906  if (bModel->GetPriorPdf()) {
907  nuisPdf = bModel->GetPriorPdf();
908  Info("StandardHypoTestInvDemo","No nuisance pdf given - try to use %s that is defined as a prior pdf in the B model",nuisPdf->GetName());
909  }
910  else {
911  Error("StandardHypoTestInvDemo","Cannot run Hybrid calculator because no prior on the nuisance parameter is specified or can be derived");
912  return 0;
913  }
914  }
915  assert(nuisPdf);
916  Info("StandardHypoTestInvDemo","Using as nuisance Pdf ... " );
917  nuisPdf->Print();
918 
919  const RooArgSet * nuisParams = (bModel->GetNuisanceParameters() ) ? bModel->GetNuisanceParameters() : sbModel->GetNuisanceParameters();
920  RooArgSet * np = nuisPdf->getObservables(*nuisParams);
921  if (np->getSize() == 0) {
922  Warning("StandardHypoTestInvDemo","Prior nuisance does not depend on nuisance parameters. They will be smeared in their full range");
923  }
924  delete np;
925 
926  hhc->ForcePriorNuisanceAlt(*nuisPdf);
927  hhc->ForcePriorNuisanceNull(*nuisPdf);
928 
929 
930  }
931  }
932  else if (type == 2 || type == 3) {
933  if (testStatType == 3) ((AsymptoticCalculator*) hc)->SetOneSided(true);
934  if (testStatType != 2 && testStatType != 3)
935  Warning("StandardHypoTestInvDemo","Only the PL test statistic can be used with AsymptoticCalculator - use by default a two-sided PL");
936  }
937  else if (type == 0 ) {
938  ((FrequentistCalculator*) hc)->SetToys(ntoys,ntoys/mNToysRatio);
939  // store also the fit information for each poi point used by calculator based on toys
940  if (mEnableDetOutput) ((FrequentistCalculator*) hc)->StoreFitInfo(true);
941  }
942  else if (type == 1 ) {
943  ((HybridCalculator*) hc)->SetToys(ntoys,ntoys/mNToysRatio);
944  // store also the fit information for each poi point used by calculator based on toys
945  //if (mEnableDetOutput) ((HybridCalculator*) hc)->StoreFitInfo(true);
946  }
947 
948  // Get the result
950 
951 
952 
953  HypoTestInverter calc(*hc);
954  calc.SetConfidenceLevel(optHTInv.confLevel);
955 
956 
957  calc.UseCLs(useCLs);
958  calc.SetVerbose(true);
959 
960  // can speed up using proof-lite
961  if (mUseProof) {
962  ProofConfig pc(*w, mNWorkers, "", kFALSE);
963  toymcs->SetProofConfig(&pc); // enable proof
964  }
965 
966 
967  if (npoints > 0) {
968  if (poimin > poimax) {
969  // if no min/max given scan between MLE and +4 sigma
970  poimin = int(poihat);
971  poimax = int(poihat + 4 * poi->getError());
972  }
973  std::cout << "Doing a fixed scan in interval : " << poimin << " , " << poimax << std::endl;
974  calc.SetFixedScan(npoints,poimin,poimax);
975  }
976  else {
977  //poi->setMax(10*int( (poihat+ 10 *poi->getError() )/10 ) );
978  std::cout << "Doing an automatic scan in interval : " << poi->getMin() << " , " << poi->getMax() << std::endl;
979  }
980 
981  tw.Start();
982  HypoTestInverterResult * r = calc.GetInterval();
983  std::cout << "Time to perform limit scan \n";
984  tw.Print();
985 
986  if (mRebuild) {
987 
988  std::cout << "\n***************************************************************\n";
989  std::cout << "Rebuild the upper limit distribution by re-generating new set of pseudo-experiment and re-compute for each of them a new upper limit\n\n";
990 
991 
992  allParams = sbModel->GetPdf()->getParameters(*data);
993 
994  // define on which value of nuisance parameters to do the rebuild
995  // default is best fit value for bmodel snapshot
996 
997 
998 
999  if (mRebuildParamValues != 0) {
1000  // set all parameters to their initial workspace values
1001  *allParams = initialParameters;
1002  }
1003  if (mRebuildParamValues == 0 || mRebuildParamValues == 1 ) {
1004  RooArgSet constrainParams;
1005  if (sbModel->GetNuisanceParameters() ) constrainParams.add(*sbModel->GetNuisanceParameters());
1006  RooStats::RemoveConstantParameters(&constrainParams);
1007 
1008  const RooArgSet * poiModel = sbModel->GetParametersOfInterest();
1009  bModel->LoadSnapshot();
1010 
1011  // do a profile using the B model snapshot
1012  if (mRebuildParamValues == 0 ) {
1013 
1014  RooStats::SetAllConstant(*poiModel,true);
1015 
1016  sbModel->GetPdf()->fitTo(*data,InitialHesse(false), Hesse(false),
1017  Minimizer(mMinimizerType.c_str(),"Migrad"), Strategy(0), PrintLevel(mPrintLevel), Constrain(constrainParams), Offset(RooStats::IsNLLOffset()) );
1018 
1019 
1020  std::cout << "rebuild using fitted parameter value for B-model snapshot" << std::endl;
1021  constrainParams.Print("v");
1022 
1023  RooStats::SetAllConstant(*poiModel,false);
1024  }
1025  }
1026  std::cout << "StandardHypoTestInvDemo: Initial parameters used for rebuilding: ";
1027  RooStats::PrintListContent(*allParams, std::cout);
1028  delete allParams;
1029 
1030  calc.SetCloseProof(1);
1031  tw.Start();
1032  SamplingDistribution * limDist = calc.GetUpperLimitDistribution(true,mNToyToRebuild);
1033  std::cout << "Time to rebuild distributions " << std::endl;
1034  tw.Print();
1035 
1036  if (limDist) {
1037  std::cout << "Expected limits after rebuild distribution " << std::endl;
1038  std::cout << "expected upper limit (median of limit distribution) " << limDist->InverseCDF(0.5) << std::endl;
1039  std::cout << "expected -1 sig limit (0.16% quantile of limit dist) " << limDist->InverseCDF(ROOT::Math::normal_cdf(-1)) << std::endl;
1040  std::cout << "expected +1 sig limit (0.84% quantile of limit dist) " << limDist->InverseCDF(ROOT::Math::normal_cdf(1)) << std::endl;
1041  std::cout << "expected -2 sig limit (.025% quantile of limit dist) " << limDist->InverseCDF(ROOT::Math::normal_cdf(-2)) << std::endl;
1042  std::cout << "expected +2 sig limit (.975% quantile of limit dist) " << limDist->InverseCDF(ROOT::Math::normal_cdf(2)) << std::endl;
1043 
1044  // Plot the upper limit distribution
1045  SamplingDistPlot limPlot( (mNToyToRebuild < 200) ? 50 : 100);
1046  limPlot.AddSamplingDistribution(limDist);
1047  limPlot.GetTH1F()->SetStats(true); // display statistics
1048  limPlot.SetLineColor(kBlue);
1049  new TCanvas("limPlot","Upper Limit Distribution");
1050  limPlot.Draw();
1051 
1052  /// save result in a file
1053  limDist->SetName("RULDist");
1054  TFile * fileOut = new TFile("RULDist.root","RECREATE");
1055  limDist->Write();
1056  fileOut->Close();
1057 
1058 
1059  //update r to a new updated result object containing the rebuilt expected p-values distributions
1060  // (it will not recompute the expected limit)
1061  if (r) delete r; // need to delete previous object since GetInterval will return a cloned copy
1062  r = calc.GetInterval();
1063 
1064  }
1065  else
1066  std::cout << "ERROR : failed to re-build distributions " << std::endl;
1067  }
1068 
1069  return r;
1070 }
1071 
1072 
1073 
1074 void ReadResult(const char * fileName, const char * resultName="", bool useCLs=true) {
1075  // read a previous stored result from a file given the result name
1076 
1077  StandardHypoTestInvDemo(fileName, resultName,"","","",0,0,useCLs);
1078 }
1079 
1080 
1081 #ifdef USE_AS_MAIN
1082 int main() {
1083  StandardHypoTestInvDemo();
1084 }
1085 #endif
1086 
1087 
1088 
1089 
virtual Double_t sumEntries() const =0
virtual Double_t getMin(const char *name=0) const
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition: TObject.cxx:785
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:1276
void UseSameAltToys()
to re-use same toys for alternate hypothesis
RooCmdArg Offset(Bool_t flag=kTRUE)
Holds configuration options for proof and proof-lite.
Definition: ProofConfig.h:46
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:30
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition: TObject.cxx:854
const RooArgSet * GetObservables() const
get RooArgSet for observables (return NULL if not existing)
Definition: ModelConfig.h:237
virtual Double_t getMax(const char *name=0) const
virtual Bool_t add(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling add() for each element in the source coll...
Definition: RooArgSet.h:86
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Definition: TStopwatch.cxx:58
void LoadSnapshot() const
load the snapshot from ws if it exists
HypoTestInverter class for performing an hypothesis test inversion by scanning the hypothesis test re...
void Print(Option_t *option="") const
Print the real and cpu time passed between the start and stop events.
Definition: TStopwatch.cxx:219
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Definition: RooAbsArg.h:194
RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name)
return c1
Definition: legend1.C:41
RooCmdArg PrintLevel(Int_t code)
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
double GetExpectedUpperLimit(double nsig=0, const char *opt="") const
get Limit value corresponding at the desired nsigma level (0) is median -1 sigma is 1 sigma ...
void convertToVectorStore()
Convert tree-based storage to vector-based storage.
Definition: RooAbsData.cxx:253
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition: TNamed.cxx:140
Double_t InverseCDF(Double_t pvalue)
get the inverse of the Cumulative distribution function
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:688
RooCmdArg Strategy(Int_t code)
virtual void SetNEventsPerToy(const Int_t nevents)
Definition: ToyMCSampler.h:147
#define gROOT
Definition: TROOT.h:410
StreamConfig & getStream(Int_t id)
static void setDefaultStorageType(StorageType s)
Definition: RooAbsData.cxx:77
void removeTopic(RooFit::MsgTopic oldTopic)
static RooMsgService & instance()
Return reference to singleton instance.
STL namespace.
void SetLogYaxis(Bool_t ly)
changes plot to log scale on y axis
virtual ModelConfig * Clone(const char *name="") const
clone
Definition: ModelConfig.h:54
void setMax(const char *name, Double_t value)
Set maximum of name range to given value.
Definition: RooRealVar.cxx:417
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=1, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3976
virtual void SetTestStatistic(TestStatistic *testStatistic, unsigned int i)
Definition: ToyMCSampler.h:184
Common base class for the Hypothesis Test Calculators.
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString...
Definition: TString.cxx:2286
virtual void SetSeed(ULong_t seed=0)
Set the random generator seed.
Definition: TRandom.cxx:589
double normal_cdf(double x, double sigma=1, double x0=0)
Cumulative distribution function of the normal (Gaussian) distribution (lower tail).
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooAbsArg.h:227
void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
int main(int argc, char **argv)
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition: RooRandom.cxx:54
friend class RooArgSet
Definition: RooAbsArg.h:471
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
Double_t LowerLimit()
lower and upper bound of the confidence interval (to get upper/lower limits, multiply the size( = 1-c...
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
RooAbsData * data(const char *name) const
Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found...
RooAbsPdf * GetPriorPdf() const
get parameters prior pdf (return NULL if not existing)
Definition: ModelConfig.h:234
static void SetDefaultStrategy(int strat)
virtual void setVal(Double_t value)
Set value of variable to &#39;value&#39;.
Definition: RooRealVar.cxx:204
Int_t getSize() const
virtual void ls(Option_t *option="") const
List TNamed name and title.
Definition: TNamed.cxx:113
RooAbsCollection * snapshot(Bool_t deepCopy=kTRUE) const
Take a snap shot of current collection contents: An owning collection is returned containing clones o...
RooCmdArg InitialHesse(Bool_t flag=kTRUE)
ROOT::R::TRInterface & r
Definition: Object.C:4
R__EXTERN TSystem * gSystem
Definition: TSystem.h:540
static const std::string & DefaultMinimizerType()
RooAbsArg * first() const
RooCmdArg Minimizer(const char *type, const char *alg=0)
void SetProofConfig(ProofConfig *pc=NULL)
Definition: ToyMCSampler.h:233
void PrintListContent(const RooArgList &l, std::ostream &os=std::cout)
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
ProfileLikelihoodTestStat is an implementation of the TestStatistic interface that calculates the pro...
ToyMCSampler is an implementation of the TestStatSampler interface.
Definition: ToyMCSampler.h:71
This class simply holds a sampling distribution of some test statistic.
TestStatistic that returns the ratio of profiled likelihoods.
bool IsTwoSided() const
query if two sided result
Bool_t canBeExtended() const
Definition: RooAbsPdf.h:216
void UseNLLOffset(bool on)
Same purpose as HybridCalculatorOriginal, but different implementation.
TObject * obj(const char *name) const
Return any type of object (RooAbsArg, RooAbsData or generic object) with given name) ...
const Bool_t kFALSE
Definition: RtypesCore.h:88
The Canvas class.
Definition: TCanvas.h:31
Namespace for the RooStats classes.
Definition: Asimov.h:20
RooAbsPdf * GetPdf() const
get model PDF (return NULL if pdf has not been specified or does not exist)
Definition: ModelConfig.h:222
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return NULL if not existing)
Definition: ModelConfig.h:225
return c2
Definition: legend2.C:14
RooCmdArg Hesse(Bool_t flag=kTRUE)
HypoTestInverterResult class holds the array of hypothesis test results and compute a confidence inte...
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
Class to plot an HypoTestInverterResult, result of the HypoTestInverter calculator.
RooArgSet * getParameters(const RooAbsData *data, Bool_t stripDisconnected=kTRUE) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don&#39;t match any of...
Definition: RooAbsArg.cxx:532
Double_t LowerLimitEstimatedError()
rough estimation of the error on the computed bound of the confidence interval Estimate of lower limi...
int type
Definition: TGX11.cxx:120
bool SetAllConstant(const RooAbsCollection &coll, bool constant=true)
Definition: RooStatsUtils.h:82
void SetToys(int toysNull, int toysAlt)
set number of toys
TestStatistic class that returns -log(L[null] / L[alt]) where L is the likelihood.
bool IsNLLOffset()
virtual Bool_t isWeighted() const
Definition: RooAbsData.h:94
This class provides simple and straightforward utilities to plot SamplingDistribution objects...
Mother of all ROOT objects.
Definition: TObject.h:37
MaxLikelihoodEstimateTestStat: TestStatistic that returns maximum likelihood estimate of a specified ...
Double_t UpperLimitEstimatedError()
Estimate of lower limit error function evaluates only a rough error on the lower limit.
RooCmdArg Save(Bool_t flag=kTRUE)
void RemoveConstantParameters(RooArgSet *set)
Definition: RooStatsUtils.h:62
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1162
const RooArgSet * GetGlobalObservables() const
get RooArgSet for global observables (return NULL if not existing)
Definition: ModelConfig.h:243
Definition: file.py:1
const RooArgSet * GetNuisanceParameters() const
get RooArgSet containing the nuisance parameters (return NULL if not existing)
Definition: ModelConfig.h:228
Does a frequentist hypothesis test.
#define gPad
Definition: TVirtualPad.h:285
virtual RooFitResult * fitTo(RooAbsData &data, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none())
Fit PDF to given dataset.
Definition: RooAbsPdf.cxx:1079
static constexpr double pc
void SetUseMultiGen(Bool_t flag)
Definition: ToyMCSampler.h:81
Definition: Rtypes.h:59
void SetGenerateBinned(bool binned=true)
Definition: ToyMCSampler.h:203
const RooArgSet * GetSnapshot() const
get RooArgSet for parameters for a particular hypothesis (return NULL if not existing) ...
Double_t Sqrt(Double_t x)
Definition: TMath.h:690
virtual void SetGlobalObservables(const RooArgSet &set)
specify the global observables
Definition: ModelConfig.h:160
static void SetDefaultMinimizer(const char *type, const char *algo=0)
Double_t getError() const
Definition: RooRealVar.h:53
void Print(Option_t *opts=0) const
Print contents of the workspace.
virtual void SetSnapshot(const RooArgSet &set)
set parameter values for a particular hypothesis if using a common PDF by saving a snapshot in the wo...
TestStatistic is an interface class to provide a facility for construction test statistics distributi...
Definition: TestStatistic.h:31
int ArraySize() const
number of entries in the results array
Hypothesis Test Calculator based on the asymptotic formulae for the profile likelihood ratio...
Int_t CeilNint(Double_t x)
Definition: TMath.h:698
char name[80]
Definition: TGX11.cxx:109
virtual void ForcePriorNuisanceAlt(RooAbsPdf &priorNuisance)
TestStatSampler * GetTestStatSampler(void) const
RooCmdArg Constrain(const RooArgSet &params)
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:866
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:43
virtual void ForcePriorNuisanceNull(RooAbsPdf &priorNuisance)
Override the distribution used for marginalizing nuisance parameters that is inferred from ModelConfi...
virtual Int_t numEntries() const
Definition: RooAbsData.cxx:285
Stopwatch class.
Definition: TStopwatch.h:28
Int_t status() const
Definition: RooFitResult.h:77
NumEventsTestStat is a simple implementation of the TestStatistic interface used for simple number co...