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