Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
FitHelpers.cxx
Go to the documentation of this file.
1/// \cond ROOFIT_INTERNAL
2
3/*
4 * Project: RooFit
5 * Authors:
6 * Jonas Rembser, CERN 2023
7 *
8 * Copyright (c) 2023, CERN
9 *
10 * Redistribution and use in source and binary forms,
11 * with or without modification, are permitted according to the terms
12 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
13 */
14
15#include "FitHelpers.h"
16
17#include <RooAbsData.h>
18#include <RooAbsPdf.h>
19#include <RooAbsReal.h>
20#include <RooAddition.h>
21#include <RooBatchCompute.h>
22#include <RooBinSamplingPdf.h>
23#include <RooCategory.h>
24#include <RooCmdConfig.h>
25#include <RooConstraintSum.h>
26#include <RooDataHist.h>
27#include <RooDataSet.h>
28#include <RooDerivative.h>
29#include <RooFit/Evaluator.h>
32#include <RooFitResult.h>
33#include <RooLinkedList.h>
34#include <RooMinimizer.h>
35#include <RooConstVar.h>
36#include <RooRealVar.h>
37#include <RooSimultaneous.h>
38#include <RooFormulaVar.h>
39
40#include <Math/CholeskyDecomp.h>
41#include <Math/Util.h>
42
43#include "ConstraintHelpers.h"
44#include "RooEvaluatorWrapper.h"
45#include "RooFitImplHelpers.h"
47
48#ifdef ROOFIT_LEGACY_EVAL_BACKEND
49#include "RooChi2Var.h"
50#include "RooNLLVar.h"
51#endif
52
53using RooFit::Detail::RooNLLVarNew;
54
55namespace {
56
57constexpr int extendedFitDefault = 2;
58
59////////////////////////////////////////////////////////////////////////////////
60/// Use the asymptotically correct approach to estimate errors in the presence of weights.
61/// This is slower but more accurate than `SumW2Error`. See also https://arxiv.org/abs/1911.01303).
62/// Applies the calculated covaraince matrix to the RooMinimizer and returns
63/// the quality of the covariance matrix.
64/// See also the documentation of RooAbsPdf::fitTo(), where this function is used.
65/// \param[in] minimizer The RooMinimizer to get the fit result from. The state
66/// of the minimizer will be altered by this function: the covariance
67/// matrix caltulated here will be applied to it via
68/// RooMinimizer::applyCovarianceMatrix().
69/// \param[in] data The dataset that was used for the fit.
71{
72 RooFormulaVar logpdf("logpdf", "log(pdf)", "log(@0)", pdf);
73 RooArgSet obs;
74 logpdf.getObservables(data.get(), obs);
75
76 // Warning if the dataset is binned. TODO: in some cases,
77 // people also use RooDataSet to encode binned data,
78 // e.g. for simultaneous fits. It would be useful to detect
79 // this in this future as well.
80 if (dynamic_cast<RooDataHist const *>(&data)) {
81 oocoutW(&pdf, InputArguments)
82 << "RooAbsPdf::fitTo(" << pdf.GetName()
83 << ") WARNING: Asymptotic error correction is requested for a binned data set. "
84 "This method is not designed to handle binned data. A standard chi2 fit will likely be more suitable.";
85 };
86
87 // Calculated corrected errors for weighted likelihood fits
88 std::unique_ptr<RooFitResult> rw(minimizer.save());
89 // Weighted inverse Hessian matrix
90 const TMatrixDSym &matV = rw->covarianceMatrix();
91 oocoutI(&pdf, Fitting)
92 << "RooAbsPdf::fitTo(" << pdf.GetName()
93 << ") Calculating covariance matrix according to the asymptotically correct approach. If you find this "
94 "method useful please consider citing https://arxiv.org/abs/1911.01303.\n";
95
96 // Initialise matrix containing first derivatives
97 int nFloatPars = rw->floatParsFinal().size();
99 for (int k = 0; k < nFloatPars; k++) {
100 for (int l = 0; l < nFloatPars; l++) {
101 num(k, l) = 0.0;
102 }
103 }
104
105 // Create derivative objects
106 std::vector<std::unique_ptr<RooDerivative>> derivatives;
107 const RooArgList &floated = rw->floatParsFinal();
109 logpdf.getParameters(data.get(), allparams);
110 std::unique_ptr<RooArgSet> floatingparams{allparams.selectByAttrib("Constant", false)};
111
112 const double eps = 1.0e-4;
113
114 // Calculate derivatives of logpdf
115 for (const auto paramresult : floated) {
116 auto paraminternal = static_cast<RooRealVar *>(floatingparams->find(*paramresult));
117 assert(floatingparams->find(*paramresult)->IsA() == RooRealVar::Class());
118 double error = static_cast<RooRealVar *>(paramresult)->getError();
119 derivatives.emplace_back(logpdf.derivative(*paraminternal, obs, 1, eps * error));
120 }
121
122 // Calculate derivatives for number of expected events, needed for extended ML fit
123 RooAbsPdf *extended_pdf = dynamic_cast<RooAbsPdf *>(&pdf);
124 std::vector<double> diffs_expected(floated.size(), 0.0);
125 if (extended_pdf && extended_pdf->expectedEvents(obs) != 0.0) {
126 for (std::size_t k = 0; k < floated.size(); k++) {
127 const auto paramresult = static_cast<RooRealVar *>(floated.at(k));
128 auto paraminternal = static_cast<RooRealVar *>(floatingparams->find(*paramresult));
129
130 *paraminternal = paramresult->getVal();
131 double error = paramresult->getError();
132 paraminternal->setVal(paramresult->getVal() + eps * error);
133 double expected_plus = log(extended_pdf->expectedEvents(obs));
134 paraminternal->setVal(paramresult->getVal() - eps * error);
135 double expected_minus = log(extended_pdf->expectedEvents(obs));
136 *paraminternal = paramresult->getVal();
137 double diff = (expected_plus - expected_minus) / (2.0 * eps * error);
138 diffs_expected[k] = diff;
139 }
140 }
141
142 // Loop over data
143 for (int j = 0; j < data.numEntries(); j++) {
144 // Sets obs to current data point, this is where the pdf will be evaluated
145 obs.assign(*data.get(j));
146 // Determine first derivatives
147 std::vector<double> diffs(floated.size(), 0.0);
148 for (std::size_t k = 0; k < floated.size(); k++) {
149 const auto paramresult = static_cast<RooRealVar *>(floated.at(k));
150 auto paraminternal = static_cast<RooRealVar *>(floatingparams->find(*paramresult));
151 // first derivative to parameter k at best estimate point for this measurement
152 double diff = derivatives[k]->getVal();
153 // need to reset to best fit point after differentiation
154 *paraminternal = paramresult->getVal();
155 diffs[k] = diff;
156 }
157
158 // Fill numerator matrix
159 for (std::size_t k = 0; k < floated.size(); k++) {
160 for (std::size_t l = 0; l < floated.size(); l++) {
161 num(k, l) += data.weightSquared() * (diffs[k] + diffs_expected[k]) * (diffs[l] + diffs_expected[l]);
162 }
163 }
164 }
165 num.Similarity(matV);
166
167 // Propagate corrected errors to parameters objects
168 minimizer.applyCovarianceMatrix(num);
169
170 // The derivatives are found in RooFit and not with the minimizer (e.g.
171 // minuit), so the quality of the corrected covariance matrix corresponds to
172 // the quality of the original covariance matrix
173 return rw->covQual();
174}
175
176////////////////////////////////////////////////////////////////////////////////
177/// Apply correction to errors and covariance matrix. This uses two covariance
178/// matrices, one with the weights, the other with squared weights, to obtain
179/// the correct errors for weighted likelihood fits.
180/// Applies the calculated covaraince matrix to the RooMinimizer and returns
181/// the quality of the covariance matrix.
182/// See also the documentation of RooAbsPdf::fitTo(), where this function is used.
183/// \param[in] minimizer The RooMinimizer to get the fit result from. The state
184/// of the minimizer will be altered by this function: the covariance
185/// matrix caltulated here will be applied to it via
186/// RooMinimizer::applyCovarianceMatrix().
187/// \param[in] nll The NLL object that was used for the fit.
188int calcSumW2CorrectedCovariance(RooAbsReal const &pdf, RooMinimizer &minimizer, RooAbsReal &nll)
189{
190 // Calculated corrected errors for weighted likelihood fits
191 std::unique_ptr<RooFitResult> rw{minimizer.save()};
192 nll.applyWeightSquared(true);
193 oocoutI(&pdf, Fitting) << "RooAbsPdf::fitTo(" << pdf.GetName()
194 << ") Calculating sum-of-weights-squared correction matrix for covariance matrix\n";
195 minimizer.hesse();
196 std::unique_ptr<RooFitResult> rw2{minimizer.save()};
197 nll.applyWeightSquared(false);
198
199 // Apply correction matrix
200 const TMatrixDSym &matV = rw->covarianceMatrix();
201 TMatrixDSym matC = rw2->covarianceMatrix();
203 if (!decomp) {
204 oocoutE(&pdf, Fitting) << "RooAbsPdf::fitTo(" << pdf.GetName()
205 << ") ERROR: Cannot apply sum-of-weights correction to covariance matrix: correction "
206 "matrix calculated with weight-squared is singular\n";
207 return -1;
208 }
209
210 // replace C by its inverse
211 decomp.Invert(matC);
212 // the class lies about the matrix being symmetric, so fill in the
213 // part above the diagonal
214 for (int i = 0; i < matC.GetNrows(); ++i) {
215 for (int j = 0; j < i; ++j) {
216 matC(j, i) = matC(i, j);
217 }
218 }
219 matC.Similarity(matV);
220 // C now contains V C^-1 V
221 // Propagate corrected errors to parameters objects
222 minimizer.applyCovarianceMatrix(matC);
223
224 return std::min(rw->covQual(), rw2->covQual());
225}
226
227/// Configuration struct for RooAbsPdf::minimizeNLL with all the default values
228/// that also should be taken as the default values for RooAbsPdf::fitTo.
229struct MinimizerConfig {
230 double recoverFromNaN = 10.;
231 int optConst = 0;
232 int verbose = 0;
233 int doSave = 0;
234 int doTimer = 0;
235 int printLevel = 1;
236 int strategy = 1;
237 int initHesse = 0;
238 int hesse = 1;
239 int minos = 0;
240 int numee = 10;
241 int doEEWall = 1;
242 int doWarn = 1;
243 int doSumW2 = -1;
244 int doAsymptotic = -1;
245 int maxCalls = -1;
246 int doOffset = -1;
247 int parallelize = 0;
248 bool enableParallelGradient = false;
249 bool enableParallelDescent = false;
250 bool timingAnalysis = false;
251 const RooArgSet *minosSet = nullptr;
252 std::string minType;
253 std::string minAlg = "minuit";
254};
255
257{
258 // Process automatic extended option
261 if (ext) {
262 oocoutI(&pdf, Minimization)
263 << "p.d.f. provides expected number of events, including extended term in likelihood." << std::endl;
264 }
265 return ext;
266 }
267 // If Extended(false) was explicitly set, but the pdf MUST be extended, then
268 // it's time to print an error. This happens when you're fitting a RooAddPdf
269 // with coefficient that represent yields, and without the additional
270 // constraint these coefficients are degenerate because the RooAddPdf
271 // normalizes itself. Nothing correct can come out of this.
272 if (extendedCmdArg == 0) {
274 std::string errMsg = "You used the Extended(false) option on a pdf where the fit MUST be extended! "
275 "The parameters are not well defined and you're getting nonsensical results.";
276 oocoutE(&pdf, InputArguments) << errMsg << std::endl;
277 }
278 }
279 return extendedCmdArg;
280}
281
282/// To set the fitrange attribute of the PDF and custom ranges for the
283/// observables so that RooPlot can automatically plot the fitting range.
284void resetFitrangeAttributes(RooAbsArg &pdf, RooAbsData const &data, std::string const &baseName, const char *rangeName,
285 bool splitRange)
286{
287 // Clear possible range attributes from previous fits.
288 pdf.removeStringAttribute("fitrange");
289
290 // No fitrange was specified, so we do nothing. Or "SplitRange" is used, and
291 // then there are no uniquely defined ranges for the observables (as they
292 // are different in each category).
293 if (!rangeName || splitRange)
294 return;
295
296 RooArgSet observables;
297 pdf.getObservables(data.get(), observables);
298
299 std::string fitrangeValue;
300 auto subranges = ROOT::Split(rangeName, ",");
301 for (auto const &subrange : subranges) {
302 if (subrange.empty())
303 continue;
304 std::string fitrangeValueSubrange = std::string("fit_") + baseName;
305 if (subranges.size() > 1) {
307 }
309 for (RooAbsArg *arg : observables) {
310
311 if (arg->isCategory())
312 continue;
313 auto &observable = static_cast<RooRealVar &>(*arg);
314
315 observable.setRange(fitrangeValueSubrange.c_str(), observable.getMin(subrange.c_str()),
316 observable.getMax(subrange.c_str()));
317 }
318 }
319 pdf.setStringAttribute("fitrange", fitrangeValue.substr(0, fitrangeValue.size() - 1).c_str());
320}
321
322std::unique_ptr<RooAbsArg> createSimultaneousNLL(RooSimultaneous const &simPdf, bool isSimPdfExtended,
323 std::string const &rangeName, RooFit::OffsetMode offset)
324{
325 RooAbsCategoryLValue const &simCat = simPdf.indexCat();
326
327 // Prepare the NLL terms for each component
329 for (auto const &catState : simCat) {
330 std::string const &catName = catState.first;
332
333 // If the channel is not in the selected range of the category variable, we
334 // won't create an for NLL this channel.
335 if (!rangeName.empty()) {
336 // Only the RooCategory supports ranges, not the other
337 // RooAbsCategoryLValue-derived classes.
338 auto simCatAsRooCategory = dynamic_cast<RooCategory const *>(&simCat);
339 if (simCatAsRooCategory && !simCatAsRooCategory->isStateInRange(rangeName.c_str(), catIndex)) {
340 continue;
341 }
342 }
343
344 if (RooAbsPdf *pdf = simPdf.getPdf(catName.c_str())) {
345 auto name = std::string("nll_") + pdf->GetName();
346 std::unique_ptr<RooArgSet> observables{
347 std::unique_ptr<RooArgSet>(pdf->getVariables())->selectByAttrib("__obs__", true)};
348 // In a simultaneous fit, it is allowed that only a subset of the pdfs
349 // are extended. Therefore, we have to make sure that we don't request
350 // extended NLL objects for channels that can't be extended.
351 const bool isPdfExtended = isSimPdfExtended && pdf->extendMode() != RooAbsPdf::CanNotBeExtended;
352 auto nll =
353 std::make_unique<RooNLLVarNew>(name.c_str(), name.c_str(), *pdf, *observables, isPdfExtended, offset);
354 // Rename the special variables
355 nll->setPrefix(std::string("_") + catName + "_");
356 nllTerms.addOwned(std::move(nll));
357 }
358 }
359
360 for (auto *nll : static_range_cast<RooNLLVarNew *>(nllTerms)) {
361 nll->setSimCount(nllTerms.size());
362 }
363
364 // Time to sum the NLLs
365 auto nll = std::make_unique<RooAddition>("mynll", "mynll", nllTerms);
366 nll->addOwnedComponents(std::move(nllTerms));
367 return nll;
368}
369
370std::unique_ptr<RooAbsReal> createNLLNew(RooAbsPdf &pdf, RooAbsData &data, std::unique_ptr<RooAbsReal> &&constraints,
371 std::string const &rangeName, RooArgSet const &projDeps, bool isExtended,
373{
374 if (constraints) {
375 // The computation graph for the constraints is very small, no need to do
376 // the tracking of clean and dirty nodes here.
377 constraints->setOperMode(RooAbsArg::ADirty);
378 }
379
380 RooArgSet observables;
381 pdf.getObservables(data.get(), observables);
382 observables.remove(projDeps, true, true);
383
384 oocxcoutI(&pdf, Fitting) << "RooAbsPdf::fitTo(" << pdf.GetName()
385 << ") fixing normalization set for coefficient determination to observables in data"
386 << "\n";
387 pdf.fixAddCoefNormalization(observables, false);
388
389 // Deal with the IntegrateBins argument
391 std::unique_ptr<RooAbsPdf> wrappedPdf = RooBinSamplingPdf::create(pdf, data, integrateOverBinsPrecision);
393 if (wrappedPdf) {
394 binSamplingPdfs.addOwned(std::move(wrappedPdf));
395 }
396 // Done dealing with the IntegrateBins option
397
399
400 auto simPdf = dynamic_cast<RooSimultaneous *>(&finalPdf);
401 if (simPdf) {
402 simPdf->wrapPdfsInBinSamplingPdfs(data, integrateOverBinsPrecision);
403 nllTerms.addOwned(createSimultaneousNLL(*simPdf, isExtended, rangeName, offset));
404 } else {
405 nllTerms.addOwned(
406 std::make_unique<RooNLLVarNew>("RooNLLVarNew", "RooNLLVarNew", finalPdf, observables, isExtended, offset));
407 }
408 if (constraints) {
409 nllTerms.addOwned(std::move(constraints));
410 }
411
412 std::string nllName = std::string("nll_") + pdf.GetName() + "_" + data.GetName();
413 auto nll = std::make_unique<RooAddition>(nllName.c_str(), nllName.c_str(), nllTerms);
414 nll->addOwnedComponents(std::move(binSamplingPdfs));
415 nll->addOwnedComponents(std::move(nllTerms));
416
417 return nll;
418}
419
420} // namespace
421
422namespace RooFit::FitHelpers {
423
425{
426 // Default-initialized instance of MinimizerConfig to get the default
427 // minimizer parameter values.
429
430 pc.defineDouble("RecoverFromUndefinedRegions", "RecoverFromUndefinedRegions", 0, minimizerDefaults.recoverFromNaN);
431 pc.defineInt("optConst", "Optimize", 0, minimizerDefaults.optConst);
432 pc.defineInt("verbose", "Verbose", 0, minimizerDefaults.verbose);
433 pc.defineInt("doSave", "Save", 0, minimizerDefaults.doSave);
434 pc.defineInt("doTimer", "Timer", 0, minimizerDefaults.doTimer);
435 pc.defineInt("printLevel", "PrintLevel", 0, minimizerDefaults.printLevel);
436 pc.defineInt("strategy", "Strategy", 0, minimizerDefaults.strategy);
437 pc.defineInt("initHesse", "InitialHesse", 0, minimizerDefaults.initHesse);
438 pc.defineInt("hesse", "Hesse", 0, minimizerDefaults.hesse);
439 pc.defineInt("minos", "Minos", 0, minimizerDefaults.minos);
440 pc.defineInt("numee", "PrintEvalErrors", 0, minimizerDefaults.numee);
441 pc.defineInt("doEEWall", "EvalErrorWall", 0, minimizerDefaults.doEEWall);
442 pc.defineInt("doWarn", "Warnings", 0, minimizerDefaults.doWarn);
443 pc.defineInt("doSumW2", "SumW2Error", 0, minimizerDefaults.doSumW2);
444 pc.defineInt("doAsymptoticError", "AsymptoticError", 0, minimizerDefaults.doAsymptotic);
445 pc.defineInt("maxCalls", "MaxCalls", 0, minimizerDefaults.maxCalls);
446 pc.defineInt("doOffset", "OffsetLikelihood", 0, minimizerDefaults.doOffset);
447 pc.defineInt("parallelize", "Parallelize", 0, minimizerDefaults.parallelize); // Three parallelize arguments
448 pc.defineInt("enableParallelGradient", "ParallelGradientOptions", 0, minimizerDefaults.enableParallelGradient);
449 pc.defineInt("enableParallelDescent", "ParallelDescentOptions", 0, minimizerDefaults.enableParallelDescent);
450 pc.defineInt("timingAnalysis", "TimingAnalysis", 0, minimizerDefaults.timingAnalysis);
451 pc.defineString("mintype", "Minimizer", 0, minimizerDefaults.minType.c_str());
452 pc.defineString("minalg", "Minimizer", 1, minimizerDefaults.minAlg.c_str());
453 pc.defineSet("minosSet", "Minos", 0, minimizerDefaults.minosSet);
454}
455
456////////////////////////////////////////////////////////////////////////////////
457/// Minimizes a given NLL variable by finding the optimal parameters with the
458/// RooMinimzer. The NLL variable can be created with RooAbsPdf::createNLL.
459/// If you are looking for a function that combines likelihood creation with
460/// fitting, see RooAbsPdf::fitTo.
461/// \param[in] nll The negative log-likelihood variable to minimize.
462/// \param[in] data The dataset that was also used for the NLL. It's a necessary
463/// parameter because it is used in the asymptotic error correction.
464/// \param[in] cfg Configuration struct with all the configuration options for
465/// the RooMinimizer. These are a subset of the options that you can
466/// also pass to RooAbsPdf::fitTo via the RooFit command arguments.
467std::unique_ptr<RooFitResult> minimize(RooAbsReal &pdf, RooAbsReal &nll, RooAbsData const &data, RooCmdConfig const &pc)
468{
469 MinimizerConfig cfg;
470 cfg.recoverFromNaN = pc.getDouble("RecoverFromUndefinedRegions");
471 cfg.optConst = pc.getInt("optConst");
472 cfg.verbose = pc.getInt("verbose");
473 cfg.doSave = pc.getInt("doSave");
474 cfg.doTimer = pc.getInt("doTimer");
475 cfg.printLevel = pc.getInt("printLevel");
476 cfg.strategy = pc.getInt("strategy");
477 cfg.initHesse = pc.getInt("initHesse");
478 cfg.hesse = pc.getInt("hesse");
479 cfg.minos = pc.getInt("minos");
480 cfg.numee = pc.getInt("numee");
481 cfg.doEEWall = pc.getInt("doEEWall");
482 cfg.doWarn = pc.getInt("doWarn");
483 cfg.doSumW2 = pc.getInt("doSumW2");
484 cfg.doAsymptotic = pc.getInt("doAsymptoticError");
485 cfg.maxCalls = pc.getInt("maxCalls");
486 cfg.minosSet = pc.getSet("minosSet");
487 cfg.minType = pc.getString("mintype", "");
488 cfg.minAlg = pc.getString("minalg", "minuit");
489 cfg.doOffset = pc.getInt("doOffset");
490 cfg.parallelize = pc.getInt("parallelize");
491 cfg.enableParallelGradient = pc.getInt("enableParallelGradient");
492 cfg.enableParallelDescent = pc.getInt("enableParallelDescent");
493 cfg.timingAnalysis = pc.getInt("timingAnalysis");
494
495 // Determine if the dataset has weights
496 bool weightedData = data.isNonPoissonWeighted();
497
498 std::string msgPrefix = std::string{"RooAbsPdf::fitTo("} + pdf.GetName() + "): ";
499
500 // Warn user that a method to determine parameter uncertainties should be provided if weighted data is offered
501 if (weightedData && cfg.doSumW2 == -1 && cfg.doAsymptotic == -1) {
502 oocoutW(&pdf, InputArguments) << msgPrefix <<
503 R"(WARNING: a likelihood fit is requested of what appears to be weighted data.
504 While the estimated values of the parameters will always be calculated taking the weights into account,
505 there are multiple ways to estimate the errors of the parameters. You are advised to make an
506 explicit choice for the error calculation:
507 - Either provide SumW2Error(true), to calculate a sum-of-weights-corrected HESSE error matrix
508 (error will be proportional to the number of events in MC).
509 - Or provide SumW2Error(false), to return errors from original HESSE error matrix
510 (which will be proportional to the sum of the weights, i.e., a dataset with <sum of weights> events).
511 - Or provide AsymptoticError(true), to use the asymptotically correct expression
512 (for details see https://arxiv.org/abs/1911.01303)."
513)";
514 }
515
516 if (cfg.minos && (cfg.doSumW2 == 1 || cfg.doAsymptotic == 1)) {
517 oocoutE(&pdf, InputArguments)
518 << msgPrefix
519 << " sum-of-weights and asymptotic error correction do not work with MINOS errors. Not fitting.\n";
520 return nullptr;
521 }
522 if (cfg.doAsymptotic == 1 && cfg.minos) {
523 oocoutW(&pdf, InputArguments) << msgPrefix << "WARNING: asymptotic correction does not apply to MINOS errors\n";
524 }
525
526 // avoid setting both SumW2 and Asymptotic for uncertainty correction
527 if (cfg.doSumW2 == 1 && cfg.doAsymptotic == 1) {
528 oocoutE(&pdf, InputArguments) << msgPrefix
529 << "ERROR: Cannot compute both asymptotically correct and SumW2 errors.\n";
530 return nullptr;
531 }
532
533 // Instantiate RooMinimizer
535 minimizerConfig.enableParallelGradient = cfg.enableParallelGradient;
536 minimizerConfig.enableParallelDescent = cfg.enableParallelDescent;
537 minimizerConfig.parallelize = cfg.parallelize;
538 minimizerConfig.timingAnalysis = cfg.timingAnalysis;
539 minimizerConfig.offsetting = cfg.doOffset;
541
542 m.setMinimizerType(cfg.minType);
543 m.setEvalErrorWall(cfg.doEEWall);
544 m.setRecoverFromNaNStrength(cfg.recoverFromNaN);
545 m.setPrintEvalErrors(cfg.numee);
546 if (cfg.maxCalls > 0)
547 m.setMaxFunctionCalls(cfg.maxCalls);
548 if (cfg.printLevel != 1)
549 m.setPrintLevel(cfg.printLevel);
550 if (cfg.optConst)
551 m.optimizeConst(cfg.optConst); // Activate constant term optimization
552 if (cfg.verbose)
553 m.setVerbose(true); // Activate verbose options
554 if (cfg.doTimer)
555 m.setProfile(true); // Activate timer options
556 if (cfg.strategy != 1)
557 m.setStrategy(cfg.strategy); // Modify fit strategy
558 if (cfg.initHesse)
559 m.hesse(); // Initialize errors with hesse
560 m.minimize(cfg.minType.c_str(), cfg.minAlg.c_str()); // Minimize using chosen algorithm
561 if (cfg.hesse)
562 m.hesse(); // Evaluate errors with Hesse
563
564 int corrCovQual = -1;
565
566 if (m.getNPar() > 0) {
567 if (cfg.doAsymptotic == 1)
568 corrCovQual = calcAsymptoticCorrectedCovariance(pdf, m, data); // Asymptotically correct
569 if (cfg.doSumW2 == 1)
571 }
572
573 if (cfg.minos)
574 cfg.minosSet ? m.minos(*cfg.minosSet) : m.minos(); // Evaluate errs with Minos
575
576 // Optionally return fit result
577 std::unique_ptr<RooFitResult> ret;
578 if (cfg.doSave) {
579 auto name = std::string("fitresult_") + pdf.GetName() + "_" + data.GetName();
580 auto title = std::string("Result of fit of p.d.f. ") + pdf.GetName() + " to dataset " + data.GetName();
581 ret = std::unique_ptr<RooFitResult>{m.save(name.c_str(), title.c_str())};
582 if ((cfg.doSumW2 == 1 || cfg.doAsymptotic == 1) && m.getNPar() > 0)
583 ret->setCovQual(corrCovQual);
584 }
585
586 if (cfg.optConst)
587 m.optimizeConst(0);
588 return ret;
589}
590
591std::unique_ptr<RooAbsReal> createNLL(RooAbsPdf &pdf, RooAbsData &data, const RooLinkedList &cmdList)
592{
593 auto timingScope = std::make_unique<ROOT::Math::Util::TimingScope>(
594 [&pdf](std::string const &msg) { oocoutI(&pdf, Fitting) << msg << std::endl; }, "Creation of NLL object took");
595
596 auto baseName = std::string("nll_") + pdf.GetName() + "_" + data.GetName();
597
598 // Select the pdf-specific commands
599 RooCmdConfig pc("RooAbsPdf::createNLL(" + std::string(pdf.GetName()) + ")");
600
601 pc.defineString("rangeName", "RangeWithName", 0, "", true);
602 pc.defineString("addCoefRange", "SumCoefRange", 0, "");
603 pc.defineString("globstag", "GlobalObservablesTag", 0, "");
604 pc.defineString("globssource", "GlobalObservablesSource", 0, "data");
605 pc.defineDouble("rangeLo", "Range", 0, -999.);
606 pc.defineDouble("rangeHi", "Range", 1, -999.);
607 pc.defineInt("splitRange", "SplitRange", 0, 0);
608 pc.defineInt("ext", "Extended", 0, extendedFitDefault);
609 pc.defineInt("numcpu", "NumCPU", 0, 1);
610 pc.defineInt("interleave", "NumCPU", 1, 0);
611 pc.defineInt("verbose", "Verbose", 0, 0);
612 pc.defineInt("optConst", "Optimize", 0, 0);
613 pc.defineInt("cloneData", "CloneData", 0, 2);
614 pc.defineSet("projDepSet", "ProjectedObservables", 0, nullptr);
615 pc.defineSet("cPars", "Constrain", 0, nullptr);
616 pc.defineSet("glObs", "GlobalObservables", 0, nullptr);
617 pc.defineInt("doOffset", "OffsetLikelihood", 0, 0);
618 pc.defineSet("extCons", "ExternalConstraints", 0, nullptr);
619 pc.defineInt("EvalBackend", "EvalBackend", 0, static_cast<int>(RooFit::EvalBackend::defaultValue()));
620 pc.defineDouble("IntegrateBins", "IntegrateBins", 0, -1.);
621 pc.defineMutex("Range", "RangeWithName");
622 pc.defineMutex("GlobalObservables", "GlobalObservablesTag");
623 pc.defineInt("ModularL", "ModularL", 0, 0);
624
625 // New style likelihoods define parallelization through Parallelize(...) on fitTo or attributes on
626 // RooMinimizer::Config.
627 pc.defineMutex("ModularL", "NumCPU");
628
629 // New style likelihoods define offsetting on minimizer, not on likelihood
630 pc.defineMutex("ModularL", "OffsetLikelihood");
631
632 // Process and check varargs
633 pc.process(cmdList);
634 if (!pc.ok(true)) {
635 return nullptr;
636 }
637
638 if (pc.getInt("ModularL")) {
639 int lut[3] = {2, 1, 0};
641 static_cast<RooFit::TestStatistics::RooAbsL::Extended>(lut[pc.getInt("ext")])};
642
646
647 if (auto tmp = pc.getSet("cPars"))
648 cParsSet.add(*tmp);
649
650 if (auto tmp = pc.getSet("extCons"))
651 extConsSet.add(*tmp);
652
653 if (auto tmp = pc.getSet("glObs"))
654 glObsSet.add(*tmp);
655
656 const std::string rangeName = pc.getString("globstag", "", false);
657
659 builder.Extended(ext)
660 .ConstrainedParameters(cParsSet)
661 .ExternalConstraints(extConsSet)
662 .GlobalObservables(glObsSet)
663 .GlobalObservablesTag(rangeName.c_str());
664
665 return std::make_unique<RooFit::TestStatistics::RooRealL>("likelihood", "", builder.build());
666 }
667
668 // Decode command line arguments
669 const char *rangeName = pc.getString("rangeName", nullptr, true);
670 const char *addCoefRangeName = pc.getString("addCoefRange", nullptr, true);
671 const bool ext = interpretExtendedCmdArg(pdf, pc.getInt("ext"));
672
673 int splitRange = pc.getInt("splitRange");
674 int optConst = pc.getInt("optConst");
675 int cloneData = pc.getInt("cloneData");
676 auto offset = static_cast<RooFit::OffsetMode>(pc.getInt("doOffset"));
677
678 // If no explicit cloneData command is specified, cloneData is set to true if optimization is activated
679 if (cloneData == 2) {
681 }
682
683 if (pc.hasProcessed("Range")) {
684 double rangeLo = pc.getDouble("rangeLo");
685 double rangeHi = pc.getDouble("rangeHi");
686
687 // Create range with name 'fit' with above limits on all observables
688 RooArgSet obs;
689 pdf.getObservables(data.get(), obs);
690 for (auto arg : obs) {
691 RooRealVar *rrv = dynamic_cast<RooRealVar *>(arg);
692 if (rrv)
693 rrv->setRange("fit", rangeLo, rangeHi);
694 }
695
696 // Set range name to be fitted to "fit"
697 rangeName = "fit";
698 }
699
700 // Set the fitrange attribute of th PDF, add observables ranges for plotting
702
703 RooArgSet projDeps;
704 auto tmp = pc.getSet("projDepSet");
705 if (tmp) {
706 projDeps.add(*tmp);
707 }
708
709 const std::string globalObservablesSource = pc.getString("globssource", "data", false);
710 if (globalObservablesSource != "data" && globalObservablesSource != "model") {
711 std::string errMsg = "RooAbsPdf::fitTo: GlobalObservablesSource can only be \"data\" or \"model\"!";
712 oocoutE(&pdf, InputArguments) << errMsg << std::endl;
713 throw std::invalid_argument(errMsg);
714 }
716
717 // Lambda function to create the correct constraint term for a PDF. In old
718 // RooFit, we use this PDF itself as the argument, for the new BatchMode
719 // we're passing a clone.
720 auto createConstr = [&]() -> std::unique_ptr<RooAbsReal> {
721 return createConstraintTerm(baseName + "_constr", // name
722 pdf, // pdf
723 data, // data
724 pc.getSet("cPars"), // Constrain RooCmdArg
725 pc.getSet("extCons"), // ExternalConstraints RooCmdArg
726 pc.getSet("glObs"), // GlobalObservables RooCmdArg
727 pc.getString("globstag", nullptr, true), // GlobalObservablesTag RooCmdArg
728 takeGlobalObservablesFromData); // From GlobalObservablesSource RooCmdArg
729 };
730
731 auto evalBackend = static_cast<RooFit::EvalBackend::Value>(pc.getInt("EvalBackend"));
732
733 // Construct BatchModeNLL if requested
735
736 // Set the normalization range. We need to do it now, because it will be
737 // considered in `compileForNormSet`.
738 std::string oldNormRange;
739 if (pdf.normRange()) {
740 oldNormRange = pdf.normRange();
741 }
743
745 pdf.getObservables(data.get(), normSet);
746
747 if (dynamic_cast<RooSimultaneous const *>(&pdf)) {
748 for (auto i : projDeps) {
749 auto res = normSet.find(i->GetName());
750 if (res != nullptr) {
751 res->setAttribute("__conditional__");
752 }
753 }
754 } else {
755 normSet.remove(projDeps);
756 }
757
758 pdf.setAttribute("SplitRange", splitRange);
759 pdf.setStringAttribute("RangeName", rangeName);
760
762 ctx.setLikelihoodMode(true);
763 std::unique_ptr<RooAbsArg> head = pdf.compileForNormSet(normSet, ctx);
764 std::unique_ptr<RooAbsPdf> pdfClone = std::unique_ptr<RooAbsPdf>{&dynamic_cast<RooAbsPdf &>(*head.release())};
765
766 // reset attributes
767 pdf.setAttribute("SplitRange", false);
768 pdf.setStringAttribute("RangeName", nullptr);
769
770 // Reset the normalization range
771 pdf.setNormRange(oldNormRange.c_str());
772
773 if (addCoefRangeName) {
774 oocxcoutI(&pdf, Fitting) << "RooAbsPdf::fitTo(" << pdf.GetName()
775 << ") fixing interpretation of coefficients of any component to range "
776 << addCoefRangeName << "\n";
777 pdfClone->fixAddCoefRange(addCoefRangeName, false);
778 }
779
780 std::unique_ptr<RooAbsReal> compiledConstr;
781 if (std::unique_ptr<RooAbsReal> constr = createConstr()) {
783 compiledConstr->addOwnedComponents(std::move(constr));
784 }
785
786 auto nll = createNLLNew(*pdfClone, data, std::move(compiledConstr), rangeName ? rangeName : "", projDeps, ext,
787 pc.getDouble("IntegrateBins"), offset);
788
789 const double correction = pdfClone->getCorrection();
790
791 if (correction > 0) {
792 oocoutI(&pdf, Fitting) << "[FitHelpers] Detected correction term from RooAbsPdf::getCorrection(). "
793 << "Adding penalty to NLL." << std::endl;
794
795 // Convert the multiplicative correction to an additive term in -log L
796 auto penaltyTerm = std::make_unique<RooConstVar>((baseName + "_Penalty").c_str(),
797 "Penalty term from getCorrection()", correction);
798
799 // add penalty and NLL
800 auto correctedNLL = std::make_unique<RooAddition>((baseName + "_corrected").c_str(), "NLL + penalty",
802
803 // transfer ownership of terms
804 correctedNLL->addOwnedComponents(std::move(nll), std::move(penaltyTerm));
805 nll = std::move(correctedNLL);
806 }
807
808 auto nllWrapper = std::make_unique<RooFit::Experimental::RooEvaluatorWrapper>(
811
812 // We destroy the timing scrope for createNLL prematurely, because we
813 // separately measure the time for jitting and gradient creation
814 // inside the RooFuncWrapper.
815 timingScope.reset();
816
818 nllWrapper->generateGradient();
819 }
821 nllWrapper->setUseGeneratedFunctionCode(true);
822 }
823
824 nllWrapper->addOwnedComponents(std::move(nll));
825 nllWrapper->addOwnedComponents(std::move(pdfClone));
826
827 return nllWrapper;
828 }
829
830 std::unique_ptr<RooAbsReal> nll;
831
832#ifdef ROOFIT_LEGACY_EVAL_BACKEND
833 bool verbose = pc.getInt("verbose");
834
835 int numcpu = pc.getInt("numcpu");
836 int numcpu_strategy = pc.getInt("interleave");
837 // strategy 3 works only for RooSimultaneous.
838 if (numcpu_strategy == 3 && !pdf.InheritsFrom("RooSimultaneous")) {
839 oocoutW(&pdf, Minimization) << "Cannot use a NumCpu Strategy = 3 when the pdf is not a RooSimultaneous, "
840 "falling back to default strategy = 0"
841 << std::endl;
842 numcpu_strategy = 0;
843 }
845
847 RooAbsPdf &actualPdf = binnedLInfo.binnedPdf ? *binnedLInfo.binnedPdf : pdf;
848
849 // Construct NLL
851 RooAbsTestStatistic::Configuration cfg;
852 cfg.addCoefRangeName = addCoefRangeName ? addCoefRangeName : "";
853 cfg.nCPU = numcpu;
854 cfg.interleave = interl;
855 cfg.verbose = verbose;
856 cfg.splitCutRange = static_cast<bool>(splitRange);
857 cfg.cloneInputData = static_cast<bool>(cloneData);
858 cfg.integrateOverBinsPrecision = pc.getDouble("IntegrateBins");
859 cfg.binnedL = binnedLInfo.isBinnedL;
860 cfg.takeGlobalObservablesFromData = takeGlobalObservablesFromData;
861 cfg.rangeName = rangeName ? rangeName : "";
862 auto nllVar = std::make_unique<RooNLLVar>(baseName.c_str(), "-log(likelihood)", actualPdf, data, projDeps, ext, cfg);
863 nllVar->enableBinOffsetting(offset == RooFit::OffsetMode::Bin);
864 nll = std::move(nllVar);
866
867 // Include constraints, if any, in likelihood
868 if (std::unique_ptr<RooAbsReal> constraintTerm = createConstr()) {
869
870 // Even though it is technically only required when the computation graph
871 // is changed because global observables are taken from data, it is safer
872 // to clone the constraint model in general to reset the normalization
873 // integral caches and avoid ASAN build failures (the PDF of the main
874 // measurement is cloned too anyway, so not much overhead). This can be
875 // reconsidered after the caching of normalization sets by pointer is changed
876 // to a more memory-safe solution.
877 constraintTerm = RooHelpers::cloneTreeWithSameParameters(*constraintTerm, data.get());
878
879 // Redirect the global observables to the ones from the dataset if applicable.
880 constraintTerm->setData(data, false);
881
882 // The computation graph for the constraints is very small, no need to do
883 // the tracking of clean and dirty nodes here.
884 constraintTerm->setOperMode(RooAbsArg::ADirty);
885
886 auto orignll = std::move(nll);
887 nll = std::make_unique<RooAddition>((baseName + "_with_constr").c_str(), "nllWithCons",
888 RooArgSet(*orignll, *constraintTerm));
889 nll->addOwnedComponents(std::move(orignll), std::move(constraintTerm));
890 }
891
892 if (optConst) {
893 nll->constOptimizeTestStatistic(RooAbsArg::Activate, optConst > 1);
894 }
895
897 nll->enableOffsetting(true);
898 }
899
900 if (const double correction = pdf.getCorrection(); correction > 0) {
901 oocoutI(&pdf, Fitting) << "[FitHelpers] Detected correction term from RooAbsPdf::getCorrection(). "
902 << "Adding penalty to NLL." << std::endl;
903
904 // Convert the multiplicative correction to an additive term in -log L
905 auto penaltyTerm = std::make_unique<RooConstVar>((baseName + "_Penalty").c_str(),
906 "Penalty term from getCorrection()", correction);
907
908 auto correctedNLL = std::make_unique<RooAddition>(
909 // add penalty and NLL
910 (baseName + "_corrected").c_str(), "NLL + penalty", RooArgSet(*nll, *penaltyTerm));
911
912 // transfer ownership of terms
913 correctedNLL->addOwnedComponents(std::move(nll), std::move(penaltyTerm));
914 nll = std::move(correctedNLL);
915 }
916#else
917 throw std::runtime_error("RooFit was not built with the legacy evaluation backend");
918#endif
919
920 return nll;
921}
922
923std::unique_ptr<RooAbsReal> createChi2(RooAbsReal &real, RooDataHist &data, const RooLinkedList &cmdList)
924{
925#ifdef ROOFIT_LEGACY_EVAL_BACKEND
926 RooCmdConfig pc("createChi2(" + std::string(real.GetName()) + ")");
927
928 pc.defineInt("numcpu", "NumCPU", 0, 1);
929 pc.defineInt("verbose", "Verbose", 0, 0);
930 pc.defineString("rangeName", "RangeWithName", 0, "", true);
931
932 RooAbsTestStatistic::Configuration cfg;
933
934 // Construct Chi2
936 std::string baseName = "chi2_" + std::string(real.GetName()) + "_" + data.GetName();
937
938 // Clear possible range attributes from previous fits.
939 real.removeStringAttribute("fitrange");
940
941 pc.defineInt("etype", "DataError", 0, (Int_t)RooDataHist::Auto);
942 pc.defineInt("extended", "Extended", 0, extendedFitDefault);
943 pc.defineInt("split_range", "SplitRange", 0, 0);
944 pc.defineDouble("integrate_bins", "IntegrateBins", 0, -1);
945 pc.defineString("addCoefRange", "SumCoefRange", 0, "");
946 pc.allowUndefined();
947
948 pc.process(cmdList);
949 if (!pc.ok(true)) {
950 return nullptr;
951 }
952
953 bool extended = false;
954 if (auto pdf = dynamic_cast<RooAbsPdf const *>(&real)) {
955 extended = interpretExtendedCmdArg(*pdf, pc.getInt("extended"));
956 }
957
958 RooDataHist::ErrorType etype = static_cast<RooDataHist::ErrorType>(pc.getInt("etype"));
959
960 const char *rangeName = pc.getString("rangeName", nullptr, true);
961 const char *addCoefRangeName = pc.getString("addCoefRange", nullptr, true);
962 int splitRange = pc.getInt("splitRange");
963
964 // Set the fitrange attribute of th PDF, add observables ranges for plotting
966
967 cfg.rangeName = rangeName ? rangeName : "";
968 cfg.nCPU = pc.getInt("numcpu");
969 cfg.interleave = RooFit::Interleave;
970 cfg.verbose = static_cast<bool>(pc.getInt("verbose"));
971 cfg.cloneInputData = false;
972 cfg.integrateOverBinsPrecision = pc.getDouble("integrate_bins");
973 cfg.addCoefRangeName = addCoefRangeName ? addCoefRangeName : "";
974 cfg.splitCutRange = static_cast<bool>(splitRange);
975 auto chi2 = std::make_unique<RooChi2Var>(baseName.c_str(), baseName.c_str(), real, static_cast<RooDataHist &>(data),
976 extended, etype, cfg);
977
979
980 return chi2;
981#else
982 throw std::runtime_error("createChi2() is not supported without the legacy evaluation backend");
983 return nullptr;
984#endif
985}
986
987std::unique_ptr<RooFitResult> fitTo(RooAbsReal &real, RooAbsData &data, const RooLinkedList &cmdList, bool chi2)
988{
989 const bool isDataHist = dynamic_cast<RooDataHist const *>(&data);
990
991 RooCmdConfig pc("fitTo(" + std::string(real.GetName()) + ")");
992
994 std::string nllCmdListString;
995 if (!chi2) {
996 nllCmdListString = "ProjectedObservables,Extended,Range,"
997 "RangeWithName,SumCoefRange,NumCPU,SplitRange,Constrained,Constrain,ExternalConstraints,"
998 "CloneData,GlobalObservables,GlobalObservablesSource,GlobalObservablesTag,"
999 "EvalBackend,IntegrateBins,ModularL";
1000
1001 if (!cmdList.FindObject("ModularL") || static_cast<RooCmdArg *>(cmdList.FindObject("ModularL"))->getInt(0) == 0) {
1002 nllCmdListString += ",OffsetLikelihood";
1003 }
1004 } else {
1005 auto createChi2DataHistCmdArgs = "Range,RangeWithName,NumCPU,Optimize,IntegrateBins,ProjectedObservables,"
1006 "AddCoefRange,SplitRange,DataError,Extended";
1007 auto createChi2DataSetCmdArgs = "YVar,Integrate,RangeWithName,NumCPU,Verbose";
1009 }
1010
1012
1013 pc.defineDouble("prefit", "Prefit", 0, 0);
1015
1016 // Process and check varargs
1017 pc.process(fitCmdList);
1018 if (!pc.ok(true)) {
1019 return nullptr;
1020 }
1021
1022 // TimingAnalysis works only for RooSimultaneous.
1023 if (pc.getInt("timingAnalysis") && !real.InheritsFrom("RooSimultaneous")) {
1024 oocoutW(&real, Minimization) << "The timingAnalysis feature was built for minimization with RooSimultaneous "
1025 "and is not implemented for other PDF's. Please create a RooSimultaneous to "
1026 "enable this feature."
1027 << std::endl;
1028 }
1029
1030 // Decode command line arguments
1031 double prefit = pc.getDouble("prefit");
1032
1033 if (prefit != 0) {
1034 size_t nEvents = static_cast<size_t>(prefit * data.numEntries());
1035 if (prefit > 0.5 || nEvents < 100) {
1036 oocoutW(&real, InputArguments) << "PrefitDataFraction should be in suitable range."
1037 << "With the current PrefitDataFraction=" << prefit
1038 << ", the number of events would be " << nEvents << " out of "
1039 << data.numEntries() << ". Skipping prefit..." << std::endl;
1040 } else {
1041 size_t step = data.numEntries() / nEvents;
1042
1043 RooDataSet tiny("tiny", "tiny", *data.get(), data.isWeighted() ? RooFit::WeightVar() : RooCmdArg());
1044
1045 for (int i = 0; i < data.numEntries(); i += step) {
1046 const RooArgSet *event = data.get(i);
1047 tiny.add(*event, data.weight());
1048 }
1050 pc.filterCmdList(tinyCmdList, "Prefit,Hesse,Minos,Verbose,Save,Timer");
1053
1056
1057 fitTo(real, tiny, tinyCmdList, chi2);
1058 }
1059 }
1060
1062 if (pc.getInt("parallelize") != 0 || pc.getInt("enableParallelGradient") || pc.getInt("enableParallelDescent")) {
1063 // Set to new style likelihood if parallelization is requested
1066 }
1067
1068 std::unique_ptr<RooAbsReal> nll;
1069 if (chi2) {
1070 if (isDataHist) {
1071 nll = std::unique_ptr<RooAbsReal>{real.createChi2(static_cast<RooDataHist &>(data), nllCmdList)};
1072 }
1073 } else {
1074 nll = std::unique_ptr<RooAbsReal>{dynamic_cast<RooAbsPdf &>(real).createNLL(data, nllCmdList)};
1075 }
1076
1077 return RooFit::FitHelpers::minimize(real, *nll, data, pc);
1078}
1079
1080} // namespace RooFit::FitHelpers
1081
1082/// \endcond
header file containing the templated implementation of matrix inversion routines for use with ROOT's ...
ROOT::RRangeCast< T, false, Range_t > static_range_cast(Range_t &&coll)
#define oocoutW(o, a)
#define oocoutE(o, a)
#define oocoutI(o, a)
#define oocxcoutI(o, a)
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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 offset
char name[80]
Definition TGX11.cxx:148
class to compute the Cholesky decomposition of a matrix
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:76
void setStringAttribute(const Text_t *key, const Text_t *value)
Associate string 'value' to this object under key 'key'.
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
void removeStringAttribute(const Text_t *key)
Delete a string attribute with a given key.
RooFit::OwningPtr< RooArgSet > getVariables(bool stripDisconnected=true) const
Return RooArgSet with all variables (tree leaf nodes of expression tree)
void setAttribute(const Text_t *name, bool value=true)
Set (default) or clear a named boolean attribute of this object.
Abstract base class for objects that represent a discrete value that can be set from the outside,...
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
void assign(const RooAbsCollection &other) const
Sets the value, cache and constant attribute of any argument in our set that also appears in the othe...
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:56
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const override
void setNormRange(const char *rangeName)
virtual double getCorrection() const
This function returns the penalty term.
@ CanBeExtended
Definition RooAbsPdf.h:208
@ MustBeExtended
Definition RooAbsPdf.h:208
@ CanNotBeExtended
Definition RooAbsPdf.h:208
const char * normRange() const
Definition RooAbsPdf.h:246
virtual ExtendMode extendMode() const
Returns ability of PDF to provide extended likelihood terms.
Definition RooAbsPdf.h:212
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
virtual void fixAddCoefNormalization(const RooArgSet &addNormSet=RooArgSet(), bool force=true)
Fix the interpretation of the coefficient of any RooAddPdf component in the expression tree headed by...
static void setEvalErrorLoggingMode(ErrorLoggingMode m)
Set evaluation error logging mode.
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
RooArgSet * selectByAttrib(const char *name, bool value) const
Use RooAbsCollection::selectByAttrib(), but return as RooArgSet.
Definition RooArgSet.h:144
static std::unique_ptr< RooAbsPdf > create(RooAbsPdf &pdf, RooAbsData const &data, double precision)
Creates a wrapping RooBinSamplingPdf if appropriate.
Object to represent discrete states.
Definition RooCategory.h:28
Named container for two doubles, two integers two object points and three string pointers that can be...
Definition RooCmdArg.h:26
Int_t getInt(Int_t idx) const
Definition RooCmdArg.h:87
Configurable parser for RooCmdArg named arguments.
void defineMutex(const char *head, Args_t &&... tail)
Define arguments where any pair is mutually exclusive.
bool process(const RooCmdArg &arg)
Process given RooCmdArg.
bool hasProcessed(const char *cmdName) const
Return true if RooCmdArg with name 'cmdName' has been processed.
double getDouble(const char *name, double defaultValue=0.0) const
Return double property registered with name 'name'.
bool defineDouble(const char *name, const char *argName, int doubleNum, double defValue=0.0)
Define double property name 'name' mapped to double in slot 'doubleNum' in RooCmdArg with name argNam...
RooArgSet * getSet(const char *name, RooArgSet *set=nullptr) const
Return RooArgSet property registered with name 'name'.
bool defineSet(const char *name, const char *argName, int setNum, const RooArgSet *set=nullptr)
Define TObject property name 'name' mapped to object in slot 'setNum' in RooCmdArg with name argName ...
bool ok(bool verbose) const
Return true of parsing was successful.
const char * getString(const char *name, const char *defaultValue="", bool convEmptyToNull=false) const
Return string property registered with name 'name'.
bool defineString(const char *name, const char *argName, int stringNum, const char *defValue="", bool appendMode=false)
Define double property name 'name' mapped to double in slot 'stringNum' in RooCmdArg with name argNam...
bool defineInt(const char *name, const char *argName, int intNum, int defValue=0)
Define integer property name 'name' mapped to integer in slot 'intNum' in RooCmdArg with name argName...
void allowUndefined(bool flag=true)
If flag is true the processing of unrecognized RooCmdArgs is not considered an error.
int getInt(const char *name, int defaultValue=0) const
Return integer property registered with name 'name'.
RooLinkedList filterCmdList(RooLinkedList &cmdInList, const char *cmdNameList, bool removeFromInList=true) const
Utility function to filter commands listed in cmdNameList from cmdInList.
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
Container class to hold unbinned data.
Definition RooDataSet.h:32
static Value & defaultValue()
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
Collection class for internal use, storing a collection of RooAbsArg pointers in a doubly linked list...
Wrapper class around ROOT::Math::Minimizer that provides a seamless interface between the minimizer f...
RooFit::OwningPtr< RooFitResult > save(const char *name=nullptr, const char *title=nullptr)
Save and return a RooFitResult snapshot of current minimizer status.
int hesse()
Execute HESSE.
void applyCovarianceMatrix(TMatrixDSym const &V)
Apply results of given external covariance matrix.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
static TClass * Class()
void setRange(const char *name, double min, double max, bool shared=true)
Set a fit or plotting range.
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:546
RooCmdArg WeightVar(const char *name="weight", bool reinterpretAsWeight=false)
RooCmdArg Hesse(bool flag=true)
RooCmdArg ModularL(bool flag=false)
RooCmdArg PrintLevel(Int_t code)
RVec< PromoteType< T > > log(const RVec< T > &v)
Definition RVec.hxx:1836
CoordSystem::Scalar get(DisplacementVector2D< CoordSystem, Tag > const &p)
std::vector< std::string > Split(std::string_view str, std::string_view delims, bool skipEmpty=false)
Splits a string at each character in delims.
double nll(double pdf, double weight, int binnedL, int doBinOffset)
Definition MathFuncs.h:404
std::unique_ptr< T > compileForNormSet(T const &arg, RooArgSet const &normSet)
OffsetMode
For setting the offset mode with the Offset() command argument to RooAbsPdf::fitTo()
std::unique_ptr< T > cloneTreeWithSameParameters(T const &arg, RooArgSet const *observables=nullptr)
Clone RooAbsArg object and reattach to original parameters.
BinnedLOutput getBinnedL(RooAbsPdf const &pdf)
Config argument to RooMinimizer constructor.
TMarker m
Definition textangle.C:8
TLine l
Definition textangle.C:4