Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
xRooHypoSpace.cxx
Go to the documentation of this file.
1/*
2 * Project: xRooFit
3 * Author:
4 * Will Buttinger, RAL 2022
5 *
6 * Copyright (c) 2022, CERN
7 *
8 * Redistribution and use in source and binary forms,
9 * with or without modification, are permitted according to the terms
10 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
11 */
12
14
15#include "RooArgSet.h"
16#include "RooArgList.h"
17#include "RooRealVar.h"
18#include "RooFitResult.h"
19#include "RooConstVar.h"
20#include "RooAbsPdf.h"
21
22#include "TCanvas.h"
23#include "TStopwatch.h"
24#include "TSystem.h"
25#include "TPRegexp.h"
26#include "TMemFile.h"
27#include "TROOT.h"
28#include "RooDataSet.h"
29#include "TKey.h"
30#include "TFile.h"
31#include "TGraphErrors.h"
32#include "TMultiGraph.h"
33#include "TH1F.h"
34#include "TStyle.h"
35#include "TLegend.h"
36#include "TLine.h"
38#include "TEnv.h"
39
41
42xRooNLLVar::xRooHypoSpace::xRooHypoSpace(const char *name, const char *title)
43 : TNamed(name, title), fPars(std::make_shared<RooArgSet>())
44{
45 if (name == nullptr || strlen(name) == 0) {
46 SetName(TUUID().AsString());
47 }
48}
49
51 : fPars(std::make_shared<RooArgSet>())
52{
53 if (!result)
54 return;
55
57
58 fPars->addClone(*std::unique_ptr<RooAbsCollection>(result->GetParameters()));
59 double spaceSize = 1;
60 for (auto p : *fPars) {
61 auto v = dynamic_cast<RooRealVar *>(p);
62 if (!v)
63 continue;
64 spaceSize *= (v->getMax() - v->getMin());
65 }
66 for (int i = 0; i < result->ArraySize(); i++) {
67 auto point = result->GetResult(i);
68 double xVal = result->GetXValue(i);
69 double ssize = spaceSize;
70 for (auto p : *fPars) {
71 auto v = dynamic_cast<RooRealVar *>(p);
72 if (!v)
73 continue;
74 ssize /= (v->getMax() - v->getMin());
75 double remain = std::fmod(xVal, ssize);
76 v->setVal((xVal - remain) / ssize);
77 xVal = remain;
78 }
79 emplace_back(xRooHypoPoint(std::make_shared<RooStats::HypoTestResult>(*point), fPars.get()));
80 }
81 // add any pars we might have missed
82 for (auto &p : *this) {
83 for (auto a : *p.coords) {
84 if (!fPars->find(a->GetName()))
85 fPars->addClone(*a);
86 }
87 }
88}
89
90std::shared_ptr<xRooNode> xRooNLLVar::xRooHypoSpace::pdf(const char *parValues) const
91{
92 return pdf(toArgs(parValues));
93}
94
95std::shared_ptr<xRooNode> xRooNLLVar::xRooHypoSpace::pdf(const RooAbsCollection &parValues) const
96{
98 rhs.add(parValues);
99 rhs.sort();
100
101 std::shared_ptr<xRooNode> out = nullptr;
102
103 for (auto &[_range, _pdf] : fPdfs) {
104 // any pars not in rhs are assumed to have infinite range in rhs
105 // and vice versa
106 bool collision = true;
107 for (auto &_lhs : *_range) {
108 auto _rhs = rhs.find(*_lhs);
109 if (!_rhs)
110 continue;
111 if (auto v = dynamic_cast<RooRealVar *>(_rhs); v) {
112 if (auto v2 = dynamic_cast<RooRealVar *>(_lhs)) {
113 if (!(v->getMin() <= v2->getMax() && v2->getMin() <= v->getMax())) {
114 collision = false;
115 break;
116 }
117 } else if (auto c2 = dynamic_cast<RooConstVar *>(_lhs)) {
118 if (!(v->getMin() <= c2->getVal() && c2->getVal() <= v->getMax())) {
119 collision = false;
120 break;
121 }
122 }
123 } else if (auto c = dynamic_cast<RooConstVar *>(_rhs); c) {
124 if (auto v2 = dynamic_cast<RooRealVar *>(_lhs)) {
125 if (!(c->getVal() <= v2->getMax() && v2->getMin() <= c->getVal())) {
126 collision = false;
127 break;
128 }
129 } else if (auto c2 = dynamic_cast<RooConstVar *>(_lhs)) {
130 if (!(c->getVal() == c2->getVal())) {
131 collision = false;
132 break;
133 }
134 }
135 }
136 }
137 if (collision) {
138 if (out) {
139 throw std::runtime_error("Multiple pdf possibilities");
140 }
141 out = _pdf;
142 }
143 }
144
145 return out;
146}
147
149{
150
151 RooArgList out;
152
153 TStringToken pattern(str, ";");
154 while (pattern.NextToken()) {
155 TString s = pattern;
156 // split by "=" sign
157 auto _idx = s.Index('=');
158 if (_idx == -1)
159 continue;
160 TString _name = s(0, _idx);
161 TString _val = s(_idx + 1, s.Length());
162
163 if (_val.IsFloat()) {
164 out.addClone(RooConstVar(_name, _name, _val.Atof()));
165 } else if (_val.BeginsWith('[')) {
166 _idx = _val.Index(',');
167 if (_idx == -1)
168 continue;
169 TString _min = _val(0, _idx);
170 TString _max = _val(_idx + 1, _val.Length() - _idx - 2);
171 out.addClone(RooRealVar(_name, _name, _min.Atof(), _max.Atof()));
172 }
173 }
174
175 return out;
176}
177
178int xRooNLLVar::xRooHypoSpace::AddPoints(const char *parName, size_t nPoints, double low, double high)
179{
180 if (nPoints == 0)
181 return nPoints;
182
183 auto _par = dynamic_cast<RooAbsRealLValue *>(fPars->find(parName));
184 if (!_par)
185 throw std::runtime_error("Unknown parameter");
186 _par->setAttribute("axis");
187
188 if (low < _par->getMin()) {
189 Warning("AddPoints", "low edge of hypoSpace %g below lower bound of parameter: %g. Changing to lower bound", low,
190 _par->getMin());
191 low = _par->getMin();
192 }
193 if (high > _par->getMax()) {
194 Warning("AddPoints", "high edge of hypoSpace %g above upper bound of parameter: %g. Changing to upper bound",
195 high, _par->getMax());
196 high = _par->getMax();
197 }
198
199 if (nPoints == 1) {
200 _par->setVal((high + low) * 0.5);
201 AddPoint();
202 return nPoints;
203 }
204
205 double step = (high - low) / (nPoints - 1);
206 if (step <= 0)
207 throw std::runtime_error("Invalid steps");
208
209 for (size_t i = 0; i < nPoints; i++) {
210 _par->setVal((i == nPoints - 1) ? high : (low + step * i));
211 AddPoint();
212 }
213 return nPoints;
214}
215
217{
218 if (axes().empty()) {
219 // set the first poi as the axis variable to scan
220 if (poi().empty()) {
221 throw std::runtime_error("No POI to scan");
222 } else {
223 poi().first()->setAttribute("axis");
224 }
225 }
226
227 if (empty()) {
228 // promote all axes to being poi and demote all non-axes to non-poi
229 poi().setAttribAll("poi", false);
230 axes().setAttribAll("poi");
231 }
232
233 return AddPoint(TString::Format("%s=%f", axes().first()->GetName(), value));
234}
235
236int xRooNLLVar::xRooHypoSpace::scan(const char *type, size_t nPoints, double low, double high,
237 const std::vector<double> &nSigmas, double relUncert)
238{
239
241 sType.ToLower();
242 if (sType.Contains("cls") && !sType.Contains("pcls"))
243 sType.ReplaceAll("cls", "pcls");
244 if (!sType.Contains("pcls") && !sType.Contains("ts") && !sType.Contains("pnull") && !sType.Contains("plr")) {
245 throw std::runtime_error("scan type must be equal to one of: plr, cls, ts, pnull");
246 }
247
248 // will scan the first axes variable ... if there is none, specify the first poi as the axis var
249 if (axes().empty()) {
250 // set the first poi as the axis variable to scan
251 if (poi().empty()) {
252 throw std::runtime_error("No POI to scan");
253 } else {
254 poi().first()->setAttribute("axis");
255 }
256 }
257
258 // promote all axes to being poi and demote all non-axes to non-poi
259 poi().setAttribAll("poi", false);
260 axes().setAttribAll("poi");
261
262 auto p = dynamic_cast<RooRealVar *>(axes().first());
263 if (!p) {
264 throw std::runtime_error(TString::Format("%s not scannable", axes().first()->GetName()));
265 }
266
267 if (sType.Contains("cls")) {
268 if (empty() && relUncert == std::numeric_limits<double>::infinity()) {
269 // use default uncertainty precision of 10%
270 ::Info("xRooHypoSpace::scan", "Using default precision of 10%% for auto-scan");
271 relUncert = 0.1;
272 }
273 for (auto a : axes()) {
274 if (!a->hasRange("physical")) {
275 ::Info("xRooHypoSpace::scan", "No physical range set for %s, setting to [0,inf]", p->GetName());
276 dynamic_cast<RooRealVar *>(a)->setRange("physical", 0, std::numeric_limits<double>::infinity());
277 }
278 if (!a->getStringAttribute("altVal") || !strlen(p->getStringAttribute("altVal"))) {
279 ::Info("xRooHypoSpace::scan", "No altVal set for %s, setting to 0", a->GetName());
280 a->setStringAttribute("altVal", "0");
281 }
282 // ensure range straddles altVal
283 double altVal = TString(a->getStringAttribute("altVal")).Atof();
284 auto v = dynamic_cast<RooRealVar *>(a);
285 if (v->getMin() >= altVal) {
286 ::Info("xRooHypoSpace::scan", "range of POI does not straddle alt value, adjusting minimum to %g",
287 altVal - 1e-5);
288 v->setMin(altVal - 1e-5);
289 }
290 if (v->getMax() <= altVal) {
291 ::Info("xRooHypoSpace::scan", "range of POI does not straddle alt value, adjusting maximum to %g",
292 altVal + 1e-5);
293 v->setMax(altVal + 1e-5);
294 }
295 for (auto &[pdf, nll] : fNlls) {
296 if (auto _v = dynamic_cast<RooRealVar *>(nll->pars()->find(*a))) {
297 _v->setRange(v->getMin(), v->getMax());
298 }
299 }
300 }
301 } else if (sType.Contains("plr")) {
302 // force use of two-sided test statistic for any new points
303 fTestStatType = xRooFit::Asymptotics::TwoSided;
304 sType.ReplaceAll("plr", "ts");
305 } else if (sType.Contains("pnull") && fTestStatType == xRooFit::Asymptotics::Unknown) {
306 // for pnull (aka discovery) "scan" (may just be a single point) default to use of
307 // uncapped test stat
308 fTestStatType = xRooFit::Asymptotics::Uncapped;
309 // and ensure altVal is set
310 for (auto a : axes()) {
311 if (!a->getStringAttribute("altVal") || !strlen(p->getStringAttribute("altVal"))) {
312 ::Info("xRooHypoSpace::scan", "No altVal set for %s, setting to 1", a->GetName());
313 a->setStringAttribute("altVal", "1");
314 }
315 }
316 }
317
318 if (/*high < low ||*/ (high == low && nPoints != 1)) {
319 // take from parameter
320 low = p->getMin("scan");
321 high = p->getMax("scan");
322 }
323 if (!std::isnan(low) && !std::isnan(high) && !(std::isinf(low) && std::isinf(high))) {
324 p->setRange("scan", std::min(low, high), std::max(low, high));
325 }
326 if (p->hasRange("scan")) {
327 ::Info("xRooHypoSpace::scan", "Using %s scan range: %g - %g", p->GetName(), p->getMin("scan"), p->getMax("scan"));
328 }
329
330 bool doObs = false;
331 for (auto nSigma : nSigmas) {
332 if (std::isnan(nSigma)) {
333 doObs = true;
334 break;
335 }
336 }
337
338 if (fNlls.empty()) {
339 // this happens when loaded hypoSpace from a hypoSpaceInverterResult
340 // set relUncert to infinity so that we don't test any new points
341 relUncert = std::numeric_limits<double>::infinity(); // no NLL available so just get whatever limit we can
342
343 // if any of the defined points are 'expected' data don't do obs
344 // for(auto& hp : *this) {
345 // if(hp.isExpected) {
346 // doObs = false; break;
347 // }
348 // }
349 } else {
350#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 26, 00)
351// bool allGen = true;
352// for (auto &[pdf, nll]: fNlls) {
353// auto _d = dynamic_cast<RooDataSet *>(nll->data());
354// if (!_d || !_d->weightVar() || !_d->weightVar()->getStringAttribute("fitResult") ||
355// !_d->weightVar()->getAttribute("expected")) {
356// allGen = false;
357// break;
358// }
359// }
360// if (allGen)
361// doObs = false;
362#endif
363 }
364
365 // create a fitDatabase if required
367 if (fFitDb) {
368 // move to the db, and unlock it if this is a TMemFile
369 fFitDb->cd();
370 if (auto myDb = dynamic_cast<TMemFile *>(fFitDb.get())) {
371 // need to unlock the database
372 myDb->SetWritable(true);
373 }
374 }
375 if (!gDirectory || !gDirectory->IsWritable()) {
376 // locate a TMemFile in the open list of files and move to that
377 // or create one if cannot find
378 /*for (auto file : *gROOT->GetListOfFiles()) {
379 if (auto f = dynamic_cast<TMemFile *>(file)) {
380 f->cd();
381 break;
382 }
383 }
384 if (!gDirectory || !gDirectory->IsWritable()) {
385 new TMemFile("fitDatabase", "RECREATE");
386 }*/
387 // now we create a TMemFile of our own, so that we don't get in the way of other hypoSpaces
388 fFitDb = std::shared_ptr<TDirectory>(
389 new TMemFile(TString::Format("fitDatabase_%s", TUUID().AsString()), "RECREATE"), [](TDirectory *) {});
390 // db can last longer than the hypoSpace, so that the fits are fully available in the browser
391 // if a scan was initiated through the browser. If user wants to cleanup they can do manually
392 // through root's GetListOfFiles()
393 // would like to clean it up ourself when the hypoSpace is destroyed, but would need way to keep alive for the
394 // browser
395 }
396
397 int out = 0;
398
399 // enable visualizing by default if scanning in non-batch mode
400 if (!gROOT->IsBatch() && !sType.Contains("visualize"))
401 sType += " visualize";
402
403 if (nPoints == 0) {
404 // automatic scan
405 if (sType.Contains("cls")) {
406 for (double nSigma : nSigmas) {
407 xValueWithError res(std::make_pair(0., 0.));
408 if (std::isnan(nSigma)) {
409 if (!doObs)
410 continue;
411 res = findlimit(TString::Format("%s obs", sType.Data()), relUncert);
412 } else {
413 res =
414 findlimit(TString::Format("%s exp%s%d", sType.Data(), nSigma > 0 ? "+" : "", int(nSigma)), relUncert);
415 }
416 if (std::isnan(res.first) || std::isnan(res.second)) {
417 out = 1;
418 } else if (std::isinf(res.second)) {
419 out = 2;
420 }
421 }
422 } else if (sType.Contains("ts")) {
423 // add a point at the current poi value
424 // use that to determine ufit value and use uncertainty to set an approx range
425 AddPoint(TString::Format("%s=%g", poi().first()->GetName(), poi().getRealValue(poi().first()->GetName())));
426 graphs(sType); // triggers computation
427 if (back().status() != 0) {
428 out += 1;
429 } else {
430 // fit was ok, so use the values to determine an appropriate range
431 low = std::max(low, back().mu_hat().getVal()-back().mu_hat().getError()*3);
432 high = std::min(high, back().mu_hat().getVal()+back().mu_hat().getError()*3);
433 nPoints = 20;
434 double step = (high - low) / (nPoints - 1);
435 for (size_t i = 0; i < nPoints; i++) {
436 AddPoint(TString::Format("%s=%g", poi().first()->GetName(), low + step * i));
437 graphs(sType); // triggers computation
438 if (back().status() != 0)
439 out += 1;
440 }
441
442 }
443 } else {
444 throw std::runtime_error(TString::Format("Automatic scanning not yet supported for %s", type));
445 }
446 } else {
447 // add the required points and then compute the required value
448 if (nPoints == 1) {
449 AddPoint(TString::Format("%s=%g", poi().first()->GetName(), (high + low) / 2.));
450 graphs(sType); // triggers computation
451 if (back().status() != 0)
452 out += 1;
453 } else {
454 double step = (high - low) / (nPoints - 1);
455 for (size_t i = 0; i < nPoints; i++) {
456 AddPoint(TString::Format("%s=%g", poi().first()->GetName(), low + step * i));
457 graphs(sType); // triggers computation
458 if (back().status() != 0)
459 out += 1;
460 }
461 }
462 }
463
464 if (origDir)
465 origDir->cd();
466
467 if (auto myDb = dynamic_cast<TMemFile *>(fFitDb.get())) {
468 // need to lock the database, because if its writable when pyroot closes it causes a crash
469 myDb->SetWritable(false);
470 }
471 return out;
472}
473
474std::map<std::string, xRooNLLVar::xValueWithError>
475xRooNLLVar::xRooHypoSpace::limits(const char *opt, const std::vector<double> &nSigmas, double relUncert)
476{
477
478 if (fNlls.empty()) {
479 // this happens when loaded hypoSpace from a hypoSpaceInverterResult
480 // set relUncert to infinity so that we don't test any new points
481 relUncert = std::numeric_limits<double>::infinity(); // no NLL available so just get whatever limit we can
482 }
483
484 scan(opt, nSigmas, relUncert);
485
486 std::map<std::string, xRooNLLVar::xValueWithError> out;
487 for (auto nSigma : nSigmas) {
488 auto lim = limit(opt, nSigma);
489 if (lim.second < 0)
490 lim.second = -lim.second; // make errors positive for this method
491 out[std::isnan(nSigma) ? "obs" : TString::Format("%d", int(nSigma)).Data()] = xRooFit::matchPrecision(lim);
492 }
493 return out;
494}
495
497{
498 // move to given coords, if any ... will mark them const too
499 std::unique_ptr<RooAbsCollection, std::function<void(RooAbsCollection *)>> _snap(fPars->snapshot(),
500 [&](RooAbsCollection *c) {
501 *fPars = *c;
502 delete c;
503 });
504 TStringToken pattern(coords, ";");
505 while (pattern.NextToken()) {
506 TString s = pattern;
507 // split by "=" sign
508 auto _idx = s.Index('=');
509 if (_idx == -1)
510 continue;
511 TString _name = s(0, _idx);
512 TString _val = s(_idx + 1, s.Length());
513 auto _v = dynamic_cast<RooRealVar *>(fPars->find(_name));
514 if (!_v)
515 continue;
516
517 if (_val.IsFloat()) {
518 _v->setConstant();
519 _v->setVal(_val.Atof());
520 }
521 }
522
523 auto _pdf = pdf();
524
525 if (!_pdf)
526 throw std::runtime_error("no model at coordinates");
527
528 // if (std::unique_ptr<RooAbsCollection>(fPars->selectByAttrib("poi", true))->size() == 0) {
529 // throw std::runtime_error(
530 // "No pars designated as POI - set with pars()->find(<parName>)->setAttribute(\"poi\",true)");
531 // }
532
533 if (fNlls.find(_pdf) == fNlls.end()) {
534 fNlls[_pdf] = std::make_shared<xRooNLLVar>(_pdf->nll("" /*TODO:allow change dataset name and nll opts*/, {}));
535 }
536
537 xRooHypoPoint out;
538
539 out.nllVar = fNlls[_pdf];
540 out.fData = fNlls[_pdf]->getData();
541 out.isExpected = dynamic_cast<RooDataSet *>(out.fData.first.get()) &&
542 dynamic_cast<RooDataSet *>(out.fData.first.get())->weightVar()->getAttribute("expected");
543 // TODO: need to access the genfit of the data and add that to the point, somehow ...
544
545 out.coords.reset(fPars->snapshot()); // should already have altVal prop on poi, and poi labelled
546 // ensure all poi are marked const ... required by xRooHypoPoint behaviour
547 out.poi().setAttribAll("Constant");
548 // and now remove anything that's marked floating
549
550 // do to bug in remove have to ensure not using the hash map otherwise will be doing an invalid read after the
551 // deletion of the owned pars
552 const_cast<RooAbsCollection *>(out.coords.get())->useHashMapForFind(false);
553 const_cast<RooAbsCollection *>(out.coords.get())
554 ->remove(*std::unique_ptr<RooAbsCollection>(out.coords->selectByAttrib("Constant", false)), true, true);
555
556 double value = out.fNullVal();
557 double alt_value = out.fAltVal();
558
559 auto _type = fTestStatType;
560 if (_type == xRooFit::Asymptotics::Unknown) {
561 // decide based on values
562 if (std::isnan(alt_value)) {
564 } else if (value >= alt_value) {
566 } else {
568 }
569 }
570
571 out.fPllType = _type;
572
573 // look for a matching point
574 for (auto &p : *this) {
575 if (p.nllVar != out.nllVar)
576 continue;
577 if (p.fData != out.fData)
578 continue;
579 if (!p.alt_poi().equals(out.alt_poi()))
580 continue;
581 bool match = true;
582 for (auto c : p.alt_poi()) {
583 if (auto v = dynamic_cast<RooAbsReal *>(c);
584 v && std::abs(v->getVal() - out.alt_poi().getRealValue(v->GetName())) > 1e-12) {
585 match = false;
586 break;
587 } else if (auto cat = dynamic_cast<RooAbsCategory *>(c);
588 cat && cat->getCurrentIndex() ==
589 out.alt_poi().getCatIndex(cat->GetName(), std::numeric_limits<int>().max())) {
590 match = false;
591 break;
592 }
593 }
594 if (!match)
595 continue;
596 if (!p.coords->equals(*out.coords))
597 continue;
598 for (auto c : *p.coords) {
599 if (c->getAttribute("poi")) {
600 continue; // check poi below
601 }
602 if (auto v = dynamic_cast<RooAbsReal *>(c);
603 v && std::abs(v->getVal() - out.coords->getRealValue(v->GetName())) > 1e-12) {
604 match = false;
605 break;
606 } else if (auto cat = dynamic_cast<RooAbsCategory *>(c);
607 cat && cat->getCurrentIndex() ==
608 out.alt_poi().getCatIndex(cat->GetName(), std::numeric_limits<int>().max())) {
609 match = false;
610 break;
611 }
612 }
613 if (!match)
614 continue;
615 // if reached here we can copy over the asimov dataset to save re-generating it
616 // first copy over cfit_alt (if its there) because that is also the same since coords and alt_poi values same
617 if (auto cfit = p.cfit_alt(true)) {
618 out.fAlt_cfit = cfit;
619 }
620 if (p.asimov(true) && p.asimov(true)->fData.first && (!out.asimov(true) || !out.asimov(true)->fData.first)) {
621 out.asimov()->fData = p.asimov(true)->fData;
622 }
623 if (!p.poi().equals(out.poi()))
624 continue;
625 for (auto c : p.poi()) {
626 if (auto v = dynamic_cast<RooAbsReal *>(c);
627 v && std::abs(v->getVal() - out.poi().getRealValue(v->GetName())) > 1e-12) {
628 match = false;
629 break;
630 }
631 }
632 if (match) {
633 // found a duplicate point, return that!
634 return p;
635 }
636 }
637
638 std::string coordString;
639 for (auto a : axes()) {
640 coordString += TString::Format("%s=%g", a->GetName(), out.coords->getRealValue(a->GetName()));
641 coordString += ",";
642 }
643 coordString.erase(coordString.end() - 1);
644
645 ::Info("xRooHypoSpace::AddPoint", "Added new point @ %s", coordString.c_str());
646 return emplace_back(out);
647}
648
650{
651
652 if (!_pdf.get<RooAbsPdf>()) {
653 throw std::runtime_error("Not a pdf");
654 }
655
656 auto pars = _pdf.pars().argList();
657
658 // replace any pars with validity pars and add new pars
659 auto vpars = toArgs(validity);
660 pars.replace(vpars);
661 vpars.remove(pars, true, true);
662 pars.add(vpars);
663
664 if (auto existing = pdf(pars)) {
665 throw std::runtime_error(std::string("Clashing model: ") + existing->GetName());
666 }
667
668 auto myPars = std::shared_ptr<RooArgList>(dynamic_cast<RooArgList *>(pars.snapshot()));
669 myPars->sort();
670
671 pars.remove(*fPars, true, true);
672
673 fPars->addClone(pars);
674
675 fPdfs.insert(std::make_pair(myPars, std::make_shared<xRooNode>(_pdf)));
676
677 return true;
678}
679
681{
682 // determine which pars are the minimal set to distinguish all points in the space
683 RooArgList out;
684 out.setName("axes");
685
686 out.add(*std::unique_ptr<RooAbsCollection>(
687 fPars->selectByAttrib("axis", true))); // start with any pars explicitly designated as axes
688
689 bool clash;
690 do {
691 clash = false;
692
693 std::set<std::vector<double>> coords;
694 for (auto &p : *this) {
695 std::vector<double> p_coords;
696 for (auto o : out) {
697 auto _v = dynamic_cast<RooRealVar *>(p.coords->find(o->GetName()));
698 p_coords.push_back(
699 (_v && _v->isConstant())
700 ? _v->getVal()
701 : std::numeric_limits<double>::infinity()); // non-const coords are treating as non-existent
702 // p_coords.push_back(p.coords->getRealValue(o->GetName(), std::numeric_limits<double>::quiet_NaN()));
703 }
704 if (coords.find(p_coords) != coords.end()) {
705 clash = true;
706 break;
707 }
708 coords.insert(p_coords);
709 }
710
711 if (clash) {
712 // add next best coordinate
713 std::map<std::string, std::unordered_set<double>> values;
714 for (auto &par : *pars()) {
715 if (out.find(*par))
716 continue;
717 for (auto p : *this) {
718 auto _v = dynamic_cast<RooRealVar *>(p.coords->find(par->GetName()));
719 values[par->GetName()].insert(
720 (_v && _v->isConstant())
721 ? _v->getVal()
722 : std::numeric_limits<double>::infinity()); // non-const coords are treating as non-existent
723 // values[par->GetName()].insert(
724 // p.coords->getRealValue(par->GetName(), std::numeric_limits<double>::quiet_NaN()));
725 }
726 }
727
728 std::string bestVar;
729 size_t maxDiff = 0;
730 bool isPOI = false;
731 for (auto &[k, v] : values) {
732 if (v.size() > maxDiff || (v.size() == maxDiff && !isPOI && pars()->find(k.c_str())->getAttribute("poi"))) {
733 bestVar = k;
734 isPOI = pars()->find(k.c_str())->getAttribute("poi");
735 maxDiff = std::max(maxDiff, v.size());
736 }
737 }
738 if (bestVar.empty()) {
739 break;
740 }
741
742 out.add(*pars()->find(bestVar.c_str()));
743 }
744 } while (clash);
745
746 // ensure poi are at the end
747 std::unique_ptr<RooAbsCollection> poi(out.selectByAttrib("poi", true));
748 out.remove(*poi);
749 out.add(*poi);
750
751 return out;
752}
753
755{
756 RooArgList out;
757 out.setName("poi");
758 out.add(*std::unique_ptr<RooAbsCollection>(pars()->selectByAttrib("poi", true)));
759 return out;
760}
761
763{
764
765 auto _axes = axes();
766
767 size_t badFits = 0;
768
769 for (size_t i = 0; i < size(); i++) {
770 std::cout << i << ") ";
771 for (auto a : _axes) {
772 if (a != _axes.first())
773 std::cout << ",";
774 std::cout << a->GetName() << "="
775 << at(i).coords->getRealValue(a->GetName(), std::numeric_limits<double>::quiet_NaN());
776 }
777 std::cout << " status=[ufit:";
778 auto ufit = const_cast<xRooHypoPoint &>(at(i)).ufit(true);
779 if (!ufit) {
780 std::cout << "-";
781 } else {
782 std::cout << ufit->status();
783 badFits += (xRooNLLVar::xRooHypoPoint::allowedStatusCodes.count(ufit->status()) == 0);
784 }
785 std::cout << ",cfit_null:";
786 auto cfit = const_cast<xRooHypoPoint &>(at(i)).cfit_null(true);
787 if (!cfit) {
788 std::cout << "-";
789 } else {
790 std::cout << cfit->status();
791 badFits += (xRooNLLVar::xRooHypoPoint::allowedStatusCodes.count(cfit->status()) == 0);
792 }
793 std::cout << ",cfit_alt:";
794 auto afit = const_cast<xRooHypoPoint &>(at(i)).cfit_alt(true);
795 if (!afit) {
796 std::cout << "-";
797 } else {
798 std::cout << afit->status();
800 }
801 if (auto asiPoint = const_cast<xRooHypoPoint &>(at(i)).asimov(true)) {
802 std::cout << ",asimov.ufit:";
803 auto asi_ufit = asiPoint->ufit(true);
804 if (!asi_ufit) {
805 std::cout << "-";
806 } else {
807 std::cout << asi_ufit->status();
809 }
810 std::cout << ",asimov.cfit_null:";
811 auto asi_cfit = asiPoint->cfit_null(true);
812 if (!asi_cfit) {
813 std::cout << "-";
814 } else {
815 std::cout << asi_cfit->status();
817 }
818 auto cfit_lbound = const_cast<xRooHypoPoint &>(at(i)).cfit_lbound(true);
819 if (!cfit_lbound) {
820 } else {
821 std::cout << ",cfit_lbound:" << cfit_lbound->status();
822 badFits += (xRooNLLVar::xRooHypoPoint::allowedStatusCodes.count(cfit_lbound->status()) == 0);
823 }
824 }
825 std::cout << "]";
826 auto sigma_mu = const_cast<xRooHypoPoint &>(at(i)).sigma_mu(true);
827 if (!std::isnan(sigma_mu.first)) {
828 std::cout << " sigma_mu=" << sigma_mu.first;
829 if (sigma_mu.second)
830 std::cout << " +/- " << sigma_mu.second;
831 }
832 std::cout << std::endl;
833 }
834 std::cout << "--------------------------" << std::endl;
835 std::cout << "Number of bad fits: " << badFits << std::endl;
836}
837
838std::shared_ptr<TGraphErrors> xRooNLLVar::xRooHypoSpace::graph(
839 const char *opt /*, const std::function<void(xRooNLLVar::xRooHypoSpace*)>& progress*/) const
840{
841
842 TString sOpt(opt);
843 sOpt.ToLower();
844
845 bool doCLs = sOpt.Contains("cls");
846 bool readOnly = sOpt.Contains("readonly");
847 bool visualize = sOpt.Contains("visualize") && !readOnly;
848
849 double nSigma =
850 (sOpt.Contains("exp"))
851 ? (TString(sOpt(sOpt.Index("exp") + 3,
852 sOpt.Index(" ", sOpt.Index("exp")) == -1 ? sOpt.Length() : sOpt.Index(" ", sOpt.Index("exp"))))
853 .Atof())
854 : std::numeric_limits<double>::quiet_NaN();
855 bool expBand =
856 !std::isnan(nSigma) && nSigma && !(sOpt(sOpt.Index("exp") + 3) == '+' || sOpt(sOpt.Index("exp") + 3) == '-');
857
858 auto _axes = axes();
859 if (_axes.size() != 1)
860 return nullptr;
861
862 auto out = std::make_shared<TGraphErrors>();
863 out->SetName(GetName());
864 out->SetEditable(false);
865 const char *sCL = (doCLs) ? "CLs" : "null";
866
867 TString title =
868 TString::Format("%s;%s;p_{%s}",
869 (std::isnan(nSigma))
870 ? "Observed"
871 : (!nSigma ? "Expected"
872 : TString::Format("%s%d#sigma Expected",
873 expBand || !nSigma ? "" : ((nSigma < 0) ? "-" : "+"), int(nSigma))
874 .Data()),
875 _axes.at(0)->GetTitle(), sCL);
876
877 if (std::isnan(nSigma)) {
878 out->SetNameTitle(TString::Format("obs_p%s", sCL), title);
879 out->SetMarkerStyle(20);
880 out->SetMarkerSize(0.5);
881 if (sOpt.Contains("ts")) {
882 out->SetNameTitle("obs_ts", TString::Format("Observed;%s;%s", _axes.at(0)->GetTitle(),
883 (empty() ? "" : front().tsTitle(true).Data())));
884 }
885 } else {
886 out->SetNameTitle(TString::Format("exp%d_p%s", int(nSigma), sCL), title);
887 out->SetMarkerStyle(0);
888 out->SetMarkerSize(0);
889 out->SetLineStyle(2 + int(nSigma));
890 if (expBand && nSigma) {
891 out->SetFillColor((nSigma == 2) ? kYellow : kGreen);
892 // out->SetFillStyle(3005);
893 out->SetLineStyle(0);
894 out->SetLineWidth(0);
895 auto x = out->Clone("up");
896 x->SetBit(kCanDelete);
897 dynamic_cast<TAttFill *>(x)->SetFillStyle(nSigma == 2 ? 3005 : 3004);
898 dynamic_cast<TAttFill *>(x)->SetFillColor(kBlack);
899 out->GetListOfFunctions()->Add(x, "F");
900 x = out->Clone("down");
901 x->SetBit(kCanDelete);
902 // dynamic_cast<TAttFill*>(x)->SetFillColor((nSigma==2) ? kYellow : kGreen);
903 // dynamic_cast<TAttFill*>(x)->SetFillStyle(1001);
904 out->GetListOfFunctions()->Add(x, "F");
905 }
906 if (sOpt.Contains("ts")) {
907 out->SetNameTitle(TString::Format("exp_ts%d", int(nSigma)),
908 TString::Format("Expected;%s;%s", _axes.at(0)->GetTitle(), front().tsTitle(true).Data()));
909 }
910 }
911
912 auto badPoints = [&]() {
913 auto badPoints2 = dynamic_cast<TGraph *>(out->GetListOfFunctions()->FindObject("badPoints"));
914 if (!badPoints2) {
915 badPoints2 = new TGraph;
916 badPoints2->SetBit(kCanDelete);
917 badPoints2->SetName("badPoints");
918 badPoints2->SetMarkerStyle(5);
919 badPoints2->SetMarkerColor(std::isnan(nSigma) ? kRed : kBlue);
920 badPoints2->SetMarkerSize(1);
921 out->GetListOfFunctions()->Add(badPoints2, "P");
922 }
923 return badPoints2;
924 };
925 int nPointsDown = 0;
926 bool above = true;
927 TStopwatch s;
928 s.Start();
929 size_t nDone = 0;
930 for (auto &p : *this) {
931 if (s.RealTime() > 5) {
932 if (visualize) {
933 // draw readonly version of the graph
934 auto gra = graph(sOpt + " readOnly");
935 if (gra && gra->GetN()) {
936 if (!gPad && gROOT->GetSelectedPad())
937 gROOT->GetSelectedPad()->cd();
938 if (gPad)
939 gPad->Clear();
940 gra->DrawClone(expBand ? "AF" : "ALP")->SetBit(kCanDelete);
941#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 30, 00)
942 if (auto pad = gROOT->GetSelectedPad()) {
943 pad->GetCanvas()->ResetUpdated(); // stops previous canvas being replaced in a jupyter notebook
944 }
945#endif
947 }
948 } else {
949 ::Info("xRooHypoSpace::graph", "Completed %d/%d points for %s", int(nDone), int(size()), sOpt.Data());
950 }
951 s.Start();
952 } else {
953 s.Continue();
954 }
955 double _x = p.coords->getRealValue(_axes.at(0)->GetName(), std::numeric_limits<double>::quiet_NaN());
956 auto pval = const_cast<xRooHypoPoint &>(p).getVal(sOpt);
957 auto idx = out->GetN() - nPointsDown;
958
959 if (std::isnan(pval.first)) {
960 if (p.status() != 0) { // if status is 0 then bad pval is really just absence of fits, not bad fits
961 badPoints()->SetPoint(badPoints()->GetN(), _x, 0);
962 }
963 } else {
964 out->InsertPointBefore(idx, _x, pval.first);
965 out->SetPointError(idx, 0, pval.second);
966 }
967
968 if (expBand && nSigma) {
970 sOpt2.ReplaceAll("exp", "exp-");
971 pval = const_cast<xRooHypoPoint &>(p).getVal(sOpt2);
972 if (std::isnan(pval.first)) {
973 if (p.status() != 0) { // if status is 0 then bad pval is really just absence of fits, not bad fits
974 badPoints()->SetPoint(badPoints()->GetN(), _x, 0);
975 }
976 } else {
977 out->InsertPointBefore(idx + 1, _x, pval.first);
978 out->SetPointError(idx + 1, 0, pval.second);
979 nPointsDown++;
980 if (out->GetPointY(idx) < pval.first)
981 above = false; // the -sigma points are actually above +sigma
982 }
983 }
984 nDone++;
985 }
986
987 if (out->GetN() == 0)
988 return out;
989
990 if (!expBand) {
991 out->Sort();
992 if (out->GetListOfFunctions()->FindObject("badPoints")) {
993 // try to interpolate the points
994 for (int i = 0; i < badPoints()->GetN(); i++) {
995 badPoints()->SetPointY(i, out->Eval(badPoints()->GetPointX(i)));
996 }
997 }
998 } else {
999 out->Sort(&TGraph::CompareX, true, 0, out->GetN() - nPointsDown - 1); // sort first half
1000 out->Sort(&TGraph::CompareX, false, out->GetN() - nPointsDown, out->GetN() - 1); // reverse sort second half
1001
1002 // now populate the up and down values
1003 auto up = dynamic_cast<TGraph *>(out->GetListOfFunctions()->FindObject("up"));
1004 auto down = dynamic_cast<TGraph *>(out->GetListOfFunctions()->FindObject("down"));
1005
1006 for (int i = 0; i < out->GetN(); i++) {
1007 if (i < out->GetN() - nPointsDown) {
1008 up->SetPoint(up->GetN(), out->GetPointX(i), out->GetPointY(i) + out->GetErrorY(i) * (above ? 1. : -1.));
1009 down->SetPoint(down->GetN(), out->GetPointX(i), out->GetPointY(i) - out->GetErrorY(i) * (above ? 1. : -1.));
1010 } else {
1011 up->SetPoint(up->GetN(), out->GetPointX(i), out->GetPointY(i) - out->GetErrorY(i) * (above ? 1. : -1.));
1012 down->SetPoint(down->GetN(), out->GetPointX(i), out->GetPointY(i) + out->GetErrorY(i) * (above ? 1. : -1.));
1013 }
1014 }
1015 }
1016
1017 if (visualize) {
1018 // draw result
1019 if (!gPad && gROOT->GetSelectedPad())
1020 gROOT->GetSelectedPad()->cd();
1021 if (gPad)
1022 gPad->Clear();
1023 out->DrawClone(expBand ? "AF" : "ALP")->SetBit(kCanDelete);
1024#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 30, 00)
1025 if (auto pad = gROOT->GetSelectedPad()) {
1026 pad->GetCanvas()->ResetUpdated(); // stops previous canvas being replaced in a jupyter notebook
1027 }
1028#endif
1030 }
1031
1032 return out;
1033}
1034
1035std::shared_ptr<TMultiGraph> xRooNLLVar::xRooHypoSpace::graphs(const char *opt)
1036{
1037 TString sOpt(opt);
1038 sOpt.ToLower();
1039 std::shared_ptr<TMultiGraph> out;
1040 if (sOpt.Contains("pcls") || sOpt.Contains("pnull") || sOpt.Contains("ts")) {
1041
1042 bool visualize = sOpt.Contains("visualize");
1043 sOpt.ReplaceAll("visualize", "");
1044
1045 auto exp2 = graph(sOpt + " exp2");
1046 auto exp1 = graph(sOpt + " exp1");
1047 auto exp = graph(sOpt + " exp");
1048 bool doObs = true;
1049 // for(auto& hp : *this) { if(hp.isExpected) {doObs=false; break;} }
1050 auto obs = (doObs) ? graph(sOpt) : nullptr;
1051
1052 out = std::make_shared<TMultiGraph>(GetName(), GetTitle());
1053 if (exp2 && exp2->GetN() > 1)
1054 out->Add(static_cast<TGraph *>(exp2->Clone()), "FP");
1055 if (exp1 && exp1->GetN() > 1)
1056 out->Add(static_cast<TGraph *>(exp1->Clone()), "FP");
1057 if (exp && exp->GetN() > 1)
1058 out->Add(static_cast<TGraph *>(exp->Clone()), "LP");
1059 if (obs && obs->GetN() > 1)
1060 out->Add(static_cast<TGraph *>(obs->Clone()), "LP");
1061
1062 if (!out->GetListOfGraphs()) {
1063 return nullptr;
1064 }
1065
1066 TGraph *testedPoints = nullptr;
1067 if (sOpt.Contains("pcls")) {
1068 TGraph *line = new TGraph;
1069 line->SetName("alpha");
1070 line->SetLineStyle(2);
1071 line->SetEditable(false);
1072 line->SetPoint(line->GetN(), out->GetHistogram()->GetXaxis()->GetXmin() - 10, 0.05);
1073 testedPoints = new TGraph;
1074 testedPoints->SetName("hypoPoints");
1075 testedPoints->SetEditable(false);
1076 testedPoints->SetMarkerStyle(24);
1077 testedPoints->SetMarkerSize(0.4); // use line to indicate tested points
1078 if (exp) {
1079 for (int i = 0; i < exp->GetN(); i++) {
1080 testedPoints->SetPoint(testedPoints->GetN(), exp->GetPointX(i), 0.05);
1081 }
1082 }
1083 line->SetPoint(line->GetN(), out->GetHistogram()->GetXaxis()->GetXmax() + 10, 0.05);
1085 out->GetListOfFunctions()->Add(line, "L");
1086 }
1087 if (exp) {
1088 out->GetHistogram()->GetXaxis()->SetTitle(exp->GetHistogram()->GetXaxis()->GetTitle());
1089 out->GetHistogram()->GetYaxis()->SetTitle(exp->GetHistogram()->GetYaxis()->GetTitle());
1090 }
1091 TLegend *leg = nullptr;
1092 if (out->GetListOfGraphs()->GetEntries() > 1) {
1093 leg = new TLegend(1. - gStyle->GetPadRightMargin() - 0.3, 1. - gStyle->GetPadTopMargin() - 0.35,
1094 1. - gStyle->GetPadRightMargin() - 0.05, 1. - gStyle->GetPadTopMargin() - 0.05);
1095 leg->SetName("legend");
1096 leg->SetBit(kCanDelete);
1097
1098 out->GetListOfFunctions()->Add(leg);
1099 // out->GetListOfFunctions()->Add(out->GetHistogram()->Clone(".axis"),"sameaxis"); // redraw axis
1100
1101 for (auto g : *out->GetListOfGraphs()) {
1102 if (auto o = dynamic_cast<TGraph *>(g)->GetListOfFunctions()->FindObject("down")) {
1103 leg->AddEntry(o, "", "F");
1104 } else {
1105 leg->AddEntry(g, "", "LPE");
1106 }
1107 }
1108 }
1109
1110 auto addToLegend = [](TLegend *l, const char *label, const std::pair<double, double> val) {
1111 if (l) {
1112 l->AddEntry((TObject *)nullptr,
1113 TString::Format("%s%s: %g #pm %g%s", std::isfinite(val.second) ? "" : "#color[2]{", label,
1114 val.first, val.second, std::isfinite(val.second) ? "" : "}"),
1115 "");
1116 }
1117 };
1118
1119 if (sOpt.Contains("pcls")) {
1120 // add current limit estimates to legend
1121 if (exp2 && exp2->GetN() > 1) {
1122 auto l = xRooFit::matchPrecision(GetLimit(*graph(sOpt + "exp-2")));
1123 addToLegend(leg, "-2#sigma", l);
1124 }
1125 if (exp1 && exp1->GetN() > 1) {
1126 auto l = xRooFit::matchPrecision(GetLimit(*graph(sOpt + "exp-1")));
1127 addToLegend(leg, "-1#sigma", l);
1128 }
1129 if (exp && exp->GetN() > 1) {
1130 auto l = xRooFit::matchPrecision(GetLimit(*exp));
1131 addToLegend(leg, "0#sigma", l);
1132 }
1133 if (exp1 && exp1->GetN() > 1) {
1134 auto l = xRooFit::matchPrecision(GetLimit(*graph(sOpt + "exp+1")));
1135 addToLegend(leg, "+1#sigma", l);
1136 }
1137 if (exp2 && exp2->GetN() > 1) {
1138 auto l = xRooFit::matchPrecision(GetLimit(*graph(sOpt + "exp+2")));
1139 addToLegend(leg, "+2#sigma", l);
1140 }
1141 if (obs && obs->GetN() > 1) {
1142 auto l = xRooFit::matchPrecision(GetLimit(*obs));
1143 addToLegend(leg, "Observed", l);
1144 }
1145 }
1146 if (testedPoints)
1147 out->Add(testedPoints, "P");
1148
1149 if (visualize) {
1150 if (!gPad && gROOT->GetSelectedPad())
1151 gROOT->GetSelectedPad()->cd();
1152 if (gPad)
1153 gPad->Clear();
1154 auto gra2 = static_cast<TMultiGraph *>(out->DrawClone("A"));
1155 gra2->SetBit(kCanDelete);
1156 if (sOpt.Contains("pcls") || sOpt.Contains("pnull")) {
1157 gra2->GetHistogram()->SetMinimum(1e-6);
1158 }
1159 if (gPad) {
1160 gPad->RedrawAxis();
1161 gPad->GetCanvas()->Paint();
1162 gPad->GetCanvas()->Update();
1163#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 30, 00)
1164 gPad->GetCanvas()->ResetUpdated(); // stops previous canvas being replaced in a jupyter notebook
1165#endif
1166 }
1168 }
1169 }
1170
1171 return out;
1172}
1173
1175{
1176
1177 if (std::isnan(target)) {
1178 target = 1. - gEnv->GetValue("xRooHypoSpace.CL", 95.) / 100.;
1179 }
1180
1181 auto gr = std::make_shared<TGraph>(pValues);
1182 // remove any nan points and duplicates
1183 int i = 0;
1184 std::set<double> existingX;
1185 while (i < gr->GetN()) {
1186 if (std::isnan(gr->GetPointY(i))) {
1187 gr->RemovePoint(i);
1188 } else if (existingX.find(gr->GetPointX(i)) != existingX.end()) {
1189 gr->RemovePoint(i);
1190 } else {
1191 existingX.insert(gr->GetPointX(i));
1192 // convert to log ....
1193 gr->SetPointY(i, log(std::max(gr->GetPointY(i), 1e-10)));
1194 i++;
1195 }
1196 }
1197
1198 gr->Sort();
1199
1200 // simple linear extrapolation to critical value ... return nan if problem
1201 if (gr->GetN() < 2) {
1202 return std::pair<double, double>(std::numeric_limits<double>::quiet_NaN(), 0);
1203 }
1204
1205 double alpha = log(target);
1206
1207 bool above = gr->GetPointY(0) > alpha;
1208 for (int ii = 1; ii < gr->GetN(); ii++) {
1209 if ((above && (gr->GetPointY(ii) <= alpha)) || (!above && (gr->GetPointY(ii) >= alpha))) {
1210 // found the limit ... return linearly extrapolated point
1211 double lim = gr->GetPointX(ii - 1) + (gr->GetPointX(ii) - gr->GetPointX(ii - 1)) *
1212 (alpha - gr->GetPointY(ii - 1)) /
1213 (gr->GetPointY(ii) - gr->GetPointY(ii - 1));
1214 // use points either side as error
1215 double err = std::max(lim - gr->GetPointX(ii - 1), gr->GetPointX(ii) - lim);
1216 // give err negative sign to indicate if error due to negative side
1217 if ((lim - gr->GetPointX(ii - 1)) > (gr->GetPointX(ii) - lim))
1218 err *= -1;
1219 return std::pair(lim, err);
1220 }
1221 }
1222 // if reach here need to extrapolate ...
1223 if ((above && gr->GetPointY(gr->GetN() - 1) <= gr->GetPointY(0)) ||
1224 (!above && gr->GetPointY(gr->GetN() - 1) >= gr->GetPointY(0))) {
1225 // extrapolating above based on last two points
1226 // in fact, if 2nd last point is a p=1 (log(p)=0) then go back
1227 int offset = 2;
1228 while (offset < gr->GetN() && gr->GetPointY(gr->GetN() - offset) == 0)
1229 offset++;
1230 double x1 = gr->GetPointX(gr->GetN() - offset);
1231 double y1 = gr->GetPointY(gr->GetN() - offset);
1232 double m = (gr->GetPointY(gr->GetN() - 1) - y1) / (gr->GetPointX(gr->GetN() - 1) - x1);
1233 if (m == 0.)
1234 return std::pair(std::numeric_limits<double>::infinity(), std::numeric_limits<double>::infinity());
1235 return std::pair((alpha - y1) / m + x1, std::numeric_limits<double>::infinity());
1236 } else {
1237 // extrapolating below based on first two points
1238 double x1 = gr->GetPointX(0);
1239 double y1 = gr->GetPointY(0);
1240 double m = (gr->GetPointY(1) - y1) / (gr->GetPointX(1) - x1);
1241 if (m == 0.)
1242 return std::pair(-std::numeric_limits<double>::infinity(), -std::numeric_limits<double>::infinity());
1243 return std::pair((alpha - y1) / m + x1, -std::numeric_limits<double>::infinity());
1244 }
1245}
1246
1248{
1249 TString sOpt = TString::Format("p%s", type);
1250 if (std::isnan(nSigma)) {
1251 sOpt += "obs";
1252 } else {
1253 sOpt += TString::Format("exp%s%d", nSigma > 0 ? "+" : "", int(nSigma));
1254 }
1255 return GetLimit(*graph(sOpt + " readonly"));
1256}
1257
1259xRooNLLVar::xRooHypoSpace::findlimit(const char *opt, double relUncert, unsigned int maxTries)
1260{
1261 TString sOpt(opt);
1262 bool visualize = sOpt.Contains("visualize");
1263 sOpt.ReplaceAll("visualize", "");
1264 std::shared_ptr<TGraphErrors> gr = graph(sOpt + " readonly");
1265 if (visualize) {
1266 auto gra = graphs(sOpt.Contains("toys") ? "pcls readonly toys" : "pcls readonly");
1267 if (gra) {
1268 if (!gPad)
1269 gra->Draw(); // in 6.28 DrawClone wont make the gPad defined :( ... so Draw then clear and Draw Clone
1270 gPad->Clear();
1271 gra->DrawClone("A")->SetBit(kCanDelete);
1272 gPad->RedrawAxis();
1273 gra->GetHistogram()->SetMinimum(1e-9);
1274 gra->GetHistogram()->GetYaxis()->SetRangeUser(1e-9, 1);
1275 gPad->Modified();
1276#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 30, 00)
1277 gPad->GetCanvas()->ResetUpdated(); // stops previous canvas being replaced in a jupyter notebook
1278#endif
1280 }
1281 }
1282
1283 // resync parameter boundaries from nlls (may have been modified by fits)
1284 for (auto p : axes()) {
1285 for (auto &[pdf, nll] : fNlls) {
1286 if (auto _v = dynamic_cast<RooRealVar *>(nll->pars()->find(*p))) {
1287 dynamic_cast<RooRealVar *>(p)->setRange(_v->getMin(), _v->getMax());
1288 }
1289 }
1290 }
1291
1292 if (!gr || gr->GetN() < 2) {
1293 auto v = (axes().empty()) ? nullptr : dynamic_cast<RooRealVar *>(*axes().rbegin());
1294 if (!v)
1295 return std::pair(std::numeric_limits<double>::quiet_NaN(), 0.);
1296 double muMax = std::min(std::min(v->getMax("physical"), v->getMax()), v->getMax("scan"));
1297 double muMin = std::max(std::max(v->getMin("physical"), v->getMin()), v->getMin("scan"));
1298 if (!gr || gr->GetN() < 1) {
1299 if (maxTries == 0 || std::isnan(AddPoint(TString::Format("%s=%g", v->GetName(), muMin)).getVal(sOpt).first)) {
1300 // first point failed ... give up
1301 ::Error("findlimit", "Problem evaluating %s @ %s=%g", sOpt.Data(), v->GetName(), muMin);
1302 return std::pair(std::numeric_limits<double>::quiet_NaN(), 0.);
1303 }
1304 gr.reset();
1305 return findlimit(opt, relUncert, maxTries - 1); // do this to resync parameter limits
1306 }
1307
1308 // can approximate expected limit using
1309 // mu_hat + sigma_mu*ROOT::Math::gaussian_quantile(1.-alpha/2.,1) for cls
1310 // or mu_hat + sigma_mu*ROOT::Math::gaussian_quantile((1.-alpha),1) for cls+b
1311 // get a very first estimate of sigma_mu from ufit to expected data, take error on mu as sigma_mu
1312 double nextPoint = muMin + (muMax - muMin) / 50;
1313
1314 // if done an expected limit, assume data is like expected and choose expected limit point as first test point
1315 if (sOpt.Contains("obs")) {
1316 TString sOpt2 = sOpt;
1317 sOpt2.ReplaceAll("obs", "exp");
1318 auto expLim = findlimit(sOpt2, std::numeric_limits<double>::infinity(), 0);
1319 if (!std::isnan(expLim.first) && expLim.first < nextPoint)
1320 nextPoint = expLim.first;
1321 }
1322
1323 auto point =
1324 (sOpt.Contains("exp")) ? back().asimov() : std::shared_ptr<xRooHypoPoint>(&back(), [](xRooHypoPoint *) {});
1325 point = nullptr;
1326 if (point && point->ufit()) {
1327 double rough_sigma_mu = point->mu_hat().getError();
1328 double another_estimate = point->mu_hat().getVal() + rough_sigma_mu * ROOT::Math::gaussian_quantile(0.95, 1);
1329 // if (another_estimate < nextPoint) {
1331 ::Info("xRooHypoSpace::findlimit", "Guessing %g based on rough sigma_mu = %g", nextPoint, rough_sigma_mu);
1332 //}
1333 }
1334
1335 if (maxTries == 0 || std::isnan(AddPoint(TString::Format("%s=%g", v->GetName(), nextPoint)).getVal(sOpt).first)) {
1336 // second point failed ... give up
1337 ::Error("xRooHypoSpace::findlimit", "Problem evaluating %s @ %s=%g", sOpt.Data(), v->GetName(), nextPoint);
1338 return std::pair(std::numeric_limits<double>::quiet_NaN(), 0.);
1339 }
1340 gr.reset();
1341 return findlimit(opt, relUncert, maxTries - 1);
1342 }
1343
1344 auto lim = GetLimit(*gr);
1345
1346 if (std::isnan(lim.first)) {
1347 return lim;
1348 }
1349
1350 auto v = dynamic_cast<RooRealVar *>(*axes().rbegin());
1351 double maxMu = std::min(std::min(v->getMax("physical"), v->getMax()), v->getMax("scan"));
1352 double minMu = std::max(std::max(v->getMin("physical"), v->getMin()), v->getMin("scan"));
1353
1354 // static double MIN_LIMIT_UNCERT = 1e-4; // stop iterating once uncert gets this small
1355 if (lim.first > -std::numeric_limits<double>::infinity() && lim.first < std::numeric_limits<double>::infinity() &&
1356 (std::abs(lim.second) <= relUncert * std::abs(lim.first) /* || std::abs(lim.second)<MIN_LIMIT_UNCERT*/))
1357 return lim;
1358
1359 double nextPoint;
1360
1361 if (lim.second == std::numeric_limits<double>::infinity()) {
1362 // limit was found by extrapolating to right
1363 nextPoint = lim.first;
1364 if (nextPoint == std::numeric_limits<double>::infinity() || nextPoint > maxMu) {
1365 nextPoint = gr->GetPointX(gr->GetN() - 1) + (maxMu - minMu) / 50;
1366 }
1367
1368 // prefer extrapolation with sigma_mu, if available, if it takes us further
1369 // as shape of p-value curve is usually
1370 auto point =
1371 (sOpt.Contains("exp")) ? back().asimov() : std::shared_ptr<xRooHypoPoint>(&back(), [](xRooHypoPoint *) {});
1372 point = nullptr;
1373 if (point && point->ufit()) {
1374 double rough_sigma_mu = point->mu_hat().getError();
1375 double another_estimate = point->mu_hat().getVal() + rough_sigma_mu * ROOT::Math::gaussian_quantile(0.95, 1);
1376 // if (another_estimate < nextPoint) {
1378 ::Info("xRooHypoSpace::findlimit", "Guessing %g based on rough sigma_mu = %g", nextPoint, rough_sigma_mu);
1379 //}
1380 }
1381 nextPoint = std::min(nextPoint + nextPoint * relUncert * 0.99, maxMu); // ensure we step over location if possible
1382
1383 if (nextPoint > maxMu)
1384 return lim;
1385 } else if (lim.second == -std::numeric_limits<double>::infinity()) {
1386 // limit from extrapolating to left
1387 nextPoint = lim.first;
1388 if (nextPoint < minMu)
1389 nextPoint = gr->GetPointX(0) - (maxMu - minMu) / 50;
1390 if (nextPoint < minMu)
1391 return lim;
1392 } else {
1393 nextPoint = lim.first + lim.second * relUncert * 0.99;
1394 }
1395
1396 // got here need a new point .... evaluate the estimated lim location +/- the relUncert (signed error takes care of
1397 // direction)
1398
1399 ::Info("xRooHypoSpace::findlimit", "%s -- Testing new point @ %s=%g (delta=%g)", sOpt.Data(), v->GetName(),
1400 nextPoint, lim.second);
1401 if (maxTries == 0 || std::isnan(AddPoint(TString::Format("%s=%g", v->GetName(), nextPoint)).getVal(sOpt).first)) {
1402 if (maxTries == 0) {
1403 ::Warning("xRooHypoSpace::findlimit", "Reached max number of point evaluations");
1404 } else {
1405 ::Error("xRooHypoSpace::findlimit", "Problem evaluating %s @ %s=%g", sOpt.Data(), v->GetName(), nextPoint);
1406 }
1407 return lim;
1408 }
1409 gr.reset();
1410 return findlimit(opt, relUncert, maxTries - 1);
1411}
1412
1414{
1415
1416 TString sOpt(opt);
1417 sOpt.ToLower();
1418
1419 if ((sOpt == "" || sOpt == "same") && !empty()) {
1420 if (front().fPllType == xRooFit::Asymptotics::OneSidedPositive) {
1421 sOpt += "pcls"; // default to showing cls p-value scan if drawing a limit
1422 for (auto &hp : *this) {
1423 if (!hp.nullToys.empty() || !hp.altToys.empty()) {
1424 sOpt += " toys";
1425 break; // default to toys if done toys
1426 }
1427 }
1428 } else if (front().fPllType == xRooFit::Asymptotics::TwoSided) {
1429 sOpt += "ts";
1430 }
1431 }
1432
1433 // split up by ; and call Draw for each (with 'same' appended)
1434 auto _axes = axes();
1435 if (_axes.empty())
1436 return;
1437
1438 if (sOpt == "status") {
1439 // draw the points in the space
1440 if (_axes.size() <= 2) {
1441 TGraphErrors *out = new TGraphErrors;
1442 out->SetBit(kCanDelete);
1443 out->SetName("points");
1444 out->SetMarkerSize(0.5);
1445 TGraph *tsAvail = new TGraph;
1446 tsAvail->SetName("ts");
1447 tsAvail->SetBit(kCanDelete);
1448 tsAvail->SetMarkerStyle(20);
1449 TGraph *expAvail = new TGraph;
1450 expAvail->SetName("exp");
1451 expAvail->SetBit(kCanDelete);
1452 expAvail->SetMarkerStyle(25);
1453 expAvail->SetMarkerSize(out->GetMarkerSize() * 1.5);
1454 TGraph *badPoints = new TGraph;
1455 badPoints->SetName("bad_ufit");
1456 badPoints->SetBit(kCanDelete);
1457 badPoints->SetMarkerStyle(5);
1458 badPoints->SetMarkerColor(kRed);
1459 badPoints->SetMarkerSize(out->GetMarkerSize());
1460 TGraph *badPoints2 = new TGraph;
1461 badPoints2->SetName("bad_cfit_null");
1462 badPoints2->SetBit(kCanDelete);
1463 badPoints2->SetMarkerStyle(2);
1464 badPoints2->SetMarkerColor(kRed);
1465 badPoints2->SetMarkerSize(out->GetMarkerSize());
1466
1467 out->SetTitle(TString::Format("%s;%s;%s", GetTitle(), _axes.at(0)->GetTitle(),
1468 (_axes.size() == 1) ? "" : _axes.at(1)->GetTitle()));
1469 for (auto &p : *this) {
1470 bool _readOnly = p.nllVar ? p.nllVar->get()->getAttribute("readOnly") : false;
1471 if (p.nllVar)
1472 p.nllVar->get()->setAttribute("readOnly", true);
1473 double x = p.coords->getRealValue(_axes.at(0)->GetName());
1474 double y = _axes.size() == 1 ? p.ts_asymp().first : p.coords->getRealValue(_axes.at(1)->GetName());
1475 out->SetPoint(out->GetN(), x, y);
1476 if (!std::isnan(p.ts_asymp().first)) {
1477 if (_axes.size() == 1)
1478 out->SetPointError(out->GetN() - 1, 0, p.ts_asymp().second);
1479 tsAvail->SetPoint(tsAvail->GetN(), x, y);
1480 } else if (p.fUfit && (std::isnan(p.fUfit->minNll()) ||
1481 xRooNLLVar::xRooHypoPoint::allowedStatusCodes.find(p.fUfit->status()) ==
1483 badPoints->SetPoint(badPoints->GetN(), x, y);
1484 } else if (p.fNull_cfit && (std::isnan(p.fNull_cfit->minNll()) ||
1485 xRooNLLVar::xRooHypoPoint::allowedStatusCodes.find(p.fNull_cfit->status()) ==
1487 badPoints2->SetPoint(badPoints2->GetN(), x, y);
1488 }
1489 if (!std::isnan(p.ts_asymp(0).first)) {
1490 expAvail->SetPoint(expAvail->GetN(), x, y);
1491 } else if (p.asimov() && p.asimov()->fUfit &&
1492 (std::isnan(p.asimov()->fUfit->minNll()) ||
1493 xRooNLLVar::xRooHypoPoint::allowedStatusCodes.find(p.asimov()->fUfit->status()) ==
1495
1496 } else if (p.asimov() && p.asimov()->fNull_cfit &&
1497 (std::isnan(p.asimov()->fNull_cfit->minNll()) ||
1498 xRooNLLVar::xRooHypoPoint::allowedStatusCodes.find(p.asimov()->fNull_cfit->status()) ==
1500 }
1501 if (p.nllVar)
1502 p.nllVar->get()->setAttribute("readOnly", _readOnly);
1503 }
1504
1505 if (_axes.size() == 1) {
1506 TGraph tmp;
1507 for (int i = 0; i < out->GetN(); i++) {
1508 if (!std::isnan(out->GetPointY(i)))
1509 tmp.SetPoint(tmp.GetN(), out->GetPointX(i), out->GetPointY(i));
1510 }
1511 auto fixPoints = [&](TGraph *g) {
1512 for (int i = 0; i < g->GetN(); i++) {
1513 if (std::isnan(g->GetPointY(i)))
1514 g->SetPointY(i, std::isnan(tmp.Eval(g->GetPointX(i))) ? 0. : tmp.Eval(g->GetPointX(i)));
1515 }
1516 };
1517 fixPoints(out);
1522 }
1523
1524 out->SetMarkerStyle(4);
1525 out->Draw("AP");
1526 auto leg = new TLegend(1. - gPad->GetRightMargin() - 0.3, 1. - gPad->GetTopMargin() - 0.35,
1527 1. - gPad->GetRightMargin() - 0.05, 1. - gPad->GetTopMargin() - 0.05);
1528 leg->SetName("legend");
1529 leg->AddEntry(out, "Uncomputed", "P");
1530
1531 if (tsAvail->GetN()) {
1532 out->GetListOfFunctions()->Add(tsAvail, "P");
1533 leg->AddEntry(tsAvail, "Computed", "P");
1534 } else {
1535 delete tsAvail;
1536 }
1537 if (expAvail->GetN()) {
1538 out->GetListOfFunctions()->Add(expAvail, "P");
1539 leg->AddEntry(expAvail, "Expected computed", "P");
1540 } else {
1541 delete expAvail;
1542 }
1543 if (badPoints->GetN()) {
1544 out->GetListOfFunctions()->Add(badPoints, "P");
1545 leg->AddEntry(badPoints, "Bad ufit", "P");
1546 } else {
1547 delete badPoints;
1548 }
1549 if (badPoints2->GetN()) {
1550 out->GetListOfFunctions()->Add(badPoints2, "P");
1551 leg->AddEntry(badPoints2, "Bad null cfit", "P");
1552 } else {
1553 delete badPoints2;
1554 }
1555 leg->SetBit(kCanDelete);
1556 leg->Draw();
1557 // if(_axes.size()==1) out->GetHistogram()->GetYaxis()->SetRangeUser(0,1);
1558 gPad->SetGrid(true, _axes.size() > 1);
1559 if (_axes.size() == 1)
1560 gPad->SetLogy(false);
1561 }
1562
1564
1565 return;
1566 }
1567
1568 if (sOpt.Contains("pcls") || sOpt.Contains("pnull") || sOpt.Contains("ts")) {
1569 auto gra = graphs(sOpt + " readonly");
1570 if (!gPad && gROOT->GetSelectedPad())
1571 gROOT->GetSelectedPad()->cd();
1572 if (!sOpt.Contains("same") && gPad) {
1573 gPad->Clear();
1574 }
1575 if (gra) {
1576 if (!gPad)
1578 auto gra2 = static_cast<TMultiGraph *>(gra->DrawClone(sOpt.Contains("same") ? "" : "A"));
1579 gra2->SetBit(kCanDelete);
1580 if (sOpt.Contains("pcls") || sOpt.Contains("pnull")) {
1581 gra2->GetHistogram()->SetMinimum(1e-6);
1582 }
1583 if (gPad) {
1584 gPad->RedrawAxis();
1585 gPad->GetCanvas()->Paint();
1586 gPad->GetCanvas()->Update();
1587#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 30, 00)
1588 gPad->GetCanvas()->ResetUpdated(); // stops previous canvas being replaced in a jupyter notebook
1589#endif
1591 }
1592 }
1593 if (!sOpt.Contains("same") && gPad) {
1594 // auto mg = static_cast<TMultiGraph*>(gPad->GetPrimitive(gra->GetName()));
1595 // mg->GetHistogram()->SetMinimum(1e-9);
1596 // mg->GetHistogram()->GetYaxis()->SetRangeUser(1e-9,1);
1597 // gPad->SetGrid(0, 0);
1598 // gPad->SetLogy(1);
1599 }
1600
1602
1603 return;
1604 }
1605
1606 // graphs("ts visualize");return;
1607
1608 TGraphErrors *out = new TGraphErrors;
1609 out->SetName(GetName());
1610
1611 TString title = (!axes().empty()) ? TString::Format(";%s", axes().first()->GetTitle()) : "";
1612
1614 if (!empty() && axes().size() == 1) {
1615 for (auto &p : *this) {
1616 if (p.fPllType != xRooFit::Asymptotics::TwoSided) {
1617 pllType = p.fPllType;
1618 }
1619 }
1620 title += ";";
1621 title += front().tsTitle(true);
1622 }
1623
1624 out->SetTitle(title);
1625 *dynamic_cast<TAttFill *>(out) = *this;
1626 *dynamic_cast<TAttLine *>(out) = *this;
1627 *dynamic_cast<TAttMarker *>(out) = *this;
1628 out->SetBit(kCanDelete);
1629
1630 if (!gPad)
1633 if (!sOpt.Contains("same"))
1634 basePad->Clear();
1635
1636 bool doFits = false;
1637 if (sOpt.Contains("fits")) {
1638 doFits = true;
1639 sOpt.ReplaceAll("fits", "");
1640 }
1641
1642 auto mainPad = gPad;
1643
1644 out->SetEditable(false);
1645
1646 if (doFits) {
1647 gPad->Divide(1, 2);
1648 gPad->cd(1);
1649 gPad->SetBottomMargin(gPad->GetBottomMargin() * 2.); // increase margin to be same as before
1650 gPad->SetGrid();
1651 out->Draw(sOpt);
1652 basePad->cd(2);
1653 mainPad = basePad->GetPad(1);
1654 } else {
1655 gPad->SetGrid();
1656 out->Draw("ALP");
1657 }
1658
1659 std::pair<double, double> minMax(std::numeric_limits<double>::infinity(), -std::numeric_limits<double>::infinity());
1660 for (auto &p : *this) {
1661 if (p.fPllType != pllType)
1662 continue; // must all have same pll type
1663 auto val = p.pll(true).first;
1664 if (std::isnan(val))
1665 continue;
1666 minMax.first = std::min(minMax.first, val);
1667 minMax.second = std::max(minMax.second, val);
1668 }
1669 if (minMax.first < std::numeric_limits<double>::infinity())
1670 out->GetHistogram()->SetMinimum(minMax.first);
1671 if (minMax.second > -std::numeric_limits<double>::infinity())
1672 out->GetHistogram()->SetMaximum(minMax.second);
1673
1674 TGraph *badPoints = nullptr;
1675
1676 TStopwatch s;
1677 s.Start();
1678 std::shared_ptr<const RooFitResult> ufr;
1679 for (auto &p : *this) {
1680 if (p.fPllType != pllType)
1681 continue; // must all have same pll type
1682 auto val = p.pll().first;
1683 if (!ufr)
1684 ufr = p.ufit();
1685 if (out->GetN() == 0 && ufr && ufr->status() == 0) {
1686 out->SetPoint(out->GetN(),
1687 ufr->floatParsFinal().getRealValue(axes().first()->GetName(),
1688 ufr->constPars().getRealValue(axes().first()->GetName())),
1689 0.);
1690 out->SetPointError(out->GetN() - 1, 0, ufr->edm());
1691 }
1692 if (auto fr = p.fNull_cfit;
1693 fr && doFits) { // access member to avoid unnecessarily creating fit result if wasnt needed
1694 // create a new subpad and draw fitResult on it
1695 auto _pad = gPad;
1696 auto pad =
1697 new TPad(fr->GetName(), TString::Format("%s = %g", poi().first()->GetTitle(), p.fNullVal()), 0, 0, 1., 1);
1698 pad->SetNumber(out->GetN() + 1); // can't use "0" for a subpad
1699 pad->cd();
1700 xRooNode(fr).Draw("goff");
1701 _pad->cd();
1702 //_pad->GetListOfPrimitives()->AddFirst(pad);
1703 pad->AppendPad();
1704 }
1705 if (std::isnan(val) && p.status() != 0) {
1706 if (!badPoints) {
1707 badPoints = new TGraph;
1708 badPoints->SetBit(kCanDelete);
1709 badPoints->SetName("badPoints");
1710 badPoints->SetMarkerStyle(5);
1711 badPoints->SetMarkerColor(kRed);
1712 badPoints->SetMarkerSize(1);
1713 out->GetListOfFunctions()->Add(badPoints, "P");
1714 }
1715 badPoints->SetPoint(badPoints->GetN(), p.fNullVal(), out->Eval(p.fNullVal()));
1716 mainPad->Modified();
1717 } else if (!std::isnan(val)) {
1718 out->SetPoint(out->GetN(), p.coords->getRealValue(axes().first()->GetName()), p.pll().first);
1719 out->SetPointError(out->GetN() - 1, 0, p.pll().second);
1720 out->Sort();
1721
1722 // reposition bad points
1723 if (badPoints) {
1724 for (int i = 0; i < badPoints->GetN(); i++)
1725 badPoints->SetPointY(i, out->Eval(badPoints->GetPointX(i)));
1726 }
1727
1728 mainPad->Modified();
1729 }
1730 if (s.RealTime() > 3) { // stops the clock
1731 basePad->GetCanvas()->Paint();
1732 basePad->GetCanvas()->Update();
1734 s.Reset();
1735 s.Start();
1736 }
1737 s.Continue();
1738 }
1739 basePad->GetCanvas()->Paint();
1740 basePad->GetCanvas()->Update();
1741#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 30, 00)
1742 basePad->GetCanvas()->ResetUpdated(); // stops previous canvas being replaced in a jupyter notebook
1743#endif
1745
1746 // finish by overlaying ufit
1747 if (ufr && doFits) {
1748 auto _pad = gPad;
1749 auto pad = new TPad(ufr->GetName(), "unconditional fit", 0, 0, 1., 1.);
1750 pad->SetNumber(-1);
1751 pad->cd();
1752 xRooNode(ufr).Draw("goff");
1753 _pad->cd();
1754 pad->AppendPad();
1755
1756 // draw one more pad to represent the selected, and draw the ufit pad onto that pad
1757 pad = new TPad("selected", "selected", 0, 0, 1, 1);
1758 pad->Draw();
1759 pad->cd();
1760 basePad->GetPad(2)->GetPad(-1)->AppendPad();
1761 pad->Modified();
1762 pad->Update();
1764
1765 basePad->cd();
1766 }
1767
1768 if (doFits) {
1769 if (!xRooNode::gIntObj) {
1771 }
1772 gPad->GetCanvas()->Connect("Highlighted(TVirtualPad*,TObject*,Int_t,Int_t)", "xRooNode::InteractiveObject",
1773 xRooNode::gIntObj, "Interactive_PLLPlot(TVirtualPad*,TObject*,Int_t,Int_t)");
1774 }
1775
1776 return;
1777}
1778
1780{
1781
1782 RooStats::HypoTestInverterResult *out = nullptr;
1783
1784 auto _axes = axes();
1785 if (_axes.empty())
1786 return out;
1787
1788 out = new RooStats::HypoTestInverterResult(GetName(), *dynamic_cast<RooRealVar *>(_axes.at(0)), 0.95);
1789 out->SetTitle(GetTitle());
1790
1791 for (auto &p : *this) {
1792 double _x = p.coords->getRealValue(_axes.at(0)->GetName(), std::numeric_limits<double>::quiet_NaN());
1793 out->Add(_x, p.result());
1794 }
1795
1796 return out;
1797}
1798
#define c(i)
Definition RSha256.hxx:101
#define g(i)
Definition RSha256.hxx:105
#define a(i)
Definition RSha256.hxx:99
#define e(i)
Definition RSha256.hxx:103
double _val
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
@ kRed
Definition Rtypes.h:67
@ kBlack
Definition Rtypes.h:66
@ kGreen
Definition Rtypes.h:67
@ kBlue
Definition Rtypes.h:67
@ kYellow
Definition Rtypes.h:67
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define gDirectory
Definition TDirectory.h:385
R__EXTERN TEnv * gEnv
Definition TEnv.h:170
void Info(const char *location, const char *msgfmt,...)
Use this function for informational messages.
Definition TError.cxx:241
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:208
void Warning(const char *location, const char *msgfmt,...)
Use this function in warning situations.
Definition TError.cxx:252
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t SetFillStyle
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h offset
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t target
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t SetFillColor
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
Option_t Option_t TPoint TPoint const char y1
char name[80]
Definition TGX11.cxx:145
@ kCanDelete
Definition TObject.h:375
#define gROOT
Definition TROOT.h:426
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
R__EXTERN TSystem * gSystem
Definition TSystem.h:582
#define gPad
static std::pair< double, double > matchPrecision(const std::pair< double, double > &in)
Definition xRooFit.cxx:2014
std::pair< std::shared_ptr< RooAbsData >, std::shared_ptr< const RooAbsCollection > > fData
Definition xRooNLLVar.h:180
std::shared_ptr< const RooAbsCollection > coords
Definition xRooNLLVar.h:263
std::shared_ptr< xRooHypoPoint > asimov(bool readOnly=false)
std::shared_ptr< const RooFitResult > fAlt_cfit
Definition xRooNLLVar.h:265
xValueWithError limit(const char *type="cls", double nSigma=std::numeric_limits< double >::quiet_NaN()) const
static xValueWithError GetLimit(const TGraph &pValues, double target=std::numeric_limits< double >::quiet_NaN())
std::map< std::string, xValueWithError > limits(const char *opt="cls", const std::vector< double > &nSigmas={0, 1, 2, -1, -2, std::numeric_limits< double >::quiet_NaN()}, double relUncert=std::numeric_limits< double >::infinity())
bool AddModel(const xRooNode &pdf, const char *validity="")
std::shared_ptr< TMultiGraph > graphs(const char *opt)
xValueWithError findlimit(const char *opt, double relUncert=std::numeric_limits< double >::infinity(), unsigned int maxTries=20)
int AddPoints(const char *parName, size_t nPoints, double low, double high)
std::shared_ptr< TGraphErrors > graph(const char *opt) const
void Print(Option_t *opt="") const override
Print TNamed name and title.
xRooHypoSpace(const char *name="", const char *title="")
int scan(const char *type, size_t nPoints, double low=std::numeric_limits< double >::quiet_NaN(), double high=std::numeric_limits< double >::quiet_NaN(), const std::vector< double > &nSigmas={0, 1, 2, -1, -2, std::numeric_limits< double >::quiet_NaN()}, double relUncert=0.1)
std::shared_ptr< xRooNode > pdf(const RooAbsCollection &parValues) const
void Draw(Option_t *opt="") override
Default Draw method for all objects.
std::shared_ptr< RooAbsPdf > pdf() const
Definition xRooNLLVar.h:452
std::shared_ptr< RooArgSet > pars(bool stripGlobalObs=true) const
The xRooNode class is designed to wrap over a TObject and provide functionality to aid with interacti...
Definition xRooNode.h:52
void Draw(Option_t *opt="") override
Default Draw method for all objects.
static InteractiveObject * gIntObj
Definition xRooNode.h:504
xRooNode pars() const
List of parameters (non-observables) of this node.
const_iterator end() const
A space to attach TBranches.
virtual value_type getCurrentIndex() const
Return index number of current state.
Abstract container object that can hold multiple RooAbsArg objects.
RooAbsCollection * selectByAttrib(const char *name, bool value) const
Create a subset of the current collection, consisting only of those elements with the specified attri...
double getRealValue(const char *name, double defVal=0.0, bool verbose=false) const
Get value of a RooAbsReal stored in set with given name.
Int_t getCatIndex(const char *name, Int_t defVal=0, bool verbose=false) const
Get index value of a RooAbsCategory stored in set with given name.
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
void setAttribAll(const Text_t *name, bool value=true)
Set given attribute in each element of the collection by calling each elements setAttribute() functio...
virtual RooAbsArg * addClone(const RooAbsArg &var, bool silent=false)
Add a clone of the specified argument to list.
void setName(const char *name)
RooAbsArg * find(const char *name) const
Find object with given name in list.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Represents a constant real-valued object.
Definition RooConstVar.h:23
Container class to hold unbinned data.
Definition RooDataSet.h:32
Variable that can be changed from the outside.
Definition RooRealVar.h:37
HypoTestInverterResult class holds the array of hypothesis test results and compute a confidence inte...
HypoTestResult * GetResult(int index) const
return a pointer to the i^th result object
int ArraySize() const
number of entries in the results array
double GetXValue(int index) const
function to return the value of the parameter of interest for the i^th entry in the results
bool Add(const HypoTestInverterResult &otherResult)
merge with the content of another HypoTestInverterResult object
RooArgSet * GetParameters() const override
return a cloned list with the parameter of interest
Fill Area Attributes class.
Definition TAttFill.h:21
Line Attributes class.
Definition TAttLine.h:21
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:46
Marker Attributes class.
Definition TAttMarker.h:21
virtual Size_t GetMarkerSize() const
Return the marker size.
Definition TAttMarker.h:35
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:43
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:48
static TCanvas * MakeDefCanvas()
Static function to build a default canvas.
Definition TCanvas.cxx:1509
Describe directory structure in memory.
Definition TDirectory.h:45
virtual Int_t GetValue(const char *name, Int_t dflt) const
Returns the integer value for a resource.
Definition TEnv.cxx:503
A TGraphErrors is a TGraph with error bars.
virtual void SetPointError(Double_t ex, Double_t ey)
Set ex and ey values for point pointed by the mouse.
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
virtual Double_t GetPointX(Int_t i) const
Get x value for point i.
Definition TGraph.cxx:1574
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
Definition TGraph.cxx:2387
static Bool_t CompareX(const TGraph *gr, Int_t left, Int_t right)
Return kTRUE if fX[left] > fX[right]. Can be used by Sort.
Definition TGraph.cxx:733
Int_t GetN() const
Definition TGraph.h:131
virtual void Sort(Bool_t(*greater)(const TGraph *, Int_t, Int_t)=&TGraph::CompareX, Bool_t ascending=kTRUE, Int_t low=0, Int_t high=-1111)
Sorts the points of this TGraph using in-place quicksort (see e.g.
Definition TGraph.cxx:2531
TList * GetListOfFunctions() const
Definition TGraph.h:125
virtual Double_t Eval(Double_t x, TSpline *spline=nullptr, Option_t *option="") const
Interpolate points in this graph at x using a TSpline.
Definition TGraph.cxx:981
virtual Int_t RemovePoint()
Delete point close to the mouse position Returns index of removed point (or -1 if nothing was changed...
Definition TGraph.cxx:2067
void SetName(const char *name="") override
Set graph name.
Definition TGraph.cxx:2426
void Draw(Option_t *chopt="") override
Draw this graph with its current attributes.
Definition TGraph.cxx:859
virtual TH1F * GetHistogram() const
Returns a pointer to the histogram used to draw the axis Takes into account the two following cases.
Definition TGraph.cxx:1458
virtual Double_t GetPointY(Int_t i) const
Get y value for point i.
Definition TGraph.cxx:1585
void SetTitle(const char *title="") override
Change (i.e.
Definition TGraph.cxx:2442
virtual void SetPointY(Int_t i, Double_t y)
Set y value for point i.
Definition TGraph.cxx:2419
virtual void SetEditable(Bool_t editable=kTRUE)
if editable=kFALSE, the graph cannot be modified with the mouse by default a TGraph is editable
Definition TGraph.cxx:2339
This class displays a legend box (TPaveText) containing several legend entries.
Definition TLegend.h:23
TObject * FindObject(const char *name) const override
Find an object in this list using its name.
Definition TList.cxx:708
void Add(TObject *obj) override
Definition TList.h:81
A TMemFile is like a normal TFile except that it reads and writes only from memory.
Definition TMemFile.h:27
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition TMultiGraph.h:34
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:173
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:50
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition TNamed.cxx:149
virtual void SetNameTitle(const char *name, const char *title)
Set all the TNamed parameters (name and title).
Definition TNamed.cxx:163
Mother of all ROOT objects.
Definition TObject.h:42
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:885
The most important graphics class in the ROOT system.
Definition TPad.h:28
Bool_t Connect(const char *signal, const char *receiver_class, void *receiver, const char *slot)
Non-static method is used to connect from the signal of this object to the receiver slot.
Definition TQObject.cxx:865
Stopwatch class.
Definition TStopwatch.h:28
Double_t RealTime()
Stop the stopwatch (if it is running) and return the realtime (in seconds) passed between the start a...
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
void Continue()
Resume a stopped stopwatch.
void Reset()
Definition TStopwatch.h:52
Provides iteration through tokens of a given string.
Definition TPRegexp.h:143
Bool_t NextToken()
Get the next token, it is stored in this TString.
Basic string class.
Definition TString.h:138
Ssiz_t Length() const
Definition TString.h:425
Double_t Atof() const
Return floating-point value contained in string.
Definition TString.cxx:2060
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:2384
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition TString.h:660
Float_t GetPadRightMargin() const
Definition TStyle.h:216
Float_t GetPadTopMargin() const
Definition TStyle.h:214
virtual Bool_t ProcessEvents()
Process pending events (GUI, timers, sockets).
Definition TSystem.cxx:418
This class defines a UUID (Universally Unique IDentifier), also known as GUIDs (Globally Unique IDent...
Definition TUUID.h:42
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition TVirtualPad.h:51
TLine * line
double gaussian_quantile(double z, double sigma)
Inverse ( ) of the cumulative distribution function of the lower tail of the normal (Gaussian) distri...
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
TGraphErrors * gr
Definition legend1.C:25
leg
Definition legend1.C:34
return c2
Definition legend2.C:14
#define BEGIN_XROOFIT_NAMESPACE
Definition Config.h:24
#define END_XROOFIT_NAMESPACE
Definition Config.h:25
TMarker m
Definition textangle.C:8
TLine l
Definition textangle.C:4