Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooAbsPdf.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * @(#)root/roofitcore:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17//////////////////////////////////////////////////////////////////////////////
18/** \class RooAbsPdf
19 \ingroup Roofitcore
20 \brief Abstract interface for all probability density functions.
21
22## RooAbsPdf, the base class of all PDFs
23
24RooAbsPdf is the base class for all probability density
25functions (PDFs). The class provides hybrid analytical/numerical
26normalization for its implementations, error tracing, and a Monte Carlo
27generator interface.
28
29### A Minimal PDF Implementation
30
31A minimal implementation of a PDF class derived from RooAbsPdf
32should override the `evaluate()` function. This function should
33return the PDF's value (which does not need to be normalised).
34
35
36#### Normalization/Integration
37
38Although the normalization of a PDF is an integral part of a
39probability density function, normalization is treated separately
40in RooAbsPdf. The reason is that a RooAbsPdf object is more than a
41PDF: it can be a building block for a more complex composite PDF
42if any of its variables are functions instead of variables. In
43such cases, the normalization of the composite PDF may not simply be
44integral over the dependents of the top-level PDF: these are
45functions with potentially non-trivial Jacobian terms themselves.
46\note Therefore, no explicit attempt should be made to normalize the
47function output in evaluate(). In particular, normalisation constants
48can be omitted to speed up the function evaluations, and included later
49in the integration of the PDF (see below), which is rarely called in
50comparison to the `evaluate()` function.
51
52In addition, RooAbsPdf objects do not have a static concept of what
53variables are parameters, and what variables are dependents (which
54need to be integrated over for a correct PDF normalization).
55Instead, the choice of normalization is always specified each time a
56normalized value is requested from the PDF via the getVal()
57method.
58
59RooAbsPdf manages the entire normalization logic of each PDF with
60the help of a RooRealIntegral object, which coordinates the integration
61of a given choice of normalization. By default, RooRealIntegral will
62perform an entirely numeric integration of all dependents. However,
63PDFs can advertise one or more (partial) analytical integrals of
64their function, and these will be used by RooRealIntegral, if it
65determines that this is safe (i.e., no hidden Jacobian terms,
66multiplication with other PDFs that have one or more dependents in
67common, etc).
68
69#### Implementing analytical integrals
70To implement analytical integrals, two functions must be implemented. First,
71
72```
73Int_t getAnalyticalIntegral(const RooArgSet& integSet, RooArgSet& anaIntSet)
74```
75should return the analytical integrals that are supported. `integSet`
76is the set of dependents for which integration is requested. The
77function should copy the subset of dependents it can analytically
78integrate to `anaIntSet`, and return a unique identification code for
79this integration configuration. If no integration can be
80performed, zero should be returned. Second,
81
82```
83double analyticalIntegral(Int_t code)
84```
85
86implements the actual analytical integral(s) advertised by
87`getAnalyticalIntegral()`. This function will only be called with
88codes returned by `getAnalyticalIntegral()`, except code zero.
89
90The integration range for each dependent to be integrated can
91be obtained from the dependent's proxy functions `min()` and
92`max()`. Never call these proxy functions for any proxy not known to
93be a dependent via the integration code. Doing so may be
94ill-defined, e.g., in case the proxy holds a function, and will
95trigger an assert. Integrated category dependents should always be
96summed over all of their states.
97
98
99
100### Direct generation of observables
101
102Distributions for any PDF can be generated with the accept/reject method,
103but for certain PDFs, more efficient methods may be implemented. To
104implement direct generation of one or more observables, two
105functions need to be implemented, similar to those for analytical
106integrals:
107
108```
109Int_t getGenerator(const RooArgSet& generateVars, RooArgSet& directVars)
110```
111and
112```
113void generateEvent(Int_t code)
114```
115
116The first function advertises observables, for which distributions can be generated,
117similar to the way analytical integrals are advertised. The second
118function implements the actual generator for the advertised observables.
119
120The generated dependent values should be stored in the proxy
121objects. For this, the assignment operator can be used (i.e. `xProxy
122= 3.0` ). Never call assign to any proxy not known to be a dependent
123via the generation code. Doing so may be ill-defined, e.g. in case
124the proxy holds a function, and will trigger an assert.
126
127### Batched function evaluations (Advanced usage)
128
129To speed up computations with large numbers of data events in unbinned fits,
130it is beneficial to override `doEval()`. Like this, large spans of
131computations can be done, without having to call `evaluate()` for each single data event.
132`doEval()` should execute the same computation as `evaluate()`, but it
133may choose an implementation that is capable of SIMD computations.
134If doEval is not implemented, the classic and slower `evaluate()` will be
135called for each data event.
136*/
137
138#include "RooAbsPdf.h"
139
140#include "FitHelpers.h"
142#include "RooMsgService.h"
143#include "RooArgSet.h"
144#include "RooArgProxy.h"
145#include "RooRealProxy.h"
146#include "RooRealVar.h"
147#include "RooGenContext.h"
148#include "RooBinnedGenContext.h"
149#include "RooPlot.h"
150#include "RooCurve.h"
151#include "RooCategory.h"
152#include "RooNameReg.h"
153#include "RooCmdConfig.h"
154#include "RooGlobalFunc.h"
155#include "RooRandom.h"
156#include "RooNumIntConfig.h"
157#include "RooProjectedPdf.h"
158#include "RooParamBinning.h"
159#include "RooNumCdf.h"
160#include "RooFitResult.h"
161#include "RooNumGenConfig.h"
162#include "RooCachedReal.h"
163#include "RooRealIntegral.h"
164#include "RooWorkspace.h"
165#include "RooNaNPacker.h"
166#include "RooFitImplHelpers.h"
167#include "RooHelpers.h"
168#include "RooFormulaVar.h"
169#include "RooDerivative.h"
170
171#include "ROOT/StringUtils.hxx"
172#include "TMath.h"
173#include "TPaveText.h"
174#include "TMatrixD.h"
175#include "TMatrixDSym.h"
176
177#include <algorithm>
178#include <iostream>
179#include <string>
180#include <cmath>
181#include <stdexcept>
182
183namespace {
184
185inline double getLog(double prob, RooAbsReal const *caller)
186{
187
188 if (prob < 0) {
189 caller->logEvalError("getLogVal() top-level p.d.f evaluates to a negative number");
191 }
192
193 if (std::isinf(prob)) {
194 oocoutW(caller, Eval) << "RooAbsPdf::getLogVal(" << caller->GetName()
195 << ") WARNING: top-level pdf has an infinite value" << std::endl;
196 }
197
198 if (prob == 0) {
199 caller->logEvalError("getLogVal() top-level p.d.f evaluates to zero");
200
201 return -std::numeric_limits<double>::infinity();
202 }
203
204 if (TMath::IsNaN(prob)) {
205 caller->logEvalError("getLogVal() top-level p.d.f evaluates to NaN");
206
207 return prob;
208 }
209
210 return std::log(prob);
211}
212
214{
215 TObject *old = lst.FindObject(obj.GetName());
216 if (old)
217 lst.Replace(old, &obj);
218 else
219 lst.Add(&obj);
220}
221
222} // namespace
223
224using std::endl, std::string, std::ostream, std::vector, std::pair, std::make_pair;
225
227
228
229
230
233
234////////////////////////////////////////////////////////////////////////////////
235/// Default constructor
236
237RooAbsPdf::RooAbsPdf() : _normMgr(this, 10) {}
238
239////////////////////////////////////////////////////////////////////////////////
240/// Constructor with name and title only
241
242RooAbsPdf::RooAbsPdf(const char *name, const char *title) :
243 RooAbsReal(name,title), _normMgr(this,10), _selectComp(true)
244{
246 setTraceCounter(0) ;
247}
248
249
250
251////////////////////////////////////////////////////////////////////////////////
252/// Constructor with name, title, and plot range
253
254RooAbsPdf::RooAbsPdf(const char *name, const char *title,
255 double plotMin, double plotMax) :
256 RooAbsReal(name,title,plotMin,plotMax), _normMgr(this,10), _selectComp(true)
257{
259 setTraceCounter(0) ;
260}
261
262
263
264////////////////////////////////////////////////////////////////////////////////
265/// Copy constructor
266
269 _normMgr(other._normMgr,this), _selectComp(other._selectComp), _normRange(other._normRange)
270{
272 setTraceCounter(other._traceCount) ;
273
274 if (other._specGeneratorConfig) {
275 _specGeneratorConfig = std::make_unique<RooNumGenConfig>(*other._specGeneratorConfig);
276 }
277}
278
279
280
281////////////////////////////////////////////////////////////////////////////////
282/// Destructor
283
287
288
290
291 if (normVal < 0. || (normVal == 0. && rawVal != 0)) {
292 //Unreasonable normalisations. A zero integral can be tolerated if the function vanishes, though.
293 const std::string msg = "p.d.f normalization integral is zero or negative: " + std::to_string(normVal);
294 logEvalError(msg.c_str());
296 return RooNaNPacker::packFloatIntoNaN(-normVal + (rawVal < 0. ? -rawVal : 0.));
297 }
298
299 if (rawVal < 0.) {
300 std::stringstream ss;
301 ss << "p.d.f value is less than zero (" << rawVal << "), trying to recover";
302 logEvalError(ss.str().c_str());
305 }
306
307 if (TMath::IsNaN(rawVal)) {
308 logEvalError("p.d.f value is Not-a-Number");
310 return rawVal;
311 }
312
313 return (rawVal == 0. && normVal == 0.) ? 0. : rawVal / normVal;
314}
315
316
317////////////////////////////////////////////////////////////////////////////////
318/// Return current value, normalized by integrating over
319/// the observables in `nset`. If `nset` is 0, the unnormalized value
320/// is returned. All elements of `nset` must be lvalues.
321///
322/// Unnormalized values are not cached.
323/// Doing so would be complicated as `_norm->getVal()` could
324/// spoil the cache and interfere with returning the cached
325/// return value. Since unnormalized calls are typically
326/// done in integration calls, there is no performance hit.
327
328double RooAbsPdf::getValV(const RooArgSet* nset) const
329{
330
331 // Special handling of case without normalization set (used in numeric integration of pdfs)
332 if (!nset) {
333 RooArgSet const* tmp = _normSet ;
334 _normSet = nullptr ;
335 double val = evaluate() ;
336 _normSet = tmp ;
337
338 return TMath::IsNaN(val) ? 0. : val;
339 }
340
341
342 // Process change in last data set used
343 bool nintChanged(false) ;
344 if (!isActiveNormSet(nset) || _norm==nullptr) {
346 }
347
348 // Return value of object. Calculated if dirty, otherwise cached value is returned.
350
351 // Evaluate numerator
352 const double rawVal = evaluate();
353
354 // Evaluate denominator
355 const double normVal = _norm->getVal();
356
358
360 }
361
362 return _value ;
363}
364
365
366////////////////////////////////////////////////////////////////////////////////
367/// Analytical integral with normalization (see RooAbsReal::analyticalIntegralWN() for further information).
368///
369/// This function applies the normalization specified by `normSet` to the integral returned
370/// by RooAbsReal::analyticalIntegral(). The passthrough scenario (code=0) is also changed
371/// to return a normalized answer.
372
373double RooAbsPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const
374{
375 cxcoutD(Eval) << "RooAbsPdf::analyticalIntegralWN(" << GetName() << ") code = " << code << " normset = " << (normSet?*normSet:RooArgSet()) << std::endl ;
376
377
378 if (code==0) return getVal(normSet) ;
379 if (normSet) {
381 } else {
382 return analyticalIntegral(code,rangeName) ;
383 }
384}
385
386
387
388////////////////////////////////////////////////////////////////////////////////
389/// Check that passed value is positive and not 'not-a-number'. If
390/// not, print an error, until the error counter reaches its set
391/// maximum.
392
394{
395 // check for a math error or negative value
396 bool error(false) ;
397 if (TMath::IsNaN(value)) {
398 logEvalError(Form("p.d.f value is Not-a-Number (%f), forcing value to zero",value)) ;
399 error=true ;
400 }
401 if (value<0) {
402 logEvalError(Form("p.d.f value is less than zero (%f), forcing value to zero",value)) ;
403 error=true ;
404 }
405
406 // do nothing if we are no longer tracing evaluations and there was no error
407 if(!error) return error ;
408
409 // otherwise, print out this evaluations input values and result
410 if(++_errorCount <= 10) {
411 cxcoutD(Tracing) << "*** Evaluation Error " << _errorCount << " ";
412 if(_errorCount == 10) cxcoutD(Tracing) << "(no more will be printed) ";
413 }
414 else {
415 return error ;
416 }
417
418 Print() ;
419 return error ;
420}
421
422
423////////////////////////////////////////////////////////////////////////////////
424/// Get normalisation term needed to normalise the raw values returned by
425/// getVal(). Note that `getVal(normalisationVariables)` will automatically
426/// apply the normalisation term returned here.
427/// \param nset Set of variables to normalise over.
428double RooAbsPdf::getNorm(const RooArgSet* nset) const
429{
430 if (!nset) return 1 ;
431
432 syncNormalization(nset,true) ;
433 if (_verboseEval>1) cxcoutD(Tracing) << ClassName() << "::getNorm(" << GetName() << "): norm(" << _norm << ") = " << _norm->getVal() << std::endl ;
434
435 double ret = _norm->getVal() ;
436 if (ret==0.) {
437 if(++_errorCount <= 10) {
438 coutW(Eval) << "RooAbsPdf::getNorm(" << GetName() << ":: WARNING normalization is zero, nset = " ; nset->Print("1") ;
439 if(_errorCount == 10) coutW(Eval) << "RooAbsPdf::getNorm(" << GetName() << ") INFO: no more messages will be printed " << std::endl ;
440 }
441 }
442
443 return ret ;
444}
445
446
447
448////////////////////////////////////////////////////////////////////////////////
449/// Return pointer to RooAbsReal object that implements calculation of integral over observables iset in range
450/// rangeName, optionally taking the integrand normalized over observables nset
451
452const RooAbsReal* RooAbsPdf::getNormObj(const RooArgSet* nset, const RooArgSet* iset, const TNamed* rangeName) const
453{
454 // Check normalization is already stored
455 CacheElem* cache = static_cast<CacheElem*>(_normMgr.getObj(nset,iset,nullptr,rangeName)) ;
456 if (cache) {
457 return cache->_norm.get();
458 }
459
460 // If not create it now
463
464 // Normalization is always over all pdf components. Overriding the global
465 // component selection temporarily makes all RooRealIntegrals created during
466 // that time always include all components.
468 RooAbsReal* norm = std::unique_ptr<RooAbsReal>{createIntegral(depList,*nset, *getIntegratorConfig(), RooNameReg::str(rangeName))}.release();
469
470 // Store it in the cache
472
473 // And return the newly created integral
474 return norm ;
475}
476
477
478
479////////////////////////////////////////////////////////////////////////////////
480/// Verify that the normalization integral cached with this PDF
481/// is valid for given set of normalization observables.
482///
483/// If not, the cached normalization integral (if any) is deleted
484/// and a new integral is constructed for use with 'nset'.
485/// Elements in 'nset' can be discrete and real, but must be lvalues.
486///
487/// For functions that declare to be self-normalized by overloading the
488/// selfNormalized() function, a unit normalization is always constructed.
489
491{
492 setActiveNormSet(nset);
493
494 // Check if data sets are identical
495 CacheElem* cache = static_cast<CacheElem*>(_normMgr.getObj(nset)) ;
496 if (cache) {
497
498 bool nintChanged = (_norm!=cache->_norm.get()) ;
499 _norm = cache->_norm.get();
500
501 // In the past, this condition read `if (nintChanged && adjustProxies)`.
502 // However, the cache checks if the nset was already cached **by content**,
503 // and not by RooArgSet instance! So it can happen that the normalization
504 // set object is different, but the integral object is the same, in which
505 // case it would be wrong to not adjust the proxies. They always have to be
506 // adjusted when the nset changed, which is always the case when
507 // `syncNormalization()` is called.
508 if (adjustProxies) {
509 // Update dataset pointers of proxies
510 const_cast<RooAbsPdf*>(this)->setProxyNormSet(nset) ;
511 }
512
513 return nintChanged ;
514 }
515
516 // Update dataset pointers of proxies
517 if (adjustProxies) {
518 const_cast<RooAbsPdf*>(this)->setProxyNormSet(nset) ;
519 }
520
522 getObservables(nset, depList);
523
524 if (_verboseEval>0) {
525 if (!selfNormalized()) {
526 cxcoutD(Tracing) << ClassName() << "::syncNormalization(" << GetName()
527 << ") recreating normalization integral " << std::endl ;
528 depList.printStream(ccoutD(Tracing),kName|kValue|kArgs,kSingleLine) ;
529 } else {
530 cxcoutD(Tracing) << ClassName() << "::syncNormalization(" << GetName() << ") selfNormalized, creating unit norm" << std::endl;
531 }
532 }
533
534 // Destroy old normalization & create new
535 if (selfNormalized() || !dependsOn(depList)) {
536 auto ntitle = std::string(GetTitle()) + " Unit Normalization";
537 auto nname = std::string(GetName()) + "_UnitNorm";
538 _norm = new RooRealVar(nname.c_str(),ntitle.c_str(),1) ;
539 } else {
540 const char* nr = (_normRangeOverride.Length()>0 ? _normRangeOverride.Data() : (_normRange.Length()>0 ? _normRange.Data() : nullptr)) ;
541
542// std::cout << "RooAbsPdf::syncNormalization(" << GetName() << ") rangeName for normalization is " << (nr?nr:"<null>") << std::endl ;
544 {
545 // Normalization is always over all pdf components. Overriding the global
546 // component selection temporarily makes all RooRealIntegrals created during
547 // that time always include all components.
549 normInt = std::unique_ptr<RooAbsReal>{createIntegral(depList,*getIntegratorConfig(),nr)}.release();
550 }
551 static_cast<RooRealIntegral*>(normInt)->setAllowComponentSelection(false);
552 normInt->getVal() ;
553// std::cout << "resulting normInt = " << normInt->GetName() << std::endl ;
554
555 const char* cacheParamsStr = getStringAttribute("CACHEPARAMINT") ;
557
558 std::unique_ptr<RooArgSet> intParams{normInt->getVariables()} ;
559
561
562 if (!cacheParams.empty()) {
563 cxcoutD(Caching) << "RooAbsReal::createIntObj(" << GetName() << ") INFO: constructing " << cacheParams.size()
564 << "-dim value cache for integral over " << depList << " as a function of " << cacheParams << " in range " << (nr?nr:"<default>") << std::endl ;
565 std::string name = normInt->GetName() + ("_CACHE_[" + cacheParams.contentsString()) + "]";
567 cachedIntegral->setInterpolationOrder(2) ;
568 cachedIntegral->addOwnedComponents(*normInt) ;
569 cachedIntegral->setCacheSource(true) ;
570 if (normInt->operMode()==ADirty) {
571 cachedIntegral->setOperMode(ADirty) ;
572 }
574 }
575
576 }
577 _norm = normInt ;
578 }
579
580 // Register new normalization with manager (takes ownership)
581 cache = new CacheElem(*_norm) ;
582 _normMgr.setObj(nset,cache) ;
583
584// std::cout << "making new object " << _norm->GetName() << std::endl ;
585
586 return true ;
587}
588
589
590
591////////////////////////////////////////////////////////////////////////////////
592/// Reset error counter to given value, limiting the number
593/// of future error messages for this pdf to 'resetValue'
594
600
601
602
603////////////////////////////////////////////////////////////////////////////////
604/// Reset trace counter to given value, limiting the
605/// number of future trace messages for this pdf to 'value'
606
608{
609 if (!allNodes) {
611 return ;
612 } else {
615 for(auto * pdf : dynamic_range_cast<RooAbsPdf*>(branchList)) {
616 if (pdf) pdf->setTraceCounter(value,false) ;
617 }
618 }
619
620}
621
622
623
624
625////////////////////////////////////////////////////////////////////////////////
626/// Return the log of the current value with given normalization
627/// An error message is printed if the argument of the log is negative.
628
629double RooAbsPdf::getLogVal(const RooArgSet* nset) const
630{
631 return getLog(getVal(nset), this);
632}
633////////////////////////////////////////////////////////////////////////////////
634/// This function returns the penalty term.
635/// Penalty terms modify the likelihood,during model parameter estimation.This penalty term is usually
636// a function of the model parameters
638{
639 return 0;
640}
641
642////////////////////////////////////////////////////////////////////////////////
643/// Check for infinity or NaN.
644/// \param[in] inputs Array to check
645/// \return True if either infinity or NaN were found.
646namespace {
647template<class T>
648bool checkInfNaNNeg(const T& inputs) {
649 // check for a math error or negative value
650 bool inf = false;
651 bool nan = false;
652 bool neg = false;
653
654 for (double val : inputs) { //CHECK_VECTORISE
655 inf |= !std::isfinite(val);
656 nan |= TMath::IsNaN(val); // Works also during fast math
657 neg |= val < 0;
658 }
659
660 return inf || nan || neg;
661}
662}
663
664
665////////////////////////////////////////////////////////////////////////////////
666/// Scan through outputs and fix+log all nans and negative values.
667/// \param[in,out] outputs Array to be scanned & fixed.
668/// \param[in] begin Begin of event range. Only needed to print the correct event number
669/// where the error occurred.
670void RooAbsPdf::logBatchComputationErrors(std::span<const double>& outputs, std::size_t begin) const {
671 for (unsigned int i=0; i<outputs.size(); ++i) {
672 const double value = outputs[i];
673 if (TMath::IsNaN(outputs[i])) {
674 logEvalError(Form("p.d.f value of (%s) is Not-a-Number (%f) for entry %zu",
675 GetName(), value, begin+i));
676 } else if (!std::isfinite(outputs[i])){
677 logEvalError(Form("p.d.f value of (%s) is (%f) for entry %zu",
678 GetName(), value, begin+i));
679 } else if (outputs[i] < 0.) {
680 logEvalError(Form("p.d.f value of (%s) is less than zero (%f) for entry %zu",
681 GetName(), value, begin+i));
682 }
683 }
684}
685
686
687void RooAbsPdf::getLogProbabilities(std::span<const double> pdfValues, double * output) const {
688 for (std::size_t i = 0; i < pdfValues.size(); ++i) {
689 output[i] = getLog(pdfValues[i], this);
690 }
691}
692
693////////////////////////////////////////////////////////////////////////////////
694/// Return the extended likelihood term (\f$ N_\mathrm{expect} - N_\mathrm{observed} \cdot \log(N_\mathrm{expect} \f$)
695/// of this PDF for the given number of observed events.
696///
697/// For successful operation, the PDF implementation must indicate that
698/// it is extendable by overloading `canBeExtended()`, and must
699/// implement the `expectedEvents()` function.
700///
701/// \param[in] sumEntries The number of observed events.
702/// \param[in] nset The normalization set when asking the pdf for the expected
703/// number of events.
704/// \param[in] sumEntriesW2 The number of observed events when weighting with
705/// squared weights. If non-zero, the weight-squared error
706/// correction is applied to the extended term.
707/// \param[in] doOffset Offset the extended term by a counterterm where the
708/// expected number of events equals the observed number of events.
709/// This constant shift results in a term closer to zero that is
710/// approximately chi-square distributed. It is useful to do this
711/// also when summing multiple NLL terms to avoid numeric precision
712/// loss that happens if you sum multiple terms of different orders
713/// of magnitude.
714///
715/// The weight-squared error correction works as follows:
716/// adjust poisson such that
717/// estimate of \f$N_\mathrm{expect}\f$ stays at the same value, but has a different variance, rescale
718/// both the observed and expected count of the Poisson with a factor \f$ \sum w_{i} / \sum w_{i}^2 \f$
719/// (the effective weight of the Poisson term),
720/// i.e., change \f$\mathrm{Poisson}(N_\mathrm{observed} = \sum w_{i} | N_\mathrm{expect} )\f$
721/// to \f$ \mathrm{Poisson}(\sum w_{i} \cdot \sum w_{i} / \sum w_{i}^2 | N_\mathrm{expect} \cdot \sum w_{i} / \sum w_{i}^2 ) \f$,
722/// weighted by the effective weight \f$ \sum w_{i}^2 / \sum w_{i} \f$ in the likelihood.
723/// Since here we compute the likelihood with the weight square, we need to multiply by the
724/// square of the effective weight:
725/// - \f$ W_\mathrm{expect} = N_\mathrm{expect} \cdot \sum w_{i} / \sum w_{i}^2 \f$ : effective expected entries
726/// - \f$ W_\mathrm{observed} = \sum w_{i} \cdot \sum w_{i} / \sum w_{i}^2 \f$ : effective observed entries
727///
728/// The extended term for the likelihood weighted by the square of the weight will be then:
729///
730/// \f$ \left(\sum w_{i}^2 / \sum w_{i}\right)^2 \cdot W_\mathrm{expect} - (\sum w_{i}^2 / \sum w_{i})^2 \cdot W_\mathrm{observed} \cdot \log{W_\mathrm{expect}} \f$
731///
732/// aund this is using the previous expressions for \f$ W_\mathrm{expect} \f$ and \f$ W_\mathrm{observed} \f$:
733///
734/// \f$ \sum w_{i}^2 / \sum w_{i} \cdot N_\mathrm{expect} - \sum w_{i}^2 \cdot \log{W_\mathrm{expect}} \f$
735///
736/// Since the weights are constants in the likelihood we can use \f$\log{N_\mathrm{expect}}\f$ instead of \f$\log{W_\mathrm{expect}}\f$.
737///
738/// See also RooAbsPdf::extendedTerm(RooAbsData const& data, bool weightSquared, bool doOffset),
739/// which takes a dataset to extract \f$N_\mathrm{observed}\f$ and the
740/// normalization set.
741double RooAbsPdf::extendedTerm(double sumEntries, RooArgSet const* nset, double sumEntriesW2, bool doOffset) const
742{
743 return extendedTerm(sumEntries, expectedEvents(nset), sumEntriesW2, doOffset);
744}
745
746double RooAbsPdf::extendedTerm(double sumEntries, double expected, double sumEntriesW2, bool doOffset) const
747{
748 // check if this PDF supports extended maximum likelihood fits
749 if(!canBeExtended()) {
750 coutE(InputArguments) << GetName() << ": this PDF does not support extended maximum likelihood"
751 << std::endl;
752 return 0.0;
753 }
754
755 if(expected < 0.0) {
756 coutE(InputArguments) << GetName() << ": calculated negative expected events: " << expected
757 << std::endl;
758 logEvalError("extendedTerm #expected events is <0 return a NaN");
759 return TMath::QuietNaN();
760 }
761
762
763 // Explicitly handle case Nobs=Nexp=0
764 if (std::abs(expected)<1e-10 && std::abs(sumEntries)<1e-10) {
765 return 0.0;
766 }
767
768 // Check for errors in Nexpected
769 if (TMath::IsNaN(expected)) {
770 logEvalError("extendedTerm #expected events is a NaN") ;
771 return TMath::QuietNaN() ;
772 }
773
774 double extra = doOffset
775 ? (expected - sumEntries) - sumEntries * (std::log(expected) - std::log(sumEntries))
776 : expected - sumEntries * std::log(expected);
777
778 if(sumEntriesW2 != 0.0) {
779 extra *= sumEntriesW2 / sumEntries;
780 }
781
782 return extra;
783}
784
785////////////////////////////////////////////////////////////////////////////////
786/// Return the extended likelihood term (\f$ N_\mathrm{expect} - N_\mathrm{observed} \cdot \log(N_\mathrm{expect} \f$)
787/// of this PDF for the given number of observed events.
788///
789/// This function is a wrapper around
790/// RooAbsPdf::extendedTerm(double, RooArgSet const *, double, bool) const,
791/// where the number of observed events and observables to be used as the
792/// normalization set for the pdf is extracted from a RooAbsData.
793///
794/// For successful operation, the PDF implementation must indicate that
795/// it is extendable by overloading `canBeExtended()`, and must
796/// implement the `expectedEvents()` function.
797///
798/// \param[in] data The RooAbsData to retrieve the set of observables and
799/// number of expected events.
800/// \param[in] weightSquared If set to `true`, the extended term will be scaled by
801/// the ratio of squared event weights over event weights:
802/// \f$ \sum w_{i}^2 / \sum w_{i} \f$.
803/// Intended to be used by fits with the `SumW2Error()` option that
804/// can be passed to RooAbsPdf::fitTo()
805/// (see the documentation of said function to learn more about the
806/// interpretation of fits with squared weights).
807/// \param[in] doOffset See RooAbsPdf::extendedTerm(double, RooArgSet const*, double, bool) const.
808
809double RooAbsPdf::extendedTerm(RooAbsData const& data, bool weightSquared, bool doOffset) const {
810 double sumW = data.sumEntries();
811 double sumW2 = 0.0;
812 if (weightSquared) {
813 sumW2 = data.sumEntriesW2();
814 }
815 return extendedTerm(sumW, data.get(), sumW2, doOffset);
816}
817
818
819/** @fn RooAbsPdf::createNLL()
820 *
821 * @brief Construct representation of -log(L) of PDF with given dataset.
822 *
823 * If dataset is unbinned, an unbinned likelihood is constructed.
824 * If the dataset is binned, a binned likelihood is constructed.
825 *
826 * @param data Reference to a RooAbsData object representing the dataset.
827 * @param cmdArgs Variadic template arguments representing optional command arguments.
828 * You can pass either an arbitrary number of RooCmdArg instances
829 * or a single RooLinkedList that points to the RooCmdArg objects.
830 * @return An owning pointer to the created RooAbsReal NLL object.
831 *
832 * @tparam CmdArgs_t Template types for optional command arguments.
833 * Can either be an arbitrary number of RooCmdArg or a single RooLinkedList.
834 *
835 * \note This front-end function should not be re-implemented in derived PDF types.
836 * If you mean to customize the NLL creation routine,
837 * you need to override the virtual RooAbsPdf::createNLLImpl() method.
838 *
839 * The following named arguments are supported:
840 *
841 * <table>
842 * <tr><th> Type of CmdArg <th> Effect on NLL
843 * <tr><td> `ConditionalObservables(Args_t &&... argsOrArgSet)` <td> Do not normalize PDF over listed observables.
844 * Arguments can either be multiple RooRealVar or a single RooArgSet containing them.
845 * <tr><td> `Extended(bool flag)` <td> Add extended likelihood term, off by default.
846 * <tr><td> `Range(const char* name)` <td> Fit only data inside range with given name. Multiple comma-separated range names can be specified.
847 * In this case, the unnormalized PDF \f$f(x)\f$ is normalized by the integral over all ranges \f$r_i\f$:
848 * \f[
849 * p(x) = \frac{f(x)}{\sum_i \int_{r_i} f(x) dx}.
850 * \f]
851 * <tr><td> `Range(double lo, double hi)` <td> Fit only data inside given range. A range named "fit" is created on the fly on all observables.
852 * <tr><td> `SumCoefRange(const char* name)` <td> Set the range in which to interpret the coefficients of RooAddPdf components
853 * <tr><td> `NumCPU(int num, int istrat)` <td> Parallelize NLL calculation on num CPUs. (Currently, this setting is ignored with the **cpu** Backend.)
854 * <table>
855 * <tr><th> Strategy <th> Effect
856 * <tr><td> 0 = RooFit::BulkPartition - *default* <td> Divide events in N equal chunks
857 * <tr><td> 1 = RooFit::Interleave <td> Process event i%N in process N. Recommended for binned data with
858 * a substantial number of zero-bins, which will be distributed across processes more equitably in this strategy
859 * <tr><td> 2 = RooFit::SimComponents <td> Process each component likelihood of a RooSimultaneous fully in a single process
860 * and distribute components over processes. This approach can be beneficial if normalization calculation time
861 * dominates the total computation time of a component (since the normalization calculation must be performed
862 * in each process in strategies 0 and 1. However beware that if the RooSimultaneous components do not share many
863 * parameters this strategy is inefficient: as most minuit-induced likelihood calculations involve changing
864 * a single parameter, only 1 of the N processes will be active most of the time if RooSimultaneous components
865 * do not share many parameters
866 * <tr><td> 3 = RooFit::Hybrid <td> Follow strategy 0 for all RooSimultaneous components, except those with less than
867 * 30 dataset entries, for which strategy 2 is followed.
868 * </table>
869 * <tr><td> `EvalBackend(std::string const&)` <td> Choose a likelihood evaluation backend:
870 * <table>
871 * <tr><th> Backend <th> Description
872 * <tr><td> **cpu** - *default* <td> New vectorized evaluation mode, using faster math functions and auto-vectorisation (currently on a single thread).
873 * Since ROOT 6.23, this is the default if `EvalBackend()` is not passed, succeeding the **legacy** backend.
874 * If all RooAbsArg objects in the model support vectorized evaluation,
875 * likelihood computations are 2 to 10 times faster than with the **legacy** backend (each on a single thread).
876 * - unless your dataset is so small that the vectorization is not worth it.
877 * The relative difference of the single log-likelihoods with respect to the legacy mode is usually better than \f$10^{-12}\f$,
878 * and for fit parameters it's usually better than \f$10^{-6}\f$. In past ROOT releases, this backend could be activated with the now deprecated `BatchMode()` option.
879 * <tr><td> **cuda** <td> Evaluate the likelihood on a GPU that supports CUDA.
880 * This backend re-uses code from the **cpu** backend, but compiled in CUDA kernels.
881 * Hence, the results are expected to be identical, modulo some numerical differences that can arise from the different order in which the GPU is summing the log probabilities.
882 * This backend can drastically speed up the fit if all RooAbsArg object in the model support it.
883 * <tr><td> **legacy** <td> The original likelihood evaluation method.
884 * Evaluates the PDF for each single data entry at a time before summing the negative log probabilities.
885 * It supports multi-threading, but you might need more than 20 threads to maybe see about 10% performance gain over the default cpu-backend (that runs currently only on a single thread).
886 * <tr><td> **codegen** <td> **Experimental** - Generates and compiles minimal C++ code for the NLL on-the-fly and wraps it in the returned RooAbsReal.
887 * Also generates and compiles the code for the gradient using Automatic Differentiation (AD) with [Clad](https://github.com/vgvassilev/clad).
888 * This analytic gradient is passed to the minimizer, which can result in significant speedups for many-parameter fits,
889 * even compared to the **cpu** backend. However, if one of the RooAbsArg objects in the model does not support the code generation,
890 * this backend can't be used.
891 * <tr><td> **codegen_no_grad** <td> **Experimental** - Same as **codegen**, but doesn't generate and compile the gradient code and use the regular numerical differentiation instead.
892 * This is expected to be slower, but useful for debugging problems with the analytic gradient.
893 * </table>
894 * <tr><td> `Optimize(bool flag)` <td> Activate constant term optimization (on by default)
895 * <tr><td> `SplitRange(bool flag)` <td> Use separate fit ranges in a simultaneous fit. Actual range name for each subsample is assumed to
896 * be `rangeName_indexState`, where `indexState` is the state of the master index category of the simultaneous fit.
897 * Using `Range("range"), SplitRange()` as switches, different ranges could be set like this:
898 * ```
899 * myVariable.setRange("range_pi0", 135, 210);
900 * myVariable.setRange("range_gamma", 50, 210);
901 * ```
902 * <tr><td> `Constrain(const RooArgSet&pars)` <td> For p.d.f.s that contain internal parameter constraint terms (that is usually product PDFs, where one
903 * term of the product depends on parameters but not on the observable(s),), only apply constraints to the given subset of parameters.
904 * <tr><td> `ExternalConstraints(const RooArgSet& )` <td> Include given external constraints to likelihood by multiplying them with the original likelihood.
905 * <tr><td> `GlobalObservables(const RooArgSet&)` <td> Define the set of normalization observables to be used for the constraint terms.
906 * If none are specified the constrained parameters are used.
907 * <tr><td> `GlobalObservablesSource(const char* sourceName)` <td> Which source to prioritize for global observable values.
908 * Can be either:
909 * - `data`: to take the values from the dataset,
910 * falling back to the pdf value if a given global observable is not available.
911 * If no `GlobalObservables` or `GlobalObservablesTag` command argument is given, the set
912 * of global observables will be automatically defined to be the set stored in the data.
913 * - `model`: to take all values from the pdf and completely ignore the set of global observables stored in the data
914 * (not even using it to automatically define the set of global observables
915 * if the `GlobalObservables` or `GlobalObservablesTag` command arguments are not given).
916 * The default option is `data`.
917 * <tr><td> `GlobalObservablesTag(const char* tagName)` <td> Define the set of normalization observables to be used for the constraint terms by
918 * a string attribute associated with pdf observables that match the given tagName.
919 * <tr><td> `Verbose(bool flag)` <td> Controls RooFit informational messages in likelihood construction
920 * <tr><td> `CloneData(bool flag)` <td> Use clone of dataset in NLL (default is true).
921 * \warning Deprecated option that is ignored. It is up to the implementation of the NLL creation method if the data is cloned or not.
922 * <tr><td> `Offset(std::string const& mode)` <td> Likelihood offsetting mode. Can be either:
923 * <table>
924 * <tr><th> Mode <th> Description
925 * <tr><td> **none** - *default* <td> No offsetting.
926 * <tr><td> **initial** <td> Offset likelihood by initial value (so that starting value of FCN in minuit is zero).
927 * This can improve numeric stability in simultaneous fits with components with large likelihood values.
928 * <tr><td> **bin** <td> Offset likelihood bin-by-bin with a template histogram model based on the obersved data.
929 * This results in per-bin values that are all in the same order of magnitude, which reduces precision loss in the sum,
930 * which can drastically improve numeric stability.
931 * Furthermore, \f$2\cdot \text{NLL}\f$ defined like this is approximately chi-square distributed, allowing for goodness-of-fit tests.
932 * </table>
933 * <tr><td> `IntegrateBins(double precision)` <td> In binned fits, integrate the PDF over the bins instead of using the probability density at the bin centre.
934 * This can reduce the bias observed when fitting functions with high curvature to binned data.
935 * - precision > 0: Activate bin integration everywhere. Use precision between 0.01 and 1.E-6, depending on binning.
936 * Note that a low precision such as 0.01 might yield identical results to 1.E-4, since the integrator might reach 1.E-4 already in its first
937 * integration step. If lower precision is desired (more speed), a RooBinSamplingPdf has to be created manually, and its integrator
938 * has to be manipulated directly.
939 * - precision = 0: Activate bin integration only for continuous PDFs fit to a RooDataHist.
940 * - precision < 0: Deactivate.
941 * \see RooBinSamplingPdf
942 * <tr><td> `ModularL(bool flag)` <td> Enable or disable modular likelihoods, which will become the default in a future release.
943 * This does not change any user-facing code, but only enables a different likelihood class in the back-end. Note that this
944 * should be set to true for parallel minimization of likelihoods!
945 * Note that it is currently not recommended to use Modular likelihoods without any parallelization enabled in the minimization, since
946 * some features such as offsetting might not yet work in this case.
947 * </table>
948 */
949
950
951/** @brief Protected implementation of the NLL creation routine.
952 *
953 * This virtual function can be overridden in case you want to change the NLL creation logic for custom PDFs.
954 *
955 * \note Never call this function directly. Instead, call RooAbsPdf::createNLL().
956 */
957
958std::unique_ptr<RooAbsReal> RooAbsPdf::createNLLImpl(RooAbsData &data, const RooLinkedList &cmdList)
959{
960 return RooFit::FitHelpers::createNLL(*this, data, cmdList);
961}
962
963
964/** @fn RooAbsPdf::fitTo()
965 *
966 * @brief Fit PDF to given dataset.
967 *
968 * If dataset is unbinned, an unbinned maximum likelihood is performed.
969 * If the dataset is binned, a binned maximum likelihood is performed.
970 * By default the fit is executed through the MINUIT commands MIGRAD, HESSE in succession.
971 *
972 * @param data Reference to a RooAbsData object representing the dataset.
973 * @param cmdArgs Variadic template arguments representing optional command arguments.
974 * You can pass either an arbitrary number of RooCmdArg instances
975 * or a single RooLinkedList that points to the RooCmdArg objects.
976 * @return An owning pointer to the created RooAbsReal NLL object.
977 * @return RooFitResult with fit status and parameters if option Save() is used, `nullptr` otherwise. The user takes ownership of the fit result.
978 *
979 * @tparam CmdArgs_t Template types for optional command arguments.
980 * Can either be an arbitrary number of RooCmdArg or a single RooLinkedList.
981 *
982 * \note This front-end function should not be re-implemented in derived PDF types.
983 * If you mean to customize the likelihood fitting routine,
984 * you need to override the virtual RooAbsPdf::fitToImpl() method.
985 *
986 * The following named arguments are supported:
987 *
988 * <table>
989 * <tr><th> Type of CmdArg <th> Options to control construction of -log(L)
990 * <tr><td> <td> All command arguments that can also be passed to the NLL creation method.
991 * \see RooAbsPdf::createNLL()
992 *
993 * <tr><th><th> Options to control flow of fit procedure
994 * <tr><td> `Minimizer("<type>", "<algo>")` <td> Choose minimization package and optionally the algorithm to use. Default is MINUIT/MIGRAD through the RooMinimizer interface,
995 * but others can be specified (through RooMinimizer interface).
996 * <table>
997 * <tr><th> Type <th> Algorithm
998 * <tr><td> Minuit <td> migrad, simplex, minimize (=migrad+simplex), migradimproved (=migrad+improve)
999 * <tr><td> Minuit2 <td> migrad, simplex, minimize, scan
1000 * <tr><td> GSLMultiMin <td> conjugatefr, conjugatepr, bfgs, bfgs2, steepestdescent
1001 * <tr><td> GSLSimAn <td> -
1002 * </table>
1003 *
1004 * <tr><td> `InitialHesse(bool flag)` <td> Flag controls if HESSE before MIGRAD as well, off by default
1005 * <tr><td> `Optimize(bool flag)` <td> Activate constant term optimization of test statistic during minimization (on by default)
1006 * <tr><td> `Hesse(bool flag)` <td> Flag controls if HESSE is run after MIGRAD, on by default
1007 * <tr><td> `Minos(bool flag)` <td> Flag controls if MINOS is run after HESSE, off by default
1008 * <tr><td> `Minos(const RooArgSet& set)` <td> Only run MINOS on given subset of arguments
1009 * <tr><td> `Save(bool flag)` <td> Flag controls if RooFitResult object is produced and returned, off by default
1010 * <tr><td> `Strategy(Int_t flag)` <td> Set Minuit strategy (0 to 2, default is 1)
1011 * <tr><td> `MaxCalls(int n)` <td> Change maximum number of likelihood function calls from MINUIT (if `n <= 0`, the default of 500 * #%parameters is used)
1012 * <tr><td> `EvalErrorWall(bool flag=true)` <td> When parameters are in disallowed regions (e.g. PDF is negative), return very high value to fitter
1013 * to force it out of that region. This can, however, mean that the fitter gets lost in this region. If
1014 * this happens, try switching it off.
1015 * <tr><td> `RecoverFromUndefinedRegions(double strength)` <td> When PDF is invalid (e.g. parameter in undefined region), try to direct minimiser away from that region.
1016 * `strength` controls the magnitude of the penalty term. Leaving out this argument defaults to 10. Switch off with `strength = 0.`.
1017 *
1018 * <tr><td> `SumW2Error(bool flag)` <td> Apply correction to errors and covariance matrix.
1019 * This uses two covariance matrices, one with the weights, the other with squared weights,
1020 * to obtain the correct errors for weighted likelihood fits. If this option is activated, the
1021 * corrected covariance matrix is calculated as \f$ V_\mathrm{corr} = V C^{-1} V \f$, where \f$ V \f$ is the original
1022 * covariance matrix and \f$ C \f$ is the inverse of the covariance matrix calculated using the
1023 * squared weights. This allows to switch between two interpretations of errors:
1024 * <table>
1025 * <tr><th> SumW2Error <th> Interpretation
1026 * <tr><td> true <td> The errors reflect the uncertainty of the Monte Carlo simulation.
1027 * Use this if you want to know how much accuracy you can get from the available Monte Carlo statistics.
1028 *
1029 * **Example**: Simulation with 1000 events, the average weight is 0.1.
1030 * The errors are as big as if one fitted to 1000 events.
1031 * <tr><td> false <td> The errors reflect the errors of a dataset, which is as big as the sum of weights.
1032 * Use this if you want to know what statistical errors you would get if you had a dataset with as many
1033 * events as the (weighted) Monte Carlo simulation represents.
1034 *
1035 * **Example** (Data as above):
1036 * The errors are as big as if one fitted to 100 events.
1037 * </table>
1038 * \note If the `SumW2Error` correction is enabled, the covariance matrix quality stored in the RooFitResult
1039 * object will be the minimum of the original covariance matrix quality and the quality of the covariance
1040 * matrix calculated with the squared weights.
1041 * <tr><td> `AsymptoticError()` <td> Use the asymptotically correct approach to estimate errors in the presence of weights.
1042 * This is slower but more accurate than `SumW2Error`. See also https://arxiv.org/abs/1911.01303).
1043 This option even correctly implements the case of extended likelihood fits
1044 (see this [writeup on extended weighted fits](https://root.cern/files/extended_weighted_fits.pdf) that complements the paper linked before).
1045 * <tr><td> `PrefitDataFraction(double fraction)`
1046 * <td> Runs a prefit on a small dataset of size fraction*(actual data size). This can speed up fits
1047 * by finding good starting values for the parameters for the actual fit.
1048 * \warning Prefitting may give bad results when used in binned analysis.
1049 *
1050 * <tr><th><th> Options to control informational output
1051 * <tr><td> `Verbose(bool flag)` <td> Flag controls if verbose output is printed (NLL, parameter changes during fit).
1052 * <tr><td> `Timer(bool flag)` <td> Time CPU and wall clock consumption of fit steps, off by default.
1053 * <tr><td> `PrintLevel(Int_t level)` <td> Set Minuit print level (-1 to 3, default is 1). At -1 all RooFit informational messages are suppressed as well.
1054 * See RooMinimizer::PrintLevel for the meaning of the levels.
1055 * <tr><td> `Warnings(bool flag)` <td> Enable or disable MINUIT warnings (enabled by default)
1056 * <tr><td> `PrintEvalErrors(Int_t numErr)` <td> Control number of p.d.f evaluation errors printed per likelihood evaluation.
1057 * A negative value suppresses output completely, a zero value will only print the error count per p.d.f component,
1058 * a positive value will print details of each error up to `numErr` messages per p.d.f component.
1059 * <tr><td> `Parallelize(Int_t nWorkers)` <td> Control global parallelization settings. Arguments 1 and above enable the use of RooFit's parallel minimization
1060 * backend and uses the number given as the number of workers to use in the parallelization. -1 also enables
1061 * RooFit's parallel minimization backend, and sets the number of workers to the number of available processes.
1062 * 0 disables this feature.
1063 * In case parallelization is requested, this option implies `ModularL(true)` in the internal call to the NLL creation method.
1064 * <tr><td> `ParallelGradientOptions(bool enable=true, int orderStrategy=0, int chainFactor=1)` <td> **Experimental** - Control gradient parallelization settings. The first argument
1065 * only disables or enables gradient parallelization, this is on by default.
1066 * The second argument determines the internal partial derivative calculation
1067 * ordering strategy. The third argument determines the number of partial
1068 * derivatives that are executed per task package on each worker.
1069 * <tr><td> `ParallelDescentOptions(bool enable=false, int splitStrategy=0, int numSplits=4)` <td> **Experimental** - Control settings related to the parallelization of likelihoods
1070 * outside of the gradient calculation but in the minimization, most prominently
1071 * in the linesearch step. The first argument this disables or enables likelihood
1072 * parallelization. The second argument determines whether to split the task batches
1073 * per event or per likelihood component. And the third argument how many events or
1074 * respectively components to include in each batch.
1075 * <tr><td> `TimingAnalysis(bool flag)` <td> **Experimental** - Log timings. This feature logs timings with NewStyle likelihoods on multiple processes simultaneously
1076 * and outputs the timings at the end of a run to json log files, which can be analyzed with the
1077 * `RooFit::MultiProcess::HeatmapAnalyzer`. Only works with simultaneous likelihoods.
1078 * </table>
1079 */
1080
1081
1082/** @brief Protected implementation of the likelihood fitting routine.
1083 *
1084 * This virtual function can be overridden in case you want to change the likelihood fitting logic for custom PDFs.
1085 *
1086 * \note Never call this function directly. Instead, call RooAbsPdf::fitTo().
1087 */
1088
1089std::unique_ptr<RooFitResult> RooAbsPdf::fitToImpl(RooAbsData& data, const RooLinkedList& cmdList)
1090{
1091 return RooFit::FitHelpers::fitTo(*this, data, cmdList, false);
1092}
1093
1094
1095////////////////////////////////////////////////////////////////////////////////
1096/// Print value of p.d.f, also print normalization integral that was last used, if any
1097
1098void RooAbsPdf::printValue(ostream& os) const
1099{
1100 // silent warning messages coming when evaluating a RooAddPdf without a normalization set
1102
1103 getVal() ;
1104
1105 if (_norm) {
1106 os << getVal() << "/" << _norm->getVal() ;
1107 } else {
1108 os << getVal();
1109 }
1110}
1111
1112
1113
1114////////////////////////////////////////////////////////////////////////////////
1115/// Print multi line detailed information of this RooAbsPdf
1116
1117void RooAbsPdf::printMultiline(ostream& os, Int_t contents, bool verbose, TString indent) const
1118{
1119 RooAbsReal::printMultiline(os,contents,verbose,indent);
1120 os << indent << "--- RooAbsPdf ---" << std::endl;
1121 os << indent << "Cached value = " << _value << std::endl ;
1122 if (_norm) {
1123 os << indent << " Normalization integral: " << std::endl ;
1124 auto moreIndent = std::string(indent.Data()) + " " ;
1126 }
1127}
1128
1129
1130
1131////////////////////////////////////////////////////////////////////////////////
1132/// Return a binned generator context
1133
1135{
1136 return new RooBinnedGenContext(*this,vars,nullptr,nullptr,verbose) ;
1137}
1138
1139
1140////////////////////////////////////////////////////////////////////////////////
1141/// Interface function to create a generator context from a p.d.f. This default
1142/// implementation returns a 'standard' context that works for any p.d.f
1143
1145 const RooArgSet* auxProto, bool verbose) const
1146{
1147 return new RooGenContext(*this,vars,prototype,auxProto,verbose) ;
1148}
1149
1150
1151////////////////////////////////////////////////////////////////////////////////
1152
1154 bool verbose, bool autoBinned, const char* binnedTag) const
1155{
1156 if (prototype || (auxProto && !auxProto->empty())) {
1157 return genContext(vars,prototype,auxProto,verbose);
1158 }
1159
1160 RooAbsGenContext *context(nullptr) ;
1161 if ( (autoBinned && isBinnedDistribution(vars)) || ( binnedTag && strlen(binnedTag) && (getAttribute(binnedTag)||string(binnedTag)=="*"))) {
1162 context = binnedGenContext(vars,verbose) ;
1163 } else {
1164 context= genContext(vars,nullptr,nullptr,verbose);
1165 }
1166 return context ;
1167}
1168
1169
1170
1171////////////////////////////////////////////////////////////////////////////////
1172/// Generate a new dataset containing the specified variables with events sampled from our distribution.
1173/// Generate the specified number of events or expectedEvents() if not specified.
1174/// \param[in] whatVars Choose variables in which to generate events. Variables not listed here will remain
1175/// constant and not be used for event generation.
1176/// \param[in] arg1,arg2,arg3,arg4,arg5,arg6 Optional RooCmdArg() to change behaviour of generate().
1177/// \return RooDataSet *, owned by caller.
1178///
1179/// Any variables of this PDF that are not in whatVars will use their
1180/// current values and be treated as fixed parameters. Returns zero
1181/// in case of an error.
1182///
1183/// <table>
1184/// <tr><th> Type of CmdArg <th> Effect on generate
1185/// <tr><td> `Name(const char* name)` <td> Name of the output dataset
1186/// <tr><td> `Verbose(bool flag)` <td> Print informational messages during event generation
1187/// <tr><td> `NumEvents(int nevt)` <td> Generate specified number of events
1188/// <tr><td> `Extended()` <td> If no number of events to be generated is given,
1189/// use expected number of events from extended likelihood term.
1190/// This evidently only works for extended PDFs.
1191/// <tr><td> `GenBinned(const char* tag)` <td> Use binned generation for all component pdfs that have 'setAttribute(tag)' set
1192/// <tr><td> `AutoBinned(bool flag)` <td> Automatically deploy binned generation for binned distributions (e.g. RooHistPdf, sums and products of
1193/// RooHistPdfs etc)
1194/// \note Datasets that are generated in binned mode are returned as weighted unbinned datasets. This means that
1195/// for each bin, there will be one event in the dataset with a weight corresponding to the (possibly randomised) bin content.
1196///
1197///
1198/// <tr><td> `AllBinned()` <td> As above, but for all components.
1199/// \note The notion of components is only meaningful for simultaneous PDFs
1200/// as binned generation is always executed at the top-level node for a regular
1201/// PDF, so for those it only mattes that the top-level node is tagged.
1202///
1203/// <tr><td> ProtoData(const RooDataSet& data, bool randOrder)
1204/// <td> Use specified dataset as prototype dataset. If randOrder in ProtoData() is set to true,
1205/// the order of the events in the dataset will be read in a random order if the requested
1206/// number of events to be generated does not match the number of events in the prototype dataset.
1207/// \note If ProtoData() is used, the specified existing dataset as a prototype: the new dataset will contain
1208/// the same number of events as the prototype (unless otherwise specified), and any prototype variables not in
1209/// whatVars will be copied into the new dataset for each generated event and also used to set our PDF parameters.
1210/// The user can specify a number of events to generate that will override the default. The result is a
1211/// copy of the prototype dataset with only variables in whatVars randomized. Variables in whatVars that
1212/// are not in the prototype will be added as new columns to the generated dataset.
1213///
1214/// </table>
1215///
1216/// #### Accessing the underlying event generator
1217/// Depending on the fit model (if it is difficult to sample), it may be necessary to change generator settings.
1218/// For the default generator (RooFoamGenerator), the number of samples or cells could be increased by e.g. using
1219/// myPdf->specialGeneratorConfig()->getConfigSection("RooFoamGenerator").setRealValue("nSample",1e4);
1220///
1221/// The foam generator e.g. has the following config options:
1222/// - nCell[123N]D
1223/// - nSample
1224/// - chatLevel
1225/// \see rf902_numgenconfig.C
1226
1228 const RooCmdArg& arg3,const RooCmdArg& arg4, const RooCmdArg& arg5,const RooCmdArg& arg6)
1229{
1230 // Select the pdf-specific commands
1231 RooCmdConfig pc("RooAbsPdf::generate(" + std::string(GetName()) + ")");
1232 pc.defineObject("proto","PrototypeData",0,nullptr) ;
1233 pc.defineString("dsetName","Name",0,"") ;
1234 pc.defineInt("randProto","PrototypeData",0,0) ;
1235 pc.defineInt("resampleProto","PrototypeData",1,0) ;
1236 pc.defineInt("verbose","Verbose",0,0) ;
1237 pc.defineInt("extended","Extended",0,0) ;
1238 pc.defineInt("nEvents","NumEvents",0,0) ;
1239 pc.defineInt("autoBinned","AutoBinned",0,1) ;
1240 pc.defineInt("expectedData","ExpectedData",0,0) ;
1241 pc.defineDouble("nEventsD","NumEventsD",0,-1.) ;
1242 pc.defineString("binnedTag","GenBinned",0,"") ;
1243 pc.defineMutex("GenBinned","ProtoData") ;
1244 pc.defineMutex("Extended", "NumEvents");
1245
1246 // Process and check varargs
1248 if (!pc.ok(true)) {
1249 return nullptr;
1250 }
1251
1252 // Decode command line arguments
1253 RooDataSet* protoData = static_cast<RooDataSet*>(pc.getObject("proto",nullptr)) ;
1254 const char* dsetName = pc.getString("dsetName") ;
1255 bool verbose = pc.getInt("verbose") ;
1256 bool randProto = pc.getInt("randProto") ;
1257 bool resampleProto = pc.getInt("resampleProto") ;
1258 bool extended = pc.getInt("extended") ;
1259 bool autoBinned = pc.getInt("autoBinned") ;
1260 const char* binnedTag = pc.getString("binnedTag") ;
1261 Int_t nEventsI = pc.getInt("nEvents") ;
1262 double nEventsD = pc.getInt("nEventsD") ;
1263 //bool verbose = pc.getInt("verbose") ;
1264 bool expectedData = pc.getInt("expectedData") ;
1265
1266 double nEvents = (nEventsD>0) ? nEventsD : double(nEventsI);
1267
1268 // Force binned mode for expected data mode
1269 if (expectedData) {
1270 binnedTag="*" ;
1271 }
1272
1273 if (extended) {
1274 if (nEvents == 0) nEvents = expectedEvents(&whatVars);
1275 } else if (nEvents==0) {
1276 cxcoutI(Generation) << "No number of events specified , number of events generated is "
1277 << GetName() << "::expectedEvents() = " << expectedEvents(&whatVars)<< std::endl ;
1278 }
1279
1280 if (extended && protoData && !randProto) {
1281 cxcoutI(Generation) << "WARNING Using generator option Extended() (Poisson distribution of #events) together "
1282 << "with a prototype dataset implies incomplete sampling or oversampling of proto data. "
1283 << "Set randomize flag in ProtoData() option to randomize prototype dataset order and thus "
1284 << "to randomize the set of over/undersampled prototype events for each generation cycle." << std::endl ;
1285 }
1286
1287
1288 // Forward to appropriate implementation
1289 std::unique_ptr<RooDataSet> data;
1290 if (protoData) {
1291 data = std::unique_ptr<RooDataSet>{generate(whatVars,*protoData,Int_t(nEvents),verbose,randProto,resampleProto)};
1292 } else {
1293 data = std::unique_ptr<RooDataSet>{generate(whatVars,nEvents,verbose,autoBinned,binnedTag,expectedData, extended)};
1294 }
1295
1296 // Rename dataset to given name if supplied
1297 if (dsetName && strlen(dsetName)>0) {
1298 data->SetName(dsetName) ;
1299 }
1300
1301 return RooFit::makeOwningPtr(std::move(data));
1302}
1303
1304
1305
1306
1307
1308
1309////////////////////////////////////////////////////////////////////////////////
1310/// \note This method does not perform any generation. To generate according to generations specification call RooAbsPdf::generate(RooAbsPdf::GenSpec&) const
1311///
1312/// Details copied from RooAbsPdf::generate():
1313/// --------------------------------------------
1314/// \copydetails RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&)
1315
1317 const RooCmdArg& arg1,const RooCmdArg& arg2,
1318 const RooCmdArg& arg3,const RooCmdArg& arg4,
1319 const RooCmdArg& arg5,const RooCmdArg& arg6)
1320{
1321
1322 // Select the pdf-specific commands
1323 RooCmdConfig pc("RooAbsPdf::generate(" + std::string(GetName()) + ")");
1324 pc.defineObject("proto","PrototypeData",0,nullptr) ;
1325 pc.defineString("dsetName","Name",0,"") ;
1326 pc.defineInt("randProto","PrototypeData",0,0) ;
1327 pc.defineInt("resampleProto","PrototypeData",1,0) ;
1328 pc.defineInt("verbose","Verbose",0,0) ;
1329 pc.defineInt("extended","Extended",0,0) ;
1330 pc.defineInt("nEvents","NumEvents",0,0) ;
1331 pc.defineInt("autoBinned","AutoBinned",0,1) ;
1332 pc.defineString("binnedTag","GenBinned",0,"") ;
1333 pc.defineMutex("GenBinned","ProtoData") ;
1334
1335
1336 // Process and check varargs
1338 if (!pc.ok(true)) {
1339 return nullptr ;
1340 }
1341
1342 // Decode command line arguments
1343 RooDataSet* protoData = static_cast<RooDataSet*>(pc.getObject("proto",nullptr)) ;
1344 const char* dsetName = pc.getString("dsetName") ;
1345 Int_t nEvents = pc.getInt("nEvents") ;
1346 bool verbose = pc.getInt("verbose") ;
1347 bool randProto = pc.getInt("randProto") ;
1348 bool resampleProto = pc.getInt("resampleProto") ;
1349 bool extended = pc.getInt("extended") ;
1350 bool autoBinned = pc.getInt("autoBinned") ;
1351 const char* binnedTag = pc.getString("binnedTag") ;
1352
1354
1355 return new GenSpec(cx,whatVars,protoData,nEvents,extended,randProto,resampleProto,dsetName) ;
1356}
1357
1358
1359////////////////////////////////////////////////////////////////////////////////
1360/// If many identical generation requests
1361/// are needed, e.g. in toy MC studies, it is more efficient to use the prepareMultiGen()/generate()
1362/// combination than calling the standard generate() multiple times as
1363/// initialization overhead is only incurred once.
1364
1366{
1367 //Int_t nEvt = spec._extended ? RooRandom::randomGenerator()->Poisson(spec._nGen) : spec._nGen ;
1368 //Int_t nEvt = spec._extended ? RooRandom::randomGenerator()->Poisson(spec._nGen==0?expectedEvents(spec._whatVars):spec._nGen) : spec._nGen ;
1369 //Int_t nEvt = spec._nGen == 0 ? RooRandom::randomGenerator()->Poisson(expectedEvents(spec._whatVars)) : spec._nGen;
1370
1371 double nEvt = spec._nGen == 0 ? expectedEvents(spec._whatVars) : spec._nGen;
1372
1373 std::unique_ptr<RooDataSet> ret{generate(*spec._genContext,spec._whatVars,spec._protoData, nEvt,false,spec._randProto,spec._resampleProto,
1374 spec._init,spec._extended)};
1375 spec._init = true ;
1376 return RooFit::makeOwningPtr(std::move(ret));
1377}
1378
1379
1380
1381
1382
1383////////////////////////////////////////////////////////////////////////////////
1384/// Generate a new dataset containing the specified variables with
1385/// events sampled from our distribution.
1386///
1387/// \param[in] whatVars Generate a dataset with the variables (and categories) in this set.
1388/// Any variables of this PDF that are not in `whatVars` will use their
1389/// current values and be treated as fixed parameters.
1390/// \param[in] nEvents Generate the specified number of events or else try to use
1391/// expectedEvents() if nEvents <= 0 (default).
1392/// \param[in] verbose Show which generator strategies are being used.
1393/// \param[in] autoBinned If original distribution is binned, return bin centers and randomise weights
1394/// instead of generating single events.
1395/// \param[in] binnedTag
1396/// \param[in] expectedData Call setExpectedData on the genContext.
1397/// \param[in] extended Randomise number of events generated according to Poisson(nEvents). Only useful
1398/// if PDF is extended.
1399/// \return New dataset. Returns zero in case of an error. The caller takes ownership of the returned
1400/// dataset.
1401
1402RooFit::OwningPtr<RooDataSet> RooAbsPdf::generate(const RooArgSet &whatVars, double nEvents, bool verbose, bool autoBinned, const char* binnedTag, bool expectedData, bool extended) const
1403{
1404 if (nEvents==0 && extendMode()==CanNotBeExtended) {
1405 return RooFit::makeOwningPtr(std::make_unique<RooDataSet>("emptyData","emptyData",whatVars));
1406 }
1407
1408 // Request for binned generation
1409 std::unique_ptr<RooAbsGenContext> context{autoGenContext(whatVars,nullptr,nullptr,verbose,autoBinned,binnedTag)};
1410 if (expectedData) {
1411 context->setExpectedData(true) ;
1412 }
1413
1414 std::unique_ptr<RooDataSet> generated;
1415 if(nullptr != context && context->isValid()) {
1416 generated = std::unique_ptr<RooDataSet>{context->generate(nEvents, false, extended)};
1417 }
1418 else {
1419 coutE(Generation) << "RooAbsPdf::generate(" << GetName() << ") cannot create a valid context" << std::endl;
1420 }
1421 return RooFit::makeOwningPtr(std::move(generated));
1422}
1423
1424
1425
1426
1427////////////////////////////////////////////////////////////////////////////////
1428/// Internal method
1429
1430std::unique_ptr<RooDataSet> RooAbsPdf::generate(RooAbsGenContext& context, const RooArgSet &whatVars, const RooDataSet *prototype,
1431 double nEvents, bool /*verbose*/, bool randProtoOrder, bool resampleProto,
1432 bool skipInit, bool extended) const
1433{
1434 if (nEvents==0 && (prototype==nullptr || prototype->numEntries()==0)) {
1435 return std::make_unique<RooDataSet>("emptyData","emptyData",whatVars);
1436 }
1437
1438 std::unique_ptr<RooDataSet> generated;
1439
1440 // Resampling implies reshuffling in the implementation
1441 if (resampleProto) {
1443 }
1444
1445 if (randProtoOrder && prototype && prototype->numEntries()!=nEvents) {
1446 coutI(Generation) << "RooAbsPdf::generate (Re)randomizing event order in prototype dataset (Nevt=" << nEvents << ")" << std::endl ;
1447 Int_t* newOrder = randomizeProtoOrder(prototype->numEntries(),Int_t(nEvents),resampleProto) ;
1448 context.setProtoDataOrder(newOrder) ;
1449 delete[] newOrder ;
1450 }
1451
1452 if(context.isValid()) {
1453 generated = std::unique_ptr<RooDataSet>{context.generate(nEvents,skipInit,extended)};
1454 }
1455 else {
1456 coutE(Generation) << "RooAbsPdf::generate(" << GetName() << ") do not have a valid generator context" << std::endl;
1457 }
1458 return generated;
1459}
1460
1461
1462
1463
1464////////////////////////////////////////////////////////////////////////////////
1465/// Generate a new dataset using a prototype dataset as a model,
1466/// with values of the variables in `whatVars` sampled from our distribution.
1467///
1468/// \param[in] whatVars Generate for these variables.
1469/// \param[in] prototype Use this dataset
1470/// as a prototype: the new dataset will contain the same number of
1471/// events as the prototype (by default), and any prototype variables not in
1472/// whatVars will be copied into the new dataset for each generated
1473/// event and also used to set our PDF parameters. The user can specify a
1474/// number of events to generate that will override the default. The result is a
1475/// copy of the prototype dataset with only variables in whatVars
1476/// randomized. Variables in whatVars that are not in the prototype
1477/// will be added as new columns to the generated dataset.
1478/// \param[in] nEvents Number of events to generate. Defaults to 0, which means number
1479/// of event in prototype dataset.
1480/// \param[in] verbose Show which generator strategies are being used.
1481/// \param[in] randProtoOrder Randomise order of retrieval of events from proto dataset.
1482/// \param[in] resampleProto Resample from the proto dataset.
1483/// \return The new dataset. Returns zero in case of an error. The caller takes ownership of the
1484/// returned dataset.
1485
1487 Int_t nEvents, bool verbose, bool randProtoOrder, bool resampleProto) const
1488{
1489 std::unique_ptr<RooAbsGenContext> context{genContext(whatVars,&prototype,nullptr,verbose)};
1490 if (context) {
1492 }
1493 coutE(Generation) << "RooAbsPdf::generate(" << GetName() << ") ERROR creating generator context" << std::endl ;
1494 return nullptr;
1495}
1496
1497
1498
1499////////////////////////////////////////////////////////////////////////////////
1500/// Return lookup table with randomized order for nProto prototype events.
1501
1503{
1504 // Make output list
1505 Int_t* lut = new Int_t[nProto] ;
1506
1507 // Randomly sample input list into output list
1508 if (!resampleProto) {
1509 // In this mode, randomization is a strict reshuffle of the order
1510 std::iota(lut, lut + nProto, 0); // fill the vector with 0 to nProto - 1
1511 // Shuffle code taken from https://en.cppreference.com/w/cpp/algorithm/random_shuffle.
1512 // The std::random_shuffle function was deprecated in C++17. We could have
1513 // used std::shuffle instead, but this is not straight-forward to use with
1514 // RooRandom::integer() and we didn't want to change the random number
1515 // generator. It might cause unwanted effects like reproducibility problems.
1516 for (int i = nProto-1; i > 0; --i) {
1517 std::swap(lut[i], lut[RooRandom::integer(i+1)]);
1518 }
1519 } else {
1520 // In this mode, we resample, i.e. events can be used more than once
1521 std::generate(lut, lut + nProto, [&]{ return RooRandom::integer(nProto); });
1522 }
1523
1524
1525 return lut ;
1526}
1527
1528
1529
1530////////////////////////////////////////////////////////////////////////////////
1531/// Load generatedVars with the subset of directVars that we can generate events for,
1532/// and return a code that specifies the generator algorithm we will use. A code of
1533/// zero indicates that we cannot generate any of the directVars (in this case, nothing
1534/// should be added to generatedVars). Any non-zero codes will be passed to our generateEvent()
1535/// implementation, but otherwise its value is arbitrary. The default implementation of
1536/// this method returns zero. Subclasses will usually implement this method using the
1537/// matchArgs() methods to advertise the algorithms they provide.
1538
1539Int_t RooAbsPdf::getGenerator(const RooArgSet &/*directVars*/, RooArgSet &/*generatedVars*/, bool /*staticInitOK*/) const
1540{
1541 return 0 ;
1542}
1543
1544
1545
1546////////////////////////////////////////////////////////////////////////////////
1547/// Interface for one-time initialization to setup the generator for the specified code.
1548
1550{
1551}
1552
1553
1554
1555////////////////////////////////////////////////////////////////////////////////
1556/// Interface for generation of an event using the algorithm
1557/// corresponding to the specified code. The meaning of each code is
1558/// defined by the getGenerator() implementation. The default
1559/// implementation does nothing.
1560
1562{
1563}
1564
1565
1566
1567////////////////////////////////////////////////////////////////////////////////
1568/// Check if given observable can be safely generated using the
1569/// pdfs internal generator mechanism (if that existsP). Observables
1570/// on which a PDF depends via more than route are not safe
1571/// for use with internal generators because they introduce
1572/// correlations not known to the internal generator
1573
1575{
1576 // Arg must be direct server of self
1577 if (!findServer(arg.GetName())) return false ;
1578
1579 // There must be no other dependency routes
1580 for (const auto server : _serverList) {
1581 if(server == &arg) continue;
1582 if(server->dependsOn(arg)) {
1583 return false ;
1584 }
1585 }
1586
1587 return true ;
1588}
1589
1590
1591////////////////////////////////////////////////////////////////////////////////
1592/// Generate a new dataset containing the specified variables with events sampled from our distribution.
1593/// \param[in] whatVars Choose variables in which to generate events. Variables not listed here will remain
1594/// constant and not be used for event generation
1595/// \param[in] arg1,arg2,arg3,arg4,arg5,arg6 Optional RooCmdArg to change behaviour of generateBinned()
1596/// \return RooDataHist *, to be managed by caller.
1597///
1598/// Generate the specified number of events or expectedEvents() if not specified.
1599///
1600/// Any variables of this PDF that are not in whatVars will use their
1601/// current values and be treated as fixed parameters. Returns zero
1602/// in case of an error. The caller takes ownership of the returned
1603/// dataset.
1604///
1605/// The following named arguments are supported
1606/// | Type of CmdArg | Effect on generation
1607/// |---------------------------|-----------------------
1608/// | `Name(const char* name)` | Name of the output dataset
1609/// | `Verbose(bool flag)` | Print informational messages during event generation
1610/// | `NumEvents(int nevt)` | Generate specified number of events
1611/// | `Extended()` | The actual number of events generated will be sampled from a Poisson distribution with mu=nevt. This can be *much* faster for peaked PDFs, but the number of events is not exactly what was requested.
1612/// | `ExpectedData()` | Return a binned dataset _without_ statistical fluctuations (also aliased as Asimov())
1613///
1614
1616 const RooCmdArg& arg3,const RooCmdArg& arg4, const RooCmdArg& arg5,const RooCmdArg& arg6) const
1617{
1618
1619 // Select the pdf-specific commands
1620 RooCmdConfig pc("RooAbsPdf::generate(" + std::string(GetName()) + ")");
1621 pc.defineString("dsetName","Name",0,"") ;
1622 pc.defineInt("verbose","Verbose",0,0) ;
1623 pc.defineInt("extended","Extended",0,0) ;
1624 pc.defineInt("nEvents","NumEvents",0,0) ;
1625 pc.defineDouble("nEventsD","NumEventsD",0,-1.) ;
1626 pc.defineInt("expectedData","ExpectedData",0,0) ;
1627
1628 // Process and check varargs
1630 if (!pc.ok(true)) {
1631 return nullptr;
1632 }
1633
1634 // Decode command line arguments
1635 double nEvents = pc.getDouble("nEventsD") ;
1636 if (nEvents<0) {
1637 nEvents = pc.getInt("nEvents") ;
1638 }
1639 //bool verbose = pc.getInt("verbose") ;
1640 bool extended = pc.getInt("extended") ;
1641 bool expectedData = pc.getInt("expectedData") ;
1642 const char* dsetName = pc.getString("dsetName") ;
1643
1644 if (extended) {
1645 //nEvents = (nEvents==0?Int_t(expectedEvents(&whatVars)+0.5):nEvents) ;
1646 nEvents = (nEvents==0 ? expectedEvents(&whatVars) :nEvents) ;
1647 cxcoutI(Generation) << " Extended mode active, number of events generated (" << nEvents << ") is Poisson fluctuation on "
1648 << GetName() << "::expectedEvents() = " << nEvents << std::endl ;
1649 // If Poisson fluctuation results in zero events, stop here
1650 if (nEvents==0) {
1651 return nullptr ;
1652 }
1653 } else if (nEvents==0) {
1654 cxcoutI(Generation) << "No number of events specified , number of events generated is "
1655 << GetName() << "::expectedEvents() = " << expectedEvents(&whatVars)<< std::endl ;
1656 }
1657
1658 // Forward to appropriate implementation
1660
1661 // Rename dataset to given name if supplied
1662 if (dsetName && strlen(dsetName)>0) {
1663 data->SetName(dsetName) ;
1664 }
1665
1666 return data;
1667}
1668
1669
1670
1671
1672////////////////////////////////////////////////////////////////////////////////
1673/// Generate a new dataset containing the specified variables with
1674/// events sampled from our distribution.
1675///
1676/// \param[in] whatVars Variables that values should be generated for.
1677/// \param[in] nEvents How many events to generate. If `nEvents <=0`, use the value returned by expectedEvents() as target.
1678/// \param[in] expectedData If set to true (false by default), the returned histogram returns the 'expected'
1679/// data sample, i.e. no statistical fluctuations are present.
1680/// \param[in] extended For each bin, generate Poisson(x, mu) events, where `mu` is chosen such that *on average*,
1681/// one would obtain `nEvents` events. This means that the true number of events will fluctuate around the desired value,
1682/// but the generation happens a lot faster.
1683/// Especially if the PDF is sharply peaked, the multinomial event generation necessary to generate *exactly* `nEvents` events can
1684/// be very slow.
1685///
1686/// The binning used for generation of events is the currently set binning for the variables.
1687/// It can e.g. be changed using
1688/// ```
1689/// x.setBins(15);
1690/// x.setRange(-5., 5.);
1691/// pdf.generateBinned(RooArgSet(x), 1000);
1692/// ```
1693///
1694/// Any variables of this PDF that are not in `whatVars` will use their
1695/// current values and be treated as fixed parameters.
1696/// \return RooDataHist* owned by the caller. Returns `nullptr` in case of an error.
1698{
1699 // Create empty RooDataHist
1700 auto hist = std::make_unique<RooDataHist>("genData","genData",whatVars);
1701
1702 // Scale to number of events and introduce Poisson fluctuations
1703 if (nEvents<=0) {
1704 if (!canBeExtended()) {
1705 coutE(InputArguments) << "RooAbsPdf::generateBinned(" << GetName() << ") ERROR: No event count provided and p.d.f does not provide expected number of events" << std::endl ;
1706 return nullptr;
1707 } else {
1708
1709 // Don't round in expectedData or extended mode
1710 if (expectedData || extended) {
1711 nEvents = expectedEvents(&whatVars) ;
1712 } else {
1713 nEvents = std::round(expectedEvents(&whatVars));
1714 }
1715 }
1716 }
1717
1718 // Sample p.d.f. distribution
1719 fillDataHist(hist.get(),&whatVars,1,true) ;
1720
1721 vector<int> histOut(hist->numEntries()) ;
1722 double histMax(-1) ;
1723 Int_t histOutSum(0) ;
1724 for (int i=0 ; i<hist->numEntries() ; i++) {
1725 hist->get(i) ;
1726 if (expectedData) {
1727
1728 // Expected data, multiply p.d.f by nEvents
1729 double w=hist->weight()*nEvents ;
1730 hist->set(i, w, sqrt(w));
1731
1732 } else if (extended) {
1733
1734 // Extended mode, set contents to Poisson(pdf*nEvents)
1735 double w = RooRandom::randomGenerator()->Poisson(hist->weight()*nEvents) ;
1736 hist->set(w,sqrt(w)) ;
1737
1738 } else {
1739
1740 // Regular mode, fill array of weights with Poisson(pdf*nEvents), but to not fill
1741 // histogram yet.
1742 if (hist->weight()>histMax) {
1743 histMax = hist->weight() ;
1744 }
1745 histOut[i] = RooRandom::randomGenerator()->Poisson(hist->weight()*nEvents) ;
1746 histOutSum += histOut[i] ;
1747 }
1748 }
1749
1750
1751 if (!expectedData && !extended) {
1752
1753 // Second pass for regular mode - Trim/Extend dataset to exact number of entries
1754
1755 // Calculate difference between what is generated so far and what is requested
1756 Int_t nEvtExtra = std::abs(Int_t(nEvents)-histOutSum) ;
1757 Int_t wgt = (histOutSum>nEvents) ? -1 : 1 ;
1758
1759 // Perform simple binned accept/reject procedure to get to exact event count
1760 std::size_t counter = 0;
1761 bool havePrintedInfo = false;
1762 while(nEvtExtra>0) {
1763
1764 Int_t ibinRand = RooRandom::randomGenerator()->Integer(hist->numEntries()) ;
1765 hist->get(ibinRand) ;
1766 double ranY = RooRandom::randomGenerator()->Uniform(histMax) ;
1767
1768 if (ranY<hist->weight()) {
1769 if (wgt==1) {
1770 histOut[ibinRand]++ ;
1771 } else {
1772 // If weight is negative, prior bin content must be at least 1
1773 if (histOut[ibinRand]>0) {
1774 histOut[ibinRand]-- ;
1775 } else {
1776 continue ;
1777 }
1778 }
1779 nEvtExtra-- ;
1780 }
1781
1782 if ((counter++ > 10*nEvents || nEvents > 1.E7) && !havePrintedInfo) {
1783 havePrintedInfo = true;
1784 coutP(Generation) << "RooAbsPdf::generateBinned(" << GetName() << ") Performing costly accept/reject sampling. If this takes too long, use "
1785 << "extended mode to speed up the process." << std::endl;
1786 }
1787 }
1788
1789 // Transfer working array to histogram
1790 for (int i=0 ; i<hist->numEntries() ; i++) {
1791 hist->get(i) ;
1792 hist->set(histOut[i],sqrt(1.0*histOut[i])) ;
1793 }
1794
1795 } else if (expectedData) {
1796
1797 // Second pass for expectedData mode -- Normalize to exact number of requested events
1798 // Minor difference may be present in first round due to difference between
1799 // bin average and bin integral in sampling bins
1800 double corr = nEvents/hist->sumEntries() ;
1801 for (int i=0 ; i<hist->numEntries() ; i++) {
1802 hist->get(i) ;
1803 hist->set(hist->weight()*corr,sqrt(hist->weight()*corr)) ;
1804 }
1805
1806 }
1807
1808 return RooFit::makeOwningPtr(std::move(hist));
1809}
1810
1811
1812
1813////////////////////////////////////////////////////////////////////////////////
1814/// Special generator interface for generation of 'global observables' -- for RooStats tools
1815
1820
1821namespace {
1822void removeRangeOverlap(std::vector<std::pair<double, double>>& ranges) {
1823 //Sort from left to right
1824 std::sort(ranges.begin(), ranges.end());
1825
1826 for (auto it = ranges.begin(); it != ranges.end(); ++it) {
1827 double& startL = it->first;
1828 double& endL = it->second;
1829
1830 for (auto innerIt = it+1; innerIt != ranges.end(); ++innerIt) {
1831 const double startR = innerIt->first;
1832 const double endR = innerIt->second;
1833
1834 if (startL <= startR && startR <= endL) {
1835 //Overlapping ranges, extend left one
1836 endL = std::max(endL, endR);
1837 *innerIt = make_pair(0., 0.);
1838 }
1839 }
1840 }
1841
1842 auto newEnd = std::remove_if(ranges.begin(), ranges.end(),
1843 [](const std::pair<double,double>& input){
1844 return input.first == input.second;
1845 });
1846 ranges.erase(newEnd, ranges.end());
1847}
1848}
1849
1850
1851////////////////////////////////////////////////////////////////////////////////
1852/// Plot (project) PDF on specified frame.
1853/// - If a PDF is plotted in an empty frame, it
1854/// will show a unit-normalized curve in the frame variable. When projecting a multi-
1855/// dimensional PDF onto the frame axis, hidden parameters are taken are taken at
1856/// their current value.
1857/// - If a PDF is plotted in a frame in which a dataset has already been plotted, it will
1858/// show a projection integrated over all variables that were present in the shown
1859/// dataset (except for the one on the x-axis). The normalization of the curve will
1860/// be adjusted to the event count of the plotted dataset. An informational message
1861/// will be printed for each projection step that is performed.
1862/// - If a PDF is plotted in a frame showing a dataset *after* a fit, the above happens,
1863/// but the PDF will be drawn and normalised only in the fit range. If this is not desired,
1864/// plotting and normalisation range can be overridden using Range() and NormRange() as
1865/// documented in the table below.
1866///
1867/// This function takes the following named arguments (for more arguments, see also
1868/// RooAbsReal::plotOn(RooPlot*,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,
1869/// const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,
1870/// const RooCmdArg&) const )
1871///
1872///
1873/// <table>
1874/// <tr><th> Type of argument <th> Controlling normalisation
1875/// <tr><td> `NormRange(const char* name)` <td> Calculate curve normalization w.r.t. specified range[s].
1876/// See the tutorial rf212_plottingInRanges_blinding.C
1877/// \note Setting a Range() by default also sets a NormRange() on the same range, meaning that the
1878/// PDF is plotted and normalised in the same range. Overriding this can be useful if the PDF was fit
1879/// in limited range[s] such as side bands, `NormRange("sidebandLeft,sidebandRight")`, but the PDF
1880/// should be drawn in the full range, `Range("")`.
1881///
1882/// <tr><td> `Normalization(double scale, ScaleType code)` <td> Adjust normalization by given scale factor.
1883/// Interpretation of number depends on code:
1884/// `RooAbsReal::Relative`: relative adjustment factor
1885/// `RooAbsReal::NumEvent`: scale to match given number of events.
1886///
1887/// <tr><th> Type of argument <th> Misc control
1888/// <tr><td> `Name(const chat* name)` <td> Give curve specified name in frame. Useful if curve is to be referenced later
1889/// <tr><td> `Asymmetry(const RooCategory& c)` <td> Show the asymmetry of the PDF in given two-state category
1890/// \f$ \frac{F(+)-F(-)}{F(+)+F(-)} \f$ rather than the PDF projection. Category must have two
1891/// states with indices -1 and +1 or three states with indices -1,0 and +1.
1892/// <tr><td> `ShiftToZero(bool flag)` <td> Shift entire curve such that lowest visible point is at exactly zero.
1893/// Mostly useful when plotting -log(L) or \f$ \chi^2 \f$ distributions
1894/// <tr><td> `AddTo(const char* name, double_t wgtSelf, double_t wgtOther)` <td> Create a projection of this PDF onto the x-axis, but
1895/// instead of plotting it directly, add it to an existing curve with given name (and relative weight factors).
1896/// <tr><td> `Components(const char* names)` <td> When plotting sums of PDFs, plot only the named components (*e.g.* only
1897/// the signal of a signal+background model).
1898/// <tr><td> `Components(const RooArgSet& compSet)` <td> As above, but pass a RooArgSet of the components themselves.
1899///
1900/// <tr><th> Type of argument <th> Projection control
1901/// <tr><td> `Slice(const RooArgSet& set)` <td> Override default projection behaviour by omitting observables listed
1902/// in set from the projection, i.e. by not integrating over these.
1903/// Slicing is usually only sensible in discrete observables, by e.g. creating a slice
1904/// of the PDF at the current value of the category observable.
1905/// <tr><td> `Slice(RooCategory& cat, const char* label)` <td> Override default projection behaviour by omitting the specified category
1906/// observable from the projection, i.e., by not integrating over all states of this category.
1907/// The slice is positioned at the given label value. Multiple Slice() commands can be given to specify slices
1908/// in multiple observables, e.g.
1909/// ```{.cpp}
1910/// pdf.plotOn(frame, Slice(tagCategory, "2tag"), Slice(jetCategory, "3jet"));
1911/// ```
1912/// <tr><td> `Project(const RooArgSet& set)` <td> Override default projection behaviour by projecting
1913/// over observables given in set, completely ignoring the default projection behavior. Advanced use only.
1914/// <tr><td> `ProjWData(const RooAbsData& d)` <td> Override default projection _technique_ (integration). For observables
1915/// present in given dataset projection of PDF is achieved by constructing an average over all observable
1916/// values in given set. Consult RooFit plotting tutorial for further explanation of meaning & use of this technique
1917/// <tr><td> `ProjWData(const RooArgSet& s, const RooAbsData& d)` <td> As above but only consider subset 's' of
1918/// observables in dataset 'd' for projection through data averaging
1919/// <tr><td> `ProjectionRange(const char* rn)` <td> When projecting the PDF onto the plot axis, it is usually integrated
1920/// over the full range of the invisible variables. The ProjectionRange overrides this.
1921/// This is useful if the PDF was fitted in a limited range in y, but it is now projected onto x. If
1922/// `ProjectionRange("<name of fit range>")` is passed, the projection is normalised correctly.
1923///
1924/// <tr><th> Type of argument <th> Plotting control
1925/// <tr><td> `LineStyle(Int_t style)` <td> Select line style by ROOT line style code, default is solid
1926/// <tr><td> `LineColor(Int_t color)` <td> Select line color by ROOT color code, default is blue
1927/// <tr><td> `LineWidth(Int_t width)` <td> Select line with in pixels, default is 3
1928/// <tr><td> `FillStyle(Int_t style)` <td> Select fill style, default is not filled. If a filled style is selected,
1929/// also use VLines() to add vertical downward lines at end of curve to ensure proper closure
1930/// <tr><td> `FillColor(Int_t color)` <td> Select fill color by ROOT color code
1931/// <tr><td> `Range(const char* name)` <td> Only draw curve in range defined by given name. Multiple comma-separated ranges can be given.
1932/// An empty string "" or `nullptr` means to use the default range of the variable.
1933/// <tr><td> `Range(double lo, double hi)` <td> Only draw curve in specified range
1934/// <tr><td> `VLines()` <td> Add vertical lines to y=0 at end points of curve
1935/// <tr><td> `Precision(double eps)` <td> Control precision of drawn curve w.r.t to scale of plot, default is 1e-3. A higher precision will
1936/// result in more and more densely spaced curve points. A negative precision value will disable
1937/// adaptive point spacing and restrict sampling to the grid point of points defined by the binning
1938/// of the plotted observable (recommended for expensive functions such as profile likelihoods)
1939/// <tr><td> `Invisible(bool flag)` <td> Add curve to frame, but do not display. Useful in combination AddTo()
1940/// <tr><td> `VisualizeError(const RooFitResult& fitres, double Z=1, bool linearMethod=true)`
1941/// <td> Visualize the uncertainty on the parameters, as given in fitres, at 'Z' sigma.
1942/// The linear method is fast but may not be accurate in the presence of strong correlations (~>0.9) and at Z>2 due to linear and Gaussian approximations made.
1943/// Intervals from the sampling method can be asymmetric, and may perform better in the presence of strong correlations, but may take (much) longer to calculate
1944/// \note To include the uncertainty from the expected number of events,
1945/// the Normalization() argument with `ScaleType` `RooAbsReal::RelativeExpected` has to be passed, e.g.
1946/// ```{.cpp}
1947/// pdf.plotOn(frame, VisualizeError(fitResult), Normalization(1.0, RooAbsReal::RelativeExpected));
1948/// ```
1949///
1950/// <tr><td> `VisualizeError(const RooFitResult& fitres, const RooArgSet& param, double Z=1, bool linearMethod=true)`
1951/// <td> Visualize the uncertainty on the subset of parameters 'param', as given in fitres, at 'Z' sigma
1952/// </table>
1953
1955{
1956
1957 // Pre-processing if p.d.f. contains a fit range and there is no command specifying one,
1958 // add a fit range as default range
1959 std::unique_ptr<RooCmdArg> plotRange;
1960 std::unique_ptr<RooCmdArg> normRange2;
1961 if (getStringAttribute("fitrange") && !cmdList.FindObject("Range") &&
1962 !cmdList.FindObject("RangeWithName")) {
1963 plotRange.reset(static_cast<RooCmdArg*>(RooFit::Range(getStringAttribute("fitrange")).Clone()));
1964 cmdList.Add(plotRange.get());
1965 }
1966
1967 if (getStringAttribute("fitrange") && !cmdList.FindObject("NormRange")) {
1968 normRange2.reset(static_cast<RooCmdArg*>(RooFit::NormRange(getStringAttribute("fitrange")).Clone()));
1969 cmdList.Add(normRange2.get());
1970 }
1971
1972 if (plotRange || normRange2) {
1973 coutI(Plotting) << "RooAbsPdf::plotOn(" << GetName() << ") p.d.f was fitted in a subrange and no explicit "
1974 << (plotRange?"Range()":"") << ((plotRange&&normRange2)?" and ":"")
1975 << (normRange2?"NormRange()":"") << " was specified. Plotting / normalising in fit range. To override, do one of the following"
1976 << "\n\t- Clear the automatic fit range attribute: <pdf>.removeStringAttribute(\"fitrange\");"
1977 << "\n\t- Explicitly specify the plotting range: Range(\"<rangeName>\")."
1978 << "\n\t- Explicitly specify where to compute the normalisation: NormRange(\"<rangeName>\")."
1979 << "\n\tThe default (full) range can be denoted with Range(\"\") / NormRange(\"\")."<< std::endl ;
1980 }
1981
1982 // Sanity checks
1983 if (plotSanityChecks(frame)) return frame ;
1984
1985 // Select the pdf-specific commands
1986 RooCmdConfig pc("RooAbsPdf::plotOn(" + std::string(GetName()) + ")");
1987 pc.defineDouble("scaleFactor","Normalization",0,1.0) ;
1988 pc.defineInt("scaleType","Normalization",0,Relative) ;
1989 pc.defineSet("compSet","SelectCompSet",0) ;
1990 pc.defineString("compSpec","SelectCompSpec",0) ;
1991 pc.defineObject("asymCat","Asymmetry",0) ;
1992 pc.defineDouble("rangeLo","Range",0,-999.) ;
1993 pc.defineDouble("rangeHi","Range",1,-999.) ;
1994 pc.defineString("rangeName","RangeWithName",0,"") ;
1995 pc.defineString("normRangeName","NormRange",0,"") ;
1996 pc.defineInt("rangeAdjustNorm","Range",0,0) ;
1997 pc.defineInt("rangeWNAdjustNorm","RangeWithName",0,0) ;
1998 pc.defineMutex("SelectCompSet","SelectCompSpec") ;
1999 pc.defineMutex("Range","RangeWithName") ;
2000 pc.allowUndefined() ; // unknowns may be handled by RooAbsReal
2001
2002 // Process and check varargs
2003 pc.process(cmdList) ;
2004 if (!pc.ok(true)) {
2005 return frame ;
2006 }
2007
2008 // Decode command line arguments
2009 ScaleType stype = (ScaleType) pc.getInt("scaleType") ;
2010 double scaleFactor = pc.getDouble("scaleFactor") ;
2011 const RooAbsCategoryLValue* asymCat = static_cast<const RooAbsCategoryLValue*>(pc.getObject("asymCat")) ;
2012 const char* compSpec = pc.getString("compSpec") ;
2013 const RooArgSet* compSet = pc.getSet("compSet");
2014 bool haveCompSel = ( (compSpec && strlen(compSpec)>0) || compSet) ;
2015
2016 // Suffix for curve name
2017 std::stringstream nameSuffix;
2018 if (compSpec && strlen(compSpec)>0) {
2019 nameSuffix << "_Comp[" << compSpec << "]";
2020 } else if (compSet) {
2021 nameSuffix << "_Comp[" << compSet->contentsString() << "]";
2022 }
2023
2024 // Remove PDF-only commands from command list
2025 RooCmdConfig::stripCmdList(cmdList,"SelectCompSet,SelectCompSpec") ;
2026
2027 // Adjust normalization, if so requested
2028 if (asymCat) {
2029 RooCmdArg cnsuffix("CurveNameSuffix", 0, 0, 0, 0, nameSuffix.str().c_str(), nullptr, nullptr, nullptr);
2030 cmdList.Add(&cnsuffix);
2031 return RooAbsReal::plotOn(frame, cmdList);
2032 }
2033
2034 // More sanity checks
2035 double nExpected(1) ;
2036 if (stype==RelativeExpected) {
2037 if (!canBeExtended()) {
2038 coutE(Plotting) << "RooAbsPdf::plotOn(" << GetName()
2039 << "): ERROR the 'Expected' scale option can only be used on extendable PDFs" << std::endl ;
2040 return frame ;
2041 }
2042 frame->updateNormVars(*frame->getPlotVar()) ;
2044 }
2045
2046 if (stype != Raw) {
2047
2048 if (frame->getFitRangeNEvt() && stype==Relative) {
2049
2050 bool hasCustomRange(false);
2051 bool adjustNorm(false);
2052
2053 std::vector<pair<double,double> > rangeLim;
2054
2055 // Retrieve plot range to be able to adjust normalization to data
2056 if (pc.hasProcessed("Range")) {
2057
2058 double rangeLo = pc.getDouble("rangeLo") ;
2059 double rangeHi = pc.getDouble("rangeHi") ;
2060 rangeLim.push_back(make_pair(rangeLo,rangeHi)) ;
2061 adjustNorm = pc.getInt("rangeAdjustNorm") ;
2063
2064 coutI(Plotting) << "RooAbsPdf::plotOn(" << GetName() << ") only plotting range ["
2065 << rangeLo << "," << rangeHi << "]" ;
2066 if (!pc.hasProcessed("NormRange")) {
2067 ccoutI(Plotting) << ", curve is normalized to data in " << (adjustNorm?"given":"full") << " range" << std::endl ;
2068 } else {
2069 ccoutI(Plotting) << std::endl ;
2070 }
2071
2072 nameSuffix << "_Range[" << rangeLo << "_" << rangeHi << "]";
2073
2074 } else if (pc.hasProcessed("RangeWithName")) {
2075
2076 for (const std::string& rangeNameToken : ROOT::Split(pc.getString("rangeName", "", false), ",")) {
2077 const char* thisRangeName = rangeNameToken.empty() ? nullptr : rangeNameToken.c_str();
2078 if (thisRangeName && !frame->getPlotVar()->hasRange(thisRangeName)) {
2079 coutE(Plotting) << "Range '" << rangeNameToken << "' not defined for variable '"
2080 << frame->getPlotVar()->GetName() << "'. Ignoring ..." << std::endl;
2081 continue;
2082 }
2083 rangeLim.push_back(frame->getPlotVar()->getRange(thisRangeName));
2084 }
2085 adjustNorm = pc.getInt("rangeWNAdjustNorm") ;
2087
2088 coutI(Plotting) << "RooAbsPdf::plotOn(" << GetName() << ") only plotting range '" << pc.getString("rangeName", "", false) << "'" ;
2089 if (!pc.hasProcessed("NormRange")) {
2090 ccoutI(Plotting) << ", curve is normalized to data in " << (adjustNorm?"given":"full") << " range" << std::endl ;
2091 } else {
2092 ccoutI(Plotting) << std::endl ;
2093 }
2094
2095 nameSuffix << "_Range[" << pc.getString("rangeName") << "]";
2096 }
2097 // Specification of a normalization range override those in a regular range
2098 if (pc.hasProcessed("NormRange")) {
2099 rangeLim.clear();
2100 for (const auto& rangeNameToken : ROOT::Split(pc.getString("normRangeName", "", false), ",")) {
2101 const char* thisRangeName = rangeNameToken.empty() ? nullptr : rangeNameToken.c_str();
2102 if (thisRangeName && !frame->getPlotVar()->hasRange(thisRangeName)) {
2103 coutE(Plotting) << "Range '" << rangeNameToken << "' not defined for variable '"
2104 << frame->getPlotVar()->GetName() << "'. Ignoring ..." << std::endl;
2105 continue;
2106 }
2107 rangeLim.push_back(frame->getPlotVar()->getRange(thisRangeName));
2108 }
2109 adjustNorm = true ;
2111 coutI(Plotting) << "RooAbsPdf::plotOn(" << GetName() << ") p.d.f. curve is normalized using explicit choice of ranges '" << pc.getString("normRangeName", "", false) << "'" << std::endl ;
2112
2113 nameSuffix << "_NormRange[" << pc.getString("rangeName") << "]";
2114 }
2115
2116 if (hasCustomRange && adjustNorm) {
2117 // If overlapping ranges were given, remove them now
2118 const std::size_t oldSize = rangeLim.size();
2120
2121 if (oldSize != rangeLim.size() && !pc.hasProcessed("NormRange")) {
2122 // User gave overlapping ranges. This leads to double-counting events and integrals, and must
2123 // therefore be avoided. If a NormRange has been given, the overlap is already gone.
2124 // It's safe to plot even with overlap now.
2125 coutE(Plotting) << "Requested plot/integration ranges overlap. For correct plotting, new ranges "
2126 "will be defined." << std::endl;
2127 auto plotVar = dynamic_cast<RooRealVar*>(frame->getPlotVar());
2128 assert(plotVar);
2129 std::string rangesNoOverlap;
2130 for (auto it = rangeLim.begin(); it != rangeLim.end(); ++it) {
2131 std::stringstream rangeName;
2132 rangeName << "Remove_overlap_range_" << it - rangeLim.begin();
2133 plotVar->setRange(rangeName.str().c_str(), it->first, it->second);
2134 if (!rangesNoOverlap.empty())
2135 rangesNoOverlap += ",";
2136 rangesNoOverlap += rangeName.str();
2137 }
2138
2139 auto rangeArg = static_cast<RooCmdArg*>(cmdList.FindObject("RangeWithName"));
2140 if (rangeArg) {
2141 rangeArg->setString(0, rangesNoOverlap.c_str());
2142 } else {
2143 plotRange = std::make_unique<RooCmdArg>(RooFit::Range(rangesNoOverlap.c_str()));
2144 cmdList.Add(plotRange.get());
2145 }
2146 }
2147
2148 double rangeNevt(0) ;
2149 for (const auto& riter : rangeLim) {
2150 double nevt= frame->getFitRangeNEvt(riter.first, riter.second);
2151 rangeNevt += nevt ;
2152 }
2153
2154 scaleFactor *= rangeNevt/nExpected ;
2155
2156 } else {
2157 // First, check if the PDF *can* be extended.
2158 if (this->canBeExtended()) {
2159 // If it can, get the expected events.
2160 const double nExp = expectedEvents(frame->getNormVars());
2161 if (nExp > 0) {
2162 // If the prediction is valid, use it for normalization.
2163 scaleFactor *= nExp / nExpected;
2164 } else {
2165 // If prediction is not valid (e.g. 0), fall back to data.
2166 scaleFactor *= frame->getFitRangeNEvt() / nExpected;
2167 }
2168 } else {
2169 // If the PDF can't be extended, just use the data.
2170 scaleFactor *= frame->getFitRangeNEvt() / nExpected;
2171 }
2172 }
2173 } else if (stype==RelativeExpected) {
2174 scaleFactor *= nExpected ;
2175 } else if (stype==NumEvent) {
2176 scaleFactor /= nExpected ;
2177 }
2178 scaleFactor *= frame->getFitRangeBinW() ;
2179 }
2180 frame->updateNormVars(*frame->getPlotVar()) ;
2181
2182 // Append overriding scale factor command at end of original command list
2183 RooCmdArg tmp = RooFit::Normalization(scaleFactor,Raw) ;
2184 tmp.setInt(1,1) ; // Flag this normalization command as created for internal use (so that VisualizeError can strip it)
2185 replaceOrAdd(cmdList, tmp);
2186
2187 // Was a component selected requested
2188 if (haveCompSel) {
2189
2190 // Get complete set of tree branch nodes
2193
2194 // Discard any non-RooAbsReal nodes
2195 for (const auto arg : branchNodeSet) {
2196 if (!dynamic_cast<RooAbsReal*>(arg)) {
2197 branchNodeSet.remove(*arg) ;
2198 }
2199 }
2200
2201 // Obtain direct selection
2202 std::unique_ptr<RooArgSet> dirSelNodes;
2203 if (compSet) {
2204 dirSelNodes = std::unique_ptr<RooArgSet>{branchNodeSet.selectCommon(*compSet)};
2205 } else {
2206 dirSelNodes = std::unique_ptr<RooArgSet>{branchNodeSet.selectByName(compSpec)};
2207 }
2208 if (!dirSelNodes->empty()) {
2209 coutI(Plotting) << "RooAbsPdf::plotOn(" << GetName() << ") directly selected PDF components: " << *dirSelNodes << std::endl ;
2210
2211 // Do indirect selection and activate both
2213 } else {
2214 if (compSet) {
2215 coutE(Plotting) << "RooAbsPdf::plotOn(" << GetName() << ") ERROR: component selection set " << *compSet << " does not match any components of p.d.f." << std::endl ;
2216 } else {
2217 coutE(Plotting) << "RooAbsPdf::plotOn(" << GetName() << ") ERROR: component selection expression '" << compSpec << "' does not select any components of p.d.f." << std::endl ;
2218 }
2219 return nullptr ;
2220 }
2221 }
2222
2223 RooCmdArg cnsuffix("CurveNameSuffix", 0, 0, 0, 0, nameSuffix.str().c_str(), nullptr, nullptr, nullptr);
2224 cmdList.Add(&cnsuffix);
2225
2227
2228 // Restore selection status ;
2229 if (haveCompSel) plotOnCompSelect(nullptr) ;
2230
2231 return ret ;
2232}
2233
2234
2235//_____________________________________________________________________________
2236/// Plot oneself on 'frame'. In addition to features detailed in RooAbsReal::plotOn(),
2237/// the scale factor for a PDF can be interpreted in three different ways. The interpretation
2238/// is controlled by ScaleType
2239/// ```
2240/// Relative - Scale factor is applied on top of PDF normalization scale factor
2241/// NumEvent - Scale factor is interpreted as a number of events. The surface area
2242/// under the PDF curve will match that of a histogram containing the specified
2243/// number of event
2244/// Raw - Scale factor is applied to the raw (projected) probability density.
2245/// Not too useful, option provided for completeness.
2246/// ```
2247// coverity[PASS_BY_VALUE]
2249{
2250
2251 // Sanity checks
2252 if (plotSanityChecks(frame)) return frame ;
2253
2254 // More sanity checks
2255 double nExpected(1) ;
2256 if (o.stype==RelativeExpected) {
2257 if (!canBeExtended()) {
2258 coutE(Plotting) << "RooAbsPdf::plotOn(" << GetName()
2259 << "): ERROR the 'Expected' scale option can only be used on extendable PDFs" << std::endl ;
2260 return frame ;
2261 }
2262 frame->updateNormVars(*frame->getPlotVar()) ;
2264 }
2265
2266 // Adjust normalization, if so requested
2267 if (o.stype != Raw) {
2268
2269 if (frame->getFitRangeNEvt() && o.stype==Relative) {
2270 // If non-default plotting range is specified, adjust number of events in fit range
2271 o.scaleFactor *= frame->getFitRangeNEvt()/nExpected ;
2272 } else if (o.stype==RelativeExpected) {
2273 o.scaleFactor *= nExpected ;
2274 } else if (o.stype==NumEvent) {
2275 o.scaleFactor /= nExpected ;
2276 }
2277 o.scaleFactor *= frame->getFitRangeBinW() ;
2278 }
2279 frame->updateNormVars(*frame->getPlotVar()) ;
2280
2281 return RooAbsReal::plotOn(frame,o) ;
2282}
2283
2284
2285
2286
2287////////////////////////////////////////////////////////////////////////////////
2288/// The following named arguments are supported
2289/// <table>
2290/// <tr><th> Type of CmdArg <th> Effect on parameter box
2291/// <tr><td> `Parameters(const RooArgSet& param)` <td> Only the specified subset of parameters will be shown. By default all non-constant parameters are shown.
2292/// <tr><td> `ShowConstants(bool flag)` <td> Also display constant parameters
2293/// <tr><td> `Format(const char* what,...)` <td> Parameter formatting options.
2294/// | Parameter | Format
2295/// | ---------------------- | --------------------------
2296/// | `const char* what` | Controls what is shown. "N" adds name (alternatively, "T" adds the title), "E" adds error, "A" shows asymmetric error, "U" shows unit, "H" hides the value
2297/// | `FixedPrecision(int n)`| Controls precision, set fixed number of digits
2298/// | `AutoPrecision(int n)` | Controls precision. Number of shown digits is calculated from error + n specified additional digits (1 is sensible default)
2299/// <tr><td> `Label(const chat* label)` <td> Add label to parameter box. Use `\n` for multi-line labels.
2300/// <tr><td> `Layout(double xmin, double xmax, double ymax)` <td> Specify relative position of left/right side of box and top of box.
2301/// Coordinates are given as position on the pad between 0 and 1.
2302/// The lower end of the box is calculated automatically from the number of lines in the box.
2303/// </table>
2304///
2305///
2306/// Example use:
2307/// ```
2308/// pdf.paramOn(frame, Label("fit result"), Format("NEU",AutoPrecision(1)) ) ;
2309/// ```
2310///
2311
2313 const RooCmdArg& arg3, const RooCmdArg& arg4, const RooCmdArg& arg5,
2314 const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8)
2315{
2316 // Stuff all arguments in a list
2318 cmdList.Add(const_cast<RooCmdArg*>(&arg1)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg2)) ;
2319 cmdList.Add(const_cast<RooCmdArg*>(&arg3)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg4)) ;
2320 cmdList.Add(const_cast<RooCmdArg*>(&arg5)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg6)) ;
2321 cmdList.Add(const_cast<RooCmdArg*>(&arg7)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg8)) ;
2322
2323 // Select the pdf-specific commands
2324 RooCmdConfig pc("RooAbsPdf::paramOn(" + std::string(GetName()) + ")");
2325 pc.defineString("label","Label",0,"") ;
2326 pc.defineDouble("xmin","Layout",0,0.65) ;
2327 pc.defineDouble("xmax","Layout",1,0.9) ;
2328 pc.defineInt("ymaxi","Layout",0,Int_t(0.9*10000)) ;
2329 pc.defineInt("showc","ShowConstants",0,0) ;
2330 pc.defineSet("params","Parameters",0,nullptr) ;
2331 pc.defineInt("dummy","FormatArgs",0,0) ;
2332
2333 // Process and check varargs
2334 pc.process(cmdList) ;
2335 if (!pc.ok(true)) {
2336 return frame ;
2337 }
2338
2339 auto formatCmd = static_cast<RooCmdArg const*>(cmdList.FindObject("FormatArgs")) ;
2340
2341 const char* label = pc.getString("label") ;
2342 double xmin = pc.getDouble("xmin") ;
2343 double xmax = pc.getDouble("xmax") ;
2344 double ymax = pc.getInt("ymaxi") / 10000. ;
2345 int showc = pc.getInt("showc") ;
2346
2347 // Decode command line arguments
2348 std::unique_ptr<RooArgSet> params{getParameters(frame->getNormVars())} ;
2349 if(RooArgSet* requestedParams = pc.getSet("params")) {
2350 params = std::unique_ptr<RooArgSet>{params->selectCommon(*requestedParams)};
2351 }
2352 paramOn(frame,*params,showc,label,xmin,xmax,ymax,formatCmd);
2353
2354 return frame ;
2355}
2356
2357
2358////////////////////////////////////////////////////////////////////////////////
2359/// Add a text box with the current parameter values and their errors to the frame.
2360/// Observables of this PDF appearing in the 'data' dataset will be omitted.
2361///
2362/// An optional label will be inserted if passed. Multi-line labels can be generated
2363/// by adding `\n` to the label string. Use 'sigDigits'
2364/// to modify the default number of significant digits printed. The 'xmin,xmax,ymax'
2365/// values specify the initial relative position of the text box in the plot frame.
2366
2367RooPlot* RooAbsPdf::paramOn(RooPlot* frame, const RooArgSet& params, bool showConstants, const char *label,
2368 double xmin, double xmax ,double ymax, const RooCmdArg* formatCmd)
2369{
2370
2371 // parse the options
2372 bool showLabel= (label != nullptr && strlen(label) > 0);
2373
2374 // calculate the box's size, adjusting for constant parameters
2375
2376 double ymin(ymax);
2377 double dy(0.06);
2378 for (const auto param : params) {
2379 auto var = static_cast<RooRealVar*>(param);
2380 if(showConstants || !var->isConstant()) ymin-= dy;
2381 }
2382
2383 std::string labelString = label;
2384 unsigned int numLines = std::count(labelString.begin(), labelString.end(), '\n') + 1;
2385 if (showLabel) ymin -= numLines * dy;
2386
2387 // create the box and set its options
2388 TPaveText *box= new TPaveText(xmin,ymax,xmax,ymin,"BRNDC");
2389 if(!box) return nullptr;
2390 box->SetName((std::string(GetName()) + "_paramBox").c_str());
2391 box->SetFillColor(0);
2392 box->SetBorderSize(0);
2393 box->SetTextAlign(12);
2394 box->SetTextSize(0.04F);
2395 box->SetFillStyle(0);
2396
2397 for (const auto param : params) {
2398 auto var = static_cast<const RooRealVar*>(param);
2399 if(var->isConstant() && !showConstants) continue;
2400
2401 box->AddText((formatCmd ? var->format(*formatCmd) : var->format(2, "NELU")).c_str());
2402 }
2403
2404 // add the optional label if specified
2405 if (showLabel) {
2406 for (const auto& line : ROOT::Split(label, "\n")) {
2407 box->AddText(line.c_str());
2408 }
2409 }
2410
2411 // Add box to frame
2412 frame->addObject(box) ;
2413
2414 return frame ;
2415}
2416
2417
2418
2419
2420////////////////////////////////////////////////////////////////////////////////
2421/// Return expected number of events from this p.d.f for use in extended
2422/// likelihood calculations. This default implementation returns zero
2423
2425{
2426 return 0 ;
2427}
2428
2429
2430
2431////////////////////////////////////////////////////////////////////////////////
2432/// Change global level of verbosity for p.d.f. evaluations
2433
2435{
2436 _verboseEval = stat ;
2437}
2438
2439
2440
2441////////////////////////////////////////////////////////////////////////////////
2442/// Return global level of verbosity for p.d.f. evaluations
2443
2445{
2446 return _verboseEval ;
2447}
2448
2449
2450
2451////////////////////////////////////////////////////////////////////////////////
2452/// Destructor of normalization cache element. If this element
2453/// provides the 'current' normalization stored in RooAbsPdf::_norm
2454/// zero _norm pointer here before object pointed to is deleted here
2455
2457{
2458 // Zero _norm pointer in RooAbsPdf if it is points to our cache payload
2459 if (_owner) {
2460 RooAbsPdf* pdfOwner = static_cast<RooAbsPdf*>(_owner) ;
2461 if (pdfOwner->_norm == _norm.get()) {
2462 pdfOwner->_norm = nullptr ;
2463 }
2464 }
2465}
2466
2467
2468
2469////////////////////////////////////////////////////////////////////////////////
2470/// Return a p.d.f that represent a projection of this p.d.f integrated over given observables
2471
2473{
2474 // Construct name for new object
2475 std::string name = std::string{GetName()} + "_Proj[" + RooHelpers::getColonSeparatedNameString(iset, ',') + "]";
2476
2477 // Return projected p.d.f.
2478 return new RooProjectedPdf(name.c_str(),name.c_str(),*this,iset) ;
2479}
2480
2481
2482
2483////////////////////////////////////////////////////////////////////////////////
2484/// Create a cumulative distribution function of this p.d.f in terms
2485/// of the observables listed in iset. If no nset argument is given
2486/// the c.d.f normalization is constructed over the integrated
2487/// observables, so that its maximum value is precisely 1. It is also
2488/// possible to choose a different normalization for
2489/// multi-dimensional p.d.f.s: eg. for a pdf f(x,y,z) one can
2490/// construct a partial cdf c(x,y) that only when integrated itself
2491/// over z results in a maximum value of 1. To construct such a cdf pass
2492/// z as argument to the optional nset argument
2493
2498
2499
2500
2501////////////////////////////////////////////////////////////////////////////////
2502/// Create an object that represents the integral of the function over one or more observables listed in `iset`.
2503/// The actual integration calculation is only performed when the return object is evaluated. The name
2504/// of the integral object is automatically constructed from the name of the input function, the variables
2505/// it integrates and the range integrates over
2506///
2507/// The following named arguments are accepted
2508/// | Type of CmdArg | Effect on CDF
2509/// | ---------------------|-------------------
2510/// | SupNormSet(const RooArgSet&) | Observables over which should be normalized _in addition_ to the integration observables
2511/// | ScanNumCdf() | Apply scanning technique if cdf integral involves numeric integration [ default ]
2512/// | ScanAllCdf() | Always apply scanning technique
2513/// | ScanNoCdf() | Never apply scanning technique
2514/// | ScanParameters(Int_t nbins, Int_t intOrder) | Parameters for scanning technique of making CDF: number of sampled bins and order of interpolation applied on numeric cdf
2515
2517 const RooCmdArg& arg3, const RooCmdArg& arg4, const RooCmdArg& arg5,
2518 const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8)
2519{
2520 // Define configuration for this method
2521 RooCmdConfig pc("RooAbsReal::createCdf(" + std::string(GetName()) + ")");
2522 pc.defineSet("supNormSet","SupNormSet",0,nullptr) ;
2523 pc.defineInt("numScanBins","ScanParameters",0,1000) ;
2524 pc.defineInt("intOrder","ScanParameters",1,2) ;
2525 pc.defineInt("doScanNum","ScanNumCdf",0,1) ;
2526 pc.defineInt("doScanAll","ScanAllCdf",0,0) ;
2527 pc.defineInt("doScanNon","ScanNoCdf",0,0) ;
2528 pc.defineMutex("ScanNumCdf","ScanAllCdf","ScanNoCdf") ;
2529
2530 // Process & check varargs
2532 if (!pc.ok(true)) {
2533 return nullptr ;
2534 }
2535
2536 // Extract values from named arguments
2537 const RooArgSet* snset = pc.getSet("supNormSet",nullptr);
2538 RooArgSet nset ;
2539 if (snset) {
2540 nset.add(*snset) ;
2541 }
2542 Int_t numScanBins = pc.getInt("numScanBins") ;
2543 Int_t intOrder = pc.getInt("intOrder") ;
2544 Int_t doScanNum = pc.getInt("doScanNum") ;
2545 Int_t doScanAll = pc.getInt("doScanAll") ;
2546 Int_t doScanNon = pc.getInt("doScanNon") ;
2547
2548 // If scanning technique is not requested make integral-based cdf and return
2549 if (doScanNon) {
2550 return createIntRI(iset,nset) ;
2551 }
2552 if (doScanAll) {
2553 return createScanCdf(iset,nset,numScanBins,intOrder) ;
2554 }
2555 if (doScanNum) {
2556 std::unique_ptr<RooAbsReal> tmp{createIntegral(iset)} ;
2557 Int_t isNum= !static_cast<RooRealIntegral&>(*tmp).numIntRealVars().empty();
2558
2559 if (isNum) {
2560 coutI(NumericIntegration) << "RooAbsPdf::createCdf(" << GetName() << ") integration over observable(s) " << iset << " involves numeric integration," << std::endl
2561 << " constructing cdf though numeric integration of sampled pdf in " << numScanBins << " bins and applying order "
2562 << intOrder << " interpolation on integrated histogram." << std::endl
2563 << " To override this choice of technique use argument ScanNone(), to change scan parameters use ScanParameters(nbins,order) argument" << std::endl ;
2564 }
2565
2567 }
2568 return nullptr ;
2569}
2570
2572{
2573 string name = string(GetName()) + "_NUMCDF_" + integralNameSuffix(iset,&nset).Data() ;
2574 RooRealVar* ivar = static_cast<RooRealVar*>(iset.first()) ;
2575 ivar->setBins(numScanBins,"numcdf") ;
2576 auto ret = std::make_unique<RooNumCdf>(name.c_str(),name.c_str(),*this,*ivar,"numcdf");
2577 ret->setInterpolationOrder(intOrder) ;
2578 return RooFit::makeOwningPtr<RooAbsReal>(std::move(ret));
2579}
2580
2581
2582
2583
2584////////////////////////////////////////////////////////////////////////////////
2585/// This helper function finds and collects all constraints terms of all component p.d.f.s
2586/// and returns a RooArgSet with all those terms.
2587
2588std::unique_ptr<RooArgSet>
2590{
2591 RooArgSet constraints;
2593
2594 std::unique_ptr<RooArgSet> comps(getComponents());
2595 for (const auto arg : *comps) {
2596 auto pdf = dynamic_cast<const RooAbsPdf*>(arg) ;
2597 if (pdf && !constraints.find(pdf->GetName())) {
2598 if (auto compRet = pdf->getConstraints(observables,constrainedParams, pdfParams)) {
2599 constraints.add(*compRet,false) ;
2600 }
2601 }
2602 }
2603
2605
2606 // Strip any constraints that are completely decoupled from the other product terms
2607 auto finalConstraints = std::make_unique<RooArgSet>("AllConstraints");
2608 for(auto * pdf : static_range_cast<RooAbsPdf*>(constraints)) {
2609
2610 RooArgSet tmp;
2611 pdf->getParameters(nullptr, tmp);
2612 conParams.add(tmp,true) ;
2613
2614 if (pdf->dependsOnValue(pdfParams) || !stripDisconnected) {
2615 finalConstraints->add(*pdf) ;
2616 } else {
2617 coutI(Minimization) << "RooAbsPdf::getAllConstraints(" << GetName() << ") omitting term " << pdf->GetName()
2618 << " as constraint term as it does not share any parameters with the other pdfs in product. "
2619 << "To force inclusion in likelihood, add an explicit Constrain() argument for the target parameter" << std::endl ;
2620 }
2621 }
2622
2623 // Now remove from constrainedParams all parameters that occur exclusively in constraint term and not in regular pdf term
2624
2626 conParams.selectCommon(constrainedParams, cexl);
2627 cexl.remove(pdfParams,true,true) ;
2628 constrainedParams.remove(cexl,true,true) ;
2629
2630 return finalConstraints ;
2631}
2632
2633
2634////////////////////////////////////////////////////////////////////////////////
2635/// Returns the default numeric MC generator configuration for all RooAbsReals
2636
2641
2642
2643////////////////////////////////////////////////////////////////////////////////
2644/// Returns the specialized integrator configuration for _this_ RooAbsReal.
2645/// If this object has no specialized configuration, a null pointer is returned
2646
2651
2652
2653
2654////////////////////////////////////////////////////////////////////////////////
2655/// Returns the specialized integrator configuration for _this_ RooAbsReal.
2656/// If this object has no specialized configuration, a null pointer is returned,
2657/// unless createOnTheFly is true in which case a clone of the default integrator
2658/// configuration is created, installed as specialized configuration, and returned
2659
2661{
2663 _specGeneratorConfig = std::make_unique<RooNumGenConfig>(*defaultGeneratorConfig()) ;
2664 }
2665 return _specGeneratorConfig.get();
2666}
2667
2668
2669
2670////////////////////////////////////////////////////////////////////////////////
2671/// Return the numeric MC generator configuration used for this object. If
2672/// a specialized configuration was associated with this object, that configuration
2673/// is returned, otherwise the default configuration for all RooAbsReals is returned
2674
2676{
2677 const RooNumGenConfig *config = specialGeneratorConfig();
2678 return config ? config : defaultGeneratorConfig();
2679}
2680
2681
2682
2683////////////////////////////////////////////////////////////////////////////////
2684/// Set the given configuration as default numeric MC generator
2685/// configuration for this object
2686
2688{
2689 _specGeneratorConfig = std::make_unique<RooNumGenConfig>(config);
2690}
2691
2692
2693
2694////////////////////////////////////////////////////////////////////////////////
2695/// Remove the specialized numeric MC generator configuration associated
2696/// with this object
2697
2702
2704
2705
2706////////////////////////////////////////////////////////////////////////////////
2707
2709 bool extended, bool randProto, bool resampleProto, TString dsetName, bool init) :
2710 _genContext(context), _whatVars(whatVars), _protoData(protoData), _nGen(nGen), _extended(extended),
2711 _randProto(randProto), _resampleProto(resampleProto), _dsetName(dsetName), _init(init)
2712{
2713}
2714
2715
2716namespace {
2717
2718void sterilizeClientCaches(RooAbsArg & arg) {
2719 auto const& clients = arg.clients();
2720 for(std::size_t iClient = 0; iClient < clients.size(); ++iClient) {
2721
2722 const std::size_t oldClientsSize = clients.size();
2723 RooAbsArg* client = clients[iClient];
2724
2725 for(int iCache = 0; iCache < client->numCaches(); ++iCache) {
2726 if(auto cacheMgr = dynamic_cast<RooObjCacheManager*>(client->getCache(iCache))) {
2727 cacheMgr->sterilize();
2728 }
2729 }
2730
2731 // It can happen that the objects cached by the client are also clients of
2732 // the arg itself! In that case, the position of the client in the client
2733 // list might have changed, and we need to find the new index.
2734 if(clients.size() != oldClientsSize) {
2735 auto clientIter = std::find(clients.begin(), clients.end(), client);
2736 if(clientIter == clients.end()) {
2737 throw std::runtime_error("After a clients caches were cleared, the client was gone! This should not happen.");
2738 }
2739 iClient = std::distance(clients.begin(), clientIter);
2740 }
2741 }
2742}
2743
2744} // namespace
2745
2746
2747////////////////////////////////////////////////////////////////////////////////
2748
2750{
2752
2753 // the stuff that the clients have cached may depend on the normalization range
2754 sterilizeClientCaches(*this);
2755
2756 if (_norm) {
2758 _norm = nullptr ;
2759 }
2760}
2761
2762
2763////////////////////////////////////////////////////////////////////////////////
2764
2766{
2768
2769 // the stuff that the clients have cached may depend on the normalization range
2770 sterilizeClientCaches(*this);
2771
2772 if (_norm) {
2774 _norm = nullptr ;
2775 }
2776}
2777
2778
2779////////////////////////////////////////////////////////////////////////////////
2780/// Hook function intercepting redirectServer calls. Discard current
2781/// normalization object if any server is redirected
2782
2784 bool nameChange, bool isRecursiveStep)
2785{
2786 // If servers are redirected, the cached normalization integrals and
2787 // normalization sets are most likely invalid.
2789
2790 // Object is own by _normCacheManager that will delete object as soon as cache
2791 // is sterilized by server redirect
2792 _norm = nullptr ;
2793
2794 // Similar to the situation with the normalization integral above: if a
2795 // server is redirected, the cached normalization set might not point to
2796 // the right observables anymore. We need to reset it.
2797 setActiveNormSet(nullptr);
2799}
2800
2801
2802std::unique_ptr<RooAbsArg>
2804{
2805 if (normSet.empty() || selfNormalized()) {
2807 }
2808 std::unique_ptr<RooAbsPdf> pdfClone(static_cast<RooAbsPdf *>(this->Clone()));
2810
2811 auto newArg = std::make_unique<RooFit::Detail::RooNormalizedPdf>(*pdfClone, normSet);
2812
2813 // The direct servers are this pdf and the normalization integral, which
2814 // don't need to be compiled further.
2815 for (RooAbsArg *server : newArg->servers()) {
2816 ctx.markAsCompiled(*server);
2817 }
2818 ctx.markAsCompiled(*newArg);
2819 newArg->addOwnedComponents(std::move(pdfClone));
2820 return newArg;
2821}
2822
2823/// Returns an object that represents the expected number of events for a given
2824/// normalization set, similar to how createIntegral() returns an object that
2825/// returns the integral. This is used to build the computation graph for the
2826/// final likelihood.
2827std::unique_ptr<RooAbsReal> RooAbsPdf::createExpectedEventsFunc(const RooArgSet * /*nset*/) const
2828{
2829 std::stringstream errMsg;
2830 errMsg << "The pdf \"" << GetName() << "\" of type " << ClassName()
2831 << " did not overload RooAbsPdf::createExpectedEventsFunc()!";
2832 coutE(InputArguments) << errMsg.str() << std::endl;
2833 return nullptr;
2834}
#define e(i)
Definition RSha256.hxx:103
#define coutI(a)
#define cxcoutI(a)
#define cxcoutD(a)
#define coutP(a)
#define oocoutW(o, a)
#define coutW(a)
#define coutE(a)
#define ccoutI(a)
#define ccoutD(a)
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
static void indent(ostringstream &buf, int indent_level)
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 input
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
char name[80]
Definition TGX11.cxx:110
float xmin
float ymin
float xmax
float ymax
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2495
const_iterator begin() const
const_iterator end() const
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:76
void clearValueAndShapeDirty() const
Definition RooAbsArg.h:535
void Print(Option_t *options=nullptr) const override
Print the object to the defaultPrintStream().
Definition RooAbsArg.h:238
bool dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=nullptr, bool valueOnly=false) const
Test whether we depend on (ie, are served by) any object in the specified collection.
RooFit::OwningPtr< RooArgSet > getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
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.
const Text_t * getStringAttribute(const Text_t *key) const
Get string attribute mapped under key 'key'.
virtual std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const
RooFit::OwningPtr< RooArgSet > getComponents() const
Create a RooArgSet with all components (branch nodes) of the expression tree headed by this object.
bool getAttribute(const Text_t *name) const
Check if a named attribute is set. By default, all attributes are unset.
RooAbsCache * getCache(Int_t index) const
Return registered cache object by index.
const RefCountList_t & clients() const
List of all clients of this object.
Definition RooAbsArg.h:137
bool isValueDirty() const
Definition RooAbsArg.h:356
void setProxyNormSet(const RooArgSet *nset)
Forward a change in the cached normalization argset to all the registered proxies.
void branchNodeServerList(RooAbsCollection *list, const RooAbsArg *arg=nullptr, bool recurseNonDerived=false) const
Fill supplied list with all branch nodes of the arg tree starting with ourself as top node.
TObject * Clone(const char *newname=nullptr) const override
Make a clone of an object using the Streamer facility.
Definition RooAbsArg.h:88
RefCountList_t _serverList
Definition RooAbsArg.h:564
Int_t numCaches() const
Return number of registered caches.
RooAbsArg * findServer(const char *name) const
Return server of this with name name. Returns nullptr if not found.
Definition RooAbsArg.h:147
RooAbsArg * _owner
! Pointer to owning RooAbsArg
Abstract base class for objects that represent a discrete value that can be set from the outside,...
Abstract container object that can hold multiple RooAbsArg objects.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
RooAbsArg * find(const char *name) const
Find object with given name in list.
void Print(Option_t *options=nullptr) const override
This method must be overridden when a class wants to print itself.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
Abstract base class for generator contexts of RooAbsPdf objects.
virtual RooDataSet * generate(double nEvents=0, bool skipInit=false, bool extendedMode=false)
Generate the specified number of events with nEvents>0 and and return a dataset containing the genera...
bool isValid() const
virtual void setProtoDataOrder(Int_t *lut)
Set the traversal order of prototype data to that in the lookup tables passed as argument.
Normalization set with for above integral.
Definition RooAbsPdf.h:318
std::unique_ptr< RooAbsReal > _norm
Definition RooAbsPdf.h:323
~CacheElem() override
Destructor of normalization cache element.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
virtual bool syncNormalization(const RooArgSet *dset, bool adjustProxies=true) const
Verify that the normalization integral cached with this PDF is valid for given set of normalization o...
double getNorm(const RooArgSet &nset) const
Get normalisation term needed to normalise the raw values returned by getVal().
Definition RooAbsPdf.h:191
std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const override
RooObjCacheManager _normMgr
Definition RooAbsPdf.h:325
std::unique_ptr< RooNumGenConfig > _specGeneratorConfig
! MC generator configuration specific for this object
Definition RooAbsPdf.h:336
double getValV(const RooArgSet *set=nullptr) const override
Return current value, normalized by integrating over the observables in nset.
virtual std::unique_ptr< RooFitResult > fitToImpl(RooAbsData &data, const RooLinkedList &cmdList)
Protected implementation of the likelihood fitting routine.
virtual void generateEvent(Int_t code)
Interface for generation of an event using the algorithm corresponding to the specified code.
RooFit::OwningPtr< RooAbsReal > createScanCdf(const RooArgSet &iset, const RooArgSet &nset, Int_t numScanBins, Int_t intOrder)
void setGeneratorConfig()
Remove the specialized numeric MC generator configuration associated with this object.
virtual void resetErrorCounters(Int_t resetValue=10)
Reset error counter to given value, limiting the number of future error messages for this pdf to 'res...
std::unique_ptr< RooArgSet > getAllConstraints(const RooArgSet &observables, RooArgSet &constrainedParams, bool stripDisconnected=true) const
This helper function finds and collects all constraints terms of all component p.d....
static int verboseEval()
Return global level of verbosity for p.d.f. evaluations.
RooFit::OwningPtr< RooAbsReal > createCdf(const RooArgSet &iset, const RooArgSet &nset=RooArgSet())
Create a cumulative distribution function of this p.d.f in terms of the observables listed in iset.
bool isActiveNormSet(RooArgSet const *normSet) const
Checks if normSet is the currently active normalization set of this PDF, meaning is exactly the same ...
Definition RooAbsPdf.h:295
virtual double expectedEvents(const RooArgSet *nset) const
Return expected number of events to be used in calculation of extended likelihood.
virtual RooAbsGenContext * binnedGenContext(const RooArgSet &vars, bool verbose=false) const
Return a binned generator context.
TString _normRange
Normalization range.
Definition RooAbsPdf.h:338
virtual bool isDirectGenSafe(const RooAbsArg &arg) const
Check if given observable can be safely generated using the pdfs internal generator mechanism (if tha...
Int_t * randomizeProtoOrder(Int_t nProto, Int_t nGen, bool resample=false) const
Return lookup table with randomized order for nProto prototype events.
void setNormRange(const char *rangeName)
~RooAbsPdf() override
Destructor.
RooArgSet const * _normSet
Normalization integral (owned by _normMgr)
Definition RooAbsPdf.h:316
RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}, const RooCmdArg &arg9={}, const RooCmdArg &arg10={}) const override
Helper calling plotOn(RooPlot*, RooLinkedList&) const.
Definition RooAbsPdf.h:116
RooNumGenConfig * specialGeneratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
virtual bool selfNormalized() const
Shows if a PDF is self-normalized, which means that no attempt is made to add a normalization term.
Definition RooAbsPdf.h:203
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Print multi line detailed information of this RooAbsPdf.
virtual double getCorrection() const
This function returns the penalty term.
Int_t _traceCount
Number of traces remaining to print.
Definition RooAbsPdf.h:331
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition RooAbsPdf.h:214
RooAbsReal * _norm
Definition RooAbsPdf.h:315
GenSpec * prepareMultiGen(const RooArgSet &whatVars, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={})
Prepare GenSpec configuration object for efficient generation of multiple datasets from identical spe...
void setTraceCounter(Int_t value, bool allNodes=false)
Reset trace counter to given value, limiting the number of future trace messages for this pdf to 'val...
Int_t _errorCount
Number of errors remaining to print.
Definition RooAbsPdf.h:330
@ CanNotBeExtended
Definition RooAbsPdf.h:208
virtual std::unique_ptr< RooAbsReal > createExpectedEventsFunc(const RooArgSet *nset) const
Returns an object that represents the expected number of events for a given normalization set,...
Int_t _negCount
Number of negative probabilities remaining to print.
Definition RooAbsPdf.h:332
virtual RooPlot * paramOn(RooPlot *frame, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={})
Add a box with parameter values (and errors) to the specified frame.
RooFit::OwningPtr< RooDataSet > generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={})
See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,...
Definition RooAbsPdf.h:49
virtual const RooAbsReal * getNormObj(const RooArgSet *set, const RooArgSet *iset, const TNamed *rangeName=nullptr) const
Return pointer to RooAbsReal object that implements calculation of integral over observables iset in ...
void setActiveNormSet(RooArgSet const *normSet) const
Setter for the _normSet member, which should never be set directly.
Definition RooAbsPdf.h:280
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Analytical integral with normalization (see RooAbsReal::analyticalIntegralWN() for further informatio...
void setNormRangeOverride(const char *rangeName)
virtual RooFit::OwningPtr< RooDataSet > generateSimGlobal(const RooArgSet &whatVars, Int_t nEvents)
Special generator interface for generation of 'global observables' – for RooStats tools.
double normalizeWithNaNPacking(double rawVal, double normVal) const
virtual RooAbsGenContext * autoGenContext(const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool verbose=false, bool autoBinned=true, const char *binnedTag="") const
const RooNumGenConfig * getGeneratorConfig() const
Return the numeric MC generator configuration used for this object.
virtual void initGenerator(Int_t code)
Interface for one-time initialization to setup the generator for the specified code.
virtual ExtendMode extendMode() const
Returns ability of PDF to provide extended likelihood terms.
Definition RooAbsPdf.h:212
RooAbsPdf()
Default constructor.
virtual RooFit::OwningPtr< RooDataHist > generateBinned(const RooArgSet &whatVars, double nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}) const
As RooAbsPdf::generateBinned(const RooArgSet&, const RooCmdArg&,const RooCmdArg&, const RooCmdArg&,...
Definition RooAbsPdf.h:102
bool traceEvalPdf(double value) const
Check that passed value is positive and not 'not-a-number'.
static RooNumGenConfig * defaultGeneratorConfig()
Returns the default numeric MC generator configuration for all RooAbsReals.
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursiveStep) override
The cache manager.
void printValue(std::ostream &os) const override
Print value of p.d.f, also print normalization integral that was last used, if any.
virtual std::unique_ptr< RooAbsReal > createNLLImpl(RooAbsData &data, const RooLinkedList &cmdList)
Protected implementation of the NLL creation routine.
void logBatchComputationErrors(std::span< const double > &outputs, std::size_t begin) const
Scan through outputs and fix+log all nans and negative values.
virtual RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool verbose=false) const
Interface function to create a generator context from a p.d.f.
void getLogProbabilities(std::span< const double > pdfValues, double *output) const
static TString _normRangeOverride
Definition RooAbsPdf.h:339
static Int_t _verboseEval
Definition RooAbsPdf.h:310
double extendedTerm(double sumEntries, double expected, double sumEntriesW2=0.0, bool doOffset=false) const
virtual Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
virtual RooAbsPdf * createProjection(const RooArgSet &iset)
Return a p.d.f that represent a projection of this p.d.f integrated over given observables.
virtual double getLogVal(const RooArgSet *set=nullptr) const
Return the log of the current value with given normalization An error message is printed if the argum...
bool hasRange(const char *name) const override
Check if variable has a binning with given name.
std::pair< double, double > getRange(const char *name=nullptr) const
Get low and high bound of the variable.
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
RooDataHist * fillDataHist(RooDataHist *hist, const RooArgSet *nset, double scaleFactor, bool correctForBinVolume=false, bool showProgress=false) const
Fill a RooDataHist with values sampled from this function at the bin centers.
void plotOnCompSelect(RooArgSet *selNodes) const
Helper function for plotting of composite p.d.fs.
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:107
bool plotSanityChecks(RooPlot *frame) const
Utility function for plotOn(), perform general sanity check on frame to ensure safe plotting operatio...
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Structure printing.
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursiveStep) override
Function that is called at the end of redirectServers().
double _value
Cache for current value of object.
Definition RooAbsReal.h:539
virtual double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
TString integralNameSuffix(const RooArgSet &iset, const RooArgSet *nset=nullptr, const char *rangeName=nullptr, bool omitEmpty=false) const
Construct string with unique suffix name to give to integral object that encodes integrated observabl...
virtual double evaluate() const =0
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
void logEvalError(const char *message, const char *serverValueString=nullptr) const
Log evaluation error message.
const RooNumIntConfig * getIntegratorConfig() const
Return the numeric integration configuration used for this object.
virtual bool isBinnedDistribution(const RooArgSet &) const
Tests if the distribution is binned. Unless overridden by derived classes, this always returns false.
Definition RooAbsReal.h:352
RooFit::OwningPtr< RooAbsReal > createIntRI(const RooArgSet &iset, const RooArgSet &nset={})
Utility function for createRunningIntegral.
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}, const RooCmdArg &arg9={}, const RooCmdArg &arg10={}) const
Plot (project) PDF on specified frame.
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}) const
Create an object that represents the integral of the function over one or more observables listed in ...
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
Efficient implementation of the generator context specific for binned pdfs.
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=nullptr)
Setter function without integration set.
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=nullptr, const TNamed *isetRangeName=nullptr)
Getter function without integration set.
Implementation of RooAbsCachedReal that can cache any external RooAbsReal input function provided in ...
Named container for two doubles, two integers two object points and three string pointers that can be...
Definition RooCmdArg.h:26
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...
static void stripCmdList(RooLinkedList &cmdList, const char *cmdsToPurge)
Utility function that strips command names listed (comma separated) in cmdsToPurge from cmdList.
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'.
TObject * getObject(const char *name, TObject *obj=nullptr) const
Return TObject property registered with name 'name'.
Container class to hold unbinned data.
Definition RooDataSet.h:32
void markAsCompiled(RooAbsArg &arg) const
void compileServers(RooAbsArg &arg, RooArgSet const &normSet)
Implements a universal generator context for all RooAbsPdf classes that do not have or need a special...
Switches the message service to a different level while the instance is alive.
Definition RooHelpers.h:37
Collection class for internal use, storing a collection of RooAbsArg pointers in a doubly linked list...
static const char * str(const TNamed *ptr)
Return C++ string corresponding to given TNamed pointer.
Definition RooNameReg.h:39
Holds the configuration parameters of the various numeric integrators used by RooRealIntegral.
static RooNumGenConfig & defaultConfig()
Return reference to instance of default numeric integrator configuration object.
Implementation of a RooCacheManager<RooAbsCacheElement> that specializes in the storage of cache elem...
void sterilize() override
Clear the cache payload but retain slot mapping w.r.t to normalization and integration sets.
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
void addObject(TObject *obj, Option_t *drawOptions="", bool invisible=false)
Add a generic object to this plot.
Definition RooPlot.cxx:325
double getFitRangeNEvt() const
Return the number of events in the fit range.
Definition RooPlot.h:139
const RooArgSet * getNormVars() const
Definition RooPlot.h:146
RooAbsRealLValue * getPlotVar() const
Definition RooPlot.h:137
void updateNormVars(const RooArgSet &vars)
Install the given set of observables are reference normalization variables for this frame.
Definition RooPlot.cxx:310
double getFitRangeBinW() const
Return the bin width that is being used to normalise the PDF.
Definition RooPlot.h:142
virtual void printStream(std::ostream &os, Int_t contents, StyleOption style, TString indent="") const
Print description of object on ostream, printing contents set by contents integer,...
A RooAbsPdf implementation that represent a projection of a given input p.d.f and the object returned...
static UInt_t integer(UInt_t max, TRandom *generator=randomGenerator())
Return an integer uniformly distributed from [0,n-1].
Definition RooRandom.cxx:95
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition RooRandom.cxx:47
Performs hybrid numerical/analytical integrals of RooAbsReal objects.
const RooArgSet & numIntRealVars() const
Variable that can be changed from the outside.
Definition RooRealVar.h:37
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:50
Mother of all ROOT objects.
Definition TObject.h:41
virtual const char * GetName() const
Returns name of object.
Definition TObject.cxx:457
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:226
A Pave (see TPave) with text, lines or/and boxes inside.
Definition TPaveText.h:21
Basic string class.
Definition TString.h:138
Ssiz_t Length() const
Definition TString.h:425
const char * Data() const
Definition TString.h:384
TLine * line
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition fillpatterns.C:1
RooCmdArg SupNormSet(const RooArgSet &nset)
RooCmdArg NormRange(const char *rangeNameList)
RooCmdArg Range(const char *rangeName, bool adjustNorm=true)
RooCmdArg Normalization(double scaleFactor)
std::vector< std::string > Split(std::string_view str, std::string_view delims, bool skipEmpty=false)
Splits a string at each character in delims.
T * OwningPtr
An alias for raw pointers for indicating that the return type of a RooFit function is an owning point...
Definition Config.h:35
OwningPtr< T > makeOwningPtr(std::unique_ptr< T > &&ptr)
Internal helper to turn a std::unique_ptr<T> into an OwningPtr.
Definition Config.h:40
RooArgSet selectFromArgSet(RooArgSet const &, std::string const &names)
std::string getColonSeparatedNameString(RooArgSet const &argSet, char delim=':')
Bool_t IsNaN(Double_t x)
Definition TMath.h:903
Double_t QuietNaN()
Returns a quiet NaN as defined by IEEE 754.
Definition TMath.h:913
static double packFloatIntoNaN(float payload)
Pack float into mantissa of a NaN.
static void output()