ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/cbrown/AnalysisFramework/limits/roostats_cl95.C
Revision: 1.1
Committed: Fri Jun 24 07:48:47 2011 UTC (13 years, 10 months ago) by fronga
Content type: text/plain
Branch: MAIN
CVS Tags: cbaf_4_98ifb_paper, cbaf_4p7ifb, Honeypot, cbaf_2p1ifb, HEAD
Error occurred while calculating annotation data.
Log Message:
Adding macros used to compute limits.

File Contents

# Content
1 static const char* desc =
2 "=====================================================================\n"
3 "| \n"
4 "|\033[1m roostats_cl95.C version 1.02 \033[0m\n"
5 "| \n"
6 "| Standard c++ routine for 95% C.L. limit calculation \n"
7 "| for cross section in a 'counting experiment' \n"
8 "| Fully backwards-compatible with the CL95 macro \n"
9 "| \n"
10 "| also known as 'CL95 with RooStats' \n"
11 "| \n"
12 "|\033[1m Gena Kukartsev, Stefan Schmitz, Gregory Schott \033[0m\n"
13 "| \n"
14 "| July 2010: first version \n"
15 "| March 2011: restructuring, interface change, expected limits \n"
16 "| May 2011: added expected limit median, \n"
17 "| 68%, 95% quantile bands and actual coverage \n"
18 "| \n"
19 "=====================================================================\n"
20 " \n"
21 "Prerequisites: \n"
22 " ROOT version 5.27/06 or higher \n"
23 " \n"
24 " \n"
25 " \n"
26 "The code should be compiled in ROOT: \n"
27 " \n"
28 "root -l \n"
29 " \n"
30 ".L roostats_cl95.C+ \n"
31 " \n"
32 "Usage: \n"
33 " limit = roostats_cl95(ilum, slum, eff, seff, bck, sbck, n, gauss = false, nuisanceModel, method, plotFileName); \n"
34 " limit_result = roostats_cla(ilum, slum, eff, seff, bck, sbck, nuisanceModel); \n"
35 " limitA = roostats_cla(ilum, slum, eff, seff, bck, sbck, nuisanceModel); \n"
36 " \n"
37 "Inputs: \n"
38 " ilum - Nominal integrated luminosity (pb^-1) \n"
39 " slum - Absolute error on the integrated luminosity \n"
40 " eff - Nominal value of the efficiency times \n"
41 " acceptance (in range 0 to 1) \n"
42 " seff - Absolute error on the efficiency times \n"
43 " acceptance \n"
44 " bck - Nominal value of the background estimate \n"
45 " sbck - Absolute error on the background \n"
46 " n - Number of observed events (not used for the \n"
47 " expected limit) \n"
48 " gauss - if true, use Gaussian statistics for signal \n"
49 " instead of Poisson; automatically false \n"
50 " for n = 0. \n"
51 " Always false for expected limit calculations \n"
52 " nuisanceModel - distribution function used in integration over\n"
53 " nuisance parameters: \n"
54 " 0 - Gaussian (default), 1 - lognormal, \n"
55 " 2 - gamma; \n"
56 " (automatically 0 when gauss == true) \n"
57 " method - method of statistical inference: \n"
58 " \"bayesian\" - Bayesian with numeric \n"
59 " integration (default), \n"
60 " \"workspace\" - only create workspace and save\n"
61 " to file, no interval calculation\n"
62 " plotFileName - file name for the control plot to be created \n"
63 " file name extension will define the format, \n"
64 " <plot_cl95.pdf> is the default value, \n"
65 " specify empty string if you do not want \n"
66 " the plot to be created (saves time) \n"
67 " \n"
68 " \n"
69 "The statistics model in this routine: the routine addresses the task \n"
70 "of a Bayesian evaluation of limits for a one-bin counting experiment \n"
71 "with systematic uncertainties on luminosity and efficiency for the \n"
72 "signal and a global uncertainty on the expected background (implying \n"
73 "no correlated error on the luminosity for signal and background, \n"
74 "which will not be suitable for all use cases!). The observable is the\n"
75 "measured number of events. Here the bayesian 90% interval (upper \n"
76 "limit corresponds to one-sided 95% upper limit) is evaluated as a \n"
77 "function of the observed number of events in a hypothetical \n"
78 "experiment. \n"
79 " \n"
80 "\033[1m Note! \033[0m\n"
81 "If you are running nonstandard ROOT environment, e.g. in CMSSW, \n"
82 "you need to make sure that the RooFit and RooStats header files \n"
83 "can be found since they might be in a nonstandard location. \n"
84 " \n"
85 "For CMSSW_4_2_0_pre8 and later, add the following line to your \n"
86 "rootlogon.C: \n"
87 " gSystem -> SetIncludePath( \"-I$ROOFITSYS/include\" ); \n";
88
89
90 #include <algorithm>
91
92 #include "TCanvas.h"
93 #include "TMath.h"
94 #include "TRandom3.h"
95 #include "TUnixSystem.h"
96 #include "TStopwatch.h"
97
98 #include "RooPlot.h"
99 #include "RooRealVar.h"
100 #include "RooProdPdf.h"
101 #include "RooWorkspace.h"
102 #include "RooDataSet.h"
103
104 #include "RooStats/ModelConfig.h"
105 #include "RooStats/SimpleInterval.h"
106 #include "RooStats/BayesianCalculator.h"
107 #include "RooRandom.h"
108
109 // FIXME: remove namespaces
110 using namespace RooFit;
111 using namespace RooStats;
112 using namespace std;
113
114 class LimitResult;
115
116 Double_t roostats_cl95(Double_t ilum, Double_t slum,
117 Double_t eff, Double_t seff,
118 Double_t bck, Double_t sbck,
119 Int_t n,
120 Bool_t gauss = kFALSE,
121 Int_t nuisanceModel = 0,
122 std::string method = "bayesian",
123 std::string plotFileName = "plot_cl95.pdf");
124
125 LimitResult roostats_clm(Double_t ilum, Double_t slum,
126 Double_t eff, Double_t seff,
127 Double_t bck, Double_t sbck,
128 Int_t nit = 200, Int_t nuisanceModel = 0,
129 std::string method = "bayesian");
130
131 // legacy support: use roostats_clm() instead
132 Double_t roostats_cla(Double_t ilum, Double_t slum,
133 Double_t eff, Double_t seff,
134 Double_t bck, Double_t sbck,
135 Int_t nuisanceModel = 0,
136 std::string method = "bayesian");
137
138
139
140
141 // ---> implementation below --------------------------------------------
142
143
144 class LimitResult{
145
146 friend class CL95Calc;
147
148 public:
149 LimitResult():
150 _expected_limit(0),
151 _low68(0),
152 _high68(0),
153 _low95(0),
154 _high95(0){};
155
156 ~LimitResult(){};
157
158 Double_t GetExpectedLimit(){return _expected_limit;};
159
160 Double_t GetOneSigmaLowRange(){return _low68;};
161 Double_t GetOneSigmaHighRange(){return _high68;};
162 Double_t GetOneSigmaCoverage(){return _cover68;};
163
164 Double_t GetTwoSigmaLowRange(){return _low95;};
165 Double_t GetTwoSigmaHighRange(){return _high95;};
166 Double_t GetTwoSigmaCoverage(){return _cover95;};
167
168 private:
169 Double_t _expected_limit;
170 Double_t _low68;
171 Double_t _high68;
172 Double_t _low95;
173 Double_t _high95;
174 Double_t _cover68;
175 Double_t _cover95;
176 };
177
178
179 class CL95Calc{
180
181 public:
182 CL95Calc();
183 ~CL95Calc();
184
185 RooWorkspace * makeWorkspace(Double_t ilum, Double_t slum,
186 Double_t eff, Double_t seff,
187 Double_t bck, Double_t sbck,
188 Bool_t gauss,
189 Int_t nuisanceModel);
190 RooWorkspace * getWorkspace(){ return ws;}
191
192 RooAbsData * makeData(Int_t n);
193
194 Double_t cl95(std::string method = "bayesian");
195
196 Double_t cla( Double_t ilum, Double_t slum,
197 Double_t eff, Double_t seff,
198 Double_t bck, Double_t sbck,
199 Int_t nuisanceModel,
200 std::string method );
201
202 LimitResult clm(Double_t ilum, Double_t slum,
203 Double_t eff, Double_t seff,
204 Double_t bck, Double_t sbck,
205 Int_t nit = 200, Int_t nuisanceModel = 0,
206 std::string method = "bayesian");
207
208 int makePlot( std::string method,
209 std::string plotFileName = "plot_cl95.pdf" );
210
211 private:
212
213 // methods
214 Double_t GetRandom( std::string pdf, std::string var );
215
216 // data members
217 RooWorkspace * ws;
218 RooStats::ModelConfig mc;
219 RooAbsData * data;
220 BayesianCalculator * bcalc;
221 RooStats::SimpleInterval * sInt;
222 double nsig_rel_err;
223 double nbkg_rel_err;
224 Int_t _nuisance_model;
225
226 // random numbers
227 TRandom3 r;
228
229 // expected limits
230 Double_t _expected_limit;
231 Double_t _low68;
232 Double_t _high68;
233 Double_t _low95;
234 Double_t _high95;
235
236 };
237
238
239 // default constructor
240 CL95Calc::CL95Calc(){
241 ws = new RooWorkspace("ws");
242 data = 0;
243
244 sInt = 0;
245 bcalc = 0;
246 mc.SetName("modelconfig");
247 mc.SetTitle("ModelConfig for roostats_cl95");
248
249 nsig_rel_err = -1.0; // default non-initialized value
250 nbkg_rel_err = -1.0; // default non-initialized value
251
252 // set random seed
253 r.SetSeed();
254 UInt_t _seed = r.GetSeed();
255 UInt_t _pid = gSystem->GetPid();
256 std::cout << "[CL95Calc]: random seed: " << _seed << std::endl;
257 std::cout << "[CL95Calc]: process ID: " << _pid << std::endl;
258 _seed = 31*_seed+_pid;
259 std::cout << "[CL95Calc]: new random seed (31*seed+pid): " << _seed << std::endl;
260 r.SetSeed(_seed);
261
262 // set RooFit random seed (it has a private copy)
263 RooRandom::randomGenerator()->SetSeed(_seed);
264
265 // default Gaussian nuisance model
266 _nuisance_model = 0;
267 }
268
269
270 CL95Calc::~CL95Calc(){
271 delete ws;
272 delete data;
273 delete sInt;
274 delete bcalc;
275 }
276
277
278 RooWorkspace * CL95Calc::makeWorkspace(Double_t ilum, Double_t slum,
279 Double_t eff, Double_t seff,
280 Double_t bck, Double_t sbck,
281 Bool_t gauss,
282 Int_t nuisanceModel){
283
284 if ( bck>0.0 && (sbck/bck)<5.0 ){
285 // check that bck is not too close to zero,
286 // so lognormal and gamma modls still make sense
287 _nuisance_model = nuisanceModel;
288 }
289 else{
290 _nuisance_model = 0;
291 std::cout << "[CL95Calc]: background expectation is too close to zero compared to its uncertainty" << std::endl;
292 std::cout << "[CL95Calc]: switching to the Gaussian nuisance model" << std::endl;
293 }
294
295 // Workspace
296 // RooWorkspace * ws = new RooWorkspace("ws",true);
297
298 // observable: number of events
299 ws->factory( "n[0]" );
300
301 // integrated luminosity
302 ws->factory( "lumi[0]" );
303
304 // cross section - parameter of interest
305 ws->factory( "xsec[0]" );
306
307 // selection efficiency * acceptance
308 ws->factory( "efficiency[0]" );
309
310 // nuisance parameter: factor 1 with combined relative uncertainty
311 ws->factory( "nsig_nuis[1.0]" ); // will adjust range below
312
313 // signal yield
314 ws->factory( "prod::nsig(lumi,xsec,efficiency, nsig_nuis)" );
315
316 // estimated background yield
317 ws->factory( "bkg_est[0]" );
318
319 // nuisance parameter: factor 1 with background relative uncertainty
320 //ws->factory( "nbkg_nuis[1.0]" ); // will adjust range below
321
322 // background yield
323 //ws->factory( "prod::nbkg(bkg_est, nbkg_nuis)" );
324 ws->factory( "nbkg[1.0]" ); // will adjust value and range below
325
326 // core model:
327 ws->factory("sum::yield(nsig,nbkg)");
328 if (gauss){
329 // Poisson probability with mean signal+bkg
330 ws->factory( "Gaussian::model_core(n,yield,expr('sqrt(yield)',yield))" );
331 }
332 else{
333 // Poisson probability with mean signal+bkg
334 ws->factory( "Poisson::model_core(n,yield)" );
335 }
336
337
338 // systematic uncertainties
339 nsig_rel_err = sqrt(slum*slum/ilum/ilum+seff*seff/eff/eff);
340 nbkg_rel_err = sbck/bck;
341 if (nuisanceModel == 0){ // gaussian model for nuisance parameters
342
343 std::cout << "[roostats_cl95]: Gaussian PDFs for nuisance parameters" << endl;
344
345 // cumulative signal uncertainty
346 ws->factory( "nsig_sigma[0.1]" );
347 ws->factory( "Gaussian::syst_nsig(nsig_nuis, 1.0, nsig_sigma)" );
348 // background uncertainty
349 ws->factory( "nbkg_sigma[0.1]" );
350 //ws->factory( "Gaussian::syst_nbkg(nbkg_nuis, 1.0, nbkg_sigma)" );
351 ws->factory( "Gaussian::syst_nbkg(nbkg, bkg_est, nbkg_sigma)" );
352
353 ws->var("nsig_sigma")->setVal(nsig_rel_err);
354 //ws->var("nbkg_sigma")->setVal(nbkg_rel_err);
355 ws->var("nbkg_sigma")->setVal(sbck);
356 ws->var("nsig_sigma")->setConstant(kTRUE);
357 ws->var("nbkg_sigma")->setConstant(kTRUE);
358 }
359 else if (nuisanceModel == 1){// Lognormal model for nuisance parameters
360
361 std::cout << "[roostats_cl95]: Lognormal PDFs for nuisance parameters" << endl;
362
363 // cumulative signal uncertainty
364 ws->factory( "nsig_kappa[1.1]" );
365 ws->factory( "Lognormal::syst_nsig(nsig_nuis, 1.0, nsig_kappa)" );
366 // background uncertainty
367 ws->factory( "nbkg_kappa[1.1]" );
368 //ws->factory( "Lognormal::syst_nbkg(nbkg_nuis, 1.0, nbkg_kappa)" );
369 ws->factory( "Lognormal::syst_nbkg(nbkg, bkg_est, nbkg_kappa)" );
370
371 ws->var("nsig_kappa")->setVal(1.0 + nsig_rel_err);
372 ws->var("nbkg_kappa")->setVal(1.0 + nbkg_rel_err);
373 ws->var("nsig_kappa")->setConstant(kTRUE);
374 ws->var("nbkg_kappa")->setConstant(kTRUE);
375 }
376 else if (nuisanceModel == 2){ // Gamma model for nuisance parameters
377
378 std::cout << "[roostats_cl95]: Gamma PDFs for nuisance parameters" << endl;
379
380 // cumulative signal uncertainty
381 ws->factory( "nsig_beta[0.01]" );
382 ws->factory( "nsig_gamma[101.0]" );
383 ws->var("nsig_beta") ->setVal(nsig_rel_err*nsig_rel_err);
384 ws->var("nsig_gamma")->setVal(1.0/nsig_rel_err/nsig_rel_err + 1.0);
385 ws->factory( "Gamma::syst_nsig(nsig_nuis, nsig_gamma, nsig_beta, 0.0)" );
386
387 // background uncertainty
388 ws->factory( "nbkg_beta[0.01]" );
389 ws->factory( "nbkg_gamma[101.0]" );
390 //ws->var("nbkg_beta") ->setVal(nbkg_rel_err*nbkg_rel_err);
391 ws->var("nbkg_beta") ->setVal(sbck*sbck/bck);
392 ws->var("nbkg_gamma")->setVal(1.0/nbkg_rel_err/nbkg_rel_err + 1.0);
393 //ws->factory( "Gamma::syst_nbkg(nbkg_nuis, nbkg_gamma, nbkg_beta, 0.0)" );
394 ws->factory( "Gamma::syst_nbkg(nbkg, nbkg_gamma, nbkg_beta, 0.0)" );
395
396 ws->var("nsig_beta") ->setConstant(kTRUE);
397 ws->var("nsig_gamma")->setConstant(kTRUE);
398 ws->var("nbkg_beta") ->setConstant(kTRUE);
399 ws->var("nbkg_gamma")->setConstant(kTRUE);
400 }
401 else{
402 std::cout <<"[roostats_cl95]: undefined nuisance parameter model specified, exiting" << std::endl;
403 }
404
405 // model with systematics
406 ws->factory( "PROD::model(model_core, syst_nsig, syst_nbkg)" );
407
408 // flat prior for the parameter of interest
409 ws->factory( "Uniform::prior(xsec)" );
410
411
412 // parameter values
413 ws->var("lumi") ->setVal(ilum);
414 ws->var("efficiency")->setVal(eff);
415 ws->var("bkg_est") ->setVal(bck);
416 ws->var("xsec") ->setVal(0.0);
417 ws->var("nsig_nuis") ->setVal(1.0);
418 //ws->var("nbkg_nuis") ->setVal(1.0);
419 ws->var("nbkg") ->setVal(bck);
420
421 // set some parameters as constants
422 ws->var("lumi") ->setConstant(kTRUE);
423 ws->var("efficiency")->setConstant(kTRUE);
424 ws->var("bkg_est") ->setConstant(kTRUE);
425 ws->var("n") ->setConstant(kFALSE); // observable
426 ws->var("xsec") ->setConstant(kFALSE); // parameter of interest
427 ws->var("nsig_nuis") ->setConstant(kFALSE); // nuisance
428 //ws->var("nbkg_nuis") ->setConstant(kFALSE); // nuisance
429 ws->var("nbkg") ->setConstant(kFALSE); // nuisance
430
431 // floating parameters ranges
432 // crude estimates! Need to know data to do better
433 ws->var("n") ->setRange( 0.0, bck+(5.0*sbck)+10.0); // ad-hoc range for obs
434 ws->var("xsec") ->setRange( 0.0, 15.0*(1.0+nsig_rel_err)/ilum/eff ); // ad-hoc range for POI
435 ws->var("nsig_nuis")->setRange( std::max(0.0, 1.0 - 5.0*nsig_rel_err), 1.0 + 5.0*nsig_rel_err);
436 //ws->var("nbkg_nuis")->setRange( std::max(0.0, 1.0 - 5.0*nbkg_rel_err), 1.0 + 5.0*nbkg_rel_err);
437 ws->var("nbkg") ->setRange( std::max(0.0, bck - 5.0*sbck), bck + 5.0*sbck);
438
439 // Definition of observables and parameters of interest
440 ws->defineSet("obsSet","n");
441 ws->defineSet("poiSet","xsec");
442 //ws->defineSet("nuisanceSet","nsig_nuis,nbkg_nuis");
443 ws->defineSet("nuisanceSet","nsig_nuis,nbkg");
444
445 // setup the ModelConfig object
446 mc.SetWorkspace(*ws);
447 mc.SetPdf(*(ws->pdf("model")));
448 mc.SetParametersOfInterest(*(ws->set("poiSet")));
449 mc.SetPriorPdf(*(ws->pdf("prior")));
450 mc.SetNuisanceParameters(*(ws->set("nuisanceSet")));
451 mc.SetObservables(*(ws->set("obsSet")));
452
453 ws->import(mc);
454
455 return ws;
456 }
457
458
459 RooAbsData * CL95Calc::makeData( Int_t n ){
460 //
461 // make the dataset owned by the class
462 // the current one is deleted
463 //
464 // set ranges as well
465 //
466
467 // floating parameters ranges
468 if (nsig_rel_err < 0.0 || nbkg_rel_err < 0.0){
469 std::cout << "[roostats_cl95]: Workspace not initialized, cannot create a dataset" << std::endl;
470 return 0;
471 }
472
473 double ilum = ws->var("lumi")->getVal();
474 double eff = ws->var("efficiency")->getVal();
475 double bck = ws->var("bkg_est")->getVal();
476 double sbck = nbkg_rel_err*bck;
477
478 ws->var("n") ->setRange( 0.0, bck+(5.0*sbck)+10.0*(n+1.0)); // ad-hoc range for obs
479 ws->var("xsec") ->setRange( 0.0, 5.0*(1.0+nsig_rel_err)*std::max(10.0,n-bck)/ilum/eff ); // ad-hoc range for POI
480 ws->var("nsig_nuis")->setRange( std::max(0.0, 1.0 - 5.0*nsig_rel_err), 1.0 + 5.0*nsig_rel_err);
481 //ws->var("nbkg_nuis")->setRange( std::max(0.0, 1.0 - 5.0*nbkg_rel_err), 1.0 + 5.0*nbkg_rel_err);
482 ws->var("nbkg") ->setRange( std::max(0.0, bck - 5.0*sbck), bck + 5.0*sbck);
483
484 // create data
485 ws->var("n") ->setVal(n);
486 delete data;
487 data = new RooDataSet("data","",*(mc.GetObservables()));
488 data->add( *(mc.GetObservables()));
489
490 return data;
491 }
492
493
494 Double_t CL95Calc::cl95( std::string method ){
495
496 // this method assumes that the workspace,
497 // data and model config are ready
498
499 Double_t upper_limit = -1.0;
500
501 // make RooFit quiet
502 RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL);
503
504 if (method.find("bayesian") != std::string::npos){
505
506 //prepare Bayesian Calulator
507 delete bcalc;
508 bcalc = new BayesianCalculator(*data, mc);
509 TString namestring = "mybc";
510 bcalc->SetName(namestring);
511 bcalc->SetConfidenceLevel(0.95);
512 bcalc->SetLeftSideTailFraction(0.0);
513
514 delete sInt;
515 sInt = bcalc->GetInterval();
516 upper_limit = sInt->UpperLimit();
517 delete sInt;
518 sInt = 0;
519
520 }
521 else{
522 std::cout << "[roostats_cl95]: method " << method
523 << "is not implemented, exiting" <<std::endl;
524 return -1.0;
525 }
526
527 return upper_limit;
528
529 }
530
531
532 Double_t CL95Calc::cla( Double_t ilum, Double_t slum,
533 Double_t eff, Double_t seff,
534 Double_t bck, Double_t sbck,
535 Int_t nuisanceModel,
536 std::string method ){
537
538 makeWorkspace( ilum, slum,
539 eff, seff,
540 bck, sbck,
541 kFALSE,
542 nuisanceModel );
543
544 Double_t CL95A = 0, precision = 1.e-4;
545
546 Int_t i;
547 for (i = bck; i >= 0; i--)
548 {
549 makeData( i );
550 //
551 Double_t s95 = cl95( method );
552 Double_t s95w =s95*TMath::Poisson( (Double_t)i, bck );
553 CL95A += s95w;
554 cout << "[roostats_cla]: n = " << i << "; 95% C.L. = " << s95 << " pb; weighted 95% C.L. = " << s95w << " pb; running <s95> = " << CL95A << " pb" << endl;
555 //
556 if (s95w < CL95A*precision) break;
557 }
558 cout << "[roostats_cla]: Lower bound on n has been found at " << i+1 << endl;
559 //
560 for (i = bck+1; ; i++)
561 {
562 makeData( i );
563 Double_t s95 = cl95( method );
564 Double_t s95w =s95*TMath::Poisson( (Double_t)i, bck );
565 CL95A += s95w;
566 cout << "[roostats_cla]: n = " << i << "; 95% C.L. = " << s95 << " pb; weighted 95% C.L. = " << s95w << " pb; running <s95> = " << CL95A << " pb" << endl;
567 //
568 if (s95w < CL95A*precision) break;
569 }
570 cout << "[roostats_cla]: Upper bound on n has been found at " << i << endl;
571 //
572 cout << "[roostats_cla]: Average upper 95% C.L. limit = " << CL95A << " pb" << endl;
573
574 return CL95A;
575 }
576
577
578
579 LimitResult CL95Calc::clm( Double_t ilum, Double_t slum,
580 Double_t eff, Double_t seff,
581 Double_t bck, Double_t sbck,
582 Int_t nit, Int_t nuisanceModel,
583 std::string method ){
584
585 makeWorkspace( ilum, slum,
586 eff, seff,
587 bck, sbck,
588 kFALSE,
589 nuisanceModel );
590
591 Double_t CLM = 0.0;
592 LimitResult _result;
593
594 Double_t b68[2] = {0.0, 0.0}; // 1-sigma expected band
595 Double_t b95[2] = {0.0, 0.0}; // 2-sigma expected band
596
597 std::vector<Double_t> pe;
598
599 // timer
600 TStopwatch t;
601 t.Start(); // start timer
602 Double_t _realtime = 0.0;
603 Double_t _cputime = 0.0;
604 Double_t _realtime_last = 0.0;
605 Double_t _cputime_last = 0.0;
606 Double_t _realtime_average = 0.0;
607 Double_t _cputime_average = 0.0;
608
609 // throw pseudoexperiments
610 if (nit <= 0)return _result;
611 for (Int_t i = 0; i < nit; i++)
612 {
613 // throw random nuisance parameter (bkg yield)
614 if (nuisanceModel == 0){ // gaussian model for nuisance parameters
615 RooRealVar * _nuis = ws->var("nbkg_sigma");
616 if (_nuis){
617 _nuis->setVal(sbck);
618 }
619 else{ // nuisance model is misconfigured - fail
620 std::cout << "[roostats_clm]: nsig_sigma missing, model misconfigured, exiting..." << std::endl;
621 exit(-1);
622 }
623 }
624 // FIXME: add non-Gaussian nuisance
625 Double_t bmean = GetRandom("syst_nbkg", "nbkg");
626
627 Int_t n = r.Poisson(bmean);
628 makeData( n );
629 std::cout << "Invoking CL95 with bmean = " << bmean << "; n = " << n << std::endl;
630
631 Double_t _pe = cl95( method );
632 pe.push_back(_pe);
633 CLM += pe[i];
634
635 _realtime_last = t.RealTime() - _realtime;
636 _cputime_last = t.CpuTime() - _cputime;
637 _realtime = t.RealTime();
638 _cputime = t.CpuTime();
639 t.Continue();
640 _realtime_average = _realtime/((Double_t)(i+1));
641 _cputime_average = _cputime/((Double_t)(i+1));
642
643 std::cout << "n = " << n << "; 95% C.L. = " << _pe << " pb; running <s95> = " << CLM/(i+1.) << std::endl;
644 std::cout << "Real time (s), this iteration: " << _realtime_last << ", average per iteration: " << _realtime_average << ", total: " << _realtime << std::endl;
645 std::cout << "CPU time (s), this iteration: " << _cputime_last << ", average per iteration: " << _cputime_average << ", total: " << _cputime << std::endl << std::endl;
646 }
647
648 CLM /= nit;
649
650 // sort the vector with limits
651 std::sort(pe.begin(), pe.end());
652
653 std::cout << pe.size() << std::endl;
654
655 // median for the expected limit
656 Double_t _median = TMath::Median(nit, &pe[0]);
657
658 // quantiles for the expected limit bands
659 Double_t _prob[4]; // array with quantile boundaries
660 _prob[0] = 0.021;
661 _prob[1] = 0.159;
662 _prob[2] = 0.841;
663 _prob[3] = 0.979;
664
665 Double_t _quantiles[4]; // array for the results
666
667 TMath::Quantiles(nit, 4, &pe[0], _quantiles, _prob); // evaluate quantiles
668
669 b68[0] = _quantiles[1];
670 b68[1] = _quantiles[2];
671 b95[0] = _quantiles[0];
672 b95[1] = _quantiles[3];
673
674 // let's get actual coverages now
675
676 // sort the vector with limits
677 std::sort(pe.begin(), pe.end());
678 Long64_t lc68 = TMath::BinarySearch(nit, &pe[0], _quantiles[1]) + 1;
679 Long64_t uc68 = nit - TMath::BinarySearch(nit, &pe[0], _quantiles[2]) - 1;
680 Long64_t lc95 = TMath::BinarySearch(nit, &pe[0], _quantiles[0]) + 1;
681 Long64_t uc95 = nit - TMath::BinarySearch(nit, &pe[0], _quantiles[3]) - 1;
682
683 Double_t _cover68 = (nit - lc68 - uc68)*100./nit;
684 Double_t _cover95 = (nit - lc95 - uc95)*100./nit;
685
686 std::cout << "[CL95Calc::clm()]: median limit: " << _median << std::endl;
687 std::cout << "[CL95Calc::clm()]: 1 sigma band: [" << b68[0] << "," << b68[1] << "]; actual coverage: " << _cover68 << "%; lower/upper percentile: " << lc68*100./nit <<"/" << uc68*100./nit << std::endl;
688 std::cout << "[CL95Calc::clm()]: 2 sigma band: [" << b95[0] << "," << b95[1] << "]; actual coverage: " << _cover95 << "%; lower/upper percentile: " << lc95*100./nit <<"/" << uc95*100./nit << std::endl;
689
690 t.Print();
691
692 _result._expected_limit = _median;
693 _result._low68 = b68[0];
694 _result._high68 = b68[1];
695 _result._low95 = b95[0];
696 _result._high95 = b95[1];
697 _result._cover68 = _cover68;
698 _result._cover95 = _cover95;
699
700 return _result;
701 }
702
703
704
705 int CL95Calc::makePlot( std::string method,
706 std::string plotFileName ){
707
708 if (method.find("bayesian") != std::string::npos){
709
710 std::cout << "[roostats_cl95]: making Bayesian posterior plot" << endl;
711
712 TCanvas c1("posterior");
713 bcalc->SetScanOfPosterior(100);
714 RooPlot * plot = bcalc->GetPosteriorPlot();
715 plot->Draw();
716 c1.SaveAs(plotFileName.c_str());
717 }
718 else{
719 std::cout << "[roostats_cl95]: method " << method
720 << "is not implemented, exiting" <<std::endl;
721 return -1;
722 }
723
724 return 0;
725 }
726
727
728
729 Double_t CL95Calc::GetRandom( std::string pdf, std::string var ){
730 //
731 // generates a random number using a pdf in the workspace
732 //
733
734 // generate a dataset with one entry
735 RooDataSet * _ds = ws->pdf(pdf.c_str())->generate(*ws->var(var.c_str()), 1);
736
737 Double_t _result = ((RooRealVar *)(_ds->get(0)->first()))->getVal();
738 delete _ds;
739
740 return _result;
741 }
742
743
744
745 Int_t banner(){
746 //#define __ROOFIT_NOBANNER // banner temporary off
747 #ifndef __EXOST_NOBANNER
748 std::cout << desc << std::endl;
749 #endif
750 return 0 ;
751 }
752
753 static Int_t dummy_ = banner() ;
754
755
756
757 Double_t roostats_cl95(Double_t ilum, Double_t slum,
758 Double_t eff, Double_t seff,
759 Double_t bck, Double_t sbck,
760 Int_t n,
761 Bool_t gauss,
762 Int_t nuisanceModel,
763 std::string method,
764 std::string plotFileName){
765
766 std::cout << "[roostats_cl95]: estimating 95% C.L. upper limit" << endl;
767 if (method.find("bayesian") != std::string::npos){
768 std::cout << "[roostats_cl95]: using Bayesian calculation via numeric integration" << endl;
769 }
770 else if (method.find("workspace") != std::string::npos){
771 std::cout << "[roostats_cl95]: no interval calculation, only create and save workspace" << endl;
772 }
773 else{
774 std::cout << "[roostats_cl95]: method " << method
775 << "is not implemented, exiting" <<std::endl;
776 return -1.0;
777 }
778
779 // some input validation
780 if (n < 0){
781 std::cout << "Negative observed number of events specified, exiting" << std::endl;
782 return -1.0;
783 }
784
785 if (n == 0) gauss = kFALSE;
786
787 if (gauss){
788 nuisanceModel = 0;
789 std::cout << "[roostats_cl95]: Gaussian statistics used" << endl;
790 }
791 else{
792 std::cout << "[roostats_cl95]: Poisson statistics used" << endl;
793 }
794
795 // limit calculation
796 CL95Calc theCalc;
797 RooWorkspace * ws = theCalc.makeWorkspace( ilum, slum,
798 eff, seff,
799 bck, sbck,
800 gauss,
801 nuisanceModel );
802 RooDataSet * data = (RooDataSet *)( theCalc.makeData( n )->Clone() );
803 data->SetName("observed_data");
804 ws->import(*data);
805
806 ws->Print();
807
808 ws->SaveAs("ws.root");
809
810 // if only workspace requested, exit here
811 if ( method.find("workspace") != std::string::npos ) return 0.0;
812
813 std::cout << "[roostats_cl95]: Range of allowed cross section values: ["
814 << ws->var("xsec")->getMin() << ", "
815 << ws->var("xsec")->getMax() << "]" << std::endl;
816 Double_t limit = theCalc.cl95( method );
817 std::cout << "[roostats_cl95]: 95% C.L. upper limit: " << limit << std::endl;
818
819 // check if the plot is requested
820 if (plotFileName.size() != 0){
821 theCalc.makePlot(method, plotFileName);
822 }
823
824 return limit;
825 }
826
827
828 Double_t roostats_cla(Double_t ilum, Double_t slum,
829 Double_t eff, Double_t seff,
830 Double_t bck, Double_t sbck,
831 Int_t nuisanceModel,
832 std::string method){
833
834 Double_t limit = -1.0;
835
836 std::cout << "[roostats_cla]: estimating average 95% C.L. upper limit" << endl;
837 if (method.find("bayesian") != std::string::npos){
838 std::cout << "[roostats_cla]: using Bayesian calculation via numeric integration" << endl;
839 }
840 else{
841 std::cout << "[roostats_cla]: method " << method
842 << "is not implemented, exiting" <<std::endl;
843 return -1.0;
844 }
845
846 std::cout << "[roostats_cla]: Poisson statistics used" << endl;
847
848 CL95Calc theCalc;
849 limit = theCalc.cla( ilum, slum,
850 eff, seff,
851 bck, sbck,
852 nuisanceModel,
853 method );
854
855 //std::cout << "[roostats_cla]: average 95% C.L. upper limit: " << limit << std::endl;
856
857 return limit;
858 }
859
860
861
862 LimitResult roostats_clm(Double_t ilum, Double_t slum,
863 Double_t eff, Double_t seff,
864 Double_t bck, Double_t sbck,
865 Int_t nit, Int_t nuisanceModel,
866 std::string method){
867
868 //Double_t limit = -1.0;
869 LimitResult limit;
870
871 std::cout << "[roostats_clm]: estimating average 95% C.L. upper limit" << endl;
872 if (method.find("bayesian") != std::string::npos){
873 std::cout << "[roostats_clm]: using Bayesian calculation via numeric integration" << endl;
874 }
875 else{
876 std::cout << "[roostats_clm]: method " << method
877 << "is not implemented, exiting" <<std::endl;
878 //return -1.0;
879 return limit;
880 }
881
882 std::cout << "[roostats_clm]: Poisson statistics used" << endl;
883
884 CL95Calc theCalc;
885 limit = theCalc.clm( ilum, slum,
886 eff, seff,
887 bck, sbck,
888 nit, nuisanceModel,
889 method );
890
891 return limit;
892 }