Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
xRooNLLVar.h
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
13#include "Config.h"
14
15#ifdef XROOFIT_USE_PRAGMA_ONCE
16#pragma once
17#endif
18#if !defined(XROOFIT_XROONLLVAR_H) || defined(XROOFIT_USE_PRAGMA_ONCE)
19#ifndef XROOFIT_USE_PRAGMA_ONCE
20#define XROOFIT_XROONLLVAR_H
21#endif
22
23#include "xRooFit.h"
24
25#include <RooFitResult.h>
26#include <RooLinkedList.h>
27
28#include <Fit/FitConfig.h>
29#include <Math/IOptions.h>
30#include <TAttFill.h>
31#include <TAttLine.h>
32#include <TAttMarker.h>
33
34#include <map>
35#include <set>
36
37class RooAbsReal;
38class RooAbsPdf;
39class RooAbsData;
41class RooNLLVar;
43class RooRealVar;
44class RooCmdArg;
45
46class TGraph;
47class TGraphErrors;
48class TMultiGraph;
49
50namespace RooStats {
51class HypoTestResult;
52class HypoTestInverterResult;
53} // namespace RooStats
54
56
57class xRooNode;
58
59class xRooNLLVar : public std::shared_ptr<RooAbsReal> {
60
61public:
62 struct xValueWithError : public std::pair<double, double> {
63 xValueWithError(const std::pair<double, double> &in = {0, 0}) : std::pair<double, double>(in) {}
64 double value() const { return std::pair<double, double>::first; }
65 double error() const { return std::pair<double, double>::second; }
66 };
67
68 void Print(Option_t *opt = "");
69
70 xRooNLLVar(RooAbsPdf &pdf, const std::pair<RooAbsData *, const RooAbsCollection *> &data,
71 const RooLinkedList &nllOpts = RooLinkedList());
72 xRooNLLVar(const std::shared_ptr<RooAbsPdf> &pdf, const std::shared_ptr<RooAbsData> &data,
73 const RooLinkedList &opts = RooLinkedList());
74 xRooNLLVar(const std::shared_ptr<RooAbsPdf> &pdf,
75 const std::pair<std::shared_ptr<RooAbsData>, std::shared_ptr<const RooAbsCollection>> &data,
76 const RooLinkedList &opts = RooLinkedList());
77
79
80 // whenever implicitly converted to a RooAbsReal we will make sure our globs are set
81 RooAbsReal *get() const { return func().get(); }
82 RooAbsReal *operator->() const { return get(); }
83
84 void reinitialize();
85
86 void AddOption(const RooCmdArg &opt);
87
88 std::pair<std::shared_ptr<RooAbsData>, std::shared_ptr<const RooAbsCollection>>
89 getData() const; // returns pointer to data and snapshot of globs
90 std::pair<std::shared_ptr<RooAbsData>, std::shared_ptr<const RooAbsCollection>>
91 generate(bool expected = false, int seed = 0);
92 // std::shared_ptr<const RooFitResult> snapshot();
93
94 class xRooFitResult : public std::shared_ptr<const RooFitResult> {
95 public:
96 xRooFitResult(const std::shared_ptr<xRooNode> &in,
97 const std::shared_ptr<xRooNLLVar> &nll = nullptr); // : fNode(in) { }
98 const RooFitResult *operator->() const;
99 // operator std::shared_ptr<const RooFitResult>() const;
100 operator const RooFitResult *() const;
101 void Draw(Option_t *opt = "");
102
103 std::shared_ptr<xRooNLLVar> nll() const { return fNll; }
104
106 {
107 return get()
108 ? RooArgList(*std::unique_ptr<RooAbsCollection>(get()->floatParsFinal().selectByAttrib("poi", true)))
109 : RooArgList();
110 }
111
112 // generate a conditional fit using the given poi set to the given values
113 // alias is used to store the fit result in the map under a different name
114 xRooFitResult cfit(const char *poiValues, const char *alias = nullptr);
115 // generate the conditional fit required for an impact calculation
116 xRooFitResult ifit(const char *np, bool up, bool prefit = false);
117 // calculate the impact on poi due to np. if approx is true, will use the covariance approximation instead
118 double impact(const char *poi, const char *np, bool up = true, bool prefit = false, bool approx = false);
119 double impact(const char *np, bool up = true, bool prefit = false, bool approx = false)
120 {
121 auto _poi = poi();
122 if (_poi.size() != 1)
123 throw std::runtime_error("xRooFitResult::impact: not one POI");
124 return impact(poi().contentsString().c_str(), np, up, prefit, approx);
125 }
126
127 // calculate error on poi conditional on the given NPs being held constant at their post-fit values
128 // The conditional error is often presented as the difference in quadrature to the total error i.e.
129 // error contribution due to conditional NPs = sqrt( pow(totError,2) - pow(condError,2) )
130 double conditionalError(const char *poi, const char *nps, bool up = true, bool approx = false);
131
132 // rank all the np based on impact ... will use the covariance approximation if full impact not available
133 // the approxThreshold sets the level below which the approximation will be returned
134 // e.g. set it to 0 to not do approximation
135 RooArgList ranknp(const char *poi, bool up = true, bool prefit = false,
136 double approxThreshold = std::numeric_limits<double>::infinity());
137 // version that assumes only one parameter is poi
139 ranknp(bool up = true, bool prefit = false, double approxThreshold = std::numeric_limits<double>::infinity())
140 {
141 auto _poi = poi();
142 if (_poi.size() != 1)
143 throw std::runtime_error("xRooFitResult::ranknp: not one POI");
144 return ranknp(poi().contentsString().c_str(), up, prefit, approxThreshold);
145 }
146
147 std::shared_ptr<xRooNode> fNode;
148 std::shared_ptr<xRooNLLVar> fNll;
149
150 std::shared_ptr<std::map<std::string, xRooFitResult>> fCfits;
151 };
152
153 xRooFitResult minimize(const std::shared_ptr<ROOT::Fit::FitConfig> & = nullptr);
154
155 void SetFitConfig(const std::shared_ptr<ROOT::Fit::FitConfig> &in) { fFitConfig = in; }
156 std::shared_ptr<ROOT::Fit::FitConfig> fitConfig(); // returns fit config, or creates a default one if not existing
157 ROOT::Math::IOptions *fitConfigOptions(); // return pointer to non-const version of the options inside the fit config
158
159 class xRooHypoPoint : public TNamed {
160 public:
161 xRooHypoPoint(std::shared_ptr<RooStats::HypoTestResult> htr = nullptr, const RooAbsCollection *_coords = nullptr);
162 static std::set<int> allowedStatusCodes;
163 void Print(Option_t *opt = "") const override;
164 void Draw(Option_t *opt = "") override;
165
166 // status bitmask of the available fit results
167 // 0 = all ok
168 int status() const;
169
170 xValueWithError pll(bool readOnly = false); // observed test statistic value
171 xValueWithError sigma_mu(bool readOnly = false); // estimate of sigma_mu parameter
172 std::shared_ptr<const RooFitResult> ufit(bool readOnly = false);
173 std::shared_ptr<const RooFitResult> cfit_null(bool readOnly = false);
174 std::shared_ptr<const RooFitResult> cfit_alt(bool readOnly = false);
175 std::shared_ptr<const RooFitResult> cfit_lbound(bool readOnly = false); // cfit @ the lower bound of mu
176 std::shared_ptr<const RooFitResult> gfit() { return fGenFit; } // non-zero if data was generated
177
178 std::pair<std::shared_ptr<RooAbsData>, std::shared_ptr<const RooAbsCollection>> fData;
179 std::pair<std::shared_ptr<RooAbsData>, std::shared_ptr<const RooAbsCollection>> data();
180
181 xValueWithError getVal(const char *what);
182
183 // leave nSigma=NaN for observed p-value
184 xValueWithError pNull_asymp(double nSigma = std::numeric_limits<double>::quiet_NaN());
185 xValueWithError pAlt_asymp(double nSigma = std::numeric_limits<double>::quiet_NaN());
186 xValueWithError pCLs_asymp(double nSigma = std::numeric_limits<double>::quiet_NaN());
187 xValueWithError ts_asymp(double nSigma = std::numeric_limits<double>::quiet_NaN()); // test statistic value
188
189 xValueWithError pNull_toys(double nSigma = std::numeric_limits<double>::quiet_NaN());
190 xValueWithError pAlt_toys(double nSigma = std::numeric_limits<double>::quiet_NaN());
191 xValueWithError pCLs_toys(double nSigma = std::numeric_limits<double>::quiet_NaN())
192 {
193 if (fNullVal() == fAltVal())
194 return std::pair<double, double>(1, 0); // by construction
195 auto null = pNull_toys(nSigma);
196 auto alt = pAlt_toys(nSigma);
197 double nom = (null.first == 0) ? 0 : null.first / alt.first;
198 // double up = (null.first + null.second == 0) ? 0 : ((alt.first-alt.second<=0) ?
199 // std::numeric_limits<double>::infinity() : (null.first + null.second)/(alt.first - alt.second)); double down
200 // = (null.first - null.second == 0) ? 0 : (null.first - null.second)/(alt.first + alt.second);
201 // old way ... now doing like in pCLs_asymp by calculating the two variations ... but this is pessimistic
202 // assumes p-values are anticorrelated!
203 // so reverting to old
204 return std::pair<double, double>(nom, (alt.first - alt.second <= 0)
205 ? std::numeric_limits<double>::infinity()
206 : (sqrt(pow(null.second, 2) + pow(alt.second * nom, 2)) / alt.first));
207 // return std::pair(nom,std::max(std::abs(up - nom), std::abs(down - nom)));
208 }
209 xValueWithError ts_toys(double nSigma = std::numeric_limits<double>::quiet_NaN()); // test statistic value
210
211 // Create a HypoTestResult representing the current state of this hypoPoint
213
214 xRooHypoPoint generateNull(int seed = 0);
215 xRooHypoPoint generateAlt(int seed = 0);
216
217 void
218 addNullToys(int nToys = 1, int seed = 0, double target = std::numeric_limits<double>::quiet_NaN(),
219 double target_nSigma = std::numeric_limits<double>::quiet_NaN()); // if seed=0 will use a random seed
220 void
221 addAltToys(int nToys = 1, int seed = 0, double target = std::numeric_limits<double>::quiet_NaN(),
222 double target_nSigma = std::numeric_limits<double>::quiet_NaN()); // if seed=0 will use a random seed
223 void
224 addCLsToys(int nToys = 1, int seed = 0, double target = std::numeric_limits<double>::quiet_NaN(),
225 double target_nSigma = std::numeric_limits<double>::quiet_NaN()); // if seed=0 will use a random seed
226
227 RooArgList poi() const;
228 RooArgList alt_poi() const; // values of the poi in the alt hypothesis (will be nans if not defined)
229 RooRealVar &mu_hat(); // throws exception if ufit not available
230
231 std::shared_ptr<xRooHypoPoint>
232 asimov(bool readOnly =
233 false); // a two-sided hypoPoint with the alt hypothesis asimov dataset (used in sigma_mu() calculation)
234
235 // std::string fPOIName;
236 const char *fPOIName();
237 xRooFit::Asymptotics::PLLType fPllType = xRooFit::Asymptotics::Unknown;
238 // double fNullVal=1; double fAltVal=0;
239 double fNullVal();
240 double fAltVal();
241
242 std::shared_ptr<const RooAbsCollection> coords; // pars of the nll that will be held const alongside POI
243
244 std::shared_ptr<const RooFitResult> fUfit, fNull_cfit, fAlt_cfit, fLbound_cfit;
245 std::shared_ptr<const RooFitResult> fGenFit; // if the data was generated, this is the fit is was generated from
246 bool isExpected = false; // if genFit, flag says is asimov or not
247
248 std::shared_ptr<xRooHypoPoint>
249 fAsimov; // same as this point but pllType is twosided and data is expected post alt-fit
250
251 // first is seed, second is ts value, third is weight
252 std::vector<std::tuple<int, double, double>> nullToys; // would have to save these vectors for specific: null_cfit
253 // (genPoint), ufit, poiName, pllType, nullVal
254 std::vector<std::tuple<int, double, double>> altToys;
255
256 std::shared_ptr<xRooNLLVar> nllVar = nullptr; // hypopoints get a copy
257 std::shared_ptr<RooStats::HypoTestResult> hypoTestResult = nullptr;
258 std::shared_ptr<const RooFitResult> retrieveFit(int type);
259
260 TString tsTitle(bool inWords = false) const;
261
262 private:
263 xValueWithError pX_toys(bool alt, double nSigma = std::numeric_limits<double>::quiet_NaN());
264 size_t addToys(bool alt, int nToys, int initialSeed = 0, double target = std::numeric_limits<double>::quiet_NaN(),
265 double target_nSigma = std::numeric_limits<double>::quiet_NaN(), bool targetCLs = false,
266 double relErrThreshold = 2., size_t maxToys = 10000);
267 };
268
269 // use alt_value = nan to skip the asimov calculations
270 xRooHypoPoint hypoPoint(const char *parName, double value,
271 double alt_value = std::numeric_limits<double>::quiet_NaN(),
272 const xRooFit::Asymptotics::PLLType &pllType = xRooFit::Asymptotics::Unknown);
273 // same as above but specify parNames and values in a string
274 xRooHypoPoint hypoPoint(const char *parValues, double alt_value, const xRooFit::Asymptotics::PLLType &pllType);
275 // this next method requires poi to be flagged in the model already (with "poi" attribute) .. must be exactly one
276 xRooHypoPoint hypoPoint(double value, double alt_value = std::numeric_limits<double>::quiet_NaN(),
277 const xRooFit::Asymptotics::PLLType &pllType = xRooFit::Asymptotics::Unknown);
278
279 class xRooHypoSpace : public TNamed,
280 public TAttFill,
281 public TAttMarker,
282 public TAttLine,
283 public std::vector<xRooHypoPoint> {
284 public:
285 friend class xRooNLLVar;
286 xRooHypoSpace(const char *name = "", const char *title = "");
288
289 bool AddModel(const xRooNode &pdf, const char *validity = "");
290
291 void LoadFits(const char *apath);
292
293 // A points over given parameter, number of points between low and high
294 int AddPoints(const char *parName, size_t nPoints, double low, double high);
295
296 void Print(Option_t *opt = "") const override;
297
298 void Draw(Option_t *opt = "") override;
299
300 RooArgList poi();
301 std::shared_ptr<RooArgSet> pars() const { return fPars; };
302 RooArgList axes() const;
303
304 xRooHypoPoint &AddPoint(double value); // adds by using the first axis var
305 xRooHypoPoint &AddPoint(const char *coords = ""); // adds a new point at given coords or returns existing
306
307 xRooHypoPoint &point(size_t i) { return at(i); }
308
309 // build a TGraphErrors of pValues over the existing points
310 // opt should include any of the following:
311 // cls: do pCLs, otherwise do pNull
312 // expX: do expected, X sigma (use +X or -X for contour, otherwise will return band unless X=0)
313 // toys: pvalues from available toys
314 // readonly: don't compute anything, just return available values
315 std::shared_ptr<TGraphErrors> graph(const char *opt) const;
316
317 // return a TMultiGraph containing the set of graphs for a particular visualization
318 std::shared_ptr<TMultiGraph> graphs(const char *opt);
319
320 // will evaluate more points until limit is below given relative uncert
321 xValueWithError findlimit(const char *opt, double relUncert = std::numeric_limits<double>::infinity(),
322 unsigned int maxTries = 20);
323
324 // get currently available limit, with error. Use nSigma = nan for observed limit
325 xValueWithError limit(const char *type = "cls", double nSigma = std::numeric_limits<double>::quiet_NaN()) const;
326 int scan(const char *type, size_t nPoints, double low = std::numeric_limits<double>::quiet_NaN(),
327 double high = std::numeric_limits<double>::quiet_NaN(),
328 const std::vector<double> &nSigmas = {0, 1, 2, -1, -2, std::numeric_limits<double>::quiet_NaN()},
329 double relUncert = 0.1);
330 int scan(const char *type = "cls",
331 const std::vector<double> &nSigmas = {0, 1, 2, -1, -2, std::numeric_limits<double>::quiet_NaN()},
332 double relUncert = 0.1)
333 {
334 return scan(type, 0, std::numeric_limits<double>::quiet_NaN(), std::numeric_limits<double>::quiet_NaN(),
335 nSigmas, relUncert);
336 }
337 int scan(const char *type, double nSigma, double relUncert = 0.1)
338 {
339 return scan(type, std::vector<double>{nSigma}, relUncert);
340 }
341
342 // key is nSigma or "obs" for observed
343 // will only do obs if "obs" dataset is not a generated dataset
344 std::map<std::string, xValueWithError>
345 limits(const char *opt = "cls",
346 const std::vector<double> &nSigmas = {0, 1, 2, -1, -2, std::numeric_limits<double>::quiet_NaN()},
347 double relUncert = std::numeric_limits<double>::infinity());
348
349 std::shared_ptr<xRooNode> pdf(const RooAbsCollection &parValues) const;
350 std::shared_ptr<xRooNode> pdf(const char *parValues = "") const;
351
352 // caller needs to take ownership of the returned object
354
355 private:
356 // estimates where corresponding pValues graph becomes equal to 0.05
357 // linearly interpolates log(pVal) when obtaining limits.
358 // returns value and error
359 static xValueWithError GetLimit(const TGraph &pValues, double target = std::numeric_limits<double>::quiet_NaN());
360
361 static RooArgList toArgs(const char *str);
362
363 xRooFit::Asymptotics::PLLType fTestStatType = xRooFit::Asymptotics::Unknown;
364 std::shared_ptr<RooArgSet> fPars;
365
366 std::map<std::shared_ptr<xRooNode>, std::shared_ptr<xRooNLLVar>> fNlls; // existing NLL functions of added pdfs;
367
368 std::set<std::pair<std::shared_ptr<RooArgList>, std::shared_ptr<xRooNode>>> fPdfs;
369 };
370
371 xRooHypoSpace hypoSpace(const char *parName, int nPoints, double low, double high,
372 double alt_value = std::numeric_limits<double>::quiet_NaN(),
373 const xRooFit::Asymptotics::PLLType &pllType = xRooFit::Asymptotics::Unknown);
374 xRooHypoSpace hypoSpace(const char *parName = "",
375 const xRooFit::Asymptotics::PLLType &pllType = xRooFit::Asymptotics::Unknown,
376 double alt_value = std::numeric_limits<double>::quiet_NaN());
377 xRooHypoSpace hypoSpace(int nPoints, double low, double high,
378 double alt_value = std::numeric_limits<double>::quiet_NaN(),
379 const xRooFit::Asymptotics::PLLType &pllType = xRooFit::Asymptotics::Unknown);
380 xRooHypoSpace hypoSpace(const char *parName, xRooFit::TestStatistic::Type tsType, int nPoints = 0)
381 {
382 return hypoSpace(parName, int(tsType), nPoints, -std::numeric_limits<double>::infinity(),
383 std::numeric_limits<double>::infinity());
384 }
385
386 std::shared_ptr<RooArgSet> pars(bool stripGlobalObs = true) const;
387
388 void Draw(Option_t *opt = "");
389
390 TObject *Scan(const RooArgList &scanPars, const std::vector<std::vector<double>> &coords,
391 const RooArgList &profilePars = RooArgList());
392 TObject *Scan(const char *scanPars, const std::vector<std::vector<double>> &coords,
393 const RooArgList &profilePars = RooArgList());
394 TObject *Scan(const char *scanPars, size_t nPoints, double low, double high, size_t nPointsY, double ylow,
395 double yhigh, const RooArgList &profilePars = RooArgList())
396 {
397 std::vector<std::vector<double>> coords;
398 if (nPoints) {
399 double step = (high - low) / (nPoints);
400 for (size_t i = 0; i < nPoints; i++) {
401 std::vector<double> coord({low + step * i});
402 if (nPointsY) {
403 double stepy = (yhigh - ylow) / (nPointsY);
404 for (size_t j = 0; j < nPointsY; j++) {
405 coord.push_back({ylow + stepy * j});
406 coords.push_back(coord);
407 coord.resize(1);
408 }
409 } else {
410 coords.push_back(coord);
411 }
412 }
413 }
414 return Scan(scanPars, coords, profilePars);
415 }
416 TObject *
417 Scan(const char *scanPars, size_t nPoints, double low, double high, const RooArgList &profilePars = RooArgList())
418 {
419 return Scan(scanPars, nPoints, low, high, 0, 0, 0, profilePars);
420 }
421
422 std::shared_ptr<RooAbsReal> func() const; // will assign globs when called
423 std::shared_ptr<RooAbsPdf> pdf() const { return fPdf; }
424 RooAbsData *data() const; // returns the data hidden inside the NLLVar if there is some
425 const RooAbsCollection *globs() const { return fGlobs.get(); }
426
427 // NLL = mainTerm + constraintTerm
428 // mainTerm = sum( entryVals ) + extendedTerm + simTerm [+ binnedDataTerm if activated binnedL option]
429 // this is what it should be, at least
430
431 // total nll should be all these values + constraint term + extended term + simTerm [+binnedDataTerm if activated
432 // binnedL option]
433 RooNLLVar *mainTerm() const;
434 RooConstraintSum *constraintTerm() const;
435
436 double getEntryVal(size_t entry) const; // get the Nll value for a specific entry
437 double extendedTerm() const;
438 double simTerm() const;
439 double binnedDataTerm() const;
440 double getEntryBinWidth(size_t entry) const;
441
442 double ndof() const;
443 double saturatedVal() const;
444 double saturatedConstraintTerm() const;
445 double saturatedMainTerm() const;
446 double pgof() const; // a goodness-of-fit pvalue based on profile likelihood of a saturated model
447 double mainTermPgof() const;
448 double mainTermNdof() const;
449
450 std::set<std::string> binnedChannels() const;
451
452 // change the dataset - will check globs are the same
453 bool setData(const std::pair<std::shared_ptr<RooAbsData>, std::shared_ptr<const RooAbsCollection>> &_data);
454 bool setData(const std::shared_ptr<RooAbsData> &data, const std::shared_ptr<const RooAbsCollection> &globs)
455 {
456 return setData(std::make_pair(data, globs));
457 }
458 bool setData(const xRooNode &data);
459
460 // using shared ptrs everywhere, even for RooLinkedList which needs custom deleter to clear itself
461 // but still work ok for assignment operations
462 std::shared_ptr<RooAbsPdf> fPdf;
463 std::shared_ptr<RooAbsData> fData;
464 std::shared_ptr<const RooAbsCollection> fGlobs;
465
466 std::shared_ptr<RooLinkedList> fOpts;
467 std::shared_ptr<ROOT::Fit::FitConfig> fFitConfig;
468
469 std::shared_ptr<RooAbsCollection> fFuncVars;
470 std::shared_ptr<RooAbsCollection> fConstVars;
471 std::shared_ptr<RooAbsCollection> fFuncGlobs;
472 std::string fFuncCreationLog; // messaging from when function was last created -- to save from printing to screen
473
474 bool kReuseNLL = true;
475};
476
477namespace cling {
478std::string printValue(const xRooNLLVar::xValueWithError *val);
479std::string printValue(const std::map<std::string, xRooNLLVar::xValueWithError> *m);
480} // namespace cling
481
483
484#endif // include guard
const char Option_t
Definition RtypesCore.h:66
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
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 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 np
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 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
char name[80]
Definition TGX11.cxx:110
void Print(GNN_Data &d, std::string txt="")
double impact(const char *poi, const char *np, bool up=true, bool prefit=false, bool approx=false)
RooArgList ranknp(bool up=true, bool prefit=false, double approxThreshold=std::numeric_limits< double >::infinity())
Definition xRooNLLVar.h:139
std::shared_ptr< xRooNLLVar > nll() const
Definition xRooNLLVar.h:103
double impact(const char *np, bool up=true, bool prefit=false, bool approx=false)
Definition xRooNLLVar.h:119
std::shared_ptr< std::map< std::string, xRooFitResult > > fCfits
Definition xRooNLLVar.h:150
std::pair< std::shared_ptr< RooAbsData >, std::shared_ptr< const RooAbsCollection > > fData
Definition xRooNLLVar.h:178
std::vector< std::tuple< int, double, double > > altToys
Definition xRooNLLVar.h:254
std::shared_ptr< const RooAbsCollection > coords
Definition xRooNLLVar.h:242
std::vector< std::tuple< int, double, double > > nullToys
Definition xRooNLLVar.h:252
std::shared_ptr< const RooFitResult > fAlt_cfit
Definition xRooNLLVar.h:244
std::shared_ptr< const RooFitResult > fGenFit
Definition xRooNLLVar.h:245
xValueWithError pCLs_toys(double nSigma=std::numeric_limits< double >::quiet_NaN())
Definition xRooNLLVar.h:191
std::shared_ptr< const RooFitResult > gfit()
Definition xRooNLLVar.h:176
std::shared_ptr< RooArgSet > pars() const
Definition xRooNLLVar.h:301
int scan(const char *type="cls", const std::vector< double > &nSigmas={0, 1, 2, -1, -2, std::numeric_limits< double >::quiet_NaN()}, double relUncert=0.1)
Definition xRooNLLVar.h:330
std::map< std::shared_ptr< xRooNode >, std::shared_ptr< xRooNLLVar > > fNlls
Definition xRooNLLVar.h:366
int scan(const char *type, double nSigma, double relUncert=0.1)
Definition xRooNLLVar.h:337
std::set< std::pair< std::shared_ptr< RooArgList >, std::shared_ptr< xRooNode > > > fPdfs
Definition xRooNLLVar.h:368
This xRooNLLVar object has several special methods, e.g.
Definition xRooNLLVar.h:59
std::shared_ptr< RooAbsCollection > fFuncGlobs
Definition xRooNLLVar.h:471
std::shared_ptr< const RooAbsCollection > fGlobs
Definition xRooNLLVar.h:464
bool setData(const std::shared_ptr< RooAbsData > &data, const std::shared_ptr< const RooAbsCollection > &globs)
Definition xRooNLLVar.h:454
xRooHypoPoint hypoPoint(const char *parValues, double alt_value, const xRooFit::Asymptotics::PLLType &pllType)
std::shared_ptr< RooLinkedList > fOpts
Definition xRooNLLVar.h:466
std::shared_ptr< ROOT::Fit::FitConfig > fFitConfig
Definition xRooNLLVar.h:467
void SetFitConfig(const std::shared_ptr< ROOT::Fit::FitConfig > &in)
Definition xRooNLLVar.h:155
xRooHypoSpace hypoSpace(const char *parName, xRooFit::TestStatistic::Type tsType, int nPoints=0)
Definition xRooNLLVar.h:380
std::shared_ptr< RooAbsCollection > fConstVars
Definition xRooNLLVar.h:470
std::shared_ptr< RooAbsPdf > pdf() const
Definition xRooNLLVar.h:423
TObject * Scan(const char *scanPars, size_t nPoints, double low, double high, size_t nPointsY, double ylow, double yhigh, const RooArgList &profilePars=RooArgList())
Definition xRooNLLVar.h:394
const RooAbsCollection * globs() const
Definition xRooNLLVar.h:425
std::shared_ptr< RooAbsCollection > fFuncVars
Definition xRooNLLVar.h:469
std::shared_ptr< RooAbsData > fData
Definition xRooNLLVar.h:463
std::shared_ptr< RooAbsPdf > fPdf
Definition xRooNLLVar.h:462
TObject * Scan(const char *scanPars, size_t nPoints, double low, double high, const RooArgList &profilePars=RooArgList())
Definition xRooNLLVar.h:417
The xRooNode class is designed to wrap over a TObject and provide functionality to aid with interacti...
Definition xRooNode.h:52
Generic interface for defining configuration options of a numerical algorithm.
Definition IOptions.h:28
Abstract container object that can hold multiple RooAbsArg objects.
Storage_t const & get() const
Const access to the underlying stl container.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
Named container for two doubles, two integers two object points and three string pointers that can be...
Definition RooCmdArg.h:26
Calculates the sum of the -(log) likelihoods of a set of RooAbsPfs that represent constraint function...
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
Collection class for internal use, storing a collection of RooAbsArg pointers in a doubly linked list...
Implements a -log(likelihood) calculation from a dataset and a PDF.
Definition RooNLLVar.h:20
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 is a base class for results from hypothesis tests.
Fill Area Attributes class.
Definition TAttFill.h:19
Line Attributes class.
Definition TAttLine.h:18
Marker Attributes class.
Definition TAttMarker.h:19
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
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
Mother of all ROOT objects.
Definition TObject.h:41
Basic string class.
Definition TString.h:139
Namespace for the RooStats classes.
Definition Asimov.h:19
Definition graph.py:1
#define BEGIN_XROOFIT_NAMESPACE
Definition Config.h:24
#define END_XROOFIT_NAMESPACE
Definition Config.h:25
static const char * what
Definition stlLoader.cc:5
xValueWithError(const std::pair< double, double > &in={0, 0})
Definition xRooNLLVar.h:63
th1 Draw()
TMarker m
Definition textangle.C:8