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