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