Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooAbsTestStatistic.cxx
Go to the documentation of this file.
1/// \cond ROOFIT_INTERNAL
2
3/*****************************************************************************
4 * Project: RooFit *
5 * Package: RooFitCore *
6 * @(#)root/roofitcore:$Id$
7 * Authors: *
8 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
9 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
10 * *
11 * Copyright (c) 2000-2005, Regents of the University of California *
12 * and Stanford University. All rights reserved. *
13 * *
14 * Redistribution and use in source and binary forms, *
15 * with or without modification, are permitted according to the terms *
16 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
17 *****************************************************************************/
18
19/**
20\file RooAbsTestStatistic.cxx
21\class RooAbsTestStatistic
22\ingroup Roofitcore
23
24Abstract base class for all test
25statistics. Test statistics that evaluate the PDF at each data
26point should inherit from the RooAbsOptTestStatistic class which
27implements several generic optimizations that can be done for such
28quantities.
29
30This test statistic base class organizes calculation of test
31statistic values for RooSimultaneous PDF as a combination of test
32statistic values for the PDF components of the simultaneous PDF and
33organizes multi-processor parallel calculation of test statistic
34values. For the latter, the test statistic value is calculated in
35partitions in parallel executing processes and a posteriori
36combined in the main thread.
37**/
38
39#include "RooAbsTestStatistic.h"
40
41#include "RooAbsPdf.h"
42#include "RooSimultaneous.h"
43#include "RooAbsData.h"
44#include "RooArgSet.h"
45#include "RooRealVar.h"
46#include "RooRealMPFE.h"
47#include "RooErrorHandler.h"
48#include "RooMsgService.h"
50#include "RooFitImplHelpers.h"
52#include "RooCategory.h"
53
54#include "TTimeStamp.h"
55#include "TClass.h"
56#include <string>
57#include <stdexcept>
58
59using std::endl, std::ostream;
60
61////////////////////////////////////////////////////////////////////////////////
62/// Create a test statistic from the given function and the data.
63/// \param[in] name Name of the test statistic
64/// \param[in] title Title (for plotting)
65/// \param[in] real Function to be used for tests
66/// \param[in] data Data to fit function to
67/// \param[in] projDeps A set of projected observables
68/// \param[in] cfg statistic configuration object
69///
70/// cfg contains:
71/// - rangeName Fit data only in range with given name
72/// - addCoefRangeName If not null, all RooAddPdf components of `real` will be instructed to fix their fraction definitions to the given named range.
73/// - nCPU If larger than one, the test statistic calculation will be parallelized over multiple processes.
74/// By default the data is split with 'bulk' partitioning (each process calculates a contiguous block of fraction 1/nCPU
75/// of the data). For binned data this approach may be suboptimal as the number of bins with >0 entries
76/// in each processing block many vary greatly thereby distributing the workload rather unevenly.
77/// - interleave is set to true, the interleave partitioning strategy is used where each partition
78/// i takes all bins for which (ibin % ncpu == i) which is more likely to result in an even workload.
79/// - verbose Be more verbose.
80/// - splitCutRange If true, a different rangeName constructed as rangeName_{catName} will be used
81/// as range definition for each index state of a RooSimultaneous. This means that a different range can be defined
82/// for each category such as
83/// ```
84/// myVariable.setRange("range_pi0", 135, 210);
85/// myVariable.setRange("range_gamma", 50, 210);
86/// ```
87/// if the categories are called "pi0" and "gamma".
88
89RooAbsTestStatistic::RooAbsTestStatistic(const char *name, const char *title, RooAbsReal& real, RooAbsData& data,
90 const RooArgSet& projDeps, RooAbsTestStatistic::Configuration const& cfg) :
91 RooAbsReal(name,title),
92 _paramSet("paramSet","Set of parameters",this),
93 _func(&real),
94 _data(&data),
95 _projDeps(static_cast<RooArgSet*>(projDeps.Clone())),
96 _rangeName(cfg.rangeName),
99 _verbose(cfg.verbose),
100 // Determine if RooAbsReal is a RooSimultaneous
101 _gofOpMode{(cfg.nCPU>1 || cfg.nCPU==-1) ? MPMaster : (dynamic_cast<RooSimultaneous*>(_func) ? SimMaster : Slave)},
102 _nEvents{data.numEntries()},
103 _nCPU(cfg.nCPU != -1 ? cfg.nCPU : 1),
104 _mpinterl(cfg.interleave),
105 _takeGlobalObservablesFromData{cfg.takeGlobalObservablesFromData}
106{
107 // Register all parameters as servers
108 _paramSet.add(*std::unique_ptr<RooArgSet>{real.getParameters(&data)});
109}
110
111
112
113////////////////////////////////////////////////////////////////////////////////
114/// Copy constructor
115
116RooAbsTestStatistic::RooAbsTestStatistic(const RooAbsTestStatistic& other, const char* name) :
118 _paramSet("paramSet","Set of parameters",this),
119 _func(other._func),
120 _data(other._data),
121 _projDeps(static_cast<RooArgSet*>(other._projDeps->Clone())),
122 _rangeName(other._rangeName),
125 _verbose(other._verbose),
126 // Determine if RooAbsReal is a RooSimultaneous
128 _nEvents{_data->numEntries()},
129 _nCPU(other._nCPU != -1 ? other._nCPU : 1),
132 _takeGlobalObservablesFromData{other._takeGlobalObservablesFromData},
133 _offset(other._offset),
135{
136 // Our parameters are those of original
137 _paramSet.add(other._paramSet) ;
138}
139
140
141
142////////////////////////////////////////////////////////////////////////////////
143/// Destructor
144
145RooAbsTestStatistic::~RooAbsTestStatistic()
146{
147 if (MPMaster == _gofOpMode && _init) {
148 for (Int_t i = 0; i < _nCPU; ++i) delete _mpfeArray[i];
149 delete[] _mpfeArray ;
150 }
151
152 delete _projDeps ;
153}
154
155
156
157////////////////////////////////////////////////////////////////////////////////
158/// Calculate and return value of test statistic. If the test statistic
159/// is calculated from a RooSimultaneous, the test statistic calculation
160/// is performed separately on each simultaneous p.d.f component and associated
161/// data, and then combined. If the test statistic calculation is parallelized,
162/// partitions are calculated in nCPU processes and combined a posteriori.
163
164double RooAbsTestStatistic::evaluate() const
165{
166 // One-time Initialization
167 if (!_init) {
168 const_cast<RooAbsTestStatistic*>(this)->initialize() ;
169 }
170
171 if (SimMaster == _gofOpMode) {
172 // Evaluate array of owned GOF objects
173 double ret = 0.;
174
176 ret = combinedValue(reinterpret_cast<RooAbsReal**>(const_cast<std::unique_ptr<RooAbsTestStatistic>*>(_gofArray.data())),_gofArray.size());
177 } else {
178 double sum = 0.;
179 double carry = 0.;
180 int i = 0;
181 for (auto& gof : _gofArray) {
182 if (i % _numSets == _setNum || (_mpinterl==RooFit::Hybrid && gof->_mpinterl != RooFit::SimComponents )) {
183 double y = gof->getValV();
184 carry += gof->getCarry();
185 y -= carry;
186 const double t = sum + y;
187 carry = (t - sum) - y;
188 sum = t;
189 }
190 ++i;
191 }
192 ret = sum ;
193 _evalCarry = carry;
194 }
195
196 // Only apply global normalization if SimMaster doesn't have MP master
197 if (numSets()==1) {
198 const double norm = globalNormalization();
199 ret /= norm;
200 _evalCarry /= norm;
201 }
202
203 return ret ;
204
205 } else if (MPMaster == _gofOpMode) {
206
207 // Start calculations in parallel
208 for (Int_t i = 0; i < _nCPU; ++i) _mpfeArray[i]->calculate();
209
210 double sum(0);
211 double carry = 0.;
212 for (Int_t i = 0; i < _nCPU; ++i) {
213 double y = _mpfeArray[i]->getValV();
214 carry += _mpfeArray[i]->getCarry();
215 y -= carry;
216 const double t = sum + y;
217 carry = (t - sum) - y;
218 sum = t;
219 }
220
221 double ret = sum ;
222 _evalCarry = carry;
223
224 const double norm = globalNormalization();
225 ret /= norm;
226 _evalCarry /= norm;
227
228 return ret ;
229
230 } else {
231
232 // Evaluate as straight FUNC
233 Int_t nFirst(0);
234 Int_t nLast(_nEvents);
235 Int_t nStep(1);
236
237 switch (_mpinterl) {
239 nFirst = _nEvents * _setNum / _numSets ;
240 nLast = _nEvents * (_setNum+1) / _numSets ;
241 nStep = 1 ;
242 break;
243
245 nFirst = _setNum ;
246 nLast = _nEvents ;
247 nStep = _numSets ;
248 break ;
249
251 nFirst = 0 ;
252 nLast = _nEvents ;
253 nStep = 1 ;
254 break ;
255
256 case RooFit::Hybrid:
257 throw std::logic_error("this should never happen");
258 break ;
259 }
260
262 double ret = evaluatePartition(nFirst,nLast,nStep);
263
264 if (numSets()==1) {
265 const double norm = globalNormalization();
266 ret /= norm;
267 _evalCarry /= norm;
268 }
269
270 return ret ;
271
272 }
273}
274
275
276
277////////////////////////////////////////////////////////////////////////////////
278/// One-time initialization of the test statistic. Setup
279/// infrastructure for simultaneous p.d.f processing and/or
280/// parallelized processing if requested
281
282bool RooAbsTestStatistic::initialize()
283{
284 if (_init) return false;
285
286 if (MPMaster == _gofOpMode) {
287 initMPMode(_func,_data,_projDeps,_rangeName,_addCoefRangeName) ;
288 } else if (SimMaster == _gofOpMode) {
289 initSimMode(static_cast<RooSimultaneous*>(_func),_data,_projDeps,_rangeName,_addCoefRangeName) ;
290 }
291 _init = true;
292 return false;
293}
294
295
296
297////////////////////////////////////////////////////////////////////////////////
298/// Forward server redirect calls to component test statistics
299
300bool RooAbsTestStatistic::redirectServersHook(const RooAbsCollection& newServerList, bool mustReplaceAll, bool nameChange, bool isRecursive)
301{
302 if (SimMaster == _gofOpMode) {
303 // Forward to slaves
304 for(auto& gof : _gofArray) {
305 gof->recursiveRedirectServers(newServerList,mustReplaceAll,nameChange);
306 }
307 } else if (MPMaster == _gofOpMode&& _mpfeArray) {
308 // Forward to slaves
309 for (Int_t i = 0; i < _nCPU; ++i) {
310 if (_mpfeArray[i]) {
311 _mpfeArray[i]->recursiveRedirectServers(newServerList,mustReplaceAll,nameChange);
312// std::cout << "redirecting servers on " << _mpfeArray[i]->GetName() << std::endl;
313 }
314 }
315 }
317}
318
319
320
321////////////////////////////////////////////////////////////////////////////////
322/// Add extra information on component test statistics when printing
323/// itself as part of a tree structure
324
325void RooAbsTestStatistic::printCompactTreeHook(ostream& os, const char* indent)
326{
327 if (SimMaster == _gofOpMode) {
328 // Forward to slaves
329 os << indent << "RooAbsTestStatistic begin GOF contents" << std::endl ;
330 for (std::size_t i = 0; i < _gofArray.size(); ++i) {
332 indent2 += "[" + std::to_string(i) + "] ";
333 _gofArray[i]->printCompactTreeHook(os,indent2);
334 }
335 os << indent << "RooAbsTestStatistic end GOF contents" << std::endl;
336 } else if (MPMaster == _gofOpMode) {
337 // WVE implement this
338 }
339}
340
341
342
343////////////////////////////////////////////////////////////////////////////////
344/// Forward constant term optimization management calls to component
345/// test statistics
346
347void RooAbsTestStatistic::constOptimizeTestStatistic(ConstOpCode opcode, bool doAlsoTrackingOpt)
348{
349 initialize();
350 if (SimMaster == _gofOpMode) {
351 // Forward to slaves
352 int i = 0;
353 for (auto& gof : _gofArray) {
354 // In SimComponents Splitting strategy only constOptimize the terms that are actually used
356 if ( (i % _numSets == _setNum) || (effSplit != RooFit::SimComponents) ) {
357 gof->constOptimizeTestStatistic(opcode,doAlsoTrackingOpt);
358 }
359 ++i;
360 }
361 } else if (MPMaster == _gofOpMode) {
362 for (Int_t i = 0; i < _nCPU; ++i) {
363 _mpfeArray[i]->constOptimizeTestStatistic(opcode,doAlsoTrackingOpt);
364 }
365 }
366}
367
368
369
370////////////////////////////////////////////////////////////////////////////////
371/// Set MultiProcessor set number identification of this instance
372
373void RooAbsTestStatistic::setMPSet(Int_t inSetNum, Int_t inNumSets)
374{
377
378 if (SimMaster == _gofOpMode) {
379 // Forward to slaves
380 initialize();
381 for(auto& gof : _gofArray) {
382 gof->setMPSet(inSetNum,inNumSets);
383 }
384 }
385}
386
387
388
389////////////////////////////////////////////////////////////////////////////////
390/// Initialize multi-processor calculation mode. Create component test statistics in separate
391/// processed that are connected to this process through a RooAbsRealMPFE front-end class.
392
393void RooAbsTestStatistic::initMPMode(RooAbsReal* real, RooAbsData* data, const RooArgSet* projDeps, std::string const& rangeName, std::string const& addCoefRangeName)
394{
396
397 // Create proto-goodness-of-fit
398 Configuration cfg;
399 cfg.rangeName = rangeName;
400 cfg.addCoefRangeName = addCoefRangeName;
401 cfg.nCPU = 1;
402 cfg.interleave = _mpinterl;
403 cfg.verbose = _verbose;
404 cfg.splitCutRange = _splitRange;
405 cfg.takeGlobalObservablesFromData = _takeGlobalObservablesFromData;
406 // This configuration parameter is stored in the RooAbsOptTestStatistic.
407 // It would have been cleaner to move the member variable into RooAbsTestStatistic,
408 // but to avoid incrementing the class version we do the dynamic_cast trick.
409 if(auto thisAsRooAbsOptTestStatistic = dynamic_cast<RooAbsOptTestStatistic const*>(this)) {
410 cfg.integrateOverBinsPrecision = thisAsRooAbsOptTestStatistic->_integrateBinsPrecision;
411 }
412 RooAbsTestStatistic* gof = create(GetName(),GetTitle(),*real,*data,*projDeps,cfg);
413 gof->recursiveRedirectServers(_paramSet);
414
415 for (Int_t i = 0; i < _nCPU; ++i) {
416 gof->setMPSet(i,_nCPU);
417 gof->SetName(Form("%s_GOF%d",GetName(),i));
418 gof->SetTitle(Form("%s_GOF%d",GetTitle(),i));
419
420 ccoutD(Eval) << "RooAbsTestStatistic::initMPMode: starting remote server process #" << i << std::endl;
421 _mpfeArray[i] = new RooRealMPFE(Form("%s_%zx_MPFE%d",GetName(),reinterpret_cast<size_t>(this),i),Form("%s_%zx_MPFE%d",GetTitle(),reinterpret_cast<size_t>(this),i),*gof,false);
422 //_mpfeArray[i]->setVerbose(true,true);
423 _mpfeArray[i]->initialize();
424 if (i > 0) {
425 _mpfeArray[i]->followAsSlave(*_mpfeArray[0]);
426 }
427 }
428 _mpfeArray[_nCPU - 1]->addOwnedComponents(*gof);
429 coutI(Eval) << "RooAbsTestStatistic::initMPMode: started " << _nCPU << " remote server process." << std::endl;
430 //cout << "initMPMode --- done" << std::endl ;
431 return ;
432}
433
434
435
436////////////////////////////////////////////////////////////////////////////////
437/// Initialize simultaneous p.d.f processing mode. Strip simultaneous
438/// p.d.f into individual components, split dataset in subset
439/// matching each component and create component test statistics for
440/// each of them.
441
442void RooAbsTestStatistic::initSimMode(RooSimultaneous* simpdf, RooAbsData* data,
443 const RooArgSet* projDeps,
444 std::string const& rangeName, std::string const& addCoefRangeName)
445{
446
447 RooAbsCategoryLValue& simCat = const_cast<RooAbsCategoryLValue&>(simpdf->indexCat());
448
449 std::unique_ptr<TList> dsetList{const_cast<RooAbsData*>(data)->split(*simpdf,processEmptyDataSets())};
450
451 // Create array of regular fit contexts, containing subset of data and single fitCat PDF
452 for (const auto& catState : simCat) {
453 const std::string& catName = catState.first;
455
456 // If the channel is not in the selected range of the category variable, we
457 // won't create a slave calculator for this channel.
458 if(!rangeName.empty()) {
459 // Only the RooCategory supports ranges, not the other
460 // RooAbsCategoryLValue-derived classes.
461 auto simCatAsRooCategory = dynamic_cast<RooCategory*>(&simCat);
462 if(simCatAsRooCategory && !simCatAsRooCategory->isStateInRange(rangeName.c_str(), catIndex)) {
463 continue;
464 }
465 }
466
467 // Retrieve the PDF for this simCat state
468 RooAbsPdf* pdf = simpdf->getPdf(catName.c_str());
469 RooAbsData* dset = static_cast<RooAbsData*>(dsetList->FindObject(catName.c_str()));
470
471 if (pdf && dset && (0. != dset->sumEntries() || processEmptyDataSets())) {
472 ccoutI(Fitting) << "RooAbsTestStatistic::initSimMode: creating slave calculator #" << _gofArray.size() << " for state " << catName
473 << " (" << dset->numEntries() << " dataset entries)" << std::endl;
474
475
476 // *** START HERE
477 // WVE HACK determine if we have a RooRealSumPdf and then treat it like a binned likelihood
479 RooAbsReal &actualPdf = binnedInfo.binnedPdf ? *binnedInfo.binnedPdf : *pdf;
480 // WVE END HACK
481 // Below here directly pass binnedPdf instead of PROD(binnedPdf,constraints) as constraints are evaluated elsewhere anyway
482 // and omitting them reduces model complexity and associated handling/cloning times
483 Configuration cfg;
484 cfg.addCoefRangeName = addCoefRangeName;
485 cfg.interleave = _mpinterl;
486 cfg.verbose = _verbose;
487 cfg.splitCutRange = _splitRange;
488 cfg.binnedL = binnedInfo.isBinnedL;
489 cfg.takeGlobalObservablesFromData = _takeGlobalObservablesFromData;
490 // This configuration parameter is stored in the RooAbsOptTestStatistic.
491 // It would have been cleaner to move the member variable into RooAbsTestStatistic,
492 // but to avoid incrementing the class version we do the dynamic_cast trick.
493 if(auto thisAsRooAbsOptTestStatistic = dynamic_cast<RooAbsOptTestStatistic const*>(this)) {
494 cfg.integrateOverBinsPrecision = thisAsRooAbsOptTestStatistic->_integrateBinsPrecision;
495 }
497 cfg.nCPU = _nCPU;
498 _gofArray.emplace_back(create(catName.c_str(),catName.c_str(),actualPdf,*dset,*projDeps,cfg));
499 // *** END HERE
500
501 // Fill per-component split mode with Bulk Partition for now so that Auto will map to bulk-splitting of all components
503 _gofArray.back()->_mpinterl = dset->numEntries()<10 ? RooFit::SimComponents : RooFit::BulkPartition;
504 }
505
506 // Servers may have been redirected between instantiation and (deferred) initialization
507
509 actualPdf.getParameters(dset->get(), actualParams);
512
513 _gofArray.back()->recursiveRedirectServers(selTargetParams);
514 }
515 }
516 for(auto& gof : _gofArray) {
517 gof->setSimCount(_gofArray.size());
518 }
519 coutI(Fitting) << "RooAbsTestStatistic::initSimMode: created " << _gofArray.size() << " slave calculators." << std::endl;
520
521 dsetList->Delete(); // delete the content.
522}
523
524
525////////////////////////////////////////////////////////////////////////////////
526/// Change dataset that is used to given one. If cloneData is true, a clone of
527/// in the input dataset is made. If the test statistic was constructed with
528/// a range specification on the data, the cloneData argument is ignored and
529/// the data is always cloned.
530bool RooAbsTestStatistic::setData(RooAbsData& indata, bool cloneData)
531{
532 // Trigger refresh of likelihood offsets
533 if (isOffsetting()) {
534 enableOffsetting(false);
535 enableOffsetting(true);
536 }
537
538 switch(operMode()) {
539 case Slave:
540 // Delegate to implementation
542 case SimMaster:
543 // Forward to slaves
544 if (indata.canSplitFast()) {
545 for(auto& gof : _gofArray) {
546 RooAbsData* compData = indata.getSimData(gof->GetName());
547 gof->setDataSlave(*compData, cloneData);
548 }
549 } else if (0 == indata.numEntries()) {
550 // For an unsplit empty dataset, simply assign empty dataset to each component
551 for(auto& gof : _gofArray) {
552 gof->setDataSlave(indata, cloneData);
553 }
554 } else {
555 std::unique_ptr<TList> dlist{indata.split(*static_cast<RooSimultaneous*>(_func), processEmptyDataSets())};
556 if (!dlist) {
557 coutE(Fitting) << "RooAbsTestStatistic::initSimMode(" << GetName() << ") ERROR: index category of simultaneous pdf is missing in dataset, aborting" << std::endl;
558 throw std::runtime_error("RooAbsTestStatistic::initSimMode() ERROR, index category of simultaneous pdf is missing in dataset, aborting");
559 }
560
561 for(auto& gof : _gofArray) {
562 if (auto compData = static_cast<RooAbsData*>(dlist->FindObject(gof->GetName()))) {
563 gof->setDataSlave(*compData,false,true);
564 } else {
565 coutE(DataHandling) << "RooAbsTestStatistic::setData(" << GetName() << ") ERROR: Cannot find component data for state " << gof->GetName() << std::endl;
566 }
567 }
568 }
569 break;
570 case MPMaster:
571 // Not supported
572 coutF(DataHandling) << "RooAbsTestStatistic::setData(" << GetName() << ") FATAL: setData() is not supported in multi-processor mode" << std::endl;
573 throw std::runtime_error("RooAbsTestStatistic::setData is not supported in MPMaster mode");
574 break;
575 }
576
577 return true;
578}
579
580
581
582void RooAbsTestStatistic::enableOffsetting(bool flag)
583{
584 // Apply internal value offsetting to control numeric precision
585 if (!_init) {
586 const_cast<RooAbsTestStatistic*>(this)->initialize() ;
587 }
588
589 switch(operMode()) {
590 case Slave:
591 _doOffset = flag ;
592 // Clear offset if feature is disabled to that it is recalculated next time it is enabled
593 if (!_doOffset) {
594 _offset = ROOT::Math::KahanSum<double>{0.} ;
595 }
596 setValueDirty() ;
597 break ;
598 case SimMaster:
599 _doOffset = flag;
600 for(auto& gof : _gofArray) {
601 gof->enableOffsetting(flag);
602 }
603 break ;
604 case MPMaster:
605 _doOffset = flag;
606 for (Int_t i = 0; i < _nCPU; ++i) {
607 _mpfeArray[i]->enableOffsetting(flag);
608 }
609 break;
610 }
611}
612
613
614double RooAbsTestStatistic::getCarry() const
615{ return _evalCarry; }
616
617/// \endcond
#define coutI(a)
#define coutF(a)
#define coutE(a)
#define ccoutI(a)
#define ccoutD(a)
int Int_t
Definition RtypesCore.h:45
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
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
The Kahan summation is a compensated summation algorithm, which significantly reduces numerical error...
Definition Util.h:136
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.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursiveStep) override
Function that is called at the end of redirectServers().
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
RooArgSet * selectCommon(const RooAbsCollection &refColl) const
Use RooAbsCollection::selecCommon(), but return as RooArgSet.
Definition RooArgSet.h:154
Object to represent discrete states.
Definition RooCategory.h:28
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
Basic string class.
Definition TString.h:139
Double_t y[n]
Definition legend1.C:17
@ SimComponents
@ BulkPartition
std::string getRangeNameForSimComponent(std::string const &rangeName, bool splitRange, std::string const &catName)
BinnedLOutput getBinnedL(RooAbsPdf const &pdf)
void initialize(typename Architecture_t::Matrix_t &A, EInitialization m)
Definition Functions.h:282
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345