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 = false;
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::FitHelpers {
424
426{
427 // Default-initialized instance of MinimizerConfig to get the default
428 // minimizer parameter values.
430
431 pc.defineDouble("RecoverFromUndefinedRegions", "RecoverFromUndefinedRegions", 0, minimizerDefaults.recoverFromNaN);
432 pc.defineInt("optConst", "Optimize", 0, minimizerDefaults.optConst);
433 pc.defineInt("verbose", "Verbose", 0, minimizerDefaults.verbose);
434 pc.defineInt("doSave", "Save", 0, minimizerDefaults.doSave);
435 pc.defineInt("doTimer", "Timer", 0, minimizerDefaults.doTimer);
436 pc.defineInt("printLevel", "PrintLevel", 0, minimizerDefaults.printLevel);
437 pc.defineInt("strategy", "Strategy", 0, minimizerDefaults.strategy);
438 pc.defineInt("initHesse", "InitialHesse", 0, minimizerDefaults.initHesse);
439 pc.defineInt("hesse", "Hesse", 0, minimizerDefaults.hesse);
440 pc.defineInt("minos", "Minos", 0, minimizerDefaults.minos);
441 pc.defineInt("numee", "PrintEvalErrors", 0, minimizerDefaults.numee);
442 pc.defineInt("doEEWall", "EvalErrorWall", 0, minimizerDefaults.doEEWall);
443 pc.defineInt("doWarn", "Warnings", 0, minimizerDefaults.doWarn);
444 pc.defineInt("doSumW2", "SumW2Error", 0, minimizerDefaults.doSumW2);
445 pc.defineInt("doAsymptoticError", "AsymptoticError", 0, minimizerDefaults.doAsymptotic);
446 pc.defineInt("maxCalls", "MaxCalls", 0, minimizerDefaults.maxCalls);
447 pc.defineInt("doOffset", "OffsetLikelihood", 0, minimizerDefaults.doOffset);
448 pc.defineInt("parallelize", "Parallelize", 0, minimizerDefaults.parallelize); // Three parallelize arguments
449 pc.defineInt("enableParallelGradient", "ParallelGradientOptions", 0, minimizerDefaults.enableParallelGradient);
450 pc.defineInt("enableParallelDescent", "ParallelDescentOptions", 0, minimizerDefaults.enableParallelDescent);
451 pc.defineInt("timingAnalysis", "TimingAnalysis", 0, minimizerDefaults.timingAnalysis);
452 pc.defineString("mintype", "Minimizer", 0, minimizerDefaults.minType.c_str());
453 pc.defineString("minalg", "Minimizer", 1, minimizerDefaults.minAlg.c_str());
454 pc.defineSet("minosSet", "Minos", 0, minimizerDefaults.minosSet);
455}
456
457////////////////////////////////////////////////////////////////////////////////
458/// Minimizes a given NLL variable by finding the optimal parameters with the
459/// RooMinimzer. The NLL variable can be created with RooAbsPdf::createNLL.
460/// If you are looking for a function that combines likelihood creation with
461/// fitting, see RooAbsPdf::fitTo.
462/// \param[in] nll The negative log-likelihood variable to minimize.
463/// \param[in] data The dataset that was also used for the NLL. It's a necessary
464/// parameter because it is used in the asymptotic error correction.
465/// \param[in] cfg Configuration struct with all the configuration options for
466/// the RooMinimizer. These are a subset of the options that you can
467/// also pass to RooAbsPdf::fitTo via the RooFit command arguments.
468std::unique_ptr<RooFitResult> minimize(RooAbsReal &pdf, RooAbsReal &nll, RooAbsData const &data, RooCmdConfig const &pc)
469{
470 MinimizerConfig cfg;
471 cfg.recoverFromNaN = pc.getDouble("RecoverFromUndefinedRegions");
472 cfg.optConst = pc.getInt("optConst");
473 cfg.verbose = pc.getInt("verbose");
474 cfg.doSave = pc.getInt("doSave");
475 cfg.doTimer = pc.getInt("doTimer");
476 cfg.printLevel = pc.getInt("printLevel");
477 cfg.strategy = pc.getInt("strategy");
478 cfg.initHesse = pc.getInt("initHesse");
479 cfg.hesse = pc.getInt("hesse");
480 cfg.minos = pc.getInt("minos");
481 cfg.numee = pc.getInt("numee");
482 cfg.doEEWall = pc.getInt("doEEWall");
483 cfg.doWarn = pc.getInt("doWarn");
484 cfg.doSumW2 = pc.getInt("doSumW2");
485 cfg.doAsymptotic = pc.getInt("doAsymptoticError");
486 cfg.maxCalls = pc.getInt("maxCalls");
487 cfg.minosSet = pc.getSet("minosSet");
488 cfg.minType = pc.getString("mintype", "");
489 cfg.minAlg = pc.getString("minalg", "minuit");
490 cfg.doOffset = pc.getInt("doOffset");
491 cfg.parallelize = pc.getInt("parallelize");
492 cfg.enableParallelGradient = pc.getInt("enableParallelGradient");
493 cfg.enableParallelDescent = pc.getInt("enableParallelDescent");
494 cfg.timingAnalysis = pc.getInt("timingAnalysis");
495
496 // Determine if the dataset has weights
497 bool weightedData = data.isNonPoissonWeighted();
498
499 std::string msgPrefix = std::string{"RooAbsPdf::fitTo("} + pdf.GetName() + "): ";
500
501 // Warn user that a method to determine parameter uncertainties should be provided if weighted data is offered
502 if (weightedData && cfg.doSumW2 == -1 && cfg.doAsymptotic == -1) {
503 oocoutW(&pdf, InputArguments) << msgPrefix <<
504 R"(WARNING: a likelihood fit is requested of what appears to be weighted data.
505 While the estimated values of the parameters will always be calculated taking the weights into account,
506 there are multiple ways to estimate the errors of the parameters. You are advised to make an
507 explicit choice for the error calculation:
508 - Either provide SumW2Error(true), to calculate a sum-of-weights-corrected HESSE error matrix
509 (error will be proportional to the number of events in MC).
510 - Or provide SumW2Error(false), to return errors from original HESSE error matrix
511 (which will be proportional to the sum of the weights, i.e., a dataset with <sum of weights> events).
512 - Or provide AsymptoticError(true), to use the asymptotically correct expression
513 (for details see https://arxiv.org/abs/1911.01303)."
514)";
515 }
516
517 if (cfg.minos && (cfg.doSumW2 == 1 || cfg.doAsymptotic == 1)) {
518 oocoutE(&pdf, InputArguments)
519 << msgPrefix
520 << " sum-of-weights and asymptotic error correction do not work with MINOS errors. Not fitting.\n";
521 return nullptr;
522 }
523 if (cfg.doAsymptotic == 1 && cfg.minos) {
524 oocoutW(&pdf, InputArguments) << msgPrefix << "WARNING: asymptotic correction does not apply to MINOS errors\n";
525 }
526
527 // avoid setting both SumW2 and Asymptotic for uncertainty correction
528 if (cfg.doSumW2 == 1 && cfg.doAsymptotic == 1) {
529 oocoutE(&pdf, InputArguments) << msgPrefix
530 << "ERROR: Cannot compute both asymptotically correct and SumW2 errors.\n";
531 return nullptr;
532 }
533
534 // Instantiate RooMinimizer
536 minimizerConfig.enableParallelGradient = cfg.enableParallelGradient;
537 minimizerConfig.enableParallelDescent = cfg.enableParallelDescent;
538 minimizerConfig.parallelize = cfg.parallelize;
539 minimizerConfig.timingAnalysis = cfg.timingAnalysis;
540 minimizerConfig.offsetting = cfg.doOffset;
542
543 m.setMinimizerType(cfg.minType);
544 m.setEvalErrorWall(cfg.doEEWall);
545 m.setRecoverFromNaNStrength(cfg.recoverFromNaN);
546 m.setPrintEvalErrors(cfg.numee);
547 if (cfg.maxCalls > 0)
548 m.setMaxFunctionCalls(cfg.maxCalls);
549 if (cfg.printLevel != 1)
550 m.setPrintLevel(cfg.printLevel);
551 if (cfg.optConst)
552 m.optimizeConst(cfg.optConst); // Activate constant term optimization
553 if (cfg.verbose)
554 m.setVerbose(true); // Activate verbose options
555 if (cfg.doTimer)
556 m.setProfile(true); // Activate timer options
557 if (cfg.strategy != 1)
558 m.setStrategy(cfg.strategy); // Modify fit strategy
559 if (cfg.initHesse)
560 m.hesse(); // Initialize errors with hesse
561 m.minimize(cfg.minType.c_str(), cfg.minAlg.c_str()); // Minimize using chosen algorithm
562 if (cfg.hesse)
563 m.hesse(); // Evaluate errors with Hesse
564
565 int corrCovQual = -1;
566
567 if (m.getNPar() > 0) {
568 if (cfg.doAsymptotic == 1)
569 corrCovQual = calcAsymptoticCorrectedCovariance(pdf, m, data); // Asymptotically correct
570 if (cfg.doSumW2 == 1)
572 }
573
574 if (cfg.minos)
575 cfg.minosSet ? m.minos(*cfg.minosSet) : m.minos(); // Evaluate errs with Minos
576
577 // Optionally return fit result
578 std::unique_ptr<RooFitResult> ret;
579 if (cfg.doSave) {
580 auto name = std::string("fitresult_") + pdf.GetName() + "_" + data.GetName();
581 auto title = std::string("Result of fit of p.d.f. ") + pdf.GetName() + " to dataset " + data.GetName();
582 ret = std::unique_ptr<RooFitResult>{m.save(name.c_str(), title.c_str())};
583 if ((cfg.doSumW2 == 1 || cfg.doAsymptotic == 1) && m.getNPar() > 0)
584 ret->setCovQual(corrCovQual);
585 }
586
587 if (cfg.optConst)
588 m.optimizeConst(0);
589 return ret;
590}
591
592std::unique_ptr<RooAbsReal> createNLL(RooAbsPdf &pdf, RooAbsData &data, const RooLinkedList &cmdList)
593{
594 auto timingScope = std::make_unique<ROOT::Math::Util::TimingScope>(
595 [&pdf](std::string const &msg) { oocoutI(&pdf, Fitting) << msg << std::endl; }, "Creation of NLL object took");
596
597 auto baseName = std::string("nll_") + pdf.GetName() + "_" + data.GetName();
598
599 // Select the pdf-specific commands
600 RooCmdConfig pc("RooAbsPdf::createNLL(" + std::string(pdf.GetName()) + ")");
601
602 pc.defineString("rangeName", "RangeWithName", 0, "", true);
603 pc.defineString("addCoefRange", "SumCoefRange", 0, "");
604 pc.defineString("globstag", "GlobalObservablesTag", 0, "");
605 pc.defineString("globssource", "GlobalObservablesSource", 0, "data");
606 pc.defineDouble("rangeLo", "Range", 0, -999.);
607 pc.defineDouble("rangeHi", "Range", 1, -999.);
608 pc.defineInt("splitRange", "SplitRange", 0, 0);
609 pc.defineInt("ext", "Extended", 0, extendedFitDefault);
610 pc.defineInt("numcpu", "NumCPU", 0, 1);
611 pc.defineInt("interleave", "NumCPU", 1, 0);
612 pc.defineInt("verbose", "Verbose", 0, 0);
613 pc.defineInt("optConst", "Optimize", 0, 0);
614 pc.defineInt("cloneData", "CloneData", 0, 2);
615 pc.defineSet("projDepSet", "ProjectedObservables", 0, nullptr);
616 pc.defineSet("cPars", "Constrain", 0, nullptr);
617 pc.defineSet("glObs", "GlobalObservables", 0, nullptr);
618 pc.defineInt("doOffset", "OffsetLikelihood", 0, 0);
619 pc.defineSet("extCons", "ExternalConstraints", 0, nullptr);
620 pc.defineInt("EvalBackend", "EvalBackend", 0, static_cast<int>(RooFit::EvalBackend::defaultValue()));
621 pc.defineDouble("IntegrateBins", "IntegrateBins", 0, -1.);
622 pc.defineMutex("Range", "RangeWithName");
623 pc.defineMutex("GlobalObservables", "GlobalObservablesTag");
624 pc.defineInt("ModularL", "ModularL", 0, 0);
625
626 // New style likelihoods define parallelization through Parallelize(...) on fitTo or attributes on
627 // RooMinimizer::Config.
628 pc.defineMutex("ModularL", "NumCPU");
629
630 // New style likelihoods define offsetting on minimizer, not on likelihood
631 pc.defineMutex("ModularL", "OffsetLikelihood");
632
633 // Process and check varargs
634 pc.process(cmdList);
635 if (!pc.ok(true)) {
636 return nullptr;
637 }
638
639 if (pc.getInt("ModularL")) {
640 int lut[3] = {2, 1, 0};
642 static_cast<RooFit::TestStatistics::RooAbsL::Extended>(lut[pc.getInt("ext")])};
643
647
648 if (auto tmp = pc.getSet("cPars"))
649 cParsSet.add(*tmp);
650
651 if (auto tmp = pc.getSet("extCons"))
652 extConsSet.add(*tmp);
653
654 if (auto tmp = pc.getSet("glObs"))
655 glObsSet.add(*tmp);
656
657 const std::string rangeName = pc.getString("globstag", "", false);
658
660 builder.Extended(ext)
661 .ConstrainedParameters(cParsSet)
662 .ExternalConstraints(extConsSet)
663 .GlobalObservables(glObsSet)
664 .GlobalObservablesTag(rangeName.c_str());
665
666 return std::make_unique<RooFit::TestStatistics::RooRealL>("likelihood", "", builder.build());
667 }
668
669 // Decode command line arguments
670 const char *rangeName = pc.getString("rangeName", nullptr, true);
671 const char *addCoefRangeName = pc.getString("addCoefRange", nullptr, true);
672 const bool ext = interpretExtendedCmdArg(pdf, pc.getInt("ext"));
673
674 int splitRange = pc.getInt("splitRange");
675 int optConst = pc.getInt("optConst");
676 int cloneData = pc.getInt("cloneData");
677 auto offset = static_cast<RooFit::OffsetMode>(pc.getInt("doOffset"));
678
679 // If no explicit cloneData command is specified, cloneData is set to true if optimization is activated
680 if (cloneData == 2) {
682 }
683
684 if (pc.hasProcessed("Range")) {
685 double rangeLo = pc.getDouble("rangeLo");
686 double rangeHi = pc.getDouble("rangeHi");
687
688 // Create range with name 'fit' with above limits on all observables
689 RooArgSet obs;
690 pdf.getObservables(data.get(), obs);
691 for (auto arg : obs) {
692 RooRealVar *rrv = dynamic_cast<RooRealVar *>(arg);
693 if (rrv)
694 rrv->setRange("fit", rangeLo, rangeHi);
695 }
696
697 // Set range name to be fitted to "fit"
698 rangeName = "fit";
699 }
700
701 // Set the fitrange attribute of th PDF, add observables ranges for plotting
703
704 RooArgSet projDeps;
705 auto getSet("projDepSet");
706 if (tmp) {
707 projDeps.add(*tmp);
708 }
709
710 const std::string globalObservablesSource = pc.getString("globssource", "data", false);
711 if (globalObservablesSource != "data" && globalObservablesSource != "model") {
712 std::string errMsg = "RooAbsPdf::fitTo: GlobalObservablesSource can only be \"data\" or \"model\"!";
713 oocoutE(&pdf, InputArguments) << errMsg << std::endl;
714 throw std::invalid_argument(errMsg);
715 }
717
718 // Lambda function to create the correct constraint term for a PDF. In old
719 // RooFit, we use this PDF itself as the argument, for the new BatchMode
720 // we're passing a clone.
721 auto createConstr = [&]() -> std::unique_ptr<RooAbsReal> {
722 return createConstraintTerm(baseName + "_constr", // name
723 pdf, // pdf
724 data, // data
725 pc.getSet("cPars"), // Constrain RooCmdArg
726 pc.getSet("extCons"), // ExternalConstraints RooCmdArg
727 pc.getSet("glObs"), // GlobalObservables RooCmdArg
728 pc.getString("globstag", nullptr, true), // GlobalObservablesTag RooCmdArg
729 takeGlobalObservablesFromData); // From GlobalObservablesSource RooCmdArg
730 };
731
732 auto evalBackend = static_cast<RooFit::EvalBackend::Value>(pc.getInt("EvalBackend"));
733
734 // Construct BatchModeNLL if requested
736
737 // Set the normalization range. We need to do it now, because it will be
738 // considered in `compileForNormSet`.
739 std::string oldNormRange;
740 if (pdf.normRange()) {
741 oldNormRange = pdf.normRange();
742 }
744
746 pdf.getObservables(data.get(), normSet);
747 normSet.remove(projDeps, true, true);
748
749 pdf.setAttribute("SplitRange", splitRange);
750 pdf.setStringAttribute("RangeName", rangeName);
751
753 ctx.setLikelihoodMode(true);
754 std::unique_ptr<RooAbsArg> head = pdf.compileForNormSet(normSet, ctx);
755 std::unique_ptr<RooAbsPdf> pdfClone = std::unique_ptr<RooAbsPdf>{static_cast<RooAbsPdf *>(head.release())};
756
757 // reset attributes
758 pdf.setAttribute("SplitRange", false);
759 pdf.setStringAttribute("RangeName", nullptr);
760
761 // Reset the normalization range
762 pdf.setNormRange(oldNormRange.c_str());
763
764 if (addCoefRangeName) {
765 oocxcoutI(&pdf, Fitting) << "RooAbsPdf::fitTo(" << pdf.GetName()
766 << ") fixing interpretation of coefficients of any component to range "
767 << addCoefRangeName << "\n";
768 pdfClone->fixAddCoefRange(addCoefRangeName, false);
769 }
770
771 std::unique_ptr<RooAbsReal> compiledConstr;
772 if (std::unique_ptr<RooAbsReal> constr = createConstr()) {
774 compiledConstr->addOwnedComponents(std::move(constr));
775 }
776
777 auto nll = createNLLNew(*pdfClone, data, std::move(compiledConstr), rangeName ? rangeName : "", projDeps, ext,
778 pc.getDouble("IntegrateBins"), offset);
779
780 std::unique_ptr<RooAbsReal> nllWrapper;
781
784 bool createGradient = evalBackend == RooFit::EvalBackend::Value::Codegen;
785 auto simPdf = dynamic_cast<RooSimultaneous const *>(pdfClone.get());
786
787 // We destroy the timing scrope for createNLL prematurely, because we
788 // separately measure the time for jitting and gradient creation
789 // inside the RooFuncWrapper.
790 timingScope.reset();
791
792 nllWrapper = std::make_unique<RooFit::Experimental::RooFuncWrapper>("nll_func_wrapper", "nll_func_wrapper",
793 *nll, &data, simPdf, createGradient);
794 if (createGradient)
795 static_cast<Experimental::RooFuncWrapper &>(*nllWrapper).createGradient();
796 } else {
797 nllWrapper = std::make_unique<RooEvaluatorWrapper>(
800 }
801
802 nllWrapper->addOwnedComponents(std::move(nll));
803 nllWrapper->addOwnedComponents(std::move(pdfClone));
804 return nllWrapper;
805 }
806
807 std::unique_ptr<RooAbsReal> nll;
808
809#ifdef ROOFIT_LEGACY_EVAL_BACKEND
810 bool verbose = pc.getInt("verbose");
811
812 int numcpu = pc.getInt("numcpu");
813 int numcpu_strategy = pc.getInt("interleave");
814 // strategy 3 works only for RooSimultaneous.
815 if (numcpu_strategy == 3 && !pdf.InheritsFrom("RooSimultaneous")) {
816 oocoutW(&pdf, Minimization) << "Cannot use a NumCpu Strategy = 3 when the pdf is not a RooSimultaneous, "
817 "falling back to default strategy = 0"
818 << std::endl;
819 numcpu_strategy = 0;
820 }
822
824 RooAbsPdf &actualPdf = binnedLInfo.binnedPdf ? *binnedLInfo.binnedPdf : pdf;
825
826 // Construct NLL
828 RooAbsTestStatistic::Configuration cfg;
829 cfg.addCoefRangeName = addCoefRangeName ? addCoefRangeName : "";
830 cfg.nCPU = numcpu;
831 cfg.interleave = interl;
832 cfg.verbose = verbose;
833 cfg.splitCutRange = static_cast<bool>(splitRange);
834 cfg.cloneInputData = static_cast<bool>(cloneData);
835 cfg.integrateOverBinsPrecision = pc.getDouble("IntegrateBins");
836 cfg.binnedL = binnedLInfo.isBinnedL;
837 cfg.takeGlobalObservablesFromData = takeGlobalObservablesFromData;
838 cfg.rangeName = rangeName ? rangeName : "";
839 auto nllVar = std::make_unique<RooNLLVar>(baseName.c_str(), "-log(likelihood)", actualPdf, data, projDeps, ext, cfg);
840 nllVar->enableBinOffsetting(offset == RooFit::OffsetMode::Bin);
841 nll = std::move(nllVar);
843
844 // Include constraints, if any, in likelihood
845 if (std::unique_ptr<RooAbsReal> constraintTerm = createConstr()) {
846
847 // Even though it is technically only required when the computation graph
848 // is changed because global observables are taken from data, it is safer
849 // to clone the constraint model in general to reset the normalization
850 // integral caches and avoid ASAN build failures (the PDF of the main
851 // measurement is cloned too anyway, so not much overhead). This can be
852 // reconsidered after the caching of normalization sets by pointer is changed
853 // to a more memory-safe solution.
854 constraintTerm = RooHelpers::cloneTreeWithSameParameters(*constraintTerm, data.get());
855
856 // Redirect the global observables to the ones from the dataset if applicable.
857 constraintTerm->setData(data, false);
858
859 // The computation graph for the constraints is very small, no need to do
860 // the tracking of clean and dirty nodes here.
861 constraintTerm->setOperMode(RooAbsArg::ADirty);
862
863 auto orignll = std::move(nll);
864 nll = std::make_unique<RooAddition>((baseName + "_with_constr").c_str(), "nllWithCons",
865 RooArgSet(*orignll, *constraintTerm));
866 nll->addOwnedComponents(std::move(orignll), std::move(constraintTerm));
867 }
868
869 if (optConst) {
870 nll->constOptimizeTestStatistic(RooAbsArg::Activate, optConst > 1);
871 }
872
874 nll->enableOffsetting(true);
875 }
876#else
877 throw std::runtime_error("RooFit was not built with the legacy evaluation backend");
878#endif
879
880 return nll;
881}
882
883std::unique_ptr<RooAbsReal> createChi2(RooAbsReal &real, RooAbsData &data, const RooLinkedList &cmdList)
884{
885#ifdef ROOFIT_LEGACY_EVAL_BACKEND
886 const bool isDataHist = dynamic_cast<RooDataHist const *>(&data);
887
888 RooCmdConfig pc("createChi2(" + std::string(real.GetName()) + ")");
889
890 pc.defineInt("numcpu", "NumCPU", 0, 1);
891 pc.defineInt("verbose", "Verbose", 0, 0);
892
893 RooAbsTestStatistic::Configuration cfg;
894
895 if (isDataHist) {
896 // Construct Chi2
898 std::string baseName = "chi2_" + std::string(real.GetName()) + "_" + data.GetName();
899
900 // Clear possible range attributes from previous fits.
901 real.removeStringAttribute("fitrange");
902
903 pc.defineInt("etype", "DataError", 0, (Int_t)RooDataHist::Auto);
904 pc.defineInt("extended", "Extended", 0, extendedFitDefault);
905 pc.defineInt("split_range", "SplitRange", 0, 0);
906 pc.defineDouble("integrate_bins", "IntegrateBins", 0, -1);
907 pc.defineString("addCoefRange", "SumCoefRange", 0, "");
908 pc.allowUndefined();
909
910 pc.process(cmdList);
911 if (!pc.ok(true)) {
912 return nullptr;
913 }
914
915 bool extended = false;
916 if (auto pdf = dynamic_cast<RooAbsPdf const *>(&real)) {
917 extended = interpretExtendedCmdArg(*pdf, pc.getInt("extended"));
918 }
919
920 RooDataHist::ErrorType etype = static_cast<RooDataHist::ErrorType>(pc.getInt("etype"));
921
922 const char *rangeName = pc.getString("rangeName", nullptr, true);
923 const char *addCoefRangeName = pc.getString("addCoefRange", nullptr, true);
924
925 cfg.rangeName = rangeName ? rangeName : "";
926 cfg.nCPU = pc.getInt("numcpu");
927 cfg.interleave = RooFit::Interleave;
928 cfg.verbose = static_cast<bool>(pc.getInt("verbose"));
929 cfg.cloneInputData = false;
930 cfg.integrateOverBinsPrecision = pc.getDouble("integrate_bins");
931 cfg.addCoefRangeName = addCoefRangeName ? addCoefRangeName : "";
932 cfg.splitCutRange = static_cast<bool>(pc.getInt("split_range"));
933 auto chi2 = std::make_unique<RooChi2Var>(baseName.c_str(), baseName.c_str(), real,
934 static_cast<RooDataHist &>(data), extended, etype, cfg);
935
937
938 return chi2;
939 } else {
940 pc.defineInt("integrate", "Integrate", 0, 0);
941 pc.defineObject("yvar", "YVar", 0, nullptr);
942 pc.defineString("rangeName", "RangeWithName", 0, "", true);
943 pc.defineInt("interleave", "NumCPU", 1, 0);
944
945 // Process and check varargs
946 pc.process(cmdList);
947 if (!pc.ok(true)) {
948 return nullptr;
949 }
950
951 // Decode command line arguments
952 bool integrate = pc.getInt("integrate");
953 RooRealVar *yvar = static_cast<RooRealVar *>(pc.getObject("yvar"));
954 const char *rangeName = pc.getString("rangeName", nullptr, true);
955 Int_t numcpu = pc.getInt("numcpu");
956 Int_t numcpu_strategy = pc.getInt("interleave");
957 // strategy 3 works only for RooSimultaneous.
958 if (numcpu_strategy == 3 && !real.InheritsFrom("RooSimultaneous")) {
959 oocoutW(&real, Minimization) << "Cannot use a NumCpu Strategy = 3 when the pdf is not a RooSimultaneous, "
960 "falling back to default strategy = 0"
961 << std::endl;
962 numcpu_strategy = 0;
963 }
965 bool verbose = pc.getInt("verbose");
966
967 cfg.rangeName = rangeName ? rangeName : "";
968 cfg.nCPU = numcpu;
969 cfg.interleave = interl;
970 cfg.verbose = verbose;
971 cfg.verbose = false;
972
973 std::string name = "chi2_" + std::string(real.GetName()) + "_" + data.GetName();
974
975 return std::make_unique<RooXYChi2Var>(name.c_str(), name.c_str(), real, static_cast<RooDataSet &>(data), yvar,
976 integrate, cfg);
977 }
978#else
979 throw std::runtime_error("createChi2() is not supported without the legacy evaluation backend");
980 return nullptr;
981#endif
982}
983
984std::unique_ptr<RooFitResult> fitTo(RooAbsReal &real, RooAbsData &data, const RooLinkedList &cmdList, bool chi2)
985{
986 const bool isDataHist = dynamic_cast<RooDataHist const *>(&data);
987
988 RooCmdConfig pc("fitTo(" + std::string(real.GetName()) + ")");
989
991 std::string nllCmdListString;
992 if (!chi2) {
993 nllCmdListString = "ProjectedObservables,Extended,Range,"
994 "RangeWithName,SumCoefRange,NumCPU,SplitRange,Constrained,Constrain,ExternalConstraints,"
995 "CloneData,GlobalObservables,GlobalObservablesSource,GlobalObservablesTag,"
996 "EvalBackend,IntegrateBins,ModularL";
997
998 if (!cmdList.FindObject("ModularL") || static_cast<RooCmdArg *>(cmdList.FindObject("ModularL"))->getInt(0) == 0) {
999 nllCmdListString += ",OffsetLikelihood";
1000 }
1001 } else {
1002 auto createChi2DataHistCmdArgs = "Range,RangeWithName,NumCPU,Optimize,IntegrateBins,ProjectedObservables,"
1003 "AddCoefRange,SplitRange,DataError,Extended";
1004 auto createChi2DataSetCmdArgs = "YVar,Integrate,RangeWithName,NumCPU,Verbose";
1006 }
1007
1009
1010 pc.defineDouble("prefit", "Prefit", 0, 0);
1012
1013 // Process and check varargs
1014 pc.process(fitCmdList);
1015 if (!pc.ok(true)) {
1016 return nullptr;
1017 }
1018
1019 // TimingAnalysis works only for RooSimultaneous.
1020 if (pc.getInt("timingAnalysis") && !real.InheritsFrom("RooSimultaneous")) {
1021 oocoutW(&real, Minimization) << "The timingAnalysis feature was built for minimization with RooSimultaneous "
1022 "and is not implemented for other PDF's. Please create a RooSimultaneous to "
1023 "enable this feature."
1024 << std::endl;
1025 }
1026
1027 // Decode command line arguments
1028 double prefit = pc.getDouble("prefit");
1029
1030 if (prefit != 0) {
1031 size_t nEvents = static_cast<size_t>(prefit * data.numEntries());
1032 if (prefit > 0.5 || nEvents < 100) {
1033 oocoutW(&real, InputArguments) << "PrefitDataFraction should be in suitable range."
1034 << "With the current PrefitDataFraction=" << prefit
1035 << ", the number of events would be " << nEvents << " out of "
1036 << data.numEntries() << ". Skipping prefit..." << std::endl;
1037 } else {
1038 size_t step = data.numEntries() / nEvents;
1039
1040 RooDataSet tiny("tiny", "tiny", *data.get(), data.isWeighted() ? RooFit::WeightVar() : RooCmdArg());
1041
1042 for (int i = 0; i < data.numEntries(); i += step) {
1043 const RooArgSet *event = data.get(i);
1044 tiny.add(*event, data.weight());
1045 }
1047 pc.filterCmdList(tinyCmdList, "Prefit,Hesse,Minos,Verbose,Save,Timer");
1050
1053
1054 fitTo(real, tiny, tinyCmdList, chi2);
1055 }
1056 }
1057
1059 if (pc.getInt("parallelize") != 0 || pc.getInt("enableParallelGradient") || pc.getInt("enableParallelDescent")) {
1060 // Set to new style likelihood if parallelization is requested
1063 }
1064
1065 std::unique_ptr<RooAbsReal> nll;
1066 if (chi2) {
1067 nll = std::unique_ptr<RooAbsReal>{isDataHist ? real.createChi2(static_cast<RooDataHist &>(data), nllCmdList)
1068 : real.createChi2(static_cast<RooDataSet &>(data), nllCmdList)};
1069 } else {
1070 nll = std::unique_ptr<RooAbsReal>{dynamic_cast<RooAbsPdf &>(real).createNLL(data, nllCmdList)};
1071 }
1072
1073 return RooFit::FitHelpers::minimize(real, *nll, data, pc);
1074}
1075
1076} // namespace RooFit::FitHelpers
1077
1078/// \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)
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