Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
FitHelpers.cxx
Go to the documentation of this file.
1/// \cond ROOFIT_INTERNAL
2
3/*
4 * Project: RooFit
5 * Authors:
6 * Jonas Rembser, CERN 2023
7 *
8 * Copyright (c) 2023, CERN
9 *
10 * Redistribution and use in source and binary forms,
11 * with or without modification, are permitted according to the terms
12 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
13 */
14
15#include "FitHelpers.h"
16
17#include <RooAbsData.h>
18#include <RooAbsPdf.h>
19#include <RooAbsReal.h>
20#include <RooAddition.h>
21#include <RooBatchCompute.h>
22#include <RooBinSamplingPdf.h>
23#include <RooCategory.h>
24#include <RooCmdConfig.h>
25#include <RooConstraintSum.h>
26#include <RooDataHist.h>
27#include <RooDataSet.h>
28#include <RooDerivative.h>
29#include <RooFit/Evaluator.h>
32#include <RooFitResult.h>
33#include <RooLinkedList.h>
34#include <RooMinimizer.h>
35#include <RooConstVar.h>
36#include <RooRealVar.h>
37#include <RooSimultaneous.h>
38#include <RooFormulaVar.h>
39
40#include <Math/CholeskyDecomp.h>
41#include <Math/Util.h>
42
43#include "ConstraintHelpers.h"
44#include "RooEvaluatorWrapper.h"
45#include "RooFitImplHelpers.h"
47
48#ifdef ROOFIT_LEGACY_EVAL_BACKEND
49#include "RooChi2Var.h"
50#include "RooNLLVar.h"
51#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 = 0;
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
323/// Iterate the simultaneous pdf's categories and build one test-statistic term
324/// per channel via `makeTerm`. Channels excluded by the (optional) category
325/// range are skipped. The per-channel term's special variables are prefixed
326/// with `_<catName>_`. The terms are summed into a RooAddition named
327/// `combinedName`.
328template <typename TermFactory>
329std::unique_ptr<RooAddition> createSimultaneousStat(RooSimultaneous const &simPdf, std::string const &rangeName,
330 std::string const &combinedName, TermFactory &&makeTerm)
331{
332 RooAbsCategoryLValue const &simCat = simPdf.indexCat();
333
334 RooArgList terms;
335 for (auto const &catState : simCat) {
336 std::string const &catName = catState.first;
338
339 // Skip channels excluded by a category range (only RooCategory supports
340 // ranges on categorical values).
341 if (!rangeName.empty()) {
342 auto simCatAsRooCategory = dynamic_cast<RooCategory const *>(&simCat);
343 if (simCatAsRooCategory && !simCatAsRooCategory->isStateInRange(rangeName.c_str(), catIndex)) {
344 continue;
345 }
346 }
347
348 RooAbsPdf *channelPdf = simPdf.getPdf(catName.c_str());
349 if (!channelPdf) {
350 continue;
351 }
352 std::unique_ptr<RooArgSet> observables{
353 std::unique_ptr<RooArgSet>(channelPdf->getVariables())->selectByAttrib("__obs__", true)};
354 std::unique_ptr<RooNLLVarNew> term = makeTerm(*channelPdf, *observables);
355 term->setPrefix(std::string("_") + catName + "_");
356 terms.addOwned(std::move(term));
357 }
358
359 auto combined = std::make_unique<RooAddition>(combinedName.c_str(), combinedName.c_str(), terms);
360 combined->addOwnedComponents(std::move(terms));
361 return combined;
362}
363
364std::unique_ptr<RooAbsArg> createSimultaneousChi2(RooSimultaneous const &simPdf, std::string const &rangeName,
366{
367 auto chi2 =
368 createSimultaneousStat(simPdf, rangeName, "simChi2", [&](RooAbsPdf &channelPdf, RooArgSet const &observables) {
369 RooNLLVarNew::Config cfg;
370 cfg.statistic = RooNLLVarNew::Statistic::Chi2;
371 cfg.extended = isSimPdfExtended && channelPdf.extendMode() != RooAbsPdf::CanNotBeExtended;
372 cfg.chi2ErrorType = etype;
373 auto name = std::string("chi2_") + channelPdf.GetName();
374 return std::make_unique<RooNLLVarNew>(name.c_str(), name.c_str(), channelPdf, observables, cfg);
375 });
376 // Flag the top node so RooEvaluatorWrapper knows not to skip zero-weight bins
377 chi2->setAttribute("Chi2EvaluationActive");
378 return chi2;
379}
380
381std::unique_ptr<RooAbsArg> createSimultaneousNLL(RooSimultaneous const &simPdf, bool isSimPdfExtended,
382 std::string const &rangeName, RooFit::OffsetMode offset)
383{
384 auto nll =
385 createSimultaneousStat(simPdf, rangeName, "mynll", [&](RooAbsPdf &channelPdf, RooArgSet const &observables) {
386 RooNLLVarNew::Config cfg;
387 // Only request extended NLLs for channels that can be extended.
388 cfg.extended = isSimPdfExtended && channelPdf.extendMode() != RooAbsPdf::CanNotBeExtended;
389 cfg.offsetMode = offset;
390 auto name = std::string("nll_") + channelPdf.GetName();
391 return std::make_unique<RooNLLVarNew>(name.c_str(), name.c_str(), channelPdf, observables, cfg);
392 });
393
394 const int simCount = nll->list().size();
395 for (auto *child : static_range_cast<RooNLLVarNew *>(nll->list())) {
396 child->setSimCount(simCount);
397 }
398 return nll;
399}
400
401/// Apply the `IntegrateBins` precision to either a single pdf or in-place to
402/// the component pdfs of a RooSimultaneous. Newly-allocated wrapper pdfs are
403/// appended to `ownedOut`. The returned reference points either at one of
404/// those wrappers or at the input pdf itself.
406{
407 if (auto *simPdf = dynamic_cast<RooSimultaneous *>(&pdf)) {
408 simPdf->wrapPdfsInBinSamplingPdfs(data, precision);
409 return pdf;
410 }
411 if (std::unique_ptr<RooAbsPdf> wrapped = RooBinSamplingPdf::create(pdf, data, precision)) {
413 ownedOut.addOwned(std::move(wrapped));
414 return ref;
415 }
416 return pdf;
417}
418
419/// RAII helper for the `setNormRange` + SplitRange/RangeName attribute setting
420/// that must be done around `compileForNormSet` so that the compiled pdf sees
421/// the restricted normalization range. All mutations are reverted on
422/// destruction.
423class NormRangeScope {
424public:
425 NormRangeScope(RooAbsPdf &pdf, const char *rangeName, bool splitRange) : _pdf{&pdf}
426 {
427 if (pdf.normRange()) {
428 _oldNormRange = pdf.normRange();
429 }
431 pdf.setAttribute("SplitRange", splitRange);
432 pdf.setStringAttribute("RangeName", rangeName);
433 }
434 NormRangeScope(NormRangeScope const &) = delete;
435 NormRangeScope &operator=(NormRangeScope const &) = delete;
437 {
438 _pdf->setAttribute("SplitRange", false);
439 _pdf->setStringAttribute("RangeName", nullptr);
440 _pdf->setNormRange(_oldNormRange.empty() ? nullptr : _oldNormRange.c_str());
441 }
442
443private:
444 RooAbsPdf *_pdf;
445 std::string _oldNormRange;
446};
447
448/// Shared `compileForNormSet` sequence used by both the NLL and chi2 CPU
449/// backends. Sets up the reduced normalization range via `NormRangeScope`,
450/// then runs the graph compilation. When `likelihoodMode` is true the
451/// CompileContext is flagged accordingly (enabling binned-likelihood
452/// optimisations in the compiled pdf). The returned compiled pdf has
453/// `fixAddCoefRange` already applied when `addCoefRangeName` is non-empty.
454std::unique_ptr<RooAbsPdf> compilePdfForFit(RooAbsPdf &pdf, RooArgSet const &normSet, const char *rangeName,
455 bool splitRange, const char *addCoefRangeName, bool likelihoodMode)
456{
458
460 ctx.setLikelihoodMode(likelihoodMode);
461 std::unique_ptr<RooAbsArg> head = pdf.compileForNormSet(normSet, ctx);
462 std::unique_ptr<RooAbsPdf> pdfClone{&dynamic_cast<RooAbsPdf &>(*head.release())};
463
465 pdfClone->fixAddCoefRange(addCoefRangeName, false);
466 }
467 return pdfClone;
468}
469
470std::unique_ptr<RooAbsReal> createNLLNew(RooAbsPdf &pdf, RooAbsData &data, std::unique_ptr<RooAbsReal> &&constraints,
471 std::string const &rangeName, RooArgSet const &projDeps, bool isExtended,
473{
474 if (constraints) {
475 // The computation graph for the constraints is very small, no need to do
476 // the tracking of clean and dirty nodes here.
477 constraints->setOperMode(RooAbsArg::ADirty);
478 }
479
480 RooArgSet observables;
481 pdf.getObservables(data.get(), observables);
482 observables.remove(projDeps, true, true);
483
484 oocxcoutI(&pdf, Fitting) << "RooAbsPdf::fitTo(" << pdf.GetName()
485 << ") fixing normalization set for coefficient determination to observables in data"
486 << "\n";
487 pdf.fixAddCoefNormalization(observables, false);
488
491
493 if (auto *simPdf = dynamic_cast<RooSimultaneous *>(&finalPdf)) {
494 nllTerms.addOwned(createSimultaneousNLL(*simPdf, isExtended, rangeName, offset));
495 } else {
496 RooNLLVarNew::Config cfg;
497 cfg.extended = isExtended;
498 cfg.offsetMode = offset;
499 nllTerms.addOwned(std::make_unique<RooNLLVarNew>("RooNLLVarNew", "RooNLLVarNew", finalPdf, observables, cfg));
500 }
501 if (constraints) {
502 nllTerms.addOwned(std::move(constraints));
503 }
504
505 std::string nllName = std::string("nll_") + pdf.GetName() + "_" + data.GetName();
506 auto nll = std::make_unique<RooAddition>(nllName.c_str(), nllName.c_str(), nllTerms);
507 nll->addOwnedComponents(std::move(binSamplingPdfs));
508 nll->addOwnedComponents(std::move(nllTerms));
509
510 return nll;
511}
512
513} // namespace
514
515namespace RooFit::FitHelpers {
516
518{
519 // Default-initialized instance of MinimizerConfig to get the default
520 // minimizer parameter values.
522
523 pc.defineDouble("RecoverFromUndefinedRegions", "RecoverFromUndefinedRegions", 0, minimizerDefaults.recoverFromNaN);
524 pc.defineInt("optConst", "Optimize", 0, minimizerDefaults.optConst);
525 pc.defineInt("verbose", "Verbose", 0, minimizerDefaults.verbose);
526 pc.defineInt("doSave", "Save", 0, minimizerDefaults.doSave);
527 pc.defineInt("doTimer", "Timer", 0, minimizerDefaults.doTimer);
528 pc.defineInt("printLevel", "PrintLevel", 0, minimizerDefaults.printLevel);
529 pc.defineInt("strategy", "Strategy", 0, minimizerDefaults.strategy);
530 pc.defineInt("initHesse", "InitialHesse", 0, minimizerDefaults.initHesse);
531 pc.defineInt("hesse", "Hesse", 0, minimizerDefaults.hesse);
532 pc.defineInt("minos", "Minos", 0, minimizerDefaults.minos);
533 pc.defineInt("numee", "PrintEvalErrors", 0, minimizerDefaults.numee);
534 pc.defineInt("doEEWall", "EvalErrorWall", 0, minimizerDefaults.doEEWall);
535 pc.defineInt("doWarn", "Warnings", 0, minimizerDefaults.doWarn);
536 pc.defineInt("doSumW2", "SumW2Error", 0, minimizerDefaults.doSumW2);
537 pc.defineInt("doAsymptoticError", "AsymptoticError", 0, minimizerDefaults.doAsymptotic);
538 pc.defineInt("maxCalls", "MaxCalls", 0, minimizerDefaults.maxCalls);
539 pc.defineInt("doOffset", "OffsetLikelihood", 0, minimizerDefaults.doOffset);
540 pc.defineInt("parallelize", "Parallelize", 0, minimizerDefaults.parallelize); // Three parallelize arguments
541 pc.defineInt("enableParallelGradient", "ParallelGradientOptions", 0, minimizerDefaults.enableParallelGradient);
542 pc.defineInt("enableParallelDescent", "ParallelDescentOptions", 0, minimizerDefaults.enableParallelDescent);
543 pc.defineInt("timingAnalysis", "TimingAnalysis", 0, minimizerDefaults.timingAnalysis);
544 pc.defineString("mintype", "Minimizer", 0, minimizerDefaults.minType.c_str());
545 pc.defineString("minalg", "Minimizer", 1, minimizerDefaults.minAlg.c_str());
546 pc.defineSet("minosSet", "Minos", 0, minimizerDefaults.minosSet);
547}
548
549////////////////////////////////////////////////////////////////////////////////
550/// Minimizes a given NLL variable by finding the optimal parameters with the
551/// RooMinimzer. The NLL variable can be created with RooAbsPdf::createNLL.
552/// If you are looking for a function that combines likelihood creation with
553/// fitting, see RooAbsPdf::fitTo.
554/// \param[in] nll The negative log-likelihood variable to minimize.
555/// \param[in] data The dataset that was also used for the NLL. It's a necessary
556/// parameter because it is used in the asymptotic error correction.
557/// \param[in] cfg Configuration struct with all the configuration options for
558/// the RooMinimizer. These are a subset of the options that you can
559/// also pass to RooAbsPdf::fitTo via the RooFit command arguments.
560std::unique_ptr<RooFitResult> minimize(RooAbsReal &pdf, RooAbsReal &nll, RooAbsData const &data, RooCmdConfig const &pc)
561{
562 MinimizerConfig cfg;
563 cfg.recoverFromNaN = pc.getDouble("RecoverFromUndefinedRegions");
564 cfg.optConst = pc.getInt("optConst");
565 cfg.verbose = pc.getInt("verbose");
566 cfg.doSave = pc.getInt("doSave");
567 cfg.doTimer = pc.getInt("doTimer");
568 cfg.printLevel = pc.getInt("printLevel");
569 cfg.strategy = pc.getInt("strategy");
570 cfg.initHesse = pc.getInt("initHesse");
571 cfg.hesse = pc.getInt("hesse");
572 cfg.minos = pc.getInt("minos");
573 cfg.numee = pc.getInt("numee");
574 cfg.doEEWall = pc.getInt("doEEWall");
575 cfg.doWarn = pc.getInt("doWarn");
576 cfg.doSumW2 = pc.getInt("doSumW2");
577 cfg.doAsymptotic = pc.getInt("doAsymptoticError");
578 cfg.maxCalls = pc.getInt("maxCalls");
579 cfg.minosSet = pc.getSet("minosSet");
580 cfg.minType = pc.getString("mintype", "");
581 cfg.minAlg = pc.getString("minalg", "minuit");
582 cfg.doOffset = pc.getInt("doOffset");
583 cfg.parallelize = pc.getInt("parallelize");
584 cfg.enableParallelGradient = pc.getInt("enableParallelGradient");
585 cfg.enableParallelDescent = pc.getInt("enableParallelDescent");
586 cfg.timingAnalysis = pc.getInt("timingAnalysis");
587
588 // Determine if the dataset has weights
589 bool weightedData = data.isNonPoissonWeighted();
590
591 // The weighted-data uncertainty correction options only apply to likelihood
592 // fits. Skip the NLL-only paths when minimizing a chi-squared test statistic.
593 const bool isChi2 = nll.getAttribute("Chi2EvaluationActive");
594
595 std::string msgPrefix = std::string{"RooAbsPdf::fitTo("} + pdf.GetName() + "): ";
596
597 // Warn user that a method to determine parameter uncertainties should be provided if weighted data is offered
598 if (!isChi2 && weightedData && cfg.doSumW2 == -1 && cfg.doAsymptotic == -1) {
599 oocoutW(&pdf, InputArguments) << msgPrefix <<
600 R"(WARNING: a likelihood fit is requested of what appears to be weighted data.
601 While the estimated values of the parameters will always be calculated taking the weights into account,
602 there are multiple ways to estimate the errors of the parameters. You are advised to make an
603 explicit choice for the error calculation:
604 - Either provide SumW2Error(true), to calculate a sum-of-weights-corrected HESSE error matrix
605 (error will be proportional to the number of events in MC).
606 - Or provide SumW2Error(false), to return errors from original HESSE error matrix
607 (which will be proportional to the sum of the weights, i.e., a dataset with <sum of weights> events).
608 - Or provide AsymptoticError(true), to use the asymptotically correct expression
609 (for details see https://arxiv.org/abs/1911.01303)."
610)";
611 }
612
613 if (cfg.minos && (cfg.doSumW2 == 1 || cfg.doAsymptotic == 1)) {
614 oocoutE(&pdf, InputArguments)
615 << msgPrefix
616 << " sum-of-weights and asymptotic error correction do not work with MINOS errors. Not fitting.\n";
617 return nullptr;
618 }
619 if (cfg.doAsymptotic == 1 && cfg.minos) {
620 oocoutW(&pdf, InputArguments) << msgPrefix << "WARNING: asymptotic correction does not apply to MINOS errors\n";
621 }
622
623 // avoid setting both SumW2 and Asymptotic for uncertainty correction
624 if (cfg.doSumW2 == 1 && cfg.doAsymptotic == 1) {
625 oocoutE(&pdf, InputArguments) << msgPrefix
626 << "ERROR: Cannot compute both asymptotically correct and SumW2 errors.\n";
627 return nullptr;
628 }
629
630 // Instantiate RooMinimizer
632 minimizerConfig.enableParallelGradient = cfg.enableParallelGradient;
633 minimizerConfig.enableParallelDescent = cfg.enableParallelDescent;
634 minimizerConfig.parallelize = cfg.parallelize;
635 minimizerConfig.timingAnalysis = cfg.timingAnalysis;
636 minimizerConfig.offsetting = cfg.doOffset;
638
639 m.setMinimizerType(cfg.minType);
640 m.setEvalErrorWall(cfg.doEEWall);
641 m.setRecoverFromNaNStrength(cfg.recoverFromNaN);
642 m.setPrintEvalErrors(cfg.numee);
643 if (cfg.maxCalls > 0)
644 m.setMaxFunctionCalls(cfg.maxCalls);
645 if (cfg.printLevel != 1)
646 m.setPrintLevel(cfg.printLevel);
647 if (cfg.optConst)
648 m.optimizeConst(cfg.optConst); // Activate constant term optimization
649 if (cfg.verbose)
650 m.setVerbose(true); // Activate verbose options
651 if (cfg.doTimer)
652 m.setProfile(true); // Activate timer options
653 if (cfg.strategy != 1)
654 m.setStrategy(cfg.strategy); // Modify fit strategy
655 if (cfg.initHesse)
656 m.hesse(); // Initialize errors with hesse
657 m.minimize(cfg.minType.c_str(), cfg.minAlg.c_str()); // Minimize using chosen algorithm
658 if (cfg.hesse)
659 m.hesse(); // Evaluate errors with Hesse
660
661 int corrCovQual = -1;
662
663 if (!isChi2 && m.getNPar() > 0) {
664 if (cfg.doAsymptotic == 1)
665 corrCovQual = calcAsymptoticCorrectedCovariance(pdf, m, data); // Asymptotically correct
666 if (cfg.doSumW2 == 1)
668 }
669
670 if (cfg.minos)
671 cfg.minosSet ? m.minos(*cfg.minosSet) : m.minos(); // Evaluate errs with Minos
672
673 // Optionally return fit result
674 std::unique_ptr<RooFitResult> ret;
675 if (cfg.doSave) {
676 auto name = std::string("fitresult_") + pdf.GetName() + "_" + data.GetName();
677 auto title = std::string("Result of fit of p.d.f. ") + pdf.GetName() + " to dataset " + data.GetName();
678 ret = std::unique_ptr<RooFitResult>{m.save(name.c_str(), title.c_str())};
679 if ((cfg.doSumW2 == 1 || cfg.doAsymptotic == 1) && m.getNPar() > 0)
680 ret->setCovQual(corrCovQual);
681 }
682
683 if (cfg.optConst)
684 m.optimizeConst(0);
685 return ret;
686}
687
688std::unique_ptr<RooAbsReal> createNLL(RooAbsPdf &pdf, RooAbsData &data, const RooLinkedList &cmdList)
689{
690 auto timingScope = std::make_unique<ROOT::Math::Util::TimingScope>(
691 [&pdf](std::string const &msg) { oocoutI(&pdf, Fitting) << msg << std::endl; }, "Creation of NLL object took");
692
693 auto baseName = std::string("nll_") + pdf.GetName() + "_" + data.GetName();
694
695 // Select the pdf-specific commands
696 RooCmdConfig pc("RooAbsPdf::createNLL(" + std::string(pdf.GetName()) + ")");
697
698 pc.defineString("rangeName", "RangeWithName", 0, "", true);
699 pc.defineString("addCoefRange", "SumCoefRange", 0, "");
700 pc.defineString("globstag", "GlobalObservablesTag", 0, "");
701 pc.defineString("globssource", "GlobalObservablesSource", 0, "data");
702 pc.defineDouble("rangeLo", "Range", 0, -999.);
703 pc.defineDouble("rangeHi", "Range", 1, -999.);
704 pc.defineInt("splitRange", "SplitRange", 0, 0);
705 pc.defineInt("ext", "Extended", 0, extendedFitDefault);
706 pc.defineInt("numcpu", "NumCPU", 0, 1);
707 pc.defineInt("interleave", "NumCPU", 1, 0);
708 pc.defineInt("verbose", "Verbose", 0, 0);
709 pc.defineInt("optConst", "Optimize", 0, 0);
710 pc.defineInt("cloneData", "CloneData", 0, 2);
711 pc.defineSet("projDepSet", "ProjectedObservables", 0, nullptr);
712 pc.defineSet("cPars", "Constrain", 0, nullptr);
713 pc.defineSet("glObs", "GlobalObservables", 0, nullptr);
714 pc.defineInt("doOffset", "OffsetLikelihood", 0, 0);
715 pc.defineSet("extCons", "ExternalConstraints", 0, nullptr);
716 pc.defineInt("EvalBackend", "EvalBackend", 0, static_cast<int>(RooFit::EvalBackend::defaultValue()));
717 pc.defineDouble("IntegrateBins", "IntegrateBins", 0, -1.);
718 pc.defineMutex("Range", "RangeWithName");
719 pc.defineMutex("GlobalObservables", "GlobalObservablesTag");
720 pc.defineInt("ModularL", "ModularL", 0, 0);
721
722 // New style likelihoods define parallelization through Parallelize(...) on fitTo or attributes on
723 // RooMinimizer::Config.
724 pc.defineMutex("ModularL", "NumCPU");
725
726 // New style likelihoods define offsetting on minimizer, not on likelihood
727 pc.defineMutex("ModularL", "OffsetLikelihood");
728
729 // Process and check varargs
730 pc.process(cmdList);
731 if (!pc.ok(true)) {
732 return nullptr;
733 }
734
735 if (pc.getInt("ModularL")) {
736 int lut[3] = {2, 1, 0};
738 static_cast<RooFit::TestStatistics::RooAbsL::Extended>(lut[pc.getInt("ext")])};
739
743
744 if (auto tmp = pc.getSet("cPars"))
745 cParsSet.add(*tmp);
746
747 if (auto tmp = pc.getSet("extCons"))
748 extConsSet.add(*tmp);
749
750 if (auto tmp = pc.getSet("glObs"))
751 glObsSet.add(*tmp);
752
753 const std::string rangeName = pc.getString("globstag", "", false);
754
756 builder.Extended(ext)
757 .ConstrainedParameters(cParsSet)
758 .ExternalConstraints(extConsSet)
759 .GlobalObservables(glObsSet)
760 .GlobalObservablesTag(rangeName.c_str());
761
762 return std::make_unique<RooFit::TestStatistics::RooRealL>("likelihood", "", builder.build());
763 }
764
765 // Decode command line arguments
766 const char *rangeName = pc.getString("rangeName", nullptr, true);
767 const char *addCoefRangeName = pc.getString("addCoefRange", nullptr, true);
768 const bool ext = interpretExtendedCmdArg(pdf, pc.getInt("ext"));
769
770 int splitRange = pc.getInt("splitRange");
771 int optConst = pc.getInt("optConst");
772 int cloneData = pc.getInt("cloneData");
773 auto offset = static_cast<RooFit::OffsetMode>(pc.getInt("doOffset"));
774
775 // If no explicit cloneData command is specified, cloneData is set to true if optimization is activated
776 if (cloneData == 2) {
778 }
779
780 if (pc.hasProcessed("Range")) {
781 double rangeLo = pc.getDouble("rangeLo");
782 double rangeHi = pc.getDouble("rangeHi");
783
784 // Create range with name 'fit' with above limits on all observables
785 RooArgSet obs;
786 pdf.getObservables(data.get(), obs);
787 for (auto arg : obs) {
788 RooRealVar *rrv = dynamic_cast<RooRealVar *>(arg);
789 if (rrv)
790 rrv->setRange("fit", rangeLo, rangeHi);
791 }
792
793 // Set range name to be fitted to "fit"
794 rangeName = "fit";
795 }
796
797 // Set the fitrange attribute of th PDF, add observables ranges for plotting
799
800 RooArgSet projDeps;
801 auto tmp = pc.getSet("projDepSet");
802 if (tmp) {
803 projDeps.add(*tmp);
804 }
805
806 const std::string globalObservablesSource = pc.getString("globssource", "data", false);
807 if (globalObservablesSource != "data" && globalObservablesSource != "model") {
808 std::string errMsg = "RooAbsPdf::fitTo: GlobalObservablesSource can only be \"data\" or \"model\"!";
809 oocoutE(&pdf, InputArguments) << errMsg << std::endl;
810 throw std::invalid_argument(errMsg);
811 }
813
814 // Lambda function to create the correct constraint term for a PDF. In old
815 // RooFit, we use this PDF itself as the argument, for the new BatchMode
816 // we're passing a clone.
817 auto createConstr = [&]() -> std::unique_ptr<RooAbsReal> {
818 return createConstraintTerm(baseName + "_constr", // name
819 pdf, // pdf
820 data, // data
821 pc.getSet("cPars"), // Constrain RooCmdArg
822 pc.getSet("extCons"), // ExternalConstraints RooCmdArg
823 pc.getSet("glObs"), // GlobalObservables RooCmdArg
824 pc.getString("globstag", nullptr, true), // GlobalObservablesTag RooCmdArg
825 takeGlobalObservablesFromData); // From GlobalObservablesSource RooCmdArg
826 };
827
828 auto evalBackend = static_cast<RooFit::EvalBackend::Value>(pc.getInt("EvalBackend"));
829
830 // Construct BatchModeNLL if requested
832
834 pdf.getObservables(data.get(), normSet);
835
836 if (dynamic_cast<RooSimultaneous const *>(&pdf)) {
837 for (auto i : projDeps) {
838 auto res = normSet.find(i->GetName());
839 if (res != nullptr) {
840 res->setAttribute("__conditional__");
841 }
842 }
843 } else {
844 normSet.remove(projDeps);
845 }
846
847 std::unique_ptr<RooAbsPdf> pdfClone =
848 compilePdfForFit(pdf, normSet, rangeName, splitRange, addCoefRangeName, /*likelihoodMode=*/true);
849
850 if (addCoefRangeName) {
851 oocxcoutI(&pdf, Fitting) << "RooAbsPdf::fitTo(" << pdf.GetName()
852 << ") fixing interpretation of coefficients of any component to range "
853 << addCoefRangeName << "\n";
854 }
855
856 std::unique_ptr<RooAbsReal> compiledConstr;
857 if (std::unique_ptr<RooAbsReal> constr = createConstr()) {
859 compiledConstr->addOwnedComponents(std::move(constr));
860 }
861
862 auto nll = createNLLNew(*pdfClone, data, std::move(compiledConstr), rangeName ? rangeName : "", projDeps, ext,
863 pc.getDouble("IntegrateBins"), offset);
864
865 const double correction = pdfClone->getCorrection();
866
867 if (correction > 0) {
868 oocoutI(&pdf, Fitting) << "[FitHelpers] Detected correction term from RooAbsPdf::getCorrection(). "
869 << "Adding penalty to NLL." << std::endl;
870
871 // Convert the multiplicative correction to an additive term in -log L
872 auto penaltyTerm = std::make_unique<RooConstVar>((baseName + "_Penalty").c_str(),
873 "Penalty term from getCorrection()", correction);
874
875 // add penalty and NLL
876 auto correctedNLL = std::make_unique<RooAddition>((baseName + "_corrected").c_str(), "NLL + penalty",
878
879 // transfer ownership of terms
880 correctedNLL->addOwnedComponents(std::move(nll), std::move(penaltyTerm));
881 nll = std::move(correctedNLL);
882 }
883
884 auto nllWrapper = std::make_unique<RooFit::Experimental::RooEvaluatorWrapper>(
887
888 // We destroy the timing scrope for createNLL prematurely, because we
889 // separately measure the time for jitting and gradient creation
890 // inside the RooFuncWrapper.
891 timingScope.reset();
892
894 nllWrapper->generateGradient();
895 }
897 nllWrapper->setUseGeneratedFunctionCode(true);
898 }
899
900 nllWrapper->addOwnedComponents(std::move(nll));
901 nllWrapper->addOwnedComponents(std::move(pdfClone));
902
903 return nllWrapper;
904 }
905
906 std::unique_ptr<RooAbsReal> nll;
907
908#ifdef ROOFIT_LEGACY_EVAL_BACKEND
909 bool verbose = pc.getInt("verbose");
910
911 int numcpu = pc.getInt("numcpu");
912 int numcpu_strategy = pc.getInt("interleave");
913 // strategy 3 works only for RooSimultaneous.
914 if (numcpu_strategy == 3 && !pdf.InheritsFrom("RooSimultaneous")) {
915 oocoutW(&pdf, Minimization) << "Cannot use a NumCpu Strategy = 3 when the pdf is not a RooSimultaneous, "
916 "falling back to default strategy = 0"
917 << std::endl;
918 numcpu_strategy = 0;
919 }
921
923 RooAbsPdf &actualPdf = binnedLInfo.binnedPdf ? *binnedLInfo.binnedPdf : pdf;
924
925 // Construct NLL
927 RooAbsTestStatistic::Configuration cfg;
928 cfg.addCoefRangeName = addCoefRangeName ? addCoefRangeName : "";
929 cfg.nCPU = numcpu;
930 cfg.interleave = interl;
931 cfg.verbose = verbose;
932 cfg.splitCutRange = static_cast<bool>(splitRange);
933 cfg.cloneInputData = static_cast<bool>(cloneData);
934 cfg.integrateOverBinsPrecision = pc.getDouble("IntegrateBins");
935 cfg.binnedL = binnedLInfo.isBinnedL;
936 cfg.takeGlobalObservablesFromData = takeGlobalObservablesFromData;
937 cfg.rangeName = rangeName ? rangeName : "";
938 auto nllVar = std::make_unique<RooNLLVar>(baseName.c_str(), "-log(likelihood)", actualPdf, data, projDeps, ext, cfg);
939 nllVar->enableBinOffsetting(offset == RooFit::OffsetMode::Bin);
940 nll = std::move(nllVar);
942
943 // Include constraints, if any, in likelihood
944 if (std::unique_ptr<RooAbsReal> constraintTerm = createConstr()) {
945
946 // Even though it is technically only required when the computation graph
947 // is changed because global observables are taken from data, it is safer
948 // to clone the constraint model in general to reset the normalization
949 // integral caches and avoid ASAN build failures (the PDF of the main
950 // measurement is cloned too anyway, so not much overhead). This can be
951 // reconsidered after the caching of normalization sets by pointer is changed
952 // to a more memory-safe solution.
953 constraintTerm = RooHelpers::cloneTreeWithSameParameters(*constraintTerm, data.get());
954
955 // Redirect the global observables to the ones from the dataset if applicable.
956 constraintTerm->setData(data, false);
957
958 // The computation graph for the constraints is very small, no need to do
959 // the tracking of clean and dirty nodes here.
960 constraintTerm->setOperMode(RooAbsArg::ADirty);
961
962 auto orignll = std::move(nll);
963 nll = std::make_unique<RooAddition>((baseName + "_with_constr").c_str(), "nllWithCons",
964 RooArgSet(*orignll, *constraintTerm));
965 nll->addOwnedComponents(std::move(orignll), std::move(constraintTerm));
966 }
967
968 if (optConst) {
969 nll->constOptimizeTestStatistic(RooAbsArg::Activate, optConst > 1);
970 }
971
973 nll->enableOffsetting(true);
974 }
975
976 if (const double correction = pdf.getCorrection(); correction > 0) {
977 oocoutI(&pdf, Fitting) << "[FitHelpers] Detected correction term from RooAbsPdf::getCorrection(). "
978 << "Adding penalty to NLL." << std::endl;
979
980 // Convert the multiplicative correction to an additive term in -log L
981 auto penaltyTerm = std::make_unique<RooConstVar>((baseName + "_Penalty").c_str(),
982 "Penalty term from getCorrection()", correction);
983
984 auto correctedNLL = std::make_unique<RooAddition>(
985 // add penalty and NLL
986 (baseName + "_corrected").c_str(), "NLL + penalty", RooArgSet(*nll, *penaltyTerm));
987
988 // transfer ownership of terms
989 correctedNLL->addOwnedComponents(std::move(nll), std::move(penaltyTerm));
990 nll = std::move(correctedNLL);
991 }
992#else
993 throw std::runtime_error("RooFit was not built with the legacy evaluation backend");
994#endif
995
996 return nll;
997}
998
999namespace {
1000
1001std::unique_ptr<RooAbsReal> createChi2FromDataHist(RooAbsReal &real, RooDataHist &data, const RooLinkedList &cmdList)
1002{
1003 RooCmdConfig pc("createChi2(" + std::string(real.GetName()) + ")");
1004
1005 pc.defineInt("EvalBackend", "EvalBackend", 0, static_cast<int>(RooFit::EvalBackend::defaultValue()));
1006 pc.defineInt("numcpu", "NumCPU", 0, 1);
1007 pc.defineInt("verbose", "Verbose", 0, 0);
1008 pc.defineString("rangeName", "RangeWithName", 0, "", true);
1009 pc.defineDouble("rangeLo", "Range", 0, -999.);
1010 pc.defineDouble("rangeHi", "Range", 1, -999.);
1011 pc.defineMutex("Range", "RangeWithName");
1012 pc.defineInt("etype", "DataError", 0, (Int_t)RooDataHist::Auto);
1013 pc.defineInt("extended", "Extended", 0, extendedFitDefault);
1014 pc.defineInt("splitRange", "SplitRange", 0, 0);
1015 pc.defineDouble("integrate_bins", "IntegrateBins", 0, -1);
1016 pc.defineString("addCoefRange", "SumCoefRange", 0, "");
1017 pc.allowUndefined();
1018
1019 pc.process(cmdList);
1020 if (!pc.ok(true)) {
1021 return nullptr;
1022 }
1023
1024 // Clear possible range attributes from previous fits.
1025 real.removeStringAttribute("fitrange");
1026
1027 std::string baseName = "chi2_" + std::string(real.GetName()) + "_" + data.GetName();
1028
1029 auto evalBackend = static_cast<RooFit::EvalBackend::Value>(pc.getInt("EvalBackend"));
1030
1031 RooDataHist::ErrorType etype = static_cast<RooDataHist::ErrorType>(pc.getInt("etype"));
1032 // Resolve Auto to a concrete mode so it's consistent across backends.
1033 if (etype == RooDataHist::Auto) {
1034 etype = data.isNonPoissonWeighted() ? RooDataHist::SumW2 : RooDataHist::Expected;
1035 }
1036
1037 auto *pdf = dynamic_cast<RooAbsPdf *>(&real);
1038 const char *rangeName = pc.getString("rangeName", nullptr, true);
1039
1040 // Translate Range(lo, hi) into a "fit" named range on all observables.
1041 if (pc.hasProcessed("Range")) {
1042 const double rangeLo = pc.getDouble("rangeLo");
1043 const double rangeHi = pc.getDouble("rangeHi");
1044 RooArgSet obs;
1045 real.getObservables(data.get(), obs);
1046 for (auto arg : obs) {
1047 if (auto *rrv = dynamic_cast<RooRealVar *>(arg)) {
1048 rrv->setRange("fit", rangeLo, rangeHi);
1049 }
1050 }
1051 rangeName = "fit";
1052 }
1053
1056
1057 const int splitRange = pc.getInt("splitRange");
1059
1060 std::unique_ptr<RooFit::Experimental::RooEvaluatorWrapper> wrapper;
1061
1062 // Function mode: the input is a non-pdf RooAbsReal. We can short-circuit
1063 // the pdf-compilation pipeline since there's no real pdf to normalize.
1064 if (!pdf) {
1065 RooArgSet observables;
1066 real.getObservables(data.get(), observables);
1067 RooNLLVarNew::Config cfg;
1068 cfg.statistic = RooNLLVarNew::Statistic::Chi2;
1069 cfg.chi2ErrorType = etype;
1070 auto chi2 = std::make_unique<RooNLLVarNew>(baseName.c_str(), baseName.c_str(), real, observables, cfg);
1071 wrapper = std::make_unique<RooFit::Experimental::RooEvaluatorWrapper>(
1073 /*simPdf=*/nullptr,
1074 /*takeGlobalObservablesFromData=*/true);
1075 wrapper->addOwnedComponents(std::move(chi2));
1076 } else {
1077 const bool extended = interpretExtendedCmdArg(*pdf, pc.getInt("extended"));
1078
1080 pdf->getObservables(data.get(), normSet);
1081
1082 oocxcoutI(pdf, Fitting) << "createChi2(" << pdf->GetName()
1083 << ") fixing normalization set for coefficient determination to observables in data\n";
1084 pdf->fixAddCoefNormalization(normSet, false);
1085
1086 std::unique_ptr<RooAbsPdf> pdfClone =
1087 compilePdfForFit(*pdf, normSet, rangeName, splitRange, pc.getString("addCoefRange", nullptr, true),
1088 /*likelihoodMode=*/false);
1089
1093
1094 std::unique_ptr<RooAbsReal> chi2;
1095 if (auto *simPdfClone = dynamic_cast<RooSimultaneous *>(&finalPdf)) {
1096 chi2 = std::unique_ptr<RooAbsReal>{dynamic_cast<RooAbsReal *>(
1097 createSimultaneousChi2(*simPdfClone, rangeName ? rangeName : "", extended, etype).release())};
1098 } else {
1099 RooArgSet observables;
1100 finalPdf.getObservables(data.get(), observables);
1101 RooNLLVarNew::Config cfg;
1102 cfg.statistic = RooNLLVarNew::Statistic::Chi2;
1103 cfg.extended = extended;
1104 cfg.chi2ErrorType = etype;
1105 chi2 = std::make_unique<RooNLLVarNew>(baseName.c_str(), baseName.c_str(), finalPdf, observables, cfg);
1106 }
1107
1108 wrapper = std::make_unique<RooFit::Experimental::RooEvaluatorWrapper>(
1110 /*takeGlobalObservablesFromData=*/true);
1111 wrapper->addOwnedComponents(std::move(binSamplingPdfs));
1112 wrapper->addOwnedComponents(std::move(chi2));
1113 wrapper->addOwnedComponents(std::move(pdfClone));
1114 }
1115
1117 wrapper->generateGradient();
1118 }
1120 wrapper->setUseGeneratedFunctionCode(true);
1121 }
1122
1124 return wrapper;
1125 }
1126
1127#ifdef ROOFIT_LEGACY_EVAL_BACKEND
1128 RooAbsTestStatistic::Configuration cfg;
1129
1131
1132 bool extended = false;
1133 if (pdf) {
1134 extended = interpretExtendedCmdArg(*pdf, pc.getInt("extended"));
1135 }
1136
1137 const char *addCoefRangeName = pc.getString("addCoefRange", nullptr, true);
1138 int splitRange = pc.getInt("splitRange");
1139
1140 // Set the fitrange attribute of th PDF, add observables ranges for plotting
1142
1143 cfg.rangeName = rangeName ? rangeName : "";
1144 cfg.nCPU = pc.getInt("numcpu");
1145 cfg.interleave = RooFit::Interleave;
1146 cfg.verbose = static_cast<bool>(pc.getInt("verbose"));
1147 cfg.cloneInputData = false;
1148 cfg.integrateOverBinsPrecision = pc.getDouble("integrate_bins");
1149 cfg.addCoefRangeName = addCoefRangeName ? addCoefRangeName : "";
1150 cfg.splitCutRange = static_cast<bool>(splitRange);
1151 auto chi2 = std::make_unique<RooChi2Var>(baseName.c_str(), baseName.c_str(), real, static_cast<RooDataHist &>(data),
1152 extended, etype, cfg);
1153
1155
1156 return chi2;
1157#else
1158 throw std::runtime_error("createChi2() is not supported without the legacy evaluation backend");
1159 return nullptr;
1160#endif
1161}
1162
1163} // namespace
1164
1165std::unique_ptr<RooAbsReal> createChi2(RooAbsReal &real, RooAbsData &data, const RooLinkedList &cmdList)
1166{
1167 if (auto dataHist = dynamic_cast<RooDataHist *>(&data))
1168 return createChi2FromDataHist(real, *dataHist, cmdList);
1169
1170#ifdef ROOFIT_LEGACY_EVAL_BACKEND
1171
1172 RooCmdConfig pc("createChi2(" + std::string(real.GetName()) + ")");
1173
1174 pc.defineInt("numcpu", "NumCPU", 0, 1);
1175 pc.defineInt("verbose", "Verbose", 0, 0);
1176 pc.defineString("rangeName", "RangeWithName", 0, "", true);
1177
1178 RooAbsTestStatistic::Configuration cfg;
1179
1180 pc.defineInt("integrate", "Integrate", 0, 0);
1181 pc.defineObject("yvar", "YVar", 0, nullptr);
1182 pc.defineInt("interleave", "NumCPU", 1, 0);
1183
1184 // Process and check varargs
1185 pc.process(cmdList);
1186 if (!pc.ok(true)) {
1187 return nullptr;
1188 }
1189
1190 // Decode command line arguments
1191 bool integrate = pc.getInt("integrate");
1192 RooRealVar *yvar = static_cast<RooRealVar *>(pc.getObject("yvar"));
1193 const char *rangeName = pc.getString("rangeName", nullptr, true);
1194 Int_t numcpu = pc.getInt("numcpu");
1195 Int_t numcpu_strategy = pc.getInt("interleave");
1196 // strategy 3 works only for RooSimultaneous.
1197 if (numcpu_strategy == 3 && !real.InheritsFrom("RooSimultaneous")) {
1198 oocoutW(&real, Minimization) << "Cannot use a NumCpu Strategy = 3 when the pdf is not a RooSimultaneous, "
1199 "falling back to default strategy = 0"
1200 << std::endl;
1201 numcpu_strategy = 0;
1202 }
1204 bool verbose = pc.getInt("verbose");
1205
1206 cfg.rangeName = rangeName ? rangeName : "";
1207 cfg.nCPU = numcpu;
1208 cfg.interleave = interl;
1209 cfg.verbose = verbose;
1210 cfg.verbose = false;
1211
1212 std::string name = "chi2_" + std::string(real.GetName()) + "_" + data.GetName();
1213
1214 return std::make_unique<RooXYChi2Var>(name.c_str(), name.c_str(), real, static_cast<RooDataSet &>(data), yvar,
1215 integrate, cfg);
1216#else
1217 throw std::runtime_error("createChi2() is not supported without the legacy evaluation backend");
1218 return nullptr;
1219#endif
1220}
1221
1222std::unique_ptr<RooFitResult> fitTo(RooAbsReal &real, RooAbsData &data, const RooLinkedList &cmdList, bool chi2)
1223{
1224 const bool isDataHist = dynamic_cast<RooDataHist const *>(&data);
1225
1226 RooCmdConfig pc("fitTo(" + std::string(real.GetName()) + ")");
1227
1229 std::string nllCmdListString;
1230 if (!chi2) {
1231 nllCmdListString = "ProjectedObservables,Extended,Range,"
1232 "RangeWithName,SumCoefRange,NumCPU,SplitRange,Constrained,Constrain,ExternalConstraints,"
1233 "CloneData,GlobalObservables,GlobalObservablesSource,GlobalObservablesTag,"
1234 "EvalBackend,IntegrateBins,ModularL";
1235
1236 if (!cmdList.FindObject("ModularL") || static_cast<RooCmdArg *>(cmdList.FindObject("ModularL"))->getInt(0) == 0) {
1237 nllCmdListString += ",OffsetLikelihood";
1238 }
1239 } else {
1240 auto createChi2DataHistCmdArgs = "Range,RangeWithName,NumCPU,Optimize,IntegrateBins,ProjectedObservables,"
1241 "AddCoefRange,SplitRange,DataError,Extended,EvalBackend";
1242 auto createChi2DataSetCmdArgs = "YVar,Integrate,RangeWithName,NumCPU,Verbose";
1244 }
1245
1247
1248 pc.defineDouble("prefit", "Prefit", 0, 0);
1250
1251 // Process and check varargs
1252 pc.process(fitCmdList);
1253 if (!pc.ok(true)) {
1254 return nullptr;
1255 }
1256
1257 // TimingAnalysis works only for RooSimultaneous.
1258 if (pc.getInt("timingAnalysis") && !real.InheritsFrom("RooSimultaneous")) {
1259 oocoutW(&real, Minimization) << "The timingAnalysis feature was built for minimization with RooSimultaneous "
1260 "and is not implemented for other PDF's. Please create a RooSimultaneous to "
1261 "enable this feature."
1262 << std::endl;
1263 }
1264
1265 // Decode command line arguments
1266 double prefit = pc.getDouble("prefit");
1267
1268 if (prefit != 0) {
1269 size_t nEvents = static_cast<size_t>(prefit * data.numEntries());
1270 if (prefit > 0.5 || nEvents < 100) {
1271 oocoutW(&real, InputArguments) << "PrefitDataFraction should be in suitable range."
1272 << "With the current PrefitDataFraction=" << prefit
1273 << ", the number of events would be " << nEvents << " out of "
1274 << data.numEntries() << ". Skipping prefit..." << std::endl;
1275 } else {
1276 size_t step = data.numEntries() / nEvents;
1277
1278 RooDataSet tiny("tiny", "tiny", *data.get(), data.isWeighted() ? RooFit::WeightVar() : RooCmdArg());
1279
1280 for (int i = 0; i < data.numEntries(); i += step) {
1281 const RooArgSet *event = data.get(i);
1282 tiny.add(*event, data.weight());
1283 }
1285 pc.filterCmdList(tinyCmdList, "Prefit,Hesse,Minos,Verbose,Save,Timer");
1288
1291
1292 fitTo(real, tiny, tinyCmdList, chi2);
1293 }
1294 }
1295
1297 if (pc.getInt("parallelize") != 0 || pc.getInt("enableParallelGradient") || pc.getInt("enableParallelDescent")) {
1298 // Set to new style likelihood if parallelization is requested
1301 }
1302
1303 std::unique_ptr<RooAbsReal> nll;
1304 if (chi2) {
1305 nll = std::unique_ptr<RooAbsReal>{isDataHist ? real.createChi2(static_cast<RooDataHist &>(data), nllCmdList)
1306 : real.createChi2(static_cast<RooDataSet &>(data), nllCmdList)};
1307 } else {
1308 nll = std::unique_ptr<RooAbsReal>{dynamic_cast<RooAbsPdf &>(real).createNLL(data, nllCmdList)};
1309 }
1310
1311 return RooFit::FitHelpers::minimize(real, *nll, data, pc);
1312}
1313
1314} // namespace RooFit::FitHelpers
1315
1316/// \endcond
header file containing the templated implementation of matrix inversion routines for use with ROOT's ...
ROOT::RRangeCast< T, false, Range_t > static_range_cast(Range_t &&coll)
#define oocoutW(o, a)
#define oocoutE(o, a)
#define oocoutI(o, a)
#define oocxcoutI(o, a)
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h offset
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 child
char name[80]
Definition TGX11.cxx:145
Binding & operator=(OUT(*fun)(void))
class to compute the Cholesky decomposition of a matrix
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:76
void setStringAttribute(const Text_t *key, const Text_t *value)
Associate string 'value' to this object under key 'key'.
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
void removeStringAttribute(const Text_t *key)
Delete a string attribute with a given key.
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...
virtual bool addOwned(RooAbsArg &var, bool silent=false)
Add an argument and transfer the ownership to the collection.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:56
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const override
void setNormRange(const char *rangeName)
virtual double getCorrection() const
This function returns the penalty term.
@ CanBeExtended
Definition RooAbsPdf.h:208
@ MustBeExtended
Definition RooAbsPdf.h:208
@ CanNotBeExtended
Definition RooAbsPdf.h:208
const char * normRange() const
Definition RooAbsPdf.h:246
virtual ExtendMode extendMode() const
Returns ability of PDF to provide extended likelihood terms.
Definition RooAbsPdf.h:212
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
virtual void fixAddCoefNormalization(const RooArgSet &addNormSet=RooArgSet(), bool force=true)
Fix the interpretation of the coefficient of any RooAddPdf component in the expression tree headed by...
static void setEvalErrorLoggingMode(ErrorLoggingMode m)
Set evaluation error logging mode.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
RooArgSet * selectByAttrib(const char *name, bool value) const
Use RooAbsCollection::selectByAttrib(), but return as RooArgSet.
Definition RooArgSet.h:144
static std::unique_ptr< RooAbsPdf > create(RooAbsPdf &pdf, RooAbsData const &data, double precision)
Creates a wrapping RooBinSamplingPdf if appropriate.
Object to represent discrete states.
Definition RooCategory.h:28
Named container for two doubles, two integers two object points and three string pointers that can be...
Definition RooCmdArg.h:26
Int_t getInt(Int_t idx) const
Definition RooCmdArg.h:87
Configurable parser for RooCmdArg named arguments.
void defineMutex(const char *head, Args_t &&... tail)
Define arguments where any pair is mutually exclusive.
bool process(const RooCmdArg &arg)
Process given RooCmdArg.
bool hasProcessed(const char *cmdName) const
Return true if RooCmdArg with name 'cmdName' has been processed.
double getDouble(const char *name, double defaultValue=0.0) const
Return double property registered with name 'name'.
bool defineDouble(const char *name, const char *argName, int doubleNum, double defValue=0.0)
Define double property name 'name' mapped to double in slot 'doubleNum' in RooCmdArg with name argNam...
RooArgSet * getSet(const char *name, RooArgSet *set=nullptr) const
Return RooArgSet property registered with name 'name'.
bool defineSet(const char *name, const char *argName, int setNum, const RooArgSet *set=nullptr)
Define TObject property name 'name' mapped to object in slot 'setNum' in RooCmdArg with name argName ...
bool ok(bool verbose) const
Return true of parsing was successful.
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:32
static Value & defaultValue()
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
Collection class for internal use, storing a collection of RooAbsArg pointers in a doubly linked list...
Wrapper class around ROOT::Math::Minimizer that provides a seamless interface between the minimizer f...
RooFit::OwningPtr< RooFitResult > save(const char *name=nullptr, const char *title=nullptr)
Save and return a RooFitResult snapshot of current minimizer status.
int hesse()
Execute HESSE.
void applyCovarianceMatrix(TMatrixDSym const &V)
Apply results of given external covariance matrix.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
static TClass * Class()
void setRange(const char *name, double min, double max, bool shared=true)
Set a fit or plotting range.
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:546
RooCmdArg WeightVar(const char *name="weight", bool reinterpretAsWeight=false)
RooCmdArg Hesse(bool flag=true)
RooCmdArg ModularL(bool flag=false)
RooCmdArg PrintLevel(Int_t code)
RVec< PromoteType< T > > log(const RVec< T > &v)
Definition RVec.hxx:1836
CoordSystem::Scalar get(DisplacementVector2D< CoordSystem, Tag > const &p)
std::vector< std::string > Split(std::string_view str, std::string_view delims, bool skipEmpty=false)
Splits a string at each character in delims.
double nll(double pdf, double weight, int binnedL, int doBinOffset)
Definition MathFuncs.h:452
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