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