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 }
306
307 if (high < low || (high == low && nPoints != 1)) {
308 // take from parameter
309 low = p->getMin("scan");
310 high = p->getMax("scan");
311 }
312 if (!std::isnan(low) && !std::isnan(high) && !(std::isinf(low) && std::isinf(high))) {
313 p->setRange("scan", low, high);
314 }
315 if (p->hasRange("scan")) {
316 ::Info("xRooHypoSpace::scan", "Using %s scan range: %g - %g", p->GetName(), p->getMin("scan"), p->getMax("scan"));
317 }
318
319 bool doObs = false;
320 for (auto nSigma : nSigmas) {
321 if (std::isnan(nSigma)) {
322 doObs = true;
323 break;
324 }
325 }
326
327 if (fNlls.empty()) {
328 // this happens when loaded hypoSpace from a hypoSpaceInverterResult
329 // set relUncert to infinity so that we don't test any new points
330 relUncert = std::numeric_limits<double>::infinity(); // no NLL available so just get whatever limit we can
331
332 // if any of the defined points are 'expected' data don't do obs
333 // for(auto& hp : *this) {
334 // if(hp.isExpected) {
335 // doObs = false; break;
336 // }
337 // }
338 } else {
339#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 26, 00)
340// bool allGen = true;
341// for (auto &[pdf, nll]: fNlls) {
342// auto _d = dynamic_cast<RooDataSet *>(nll->data());
343// if (!_d || !_d->weightVar() || !_d->weightVar()->getStringAttribute("fitResult") ||
344// !_d->weightVar()->getAttribute("expected")) {
345// allGen = false;
346// break;
347// }
348// }
349// if (allGen)
350// doObs = false;
351#endif
352 }
353
354 // create a fitDatabase if required
356 if (!gDirectory || !gDirectory->IsWritable()) {
357 // locate a TMemFile in the open list of files and move to that
358 // or create one if cannot find
359 for (auto file : *gROOT->GetListOfFiles()) {
360 if (auto f = dynamic_cast<TMemFile *>(file)) {
361 f->cd();
362 break;
363 }
364 }
365 if (!gDirectory || !gDirectory->IsWritable()) {
366 new TMemFile("fitDatabase", "RECREATE");
367 }
368 }
369
370 int out = 0;
371
372 if (nPoints == 0) {
373 // automatic scan
374 if (sType.Contains("cls")) {
375 for (double nSigma : nSigmas) {
376 xValueWithError res(std::make_pair(0., 0.));
377 if (std::isnan(nSigma)) {
378 if (!doObs)
379 continue;
380 res = findlimit(TString::Format("%s obs", sType.Data()), relUncert);
381 } else {
382 res =
383 findlimit(TString::Format("%s exp%s%d", sType.Data(), nSigma > 0 ? "+" : "", int(nSigma)), relUncert);
384 }
385 if (std::isnan(res.first) || std::isnan(res.second)) {
386 out = 1;
387 } else if (std::isinf(res.second)) {
388 out = 2;
389 }
390 }
391 } else {
392 throw std::runtime_error(TString::Format("Automatic scanning not yet supported for %s", type));
393 }
394 } else {
395 // add the required points and then compute the required value
396 if (nPoints == 1) {
397 AddPoint(TString::Format("%s=%g", poi().first()->GetName(), (high + low) / 2.));
398 graphs(sType); // triggers computation
399 } else {
400 double step = (high - low) / (nPoints - 1);
401 for (size_t i = 0; i < nPoints; i++) {
402 AddPoint(TString::Format("%s=%g", poi().first()->GetName(), low + step * i));
403 graphs(sType); // triggers computation
404 }
405 }
406 }
407
408 if (origDir)
409 origDir->cd();
410
411 return out;
412}
413
414std::map<std::string, xRooNLLVar::xValueWithError>
415xRooNLLVar::xRooHypoSpace::limits(const char *opt, const std::vector<double> &nSigmas, double relUncert)
416{
417
418 if (fNlls.empty()) {
419 // this happens when loaded hypoSpace from a hypoSpaceInverterResult
420 // set relUncert to infinity so that we don't test any new points
421 relUncert = std::numeric_limits<double>::infinity(); // no NLL available so just get whatever limit we can
422 }
423
424 scan(opt, nSigmas, relUncert);
425
426 std::map<std::string, xRooNLLVar::xValueWithError> out;
427 for (auto nSigma : nSigmas) {
428 auto lim = limit(opt, nSigma);
429 if (lim.second < 0)
430 lim.second = -lim.second; // make errors positive for this method
431 out[std::isnan(nSigma) ? "obs" : TString::Format("%d", int(nSigma)).Data()] = xRooFit::matchPrecision(lim);
432 }
433 return out;
434}
435
437{
438 // move to given coords, if any ... will mark them const too
439 std::unique_ptr<RooAbsCollection, std::function<void(RooAbsCollection *)>> _snap(fPars->snapshot(),
440 [&](RooAbsCollection *c) {
441 *fPars = *c;
442 delete c;
443 });
444 TStringToken pattern(coords, ";");
445 while (pattern.NextToken()) {
446 TString s = pattern;
447 // split by "=" sign
448 auto _idx = s.Index('=');
449 if (_idx == -1)
450 continue;
451 TString _name = s(0, _idx);
452 TString _val = s(_idx + 1, s.Length());
453 auto _v = dynamic_cast<RooRealVar *>(fPars->find(_name));
454 if (!_v)
455 continue;
456
457 if (_val.IsFloat()) {
458 _v->setConstant();
459 _v->setVal(_val.Atof());
460 }
461 }
462
463 auto _pdf = pdf();
464
465 if (!_pdf)
466 throw std::runtime_error("no model at coordinates");
467
468 // if (std::unique_ptr<RooAbsCollection>(fPars->selectByAttrib("poi", true))->size() == 0) {
469 // throw std::runtime_error(
470 // "No pars designated as POI - set with pars()->find(<parName>)->setAttribute(\"poi\",true)");
471 // }
472
473 if (fNlls.find(_pdf) == fNlls.end()) {
474 fNlls[_pdf] = std::make_shared<xRooNLLVar>(_pdf->nll("" /*TODO:allow change dataset name and nll opts*/, {}));
475 }
476
477 xRooHypoPoint out;
478
479 out.nllVar = fNlls[_pdf];
480 out.fData = fNlls[_pdf]->getData();
481 out.isExpected = dynamic_cast<RooDataSet *>(out.fData.first.get()) &&
482 dynamic_cast<RooDataSet *>(out.fData.first.get())->weightVar()->getAttribute("expected");
483 // TODO: need to access the genfit of the data and add that to the point, somehow ...
484
485 out.coords.reset(fPars->snapshot()); // should already have altVal prop on poi, and poi labelled
486 // ensure all poi are marked const ... required by xRooHypoPoint behaviour
487 out.poi().setAttribAll("Constant");
488 // and now remove anything that's marked floating
489
490 // do to bug in remove have to ensure not using the hash map otherwise will be doing an invalid read after the
491 // deletion of the owned pars
492 const_cast<RooAbsCollection *>(out.coords.get())->useHashMapForFind(false);
493 const_cast<RooAbsCollection *>(out.coords.get())
494 ->remove(*std::unique_ptr<RooAbsCollection>(out.coords->selectByAttrib("Constant", false)), true, true);
495
496 double value = out.fNullVal();
497 double alt_value = out.fAltVal();
498
499 auto _type = fTestStatType;
500 if (_type == xRooFit::Asymptotics::Unknown) {
501 // decide based on values
502 if (std::isnan(alt_value)) {
504 } else if (value >= alt_value) {
506 } else {
508 }
509 }
510
511 out.fPllType = _type;
512
513 // look for a matching point
514 for (auto &p : *this) {
515 if (p.nllVar != out.nllVar)
516 continue;
517 if (p.fData != out.fData)
518 continue;
519 if (!p.alt_poi().equals(out.alt_poi()))
520 continue;
521 bool match = true;
522 for (auto c : p.alt_poi()) {
523 if (auto v = dynamic_cast<RooAbsReal *>(c);
524 v && std::abs(v->getVal() - out.alt_poi().getRealValue(v->GetName())) > 1e-12) {
525 match = false;
526 break;
527 } else if (auto cat = dynamic_cast<RooAbsCategory *>(c);
528 cat && cat->getCurrentIndex() ==
529 out.alt_poi().getCatIndex(cat->GetName(), std::numeric_limits<int>().max())) {
530 match = false;
531 break;
532 }
533 }
534 if (!match)
535 continue;
536 if (!p.coords->equals(*out.coords))
537 continue;
538 for (auto c : *p.coords) {
539 if (c->getAttribute("poi")) {
540 continue; // check poi below
541 }
542 if (auto v = dynamic_cast<RooAbsReal *>(c);
543 v && std::abs(v->getVal() - out.coords->getRealValue(v->GetName())) > 1e-12) {
544 match = false;
545 break;
546 } else if (auto cat = dynamic_cast<RooAbsCategory *>(c);
547 cat && cat->getCurrentIndex() ==
548 out.alt_poi().getCatIndex(cat->GetName(), std::numeric_limits<int>().max())) {
549 match = false;
550 break;
551 }
552 }
553 if (!match)
554 continue;
555 // if reached here we can copy over the asimov dataset to save re-generating it
556 // first copy over cfit_alt (if its there) because that is also the same since coords and alt_poi values same
557 if (auto cfit = p.cfit_alt(true)) {
558 out.fAlt_cfit = cfit;
559 }
560 if (p.asimov(true) && p.asimov(true)->fData.first && (!out.asimov(true) || !out.asimov(true)->fData.first)) {
561 out.asimov()->fData = p.asimov(true)->fData;
562 }
563 if (!p.poi().equals(out.poi()))
564 continue;
565 for (auto c : p.poi()) {
566 if (auto v = dynamic_cast<RooAbsReal *>(c);
567 v && std::abs(v->getVal() - out.poi().getRealValue(v->GetName())) > 1e-12) {
568 match = false;
569 break;
570 }
571 }
572 if (match) {
573 // found a duplicate point, return that!
574 return p;
575 }
576 }
577
578 std::string coordString;
579 for (auto a : axes()) {
580 coordString += TString::Format("%s=%g", a->GetName(), out.coords->getRealValue(a->GetName()));
581 coordString += ",";
582 }
583 coordString.erase(coordString.end() - 1);
584
585 ::Info("xRooHypoSpace::AddPoint", "Added new point @ %s", coordString.c_str());
586 return emplace_back(out);
587}
588
590{
591
592 if (!_pdf.get<RooAbsPdf>()) {
593 throw std::runtime_error("Not a pdf");
594 }
595
596 auto pars = _pdf.pars().argList();
597
598 // replace any pars with validity pars and add new pars
599 auto vpars = toArgs(validity);
600 pars.replace(vpars);
601 vpars.remove(pars, true, true);
602 pars.add(vpars);
603
604 if (auto existing = pdf(pars)) {
605 throw std::runtime_error(std::string("Clashing model: ") + existing->GetName());
606 }
607
608 auto myPars = std::shared_ptr<RooArgList>(dynamic_cast<RooArgList *>(pars.snapshot()));
609 myPars->sort();
610
611 pars.remove(*fPars, true, true);
612
613 fPars->addClone(pars);
614
615 fPdfs.insert(std::make_pair(myPars, std::make_shared<xRooNode>(_pdf)));
616
617 return true;
618}
619
621{
622 // determine which pars are the minimal set to distinguish all points in the space
623 RooArgList out;
624 out.setName("axes");
625
626 out.add(*std::unique_ptr<RooAbsCollection>(
627 fPars->selectByAttrib("axis", true))); // start with any pars explicitly designated as axes
628
629 bool clash;
630 do {
631 clash = false;
632
633 std::set<std::vector<double>> coords;
634 for (auto &p : *this) {
635 std::vector<double> p_coords;
636 for (auto o : out) {
637 auto _v = dynamic_cast<RooRealVar *>(p.coords->find(o->GetName()));
638 p_coords.push_back(
639 (_v && _v->isConstant())
640 ? _v->getVal()
641 : std::numeric_limits<double>::infinity()); // non-const coords are treating as non-existent
642 // p_coords.push_back(p.coords->getRealValue(o->GetName(), std::numeric_limits<double>::quiet_NaN()));
643 }
644 if (coords.find(p_coords) != coords.end()) {
645 clash = true;
646 break;
647 }
648 coords.insert(p_coords);
649 }
650
651 if (clash) {
652 // add next best coordinate
653 std::map<std::string, std::unordered_set<double>> values;
654 for (auto &par : *pars()) {
655 if (out.find(*par))
656 continue;
657 for (auto p : *this) {
658 auto _v = dynamic_cast<RooRealVar *>(p.coords->find(par->GetName()));
659 values[par->GetName()].insert(
660 (_v && _v->isConstant())
661 ? _v->getVal()
662 : std::numeric_limits<double>::infinity()); // non-const coords are treating as non-existent
663 // values[par->GetName()].insert(
664 // p.coords->getRealValue(par->GetName(), std::numeric_limits<double>::quiet_NaN()));
665 }
666 }
667
668 std::string bestVar;
669 size_t maxDiff = 0;
670 bool isPOI = false;
671 for (auto &[k, v] : values) {
672 if (v.size() > maxDiff || (v.size() == maxDiff && !isPOI && pars()->find(k.c_str())->getAttribute("poi"))) {
673 bestVar = k;
674 isPOI = pars()->find(k.c_str())->getAttribute("poi");
675 maxDiff = std::max(maxDiff, v.size());
676 }
677 }
678 if (bestVar.empty()) {
679 break;
680 }
681
682 out.add(*pars()->find(bestVar.c_str()));
683 }
684 } while (clash);
685
686 // ensure poi are at the end
687 std::unique_ptr<RooAbsCollection> poi(out.selectByAttrib("poi", true));
688 out.remove(*poi);
689 out.add(*poi);
690
691 return out;
692}
693
695{
696 RooArgList out;
697 out.setName("poi");
698 out.add(*std::unique_ptr<RooAbsCollection>(pars()->selectByAttrib("poi", true)));
699 return out;
700}
701
703{
704
705 if (!gDirectory)
706 return;
707 auto dir = gDirectory->GetDirectory(apath);
708 if (!dir) {
709 // try open file first
710 TString s(apath);
711 auto f = TFile::Open(s.Contains(":") ? TString(s(0, s.Index(":"))) : s);
712 if (f) {
713 if (!s.Contains(":"))
714 s += ":";
715 dir = gDirectory->GetDirectory(s);
716 if (dir) {
717 LoadFits(s);
718 return;
719 }
720 }
721 if (!dir) {
722 Error("LoadFits", "Path not found %s", apath);
723 return;
724 }
725 }
726
727 // assume for now all fits in given path will have the same pars
728 // so can just look at the float and const pars of first fit result to get all of them
729 // tuple is: parName, parValue, parAltValue (blank if nan)
730 // key represents the ufit values, value represents the sets of poi for the available cfits (subfits of the ufit)
731
732 std::map<std::set<std::tuple<std::string, double, std::string>>, std::set<std::set<std::string>>> cfits;
733 std::set<std::string> allpois;
734
735 int nFits = 0;
736 std::function<void(TDirectory *)> processDir;
737 processDir = [&](TDirectory *_dir) {
738 std::cout << "Processing " << _dir->GetName() << std::endl;
739 if (auto keys = _dir->GetListOfKeys(); keys) {
740 // first check if dir doesn't contain any RooLinkedList ... this identifies it as not an nll dir
741 // so treat any sub-dirs as new nll
742 bool isNllDir = false;
743 for (auto &&k : *keys) {
744 TKey *key = dynamic_cast<TKey *>(k);
745 if (strcmp(key->GetClassName(), "RooLinkedList") == 0) {
746 isNllDir = true;
747 break;
748 }
749 }
750
751 for (auto &&k : *keys) {
752 if (auto subdir = _dir->GetDirectory(k->GetName()); subdir) {
753 if (!isNllDir) {
754 LoadFits(subdir->GetPath());
755 } else {
757 }
758 continue;
759 }
760 auto cl = TClass::GetClass((static_cast<TKey *>(k))->GetClassName());
761 if (cl->InheritsFrom("RooFitResult")) {
762 if (auto cachedFit = _dir->Get<RooFitResult>(k->GetName()); cachedFit) {
763 nFits++;
764 if (nFits == 1) {
765 // for first fit add any missing float pars
766 std::unique_ptr<RooAbsCollection> snap(cachedFit->floatParsFinal().snapshot());
767 snap->remove(*fPars, true, true);
768 fPars->addClone(*snap);
769 // add also the non-string const pars
770 for (auto &p : cachedFit->constPars()) {
771 if (p->getAttribute("global"))
772 continue; // don't consider globals
773 auto v = dynamic_cast<RooAbsReal *>(p);
774 if (!v) {
775 continue;
776 };
777 if (!fPars->contains(*v))
778 fPars->addClone(*v);
779 }
780 }
781 // get names of all the floats
782 std::set<std::string> floatPars;
783 for (auto &p : cachedFit->floatParsFinal())
784 floatPars.insert(p->GetName());
785 // see if
786
787 // build a set of the const par values
788 std::set<std::tuple<std::string, double, std::string>> constPars;
789 for (auto &p : cachedFit->constPars()) {
790 if (p->getAttribute("global"))
791 continue; // don't consider globals when looking for cfits
792 auto v = dynamic_cast<RooAbsReal *>(p);
793 if (!v) {
794 continue;
795 };
796 constPars.insert(
797 std::make_tuple(v->GetName(), v->getVal(),
798 v->getStringAttribute("altVal") ? v->getStringAttribute("altVal") : ""));
799 }
800 // now see if this is a subset of any existing cfit ...
801 for (auto &&[key, value] : cfits) {
802 if (constPars == key)
803 continue; // ignore cases where we already recorded this list of constPars
804 if (std::includes(constPars.begin(), constPars.end(), key.begin(), key.end())) {
805 // usual case ... cachedFit has more constPars than one of the fits we have already encountered
806 // (the ufit)
807 // => cachedFit is a cfit of key fr ...
808 std::set<std::string> pois;
809 for (auto &&par : constPars) {
810 if (key.find(par) == key.end()) {
811 pois.insert(std::get<0>(par));
812 allpois.insert(std::get<0>(par));
813 }
814 }
815 if (!pois.empty()) {
816 cfits[constPars].insert(pois);
817 // std::cout << cachedFit->GetName() << " ";
818 // for(auto ff: constPars) std::cout << ff.first << "=" <<
819 // ff.second << " "; std::cout << std::endl;
820 }
821 }
822 /* FOR NOW we will skip cases where we encounter the cfit before the ufit - usually should eval the
823 ufit first
824 * else if (std::includes(key.begin(), key.end(), constPars.begin(), constPars.end())) {
825 // constPars are subset of key
826 // => key is a ufit of the cachedFit
827 // add all par names of key that aren't in constPars ... these are the poi
828 std::set<std::string> pois;
829 for (auto &&par: key) {
830 if (constPars.find(par) == constPars.end()) {
831 pois.insert(std::get<0>(par));
832 allpois.insert(std::get<0>(par));
833 }
834 }
835 if (!pois.empty()) {
836 std::cout << "found cfit BEFORE ufit??" << std::endl;
837 value.insert(pois);
838 }
839 } */
840 }
841 // ensure that this combination of constPars has entry in map,
842 // even if it doesn't end up with any poi identified from cfits to it
843 cfits[constPars];
844 delete cachedFit;
845 }
846 }
847 }
848 }
849 };
850 processDir(dir);
851 ::Info("xRooHypoSpace::xRooHypoSpace", "%s - Loaded %d fits", apath, nFits);
852
853 if (allpois.size() == 1) {
854 ::Info("xRooHypoSpace::xRooHypoSpace", "Detected POI: %s", allpois.begin()->c_str());
855
856 auto nll = std::make_shared<xRooNLLVar>(nullptr, nullptr);
857 auto dummyNll = std::make_shared<RooRealVar>(apath, "Dummy NLL", std::numeric_limits<double>::quiet_NaN());
858 nll->std::shared_ptr<RooAbsReal>::operator=(dummyNll);
859 dummyNll->setAttribute("readOnly");
860 // add pars as 'servers' on the dummy NLL
861 if (fPars) {
862 for (auto &&p : *fPars) {
863 dummyNll->addServer(
864 *p); // this is ok provided fPars (i.e. hypoSpace) stays alive as long as the hypoPoint ...
865 }
866 // flag poi
867 for (auto &p : allpois) {
868 fPars->find(p.c_str())->setAttribute("poi", true);
869 }
870 }
871 nll->reinitialize(); // triggers filling of par lists etc
872
873 for (auto &&[key, value] : cfits) {
874 if (value.find(allpois) != value.end()) {
875 // get the value of the poi in the key set
876 auto _coords = std::make_shared<RooArgSet>();
877 for (auto &k : key) {
878 auto v = _coords->addClone(RooRealVar(std::get<0>(k).c_str(), std::get<0>(k).c_str(), std::get<1>(k)));
879 v->setAttribute("poi", allpois.find(std::get<0>(k)) != allpois.end());
880 if (!std::get<2>(k).empty()) {
881 v->setStringAttribute("altVal", std::get<2>(k).c_str());
882 }
883 }
885 // hp.fPOIName = allpois.begin()->c_str();
886 // hp.fNullVal = _coords->getRealValue(hp.fPOIName.c_str());
887 hp.coords = _coords;
888 hp.nllVar = nll;
889
890 // auto altVal =
891 // hp.null_cfit()->constPars().find(hp.fPOIName.c_str())->getStringAttribute("altVal");
892 // if(altVal) hp.fAltVal = TString(altVal).Atof();
893 // else hp.fAltVal = std::numeric_limits<double>::quiet_NaN();
894
895 // decide based on values
896 if (std::isnan(hp.fAltVal())) {
898 } else if (hp.fNullVal() >= hp.fAltVal()) {
900 } else {
902 }
903
904 emplace_back(hp);
905 }
906 }
907 } else if (nFits > 0) {
908 std::cout << "possible POI: ";
909 for (auto p : allpois)
910 std::cout << p << ",";
911 std::cout << std::endl;
912 }
913}
914
916{
917
918 auto _axes = axes();
919
920 size_t badFits = 0;
921
922 for (size_t i = 0; i < size(); i++) {
923 std::cout << i << ") ";
924 for (auto a : _axes) {
925 if (a != _axes.first())
926 std::cout << ",";
927 std::cout << a->GetName() << "="
928 << at(i).coords->getRealValue(a->GetName(), std::numeric_limits<double>::quiet_NaN());
929 }
930 std::cout << " status=[ufit:";
931 auto ufit = const_cast<xRooHypoPoint &>(at(i)).ufit(true);
932 if (!ufit) {
933 std::cout << "-";
934 } else {
935 std::cout << ufit->status();
936 badFits += (xRooNLLVar::xRooHypoPoint::allowedStatusCodes.count(ufit->status()) == 0);
937 }
938 std::cout << ",cfit_null:";
939 auto cfit = const_cast<xRooHypoPoint &>(at(i)).cfit_null(true);
940 if (!cfit) {
941 std::cout << "-";
942 } else {
943 std::cout << cfit->status();
944 badFits += (xRooNLLVar::xRooHypoPoint::allowedStatusCodes.count(cfit->status()) == 0);
945 }
946 std::cout << ",cfit_alt:";
947 auto afit = const_cast<xRooHypoPoint &>(at(i)).cfit_alt(true);
948 if (!afit) {
949 std::cout << "-";
950 } else {
951 std::cout << afit->status();
953 }
954 if (auto asiPoint = const_cast<xRooHypoPoint &>(at(i)).asimov(true)) {
955 std::cout << ",asimov.ufit:";
956 auto asi_ufit = asiPoint->ufit(true);
957 if (!asi_ufit) {
958 std::cout << "-";
959 } else {
960 std::cout << asi_ufit->status();
962 }
963 std::cout << ",asimov.cfit_null:";
964 auto asi_cfit = asiPoint->cfit_null(true);
965 if (!asi_cfit) {
966 std::cout << "-";
967 } else {
968 std::cout << asi_cfit->status();
970 }
971 }
972 std::cout << "]";
973 auto sigma_mu = const_cast<xRooHypoPoint &>(at(i)).sigma_mu(true);
974 if (!std::isnan(sigma_mu.first)) {
975 std::cout << " sigma_mu=" << sigma_mu.first;
976 if (sigma_mu.second)
977 std::cout << " +/- " << sigma_mu.second;
978 }
979 std::cout << std::endl;
980 }
981 std::cout << "--------------------------" << std::endl;
982 std::cout << "Number of bad fits: " << badFits << std::endl;
983}
984
985std::shared_ptr<TGraphErrors> xRooNLLVar::xRooHypoSpace::graph(
986 const char *opt /*, const std::function<void(xRooNLLVar::xRooHypoSpace*)>& progress*/) const
987{
988
989 TString sOpt(opt);
990 sOpt.ToLower();
991
992 bool doCLs = sOpt.Contains("cls");
993 bool readOnly = sOpt.Contains("readonly");
994 bool visualize = sOpt.Contains("visualize") && !readOnly;
995
996 double nSigma =
997 (sOpt.Contains("exp"))
998 ? (TString(sOpt(sOpt.Index("exp") + 3,
999 sOpt.Index(" ", sOpt.Index("exp")) == -1 ? sOpt.Length() : sOpt.Index(" ", sOpt.Index("exp"))))
1000 .Atof())
1001 : std::numeric_limits<double>::quiet_NaN();
1002 bool expBand =
1003 !std::isnan(nSigma) && nSigma && !(sOpt(sOpt.Index("exp") + 3) == '+' || sOpt(sOpt.Index("exp") + 3) == '-');
1004
1005 auto _axes = axes();
1006 if (_axes.size() != 1)
1007 return nullptr;
1008
1009 auto out = std::make_shared<TGraphErrors>();
1010 out->SetName(GetName());
1011 out->SetEditable(false);
1012 const char *sCL = (doCLs) ? "CLs" : "null";
1013
1014 TString title =
1015 TString::Format("%s;%s;p_{%s}",
1016 (std::isnan(nSigma))
1017 ? "Observed"
1018 : (!nSigma ? "Expected"
1019 : TString::Format("%s%d#sigma Expected",
1020 expBand || !nSigma ? "" : ((nSigma < 0) ? "-" : "+"), int(nSigma))
1021 .Data()),
1022 _axes.at(0)->GetTitle(), sCL);
1023
1024 if (std::isnan(nSigma)) {
1025 out->SetNameTitle(TString::Format("obs_p%s", sCL), title);
1026 out->SetMarkerStyle(20);
1027 out->SetMarkerSize(0.5);
1028 if (sOpt.Contains("ts")) {
1029 out->SetNameTitle("obs_ts", TString::Format("Observed;%s;%s", _axes.at(0)->GetTitle(),
1030 (empty() ? "" : front().tsTitle(true).Data())));
1031 }
1032 } else {
1033 out->SetNameTitle(TString::Format("exp%d_p%s", int(nSigma), sCL), title);
1034 out->SetMarkerStyle(0);
1035 out->SetMarkerSize(0);
1036 out->SetLineStyle(2 + int(nSigma));
1037 if (expBand && nSigma) {
1038 out->SetFillColor((nSigma == 2) ? kYellow : kGreen);
1039 // out->SetFillStyle(3005);
1040 out->SetLineStyle(0);
1041 out->SetLineWidth(0);
1042 auto x = out->Clone("up");
1043 x->SetBit(kCanDelete);
1044 dynamic_cast<TAttFill *>(x)->SetFillStyle(nSigma == 2 ? 3005 : 3004);
1045 dynamic_cast<TAttFill *>(x)->SetFillColor(kBlack);
1046 out->GetListOfFunctions()->Add(x, "F");
1047 x = out->Clone("down");
1048 x->SetBit(kCanDelete);
1049 // dynamic_cast<TAttFill*>(x)->SetFillColor((nSigma==2) ? kYellow : kGreen);
1050 // dynamic_cast<TAttFill*>(x)->SetFillStyle(1001);
1051 out->GetListOfFunctions()->Add(x, "F");
1052 }
1053 if (sOpt.Contains("ts")) {
1054 out->SetNameTitle(TString::Format("exp_ts%d", int(nSigma)),
1055 TString::Format("Expected;%s;%s", _axes.at(0)->GetTitle(), front().tsTitle(true).Data()));
1056 }
1057 }
1058
1059 auto badPoints = [&]() {
1060 auto badPoints2 = dynamic_cast<TGraph *>(out->GetListOfFunctions()->FindObject("badPoints"));
1061 if (!badPoints2) {
1062 badPoints2 = new TGraph;
1063 badPoints2->SetBit(kCanDelete);
1064 badPoints2->SetName("badPoints");
1065 badPoints2->SetMarkerStyle(5);
1066 badPoints2->SetMarkerColor(std::isnan(nSigma) ? kRed : kBlue);
1067 badPoints2->SetMarkerSize(1);
1068 out->GetListOfFunctions()->Add(badPoints2, "P");
1069 }
1070 return badPoints2;
1071 };
1072 int nPointsDown = 0;
1073 bool above = true;
1074 TStopwatch s;
1075 s.Start();
1076 size_t nDone = 0;
1077 for (auto &p : *this) {
1078 if (s.RealTime() > 5) {
1079 if (visualize) {
1080 // draw readonly version of the graph
1081 auto gra = graph(sOpt + " readOnly");
1082 if (gra && gra->GetN()) {
1083 if (!gPad && gROOT->GetSelectedPad())
1084 gROOT->GetSelectedPad()->cd();
1085 if (gPad)
1086 gPad->Clear();
1087 gra->DrawClone(expBand ? "AF" : "ALP")->SetBit(kCanDelete);
1088#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 30, 00)
1089 if (auto pad = gROOT->GetSelectedPad()) {
1090 pad->GetCanvas()->ResetUpdated(); // stops previous canvas being replaced in a jupyter notebook
1091 }
1092#endif
1094 }
1095 } else {
1096 ::Info("xRooHypoSpace::graph", "Completed %d/%d points for %s", int(nDone), int(size()), sOpt.Data());
1097 }
1098 s.Start();
1099 } else {
1100 s.Continue();
1101 }
1102 double _x = p.coords->getRealValue(_axes.at(0)->GetName(), std::numeric_limits<double>::quiet_NaN());
1103 auto pval = const_cast<xRooHypoPoint &>(p).getVal(sOpt);
1104 auto idx = out->GetN() - nPointsDown;
1105
1106 if (std::isnan(pval.first)) {
1107 if (p.status() != 0) { // if status is 0 then bad pval is really just absence of fits, not bad fits
1108 badPoints()->SetPoint(badPoints()->GetN(), _x, 0);
1109 }
1110 } else {
1111 out->InsertPointBefore(idx, _x, pval.first);
1112 out->SetPointError(idx, 0, pval.second);
1113 }
1114
1115 if (expBand && nSigma) {
1116 TString sOpt2 = sOpt;
1117 sOpt2.ReplaceAll("exp", "exp-");
1118 pval = const_cast<xRooHypoPoint &>(p).getVal(sOpt2);
1119 if (std::isnan(pval.first)) {
1120 if (p.status() != 0) { // if status is 0 then bad pval is really just absence of fits, not bad fits
1121 badPoints()->SetPoint(badPoints()->GetN(), _x, 0);
1122 }
1123 } else {
1124 out->InsertPointBefore(idx + 1, _x, pval.first);
1125 out->SetPointError(idx + 1, 0, pval.second);
1126 nPointsDown++;
1127 if (out->GetPointY(idx) < pval.first)
1128 above = false; // the -sigma points are actually above +sigma
1129 }
1130 }
1131 nDone++;
1132 }
1133
1134 if (out->GetN() == 0)
1135 return out;
1136
1137 if (!expBand) {
1138 out->Sort();
1139 if (out->GetListOfFunctions()->FindObject("badPoints")) {
1140 // try to interpolate the points
1141 for (int i = 0; i < badPoints()->GetN(); i++) {
1142 badPoints()->SetPointY(i, out->Eval(badPoints()->GetPointX(i)));
1143 }
1144 }
1145 } else {
1146 out->Sort(&TGraph::CompareX, true, 0, out->GetN() - nPointsDown - 1); // sort first half
1147 out->Sort(&TGraph::CompareX, false, out->GetN() - nPointsDown, out->GetN() - 1); // reverse sort second half
1148
1149 // now populate the up and down values
1150 auto up = dynamic_cast<TGraph *>(out->GetListOfFunctions()->FindObject("up"));
1151 auto down = dynamic_cast<TGraph *>(out->GetListOfFunctions()->FindObject("down"));
1152
1153 for (int i = 0; i < out->GetN(); i++) {
1154 if (i < out->GetN() - nPointsDown) {
1155 up->SetPoint(up->GetN(), out->GetPointX(i), out->GetPointY(i) + out->GetErrorY(i) * (above ? 1. : -1.));
1156 down->SetPoint(down->GetN(), out->GetPointX(i), out->GetPointY(i) - out->GetErrorY(i) * (above ? 1. : -1.));
1157 } else {
1158 up->SetPoint(up->GetN(), out->GetPointX(i), out->GetPointY(i) - out->GetErrorY(i) * (above ? 1. : -1.));
1159 down->SetPoint(down->GetN(), out->GetPointX(i), out->GetPointY(i) + out->GetErrorY(i) * (above ? 1. : -1.));
1160 }
1161 }
1162 }
1163
1164 if (visualize) {
1165 // draw result
1166 if (!gPad && gROOT->GetSelectedPad())
1167 gROOT->GetSelectedPad()->cd();
1168 if (gPad)
1169 gPad->Clear();
1170 out->DrawClone(expBand ? "AF" : "ALP")->SetBit(kCanDelete);
1171#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 30, 00)
1172 if (auto pad = gROOT->GetSelectedPad()) {
1173 pad->GetCanvas()->ResetUpdated(); // stops previous canvas being replaced in a jupyter notebook
1174 }
1175#endif
1177 }
1178
1179 return out;
1180}
1181
1182std::shared_ptr<TMultiGraph> xRooNLLVar::xRooHypoSpace::graphs(const char *opt)
1183{
1184 TString sOpt(opt);
1185 sOpt.ToLower();
1186 std::shared_ptr<TMultiGraph> out;
1187 if (sOpt.Contains("pcls") || sOpt.Contains("pnull") || sOpt.Contains("ts")) {
1188
1189 bool visualize = sOpt.Contains("visualize");
1190 sOpt.ReplaceAll("visualize", "");
1191
1192 auto exp2 = graph(sOpt + " exp2");
1193 auto exp1 = graph(sOpt + " exp1");
1194 auto exp = graph(sOpt + " exp");
1195 bool doObs = true;
1196 // for(auto& hp : *this) { if(hp.isExpected) {doObs=false; break;} }
1197 auto obs = (doObs) ? graph(sOpt) : nullptr;
1198
1199 out = std::make_shared<TMultiGraph>(GetName(), GetTitle());
1200 if (exp2 && exp2->GetN() > 1)
1201 out->Add(static_cast<TGraph *>(exp2->Clone()), "FP");
1202 if (exp1 && exp1->GetN() > 1)
1203 out->Add(static_cast<TGraph *>(exp1->Clone()), "FP");
1204 if (exp && exp->GetN() > 1)
1205 out->Add(static_cast<TGraph *>(exp->Clone()), "LP");
1206 if (obs && obs->GetN() > 1)
1207 out->Add(static_cast<TGraph *>(obs->Clone()), "LP");
1208
1209 if (!out->GetListOfGraphs()) {
1210 return nullptr;
1211 }
1212
1213 TGraph *testedPoints = nullptr;
1214 if (sOpt.Contains("pcls")) {
1215 TGraph *line = new TGraph;
1216 line->SetName("alpha");
1217 line->SetLineStyle(2);
1218 line->SetEditable(false);
1219 line->SetPoint(line->GetN(), out->GetHistogram()->GetXaxis()->GetXmin() - 10, 0.05);
1220 testedPoints = new TGraph;
1221 testedPoints->SetName("hypoPoints");
1222 testedPoints->SetEditable(false);
1223 testedPoints->SetMarkerStyle(24);
1224 testedPoints->SetMarkerSize(0.4); // use line to indicate tested points
1225 if (exp) {
1226 for (int i = 0; i < exp->GetN(); i++) {
1227 testedPoints->SetPoint(testedPoints->GetN(), exp->GetPointX(i), 0.05);
1228 }
1229 }
1230 line->SetPoint(line->GetN(), out->GetHistogram()->GetXaxis()->GetXmax() + 10, 0.05);
1232 out->GetListOfFunctions()->Add(line, "L");
1233 }
1234 if (exp) {
1235 out->GetHistogram()->GetXaxis()->SetTitle(exp->GetHistogram()->GetXaxis()->GetTitle());
1236 out->GetHistogram()->GetYaxis()->SetTitle(exp->GetHistogram()->GetYaxis()->GetTitle());
1237 }
1238 auto leg = new TLegend(1. - gStyle->GetPadRightMargin() - 0.3, 1. - gStyle->GetPadTopMargin() - 0.35,
1239 1. - gStyle->GetPadRightMargin() - 0.05, 1. - gStyle->GetPadTopMargin() - 0.05);
1240 leg->SetName("legend");
1241 leg->SetBit(kCanDelete);
1242
1243 out->GetListOfFunctions()->Add(leg);
1244 // out->GetListOfFunctions()->Add(out->GetHistogram()->Clone(".axis"),"sameaxis"); // redraw axis
1245
1246 for (auto g : *out->GetListOfGraphs()) {
1247 if (auto o = dynamic_cast<TGraph *>(g)->GetListOfFunctions()->FindObject("down")) {
1248 leg->AddEntry(o, "", "F");
1249 } else {
1250 leg->AddEntry(g, "", "LPE");
1251 }
1252 }
1253
1254 if (sOpt.Contains("pcls")) {
1255 // add current limit estimates to legend
1256 if (exp2 && exp2->GetN() > 1) {
1257 auto l = xRooFit::matchPrecision(GetLimit(*graph(sOpt + "exp-2")));
1258 leg->AddEntry((TObject *)nullptr, TString::Format("-2#sigma: %g +/- %g", l.first, l.second), "");
1259 }
1260 if (exp1 && exp1->GetN() > 1) {
1261 auto l = xRooFit::matchPrecision(GetLimit(*graph(sOpt + "exp-1")));
1262 leg->AddEntry((TObject *)nullptr, TString::Format("-1#sigma: %g +/- %g", l.first, l.second), "");
1263 }
1264 if (exp && exp->GetN() > 1) {
1265 auto l = xRooFit::matchPrecision(GetLimit(*exp));
1266 leg->AddEntry((TObject *)nullptr, TString::Format("0#sigma: %g +/- %g", l.first, l.second), "");
1267 }
1268 if (exp1 && exp1->GetN() > 1) {
1269 auto l = xRooFit::matchPrecision(GetLimit(*graph(sOpt + "exp+1")));
1270 leg->AddEntry((TObject *)nullptr, TString::Format("+1#sigma: %g +/- %g", l.first, l.second), "");
1271 }
1272 if (exp2 && exp2->GetN() > 1) {
1273 auto l = xRooFit::matchPrecision(GetLimit(*graph(sOpt + "exp+2")));
1274 leg->AddEntry((TObject *)nullptr, TString::Format("+2#sigma: %g +/- %g", l.first, l.second), "");
1275 }
1276 if (obs && obs->GetN() > 1) {
1277 auto l = xRooFit::matchPrecision(GetLimit(*obs));
1278 leg->AddEntry((TObject *)nullptr, TString::Format("Observed: %g +/- %g", l.first, l.second), "");
1279 }
1280 }
1281 if (testedPoints)
1282 out->Add(testedPoints, "P");
1283
1284 if (visualize) {
1285 if (!gPad && gROOT->GetSelectedPad())
1286 gROOT->GetSelectedPad()->cd();
1287 if (gPad)
1288 gPad->Clear();
1289 auto gra2 = static_cast<TMultiGraph *>(out->DrawClone("A"));
1290 gra2->SetBit(kCanDelete);
1291 if (sOpt.Contains("pcls") || sOpt.Contains("pnull")) {
1292 gra2->GetHistogram()->SetMinimum(1e-6);
1293 }
1294 if (gPad) {
1295 gPad->RedrawAxis();
1296 gPad->GetCanvas()->Paint();
1297 gPad->GetCanvas()->Update();
1298#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 30, 00)
1299 gPad->GetCanvas()->ResetUpdated(); // stops previous canvas being replaced in a jupyter notebook
1300#endif
1301 }
1303 }
1304 }
1305
1306 return out;
1307}
1308
1310{
1311
1312 if (std::isnan(target)) {
1313 target = 1. - gEnv->GetValue("xRooHypoSpace.CL", 95.) / 100.;
1314 }
1315
1316 auto gr = std::make_shared<TGraph>(pValues);
1317 // remove any nan points and duplicates
1318 int i = 0;
1319 std::set<double> existingX;
1320 while (i < gr->GetN()) {
1321 if (std::isnan(gr->GetPointY(i))) {
1322 gr->RemovePoint(i);
1323 } else if (existingX.find(gr->GetPointX(i)) != existingX.end()) {
1324 gr->RemovePoint(i);
1325 } else {
1326 existingX.insert(gr->GetPointX(i));
1327 // convert to log ....
1328 gr->SetPointY(i, log(std::max(gr->GetPointY(i), 1e-10)));
1329 i++;
1330 }
1331 }
1332
1333 gr->Sort();
1334
1335 // simple linear extrapolation to critical value ... return nan if problem
1336 if (gr->GetN() < 2) {
1337 return std::pair<double, double>(std::numeric_limits<double>::quiet_NaN(), 0);
1338 }
1339
1340 double alpha = log(target);
1341
1342 bool above = gr->GetPointY(0) > alpha;
1343 for (int ii = 1; ii < gr->GetN(); ii++) {
1344 if ((above && (gr->GetPointY(ii) <= alpha)) || (!above && (gr->GetPointY(ii) >= alpha))) {
1345 // found the limit ... return linearly extrapolated point
1346 double lim = gr->GetPointX(ii - 1) + (gr->GetPointX(ii) - gr->GetPointX(ii - 1)) *
1347 (alpha - gr->GetPointY(ii - 1)) /
1348 (gr->GetPointY(ii) - gr->GetPointY(ii - 1));
1349 // use points either side as error
1350 double err = std::max(lim - gr->GetPointX(ii - 1), gr->GetPointX(ii) - lim);
1351 // give err negative sign to indicate if error due to negative side
1352 if ((lim - gr->GetPointX(ii - 1)) > (gr->GetPointX(ii) - lim))
1353 err *= -1;
1354 return std::pair(lim, err);
1355 }
1356 }
1357 // if reach here need to extrapolate ...
1358 if ((above && gr->GetPointY(gr->GetN() - 1) <= gr->GetPointY(0)) ||
1359 (!above && gr->GetPointY(gr->GetN() - 1) >= gr->GetPointY(0))) {
1360 // extrapolating above based on last two points
1361 // in fact, if 2nd last point is a p=1 (log(p)=0) then go back
1362 int offset = 2;
1363 while (offset < gr->GetN() && gr->GetPointY(gr->GetN() - offset) == 0)
1364 offset++;
1365 double x1 = gr->GetPointX(gr->GetN() - offset);
1366 double y1 = gr->GetPointY(gr->GetN() - offset);
1367 double m = (gr->GetPointY(gr->GetN() - 1) - y1) / (gr->GetPointX(gr->GetN() - 1) - x1);
1368 if (m == 0.)
1369 return std::pair(std::numeric_limits<double>::infinity(), std::numeric_limits<double>::infinity());
1370 return std::pair((alpha - y1) / m + x1, std::numeric_limits<double>::infinity());
1371 } else {
1372 // extrapolating below based on first two points
1373 double x1 = gr->GetPointX(0);
1374 double y1 = gr->GetPointY(0);
1375 double m = (gr->GetPointY(1) - y1) / (gr->GetPointX(1) - x1);
1376 if (m == 0.)
1377 return std::pair(-std::numeric_limits<double>::infinity(), -std::numeric_limits<double>::infinity());
1378 return std::pair((alpha - y1) / m + x1, -std::numeric_limits<double>::infinity());
1379 }
1380}
1381
1383{
1384 TString sOpt = TString::Format("p%s", type);
1385 if (std::isnan(nSigma)) {
1386 sOpt += "obs";
1387 } else {
1388 sOpt += TString::Format("exp%s%d", nSigma > 0 ? "+" : "", int(nSigma));
1389 }
1390 return GetLimit(*graph(sOpt + " readonly"));
1391}
1392
1394xRooNLLVar::xRooHypoSpace::findlimit(const char *opt, double relUncert, unsigned int maxTries)
1395{
1396 TString sOpt(opt);
1397 bool visualize = sOpt.Contains("visualize");
1398 sOpt.ReplaceAll("visualize", "");
1399 std::shared_ptr<TGraphErrors> gr = graph(sOpt + " readonly");
1400 if (visualize) {
1401 auto gra = graphs(sOpt.Contains("toys") ? "pcls readonly toys" : "pcls readonly");
1402 if (gra) {
1403 if (!gPad)
1404 gra->Draw(); // in 6.28 DrawClone wont make the gPad defined :( ... so Draw then clear and Draw Clone
1405 gPad->Clear();
1406 gra->DrawClone("A")->SetBit(kCanDelete);
1407 gPad->RedrawAxis();
1408 gra->GetHistogram()->SetMinimum(1e-9);
1409 gra->GetHistogram()->GetYaxis()->SetRangeUser(1e-9, 1);
1410 gPad->Modified();
1411#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 30, 00)
1412 gPad->GetCanvas()->ResetUpdated(); // stops previous canvas being replaced in a jupyter notebook
1413#endif
1415 }
1416 }
1417
1418 // resync parameter boundaries from nlls (may have been modified by fits)
1419 for (auto p : axes()) {
1420 for (auto &[pdf, nll] : fNlls) {
1421 if (auto _v = dynamic_cast<RooRealVar *>(nll->pars()->find(*p))) {
1422 dynamic_cast<RooRealVar *>(p)->setRange(_v->getMin(), _v->getMax());
1423 }
1424 }
1425 }
1426
1427 if (!gr || gr->GetN() < 2) {
1428 auto v = (axes().empty()) ? nullptr : dynamic_cast<RooRealVar *>(*axes().rbegin());
1429 if (!v)
1430 return std::pair(std::numeric_limits<double>::quiet_NaN(), 0.);
1431 double muMax = std::min(std::min(v->getMax("physical"), v->getMax()), v->getMax("scan"));
1432 double muMin = std::max(std::max(v->getMin("physical"), v->getMin()), v->getMin("scan"));
1433 if (!gr || gr->GetN() < 1) {
1434 if (maxTries == 0 || std::isnan(AddPoint(TString::Format("%s=%g", v->GetName(), muMin)).getVal(sOpt).first)) {
1435 // first point failed ... give up
1436 Error("findlimit", "Problem evaluating %s @ %s=%g", sOpt.Data(), v->GetName(), muMin);
1437 return std::pair(std::numeric_limits<double>::quiet_NaN(), 0.);
1438 }
1439 gr.reset();
1440 return findlimit(opt, relUncert, maxTries - 1); // do this to resync parameter limits
1441 }
1442
1443 // can approximate expected limit using
1444 // mu_hat + sigma_mu*ROOT::Math::gaussian_quantile(1.-alpha/2.,1) for cls
1445 // or mu_hat + sigma_mu*ROOT::Math::gaussian_quantile((1.-alpha),1) for cls+b
1446 // get a very first estimate of sigma_mu from ufit to expected data, take error on mu as sigma_mu
1447 double nextPoint = muMin + (muMax - muMin) / 50;
1448
1449 // if done an expected limit, assume data is like expected and choose expected limit point as first test point
1450 if (sOpt.Contains("obs")) {
1451 TString sOpt2 = sOpt;
1452 sOpt2.ReplaceAll("obs", "exp");
1453 auto expLim = findlimit(sOpt2, std::numeric_limits<double>::infinity(), 0);
1454 if (!std::isnan(expLim.first) && expLim.first < nextPoint)
1455 nextPoint = expLim.first;
1456 }
1457
1458 auto point =
1459 (sOpt.Contains("exp")) ? back().asimov() : std::shared_ptr<xRooHypoPoint>(&back(), [](xRooHypoPoint *) {});
1460 point = nullptr;
1461 if (point && point->ufit()) {
1462 double rough_sigma_mu = point->mu_hat().getError();
1463 double another_estimate = point->mu_hat().getVal() + rough_sigma_mu * ROOT::Math::gaussian_quantile(0.95, 1);
1464 // if (another_estimate < nextPoint) {
1466 ::Info("xRooHypoSpace::findlimit", "Guessing %g based on rough sigma_mu = %g", nextPoint, rough_sigma_mu);
1467 //}
1468 }
1469
1470 if (maxTries == 0 || std::isnan(AddPoint(TString::Format("%s=%g", v->GetName(), nextPoint)).getVal(sOpt).first)) {
1471 // second point failed ... give up
1472 Error("findlimit", "Problem evaluating %s @ %s=%g", sOpt.Data(), v->GetName(), nextPoint);
1473 return std::pair(std::numeric_limits<double>::quiet_NaN(), 0.);
1474 }
1475 gr.reset();
1476 return findlimit(opt, relUncert, maxTries - 1);
1477 }
1478
1479 auto lim = GetLimit(*gr);
1480
1481 if (std::isnan(lim.first)) {
1482 return lim;
1483 }
1484
1485 auto v = dynamic_cast<RooRealVar *>(*axes().rbegin());
1486 double maxMu = std::min(std::min(v->getMax("physical"), v->getMax()), v->getMax("scan"));
1487 double minMu = std::max(std::max(v->getMin("physical"), v->getMin()), v->getMin("scan"));
1488
1489 // static double MIN_LIMIT_UNCERT = 1e-4; // stop iterating once uncert gets this small
1490 if (lim.first > -std::numeric_limits<double>::infinity() && lim.first < std::numeric_limits<double>::infinity() &&
1491 (std::abs(lim.second) <= relUncert * std::abs(lim.first) /* || std::abs(lim.second)<MIN_LIMIT_UNCERT*/))
1492 return lim;
1493
1494 double nextPoint;
1495
1496 if (lim.second == std::numeric_limits<double>::infinity()) {
1497 // limit was found by extrapolating to right
1498 nextPoint = lim.first;
1499 if (nextPoint == std::numeric_limits<double>::infinity() || nextPoint > maxMu) {
1500 nextPoint = gr->GetPointX(gr->GetN() - 1) + (maxMu - minMu) / 50;
1501 }
1502
1503 // prefer extrapolation with sigma_mu, if available, if it takes us further
1504 // as shape of p-value curve is usually
1505 auto point =
1506 (sOpt.Contains("exp")) ? back().asimov() : std::shared_ptr<xRooHypoPoint>(&back(), [](xRooHypoPoint *) {});
1507 point = nullptr;
1508 if (point && point->ufit()) {
1509 double rough_sigma_mu = point->mu_hat().getError();
1510 double another_estimate = point->mu_hat().getVal() + rough_sigma_mu * ROOT::Math::gaussian_quantile(0.95, 1);
1511 // if (another_estimate < nextPoint) {
1513 ::Info("xRooHypoSpace::findlimit", "Guessing %g based on rough sigma_mu = %g", nextPoint, rough_sigma_mu);
1514 //}
1515 }
1516 nextPoint = std::min(nextPoint + nextPoint * relUncert * 0.99, maxMu); // ensure we step over location if possible
1517
1518 if (nextPoint > maxMu)
1519 return lim;
1520 } else if (lim.second == -std::numeric_limits<double>::infinity()) {
1521 // limit from extrapolating to left
1522 nextPoint = lim.first;
1523 if (nextPoint < minMu)
1524 nextPoint = gr->GetPointX(0) - (maxMu - minMu) / 50;
1525 if (nextPoint < minMu)
1526 return lim;
1527 } else {
1528 nextPoint = lim.first + lim.second * relUncert * 0.99;
1529 }
1530
1531 // got here need a new point .... evaluate the estimated lim location +/- the relUncert (signed error takes care of
1532 // direction)
1533
1534 ::Info("xRooHypoSpace::findlimit", "%s -- Testing new point @ %s=%g (delta=%g)", sOpt.Data(), v->GetName(),
1535 nextPoint, lim.second);
1536 if (maxTries == 0 || std::isnan(AddPoint(TString::Format("%s=%g", v->GetName(), nextPoint)).getVal(sOpt).first)) {
1537 if (maxTries == 0) {
1538 Warning("findlimit", "Reached max number of point evaluations");
1539 } else {
1540 Error("findlimit", "Problem evaluating %s @ %s=%g", sOpt.Data(), v->GetName(), nextPoint);
1541 }
1542 return lim;
1543 }
1544 gr.reset();
1545 return findlimit(opt, relUncert, maxTries - 1);
1546}
1547
1549{
1550
1551 TString sOpt(opt);
1552 sOpt.ToLower();
1553
1554 if ((sOpt == "" || sOpt == "same") && !empty()) {
1555 if (front().fPllType == xRooFit::Asymptotics::OneSidedPositive) {
1556 sOpt += "pcls"; // default to showing cls p-value scan if drawing a limit
1557 for (auto &hp : *this) {
1558 if (!hp.nullToys.empty() || !hp.altToys.empty()) {
1559 sOpt += " toys";
1560 break; // default to toys if done toys
1561 }
1562 }
1563 } else if (front().fPllType == xRooFit::Asymptotics::TwoSided) {
1564 sOpt += "ts";
1565 }
1566 }
1567
1568 // split up by ; and call Draw for each (with 'same' appended)
1569 auto _axes = axes();
1570 if (_axes.empty())
1571 return;
1572
1573 if (sOpt == "status") {
1574 // draw the points in the space
1575 if (_axes.size() <= 2) {
1576 TGraphErrors *out = new TGraphErrors;
1577 out->SetBit(kCanDelete);
1578 out->SetName("points");
1579 out->SetMarkerSize(0.5);
1580 TGraph *tsAvail = new TGraph;
1581 tsAvail->SetName("ts");
1582 tsAvail->SetBit(kCanDelete);
1583 tsAvail->SetMarkerStyle(20);
1584 TGraph *expAvail = new TGraph;
1585 expAvail->SetName("exp");
1586 expAvail->SetBit(kCanDelete);
1587 expAvail->SetMarkerStyle(25);
1588 expAvail->SetMarkerSize(out->GetMarkerSize() * 1.5);
1589 TGraph *badPoints = new TGraph;
1590 badPoints->SetName("bad_ufit");
1591 badPoints->SetBit(kCanDelete);
1592 badPoints->SetMarkerStyle(5);
1593 badPoints->SetMarkerColor(kRed);
1594 badPoints->SetMarkerSize(out->GetMarkerSize());
1595 TGraph *badPoints2 = new TGraph;
1596 badPoints2->SetName("bad_cfit_null");
1597 badPoints2->SetBit(kCanDelete);
1598 badPoints2->SetMarkerStyle(2);
1599 badPoints2->SetMarkerColor(kRed);
1600 badPoints2->SetMarkerSize(out->GetMarkerSize());
1601
1602 out->SetTitle(TString::Format("%s;%s;%s", GetTitle(), _axes.at(0)->GetTitle(),
1603 (_axes.size() == 1) ? "" : _axes.at(1)->GetTitle()));
1604 for (auto &p : *this) {
1605 bool _readOnly = p.nllVar ? p.nllVar->get()->getAttribute("readOnly") : false;
1606 if (p.nllVar)
1607 p.nllVar->get()->setAttribute("readOnly", true);
1608 double x = p.coords->getRealValue(_axes.at(0)->GetName());
1609 double y = _axes.size() == 1 ? p.ts_asymp().first : p.coords->getRealValue(_axes.at(1)->GetName());
1610 out->SetPoint(out->GetN(), x, y);
1611 if (!std::isnan(p.ts_asymp().first)) {
1612 if (_axes.size() == 1)
1613 out->SetPointError(out->GetN() - 1, 0, p.ts_asymp().second);
1614 tsAvail->SetPoint(tsAvail->GetN(), x, y);
1615 } else if (p.fUfit && (std::isnan(p.fUfit->minNll()) ||
1616 xRooNLLVar::xRooHypoPoint::allowedStatusCodes.find(p.fUfit->status()) ==
1618 badPoints->SetPoint(badPoints->GetN(), x, y);
1619 } else if (p.fNull_cfit && (std::isnan(p.fNull_cfit->minNll()) ||
1620 xRooNLLVar::xRooHypoPoint::allowedStatusCodes.find(p.fNull_cfit->status()) ==
1622 badPoints2->SetPoint(badPoints2->GetN(), x, y);
1623 }
1624 if (!std::isnan(p.ts_asymp(0).first)) {
1625 expAvail->SetPoint(expAvail->GetN(), x, y);
1626 } else if (p.asimov() && p.asimov()->fUfit &&
1627 (std::isnan(p.asimov()->fUfit->minNll()) ||
1628 xRooNLLVar::xRooHypoPoint::allowedStatusCodes.find(p.asimov()->fUfit->status()) ==
1630
1631 } else if (p.asimov() && p.asimov()->fNull_cfit &&
1632 (std::isnan(p.asimov()->fNull_cfit->minNll()) ||
1633 xRooNLLVar::xRooHypoPoint::allowedStatusCodes.find(p.asimov()->fNull_cfit->status()) ==
1635 }
1636 if (p.nllVar)
1637 p.nllVar->get()->setAttribute("readOnly", _readOnly);
1638 }
1639
1640 if (_axes.size() == 1) {
1641 TGraph tmp;
1642 for (int i = 0; i < out->GetN(); i++) {
1643 if (!std::isnan(out->GetPointY(i)))
1644 tmp.SetPoint(tmp.GetN(), out->GetPointX(i), out->GetPointY(i));
1645 }
1646 auto fixPoints = [&](TGraph *g) {
1647 for (int i = 0; i < g->GetN(); i++) {
1648 if (std::isnan(g->GetPointY(i)))
1649 g->SetPointY(i, std::isnan(tmp.Eval(g->GetPointX(i))) ? 0. : tmp.Eval(g->GetPointX(i)));
1650 }
1651 };
1652 fixPoints(out);
1657 }
1658
1659 out->SetMarkerStyle(4);
1660 out->Draw("AP");
1661 auto leg = new TLegend(1. - gPad->GetRightMargin() - 0.3, 1. - gPad->GetTopMargin() - 0.35,
1662 1. - gPad->GetRightMargin() - 0.05, 1. - gPad->GetTopMargin() - 0.05);
1663 leg->SetName("legend");
1664 leg->AddEntry(out, "Uncomputed", "P");
1665
1666 if (tsAvail->GetN()) {
1667 out->GetListOfFunctions()->Add(tsAvail, "P");
1668 leg->AddEntry(tsAvail, "Computed", "P");
1669 } else {
1670 delete tsAvail;
1671 }
1672 if (expAvail->GetN()) {
1673 out->GetListOfFunctions()->Add(expAvail, "P");
1674 leg->AddEntry(expAvail, "Expected computed", "P");
1675 } else {
1676 delete expAvail;
1677 }
1678 if (badPoints->GetN()) {
1679 out->GetListOfFunctions()->Add(badPoints, "P");
1680 leg->AddEntry(badPoints, "Bad ufit", "P");
1681 } else {
1682 delete badPoints;
1683 }
1684 if (badPoints2->GetN()) {
1685 out->GetListOfFunctions()->Add(badPoints2, "P");
1686 leg->AddEntry(badPoints2, "Bad null cfit", "P");
1687 } else {
1688 delete badPoints2;
1689 }
1690 leg->SetBit(kCanDelete);
1691 leg->Draw();
1692 // if(_axes.size()==1) out->GetHistogram()->GetYaxis()->SetRangeUser(0,1);
1693 gPad->SetGrid(true, _axes.size() > 1);
1694 if (_axes.size() == 1)
1695 gPad->SetLogy(false);
1696 }
1697
1699
1700 return;
1701 }
1702
1703 if (sOpt.Contains("pcls") || sOpt.Contains("pnull") || sOpt.Contains("ts")) {
1704 auto gra = graphs(sOpt + " readonly");
1705 if (!gPad && gROOT->GetSelectedPad())
1706 gROOT->GetSelectedPad()->cd();
1707 if (!sOpt.Contains("same") && gPad) {
1708 gPad->Clear();
1709 }
1710 if (gra) {
1711 auto gra2 = static_cast<TMultiGraph *>(gra->DrawClone(sOpt.Contains("same") ? "" : "A"));
1712 gra2->SetBit(kCanDelete);
1713 if (sOpt.Contains("pcls") || sOpt.Contains("pnull")) {
1714 gra2->GetHistogram()->SetMinimum(1e-6);
1715 }
1716 if (gPad) {
1717 gPad->RedrawAxis();
1718 gPad->GetCanvas()->Paint();
1719 gPad->GetCanvas()->Update();
1720#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 30, 00)
1721 gPad->GetCanvas()->ResetUpdated(); // stops previous canvas being replaced in a jupyter notebook
1722#endif
1724 }
1725 }
1726 if (!sOpt.Contains("same") && gPad) {
1727 // auto mg = static_cast<TMultiGraph*>(gPad->GetPrimitive(gra->GetName()));
1728 // mg->GetHistogram()->SetMinimum(1e-9);
1729 // mg->GetHistogram()->GetYaxis()->SetRangeUser(1e-9,1);
1730 // gPad->SetGrid(0, 0);
1731 // gPad->SetLogy(1);
1732 }
1733
1735
1736 return;
1737 }
1738
1739 // graphs("ts visualize");return;
1740
1741 TGraphErrors *out = new TGraphErrors;
1742 out->SetName(GetName());
1743
1744 TString title = (!axes().empty()) ? TString::Format(";%s", axes().first()->GetTitle()) : "";
1745
1747 if (!empty() && axes().size() == 1) {
1748 for (auto &p : *this) {
1749 if (p.fPllType != xRooFit::Asymptotics::TwoSided) {
1750 pllType = p.fPllType;
1751 }
1752 }
1753 title += ";";
1754 title += front().tsTitle(true);
1755 }
1756
1757 out->SetTitle(title);
1758 *dynamic_cast<TAttFill *>(out) = *this;
1759 *dynamic_cast<TAttLine *>(out) = *this;
1760 *dynamic_cast<TAttMarker *>(out) = *this;
1761 out->SetBit(kCanDelete);
1762
1763 if (!gPad)
1766 if (!sOpt.Contains("same"))
1767 basePad->Clear();
1768
1769 bool doFits = false;
1770 if (sOpt.Contains("fits")) {
1771 doFits = true;
1772 sOpt.ReplaceAll("fits", "");
1773 }
1774
1775 auto mainPad = gPad;
1776
1777 out->SetEditable(false);
1778
1779 if (doFits) {
1780 gPad->Divide(1, 2);
1781 gPad->cd(1);
1782 gPad->SetBottomMargin(gPad->GetBottomMargin() * 2.); // increase margin to be same as before
1783 gPad->SetGrid();
1784 out->Draw(sOpt);
1785 basePad->cd(2);
1786 mainPad = basePad->GetPad(1);
1787 } else {
1788 gPad->SetGrid();
1789 out->Draw("ALP");
1790 }
1791
1792 std::pair<double, double> minMax(std::numeric_limits<double>::infinity(), -std::numeric_limits<double>::infinity());
1793 for (auto &p : *this) {
1794 if (p.fPllType != pllType)
1795 continue; // must all have same pll type
1796 auto val = p.pll(true).first;
1797 if (std::isnan(val))
1798 continue;
1799 minMax.first = std::min(minMax.first, val);
1800 minMax.second = std::max(minMax.second, val);
1801 }
1802 if (minMax.first < std::numeric_limits<double>::infinity())
1803 out->GetHistogram()->SetMinimum(minMax.first);
1804 if (minMax.second > -std::numeric_limits<double>::infinity())
1805 out->GetHistogram()->SetMaximum(minMax.second);
1806
1807 TGraph *badPoints = nullptr;
1808
1809 TStopwatch s;
1810 s.Start();
1811 std::shared_ptr<const RooFitResult> ufr;
1812 for (auto &p : *this) {
1813 if (p.fPllType != pllType)
1814 continue; // must all have same pll type
1815 auto val = p.pll().first;
1816 if (!ufr)
1817 ufr = p.ufit();
1818 if (out->GetN() == 0 && ufr && ufr->status() == 0) {
1819 out->SetPoint(out->GetN(),
1820 ufr->floatParsFinal().getRealValue(axes().first()->GetName(),
1821 ufr->constPars().getRealValue(axes().first()->GetName())),
1822 0.);
1823 out->SetPointError(out->GetN() - 1, 0, ufr->edm());
1824 }
1825 if (auto fr = p.fNull_cfit;
1826 fr && doFits) { // access member to avoid unnecessarily creating fit result if wasnt needed
1827 // create a new subpad and draw fitResult on it
1828 auto _pad = gPad;
1829 auto pad =
1830 new TPad(fr->GetName(), TString::Format("%s = %g", poi().first()->GetTitle(), p.fNullVal()), 0, 0, 1., 1);
1831 pad->SetNumber(out->GetN() + 1); // can't use "0" for a subpad
1832 pad->cd();
1833 xRooNode(fr).Draw("goff");
1834 _pad->cd();
1835 //_pad->GetListOfPrimitives()->AddFirst(pad);
1836 pad->AppendPad();
1837 }
1838 if (std::isnan(val) && p.status() != 0) {
1839 if (!badPoints) {
1840 badPoints = new TGraph;
1841 badPoints->SetBit(kCanDelete);
1842 badPoints->SetName("badPoints");
1843 badPoints->SetMarkerStyle(5);
1844 badPoints->SetMarkerColor(kRed);
1845 badPoints->SetMarkerSize(1);
1846 out->GetListOfFunctions()->Add(badPoints, "P");
1847 }
1848 badPoints->SetPoint(badPoints->GetN(), p.fNullVal(), out->Eval(p.fNullVal()));
1849 mainPad->Modified();
1850 } else if (!std::isnan(val)) {
1851 out->SetPoint(out->GetN(), p.coords->getRealValue(axes().first()->GetName()), p.pll().first);
1852 out->SetPointError(out->GetN() - 1, 0, p.pll().second);
1853 out->Sort();
1854
1855 // reposition bad points
1856 if (badPoints) {
1857 for (int i = 0; i < badPoints->GetN(); i++)
1858 badPoints->SetPointY(i, out->Eval(badPoints->GetPointX(i)));
1859 }
1860
1861 mainPad->Modified();
1862 }
1863 if (s.RealTime() > 3) { // stops the clock
1864 basePad->GetCanvas()->Paint();
1865 basePad->GetCanvas()->Update();
1867 s.Reset();
1868 s.Start();
1869 }
1870 s.Continue();
1871 }
1872 basePad->GetCanvas()->Paint();
1873 basePad->GetCanvas()->Update();
1874#if ROOT_VERSION_CODE >= ROOT_VERSION(6, 30, 00)
1875 basePad->GetCanvas()->ResetUpdated(); // stops previous canvas being replaced in a jupyter notebook
1876#endif
1878
1879 // finish by overlaying ufit
1880 if (ufr && doFits) {
1881 auto _pad = gPad;
1882 auto pad = new TPad(ufr->GetName(), "unconditional fit", 0, 0, 1., 1.);
1883 pad->SetNumber(-1);
1884 pad->cd();
1885 xRooNode(ufr).Draw("goff");
1886 _pad->cd();
1887 pad->AppendPad();
1888
1889 // draw one more pad to represent the selected, and draw the ufit pad onto that pad
1890 pad = new TPad("selected", "selected", 0, 0, 1, 1);
1891 pad->Draw();
1892 pad->cd();
1893 basePad->GetPad(2)->GetPad(-1)->AppendPad();
1894 pad->Modified();
1895 pad->Update();
1897
1898 basePad->cd();
1899 }
1900
1901 if (doFits) {
1902 if (!xRooNode::gIntObj) {
1904 }
1905 gPad->GetCanvas()->Connect("Highlighted(TVirtualPad*,TObject*,Int_t,Int_t)", "xRooNode::InteractiveObject",
1906 xRooNode::gIntObj, "Interactive_PLLPlot(TVirtualPad*,TObject*,Int_t,Int_t)");
1907 }
1908
1909 return;
1910}
1911
1913{
1914
1915 RooStats::HypoTestInverterResult *out = nullptr;
1916
1917 auto _axes = axes();
1918 if (_axes.empty())
1919 return out;
1920
1921 out = new RooStats::HypoTestInverterResult(GetName(), *dynamic_cast<RooRealVar *>(_axes.at(0)), 0.95);
1922 out->SetTitle(GetTitle());
1923
1924 for (auto &p : *this) {
1925 double _x = p.coords->getRealValue(_axes.at(0)->GetName(), std::numeric_limits<double>::quiet_NaN());
1926 out->Add(_x, p.result());
1927 }
1928
1929 return out;
1930}
1931
#define f(i)
Definition RSha256.hxx:104
#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
Definition RtypesCore.h:66
@ kRed
Definition Rtypes.h:66
@ kBlack
Definition Rtypes.h:65
@ kGreen
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
@ kYellow
Definition Rtypes.h:66
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define gDirectory
Definition TDirectory.h:384
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:218
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:185
void Warning(const char *location, const char *msgfmt,...)
Use this function in warning situations.
Definition TError.cxx:229
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:110
@ kCanDelete
Definition TObject.h:367
#define gROOT
Definition TROOT.h:414
R__EXTERN TStyle * gStyle
Definition TStyle.h:436
R__EXTERN TSystem * gSystem
Definition TSystem.h:572
#define gPad
static std::pair< double, double > matchPrecision(const std::pair< double, double > &in)
Definition xRooFit.cxx:2000
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:423
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:498
xRooNode pars() const
List of parameters (non-observables) of this node.
const_iterator begin() const
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.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
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:59
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:33
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
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
RooArgSet * GetParameters() const override
return a cloned list with the parameter of interest
Fill Area Attributes class.
Definition TAttFill.h:19
Line Attributes class.
Definition TAttLine.h:18
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:42
Marker Attributes class.
Definition TAttMarker.h:19
static TCanvas * MakeDefCanvas()
Static function to build a default canvas.
Definition TCanvas.cxx:1516
static TClass * GetClass(const char *name, Bool_t load=kTRUE, Bool_t silent=kFALSE)
Static method returning pointer to TClass of the specified class name.
Definition TClass.cxx:3038
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:491
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition TFile.cxx:4094
A TGraphErrors is a TGraph with error bars.
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:1544
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:705
Int_t GetN() const
Definition TGraph.h:132
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:2486
TList * GetListOfFunctions() const
Definition TGraph.h:126
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:2037
virtual Double_t GetPointY(Int_t i) const
Get y value for point i.
Definition TGraph.cxx:1555
virtual void SetPointY(Int_t i, Double_t y)
Set y value for point i.
Definition TGraph.cxx:2374
Book space in a file, create I/O buffers, to fill them, (un)compress them.
Definition TKey.h:28
virtual const char * GetClassName() const
Definition TKey.h:75
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:576
A TMemFile is like a normal TFile except that it reads and writes only from memory.
Definition TMemFile.h:19
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
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition TNamed.cxx:140
virtual void SetNameTitle(const char *name, const char *title)
Set all the TNamed parameters (name and title).
Definition TNamed.cxx:154
Mother of all ROOT objects.
Definition TObject.h:41
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:786
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:869
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:139
Ssiz_t Length() const
Definition TString.h:417
Double_t Atof() const
Return floating-point value contained in string.
Definition TString.cxx:2054
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:2378
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:632
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition TString.h:651
Float_t GetPadRightMargin() const
Definition TStyle.h:214
Float_t GetPadTopMargin() const
Definition TStyle.h:212
virtual Bool_t ProcessEvents()
Process pending events (GUI, timers, sockets).
Definition TSystem.cxx:416
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
Definition graph.py:1
#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