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