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