Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooSimultaneous.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\file RooSimultaneous.cxx
19\class RooSimultaneous
20\ingroup Roofitcore
21
22Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
23The class takes an index category, which is used as a selector
24for PDFs, and a list of PDFs, each associated
25with a state of the index category. RooSimultaneous always returns
26the value of the PDF that is associated with the current value
27of the index category.
28
29Extended likelihood fitting is supported if all components support
30extended likelihood mode. The expected number of events by a RooSimultaneous
31is that of the component p.d.f. selected by the index category.
32
33The index category can be accessed using indexCategory().
34
35###Generating events
36When generating events from a RooSimultaneous, the index category has to be added to
37the dataset. Further, the PDF needs to know the relative probabilities of each category, i.e.,
38how many events are in which category. This can be achieved in two ways:
39- Generating with proto data that have category entries: An event from the same category as
40in the proto data is created for each event in the proto data.
41See RooAbsPdf::generate(const RooArgSet&,const RooDataSet&,Int_t,bool,bool,bool) const.
42- No proto data: A category is chosen randomly.
43\note This requires that the PDFs building the simultaneous are extended. In this way,
44the relative probability of each category can be calculated from the number of events
45in each category.
46**/
47
48#include "RooSimultaneous.h"
49
50#include "Roo1DTable.h"
52#include "RooAbsData.h"
53#include "RooAddPdf.h"
54#include "RooArgSet.h"
55#include "RooBinSamplingPdf.h"
56#include "RooCategory.h"
57#include "RooCmdConfig.h"
58#include "RooDataHist.h"
59#include "RooDataSet.h"
60#include "RooGlobalFunc.h"
61#include "RooMsgService.h"
62#include "RooNameReg.h"
63#include "RooPlot.h"
64#include "RooRandom.h"
65#include "RooRealVar.h"
66#include "RooSimGenContext.h"
68#include "RooSuperCategory.h"
69
70#include "RooFitImplHelpers.h"
71
72#include <ROOT/StringUtils.hxx>
73
74#include <iostream>
75
76namespace {
77
78std::map<std::string, RooAbsPdf *> createPdfMap(const RooArgList &inPdfList, RooAbsCategoryLValue &inIndexCat)
79{
80 std::map<std::string, RooAbsPdf *> pdfMap;
82 for (unsigned int i = 0; i < inPdfList.size(); ++i) {
83 auto pdf = static_cast<RooAbsPdf *>(&inPdfList[i]);
84 const auto &nameIdx = (*indexCatIt++);
85 pdfMap[nameIdx.first] = pdf;
86 }
87 return pdfMap;
88}
89
91{
92 TObject *old = lst.FindObject(obj.GetName());
93 if (old)
94 lst.Replace(old, &obj);
95 else
96 lst.Add(&obj);
97}
98
99} // namespace
100
102
104{
105 finalPdfs.push_back(&pdf);
106 finalCatLabels.emplace_back(catLabel);
107}
108
109using std::string;
110
111
112
113
114////////////////////////////////////////////////////////////////////////////////
115/// Constructor with index category. PDFs associated with indexCat
116/// states can be added after construction with the addPdf() function.
117///
118/// RooSimultaneous can function without having a PDF associated
119/// with every single state. The normalization in such cases is taken
120/// from the number of registered PDFs, but getVal() will assert if
121/// when called for an unregistered index state.
122
123RooSimultaneous::RooSimultaneous(const char *name, const char *title,
125 RooSimultaneous{name, title, std::map<std::string, RooAbsPdf*>{}, inIndexCat}
126{
127}
128
129
130////////////////////////////////////////////////////////////////////////////////
131/// Constructor from index category and full list of PDFs.
132/// In this constructor form, a PDF must be supplied for each indexCat state
133/// to avoid ambiguities. The PDFs are associated with the states of the
134/// index category as they appear when iterating through the category states
135/// with RooAbsCategory::begin() and RooAbsCategory::end(). This usually means
136/// they are associated by ascending index numbers.
137///
138/// PDFs may not overlap (i.e. share any variables) with the index category (function)
139
140RooSimultaneous::RooSimultaneous(const char *name, const char *title,
143{
144 if (inPdfList.size() != inIndexCat.size()) {
145 std::stringstream errMsg;
146 errMsg << "RooSimultaneous::ctor(" << GetName()
147 << " ERROR: Number PDF list entries must match number of index category states, no PDFs added";
148 coutE(InputArguments) << errMsg.str() << std::endl;
149 throw std::invalid_argument(errMsg.str());
150 }
151}
152
153
154////////////////////////////////////////////////////////////////////////////////
155
156RooSimultaneous::RooSimultaneous(const char *name, const char *title, std::map<string, RooAbsPdf *> pdfMap,
158 : RooSimultaneous(name, title, std::move(*initialize(name ? name : "", inIndexCat, pdfMap)))
159{
160}
161
162/// For internal use in RooFit.
163RooSimultaneous::RooSimultaneous(const char *name, const char *title,
166 : RooSimultaneous(name, title, RooFit::Detail::flatMapToStdMap(pdfMap), inIndexCat)
167{
168}
169
171 : RooAbsPdf(name, title),
172 _plotCoefNormSet("!plotCoefNormSet", "plotCoefNormSet", this, false, false),
173 _partIntMgr(this, 10),
174 _indexCat("indexCat", "Index category", this, *initInfo.indexCat)
175{
176 for (std::size_t i = 0; i < initInfo.finalPdfs.size(); ++i) {
177 addPdf(*initInfo.finalPdfs[i], initInfo.finalCatLabels[i].c_str());
178 }
179
180 // Take ownership of eventual super category
181 if (initInfo.superIndex) {
182 addOwnedComponents(std::move(initInfo.superIndex));
183 }
184}
185
186/// \cond ROOFIT_INTERNAL
187
188// This class cannot be locally defined in initialize as it cannot be
189// used as a template argument in that case
190namespace RooSimultaneousAux {
191 struct CompInfo {
192 RooAbsPdf* pdf ;
195 std::unique_ptr<RooArgSet> subIndexComps;
196 } ;
197}
198
199/// \endcond
200
201std::unique_ptr<RooSimultaneous::InitializationOutput>
203 std::map<std::string, RooAbsPdf *> const& pdfMap)
204
205{
206 auto out = std::make_unique<RooSimultaneous::InitializationOutput>();
207 out->indexCat = &inIndexCat;
208
209 // First see if there are any RooSimultaneous input components
210 bool simComps(false) ;
211 for (auto const& item : pdfMap) {
212 if (dynamic_cast<RooSimultaneous*>(item.second)) {
213 simComps = true ;
214 break ;
215 }
216 }
217
218 // If there are no simultaneous component p.d.f. do simple processing through addPdf()
219 if (!simComps) {
220 for (auto const& item : pdfMap) {
221 out->addPdf(*item.second,item.first);
222 }
223 return out;
224 }
225
226 std::string msgPrefix = "RooSimultaneous::initialize(" + name + ") ";
227
228 // Issue info message that we are about to do some rearranging
229 oocoutI(nullptr, InputArguments) << msgPrefix << "INFO: one or more input component of simultaneous p.d.f.s are"
230 << " simultaneous p.d.f.s themselves, rewriting composite expressions as one-level simultaneous p.d.f. in terms of"
231 << " final constituents and extended index category" << std::endl;
232
233
235 std::map<string,RooSimultaneousAux::CompInfo> compMap ;
236 for (auto const& item : pdfMap) {
237 RooSimultaneousAux::CompInfo ci ;
238 ci.pdf = item.second ;
239 RooSimultaneous* simComp = dynamic_cast<RooSimultaneous*>(item.second) ;
240 if (simComp) {
241 ci.simPdf = simComp ;
242 ci.subIndex = &simComp->indexCat() ;
243 ci.subIndexComps = simComp->indexCat().isFundamental()
244 ? std::make_unique<RooArgSet>(simComp->indexCat())
245 : std::unique_ptr<RooArgSet>(simComp->indexCat().getVariables());
246 allAuxCats.add(*ci.subIndexComps,true) ;
247 } else {
248 ci.simPdf = nullptr;
249 ci.subIndex = nullptr;
250 }
251 compMap[item.first] = std::move(ci);
252 }
253
254 // Construct the 'superIndex' from the nominal index category and all auxiliary components
255 RooArgSet allCats(inIndexCat) ;
256 allCats.add(allAuxCats) ;
257 std::string siname = name + "_index";
258 out->superIndex = std::make_unique<RooSuperCategory>(siname.c_str(),siname.c_str(),allCats) ;
259 auto *superIndex = out->superIndex.get();
260 out->indexCat = superIndex;
261
262 // Now process each of original pdf/state map entries
263 for (auto const& citem : compMap) {
264
266 if (citem.second.subIndexComps) {
267 repliCats.remove(*citem.second.subIndexComps) ;
268 }
269 inIndexCat.setLabel(citem.first.c_str()) ;
270
271 if (!citem.second.simPdf) {
272
273 // Entry is a plain p.d.f. assign it to every state permutation of the repliCats set
275
276 // Iterator over all states of repliSuperCat
277 for (const auto& nameIdx : repliSuperCat) {
278 // Set value
279 repliSuperCat.setLabel(nameIdx.first) ;
280 // Retrieve corresponding label of superIndex
281 string superLabel = superIndex->getCurrentLabel() ;
282 out->addPdf(*citem.second.pdf,superLabel);
283 oocxcoutD(static_cast<RooAbsArg*>(nullptr), InputArguments) << msgPrefix
284 << "assigning pdf " << citem.second.pdf->GetName() << " to super label " << superLabel << std::endl ;
285 }
286 } else {
287
288 // Entry is a simultaneous p.d.f
289
290 if (repliCats.empty()) {
291
292 // Case 1 -- No replication of components of RooSim component are required
293
294 for (const auto& type : *citem.second.subIndex) {
295 const_cast<RooAbsCategoryLValue*>(citem.second.subIndex)->setLabel(type.first.c_str());
296 string superLabel = superIndex->getCurrentLabel() ;
297 RooAbsPdf* compPdf = citem.second.simPdf->getPdf(type.first);
298 if (compPdf) {
299 out->addPdf(*compPdf,superLabel);
300 oocxcoutD(static_cast<RooAbsArg*>(nullptr), InputArguments) << msgPrefix
301 << "assigning pdf " << compPdf->GetName() << "(member of " << citem.second.pdf->GetName()
302 << ") to super label " << superLabel << std::endl ;
303 } else {
304 oocoutW(nullptr, InputArguments) << msgPrefix << "WARNING: No p.d.f. associated with label "
305 << type.second << " for component RooSimultaneous p.d.f " << citem.second.pdf->GetName()
306 << "which is associated with master index label " << citem.first << std::endl ;
307 }
308 }
309
310 } else {
311
312 // Case 2 -- Replication of components of RooSim component are required
313
314 // Make replication supercat
316
317 for (const auto& stype : *citem.second.subIndex) {
318 const_cast<RooAbsCategoryLValue*>(citem.second.subIndex)->setLabel(stype.first.c_str());
319
320 for (const auto& nameIdx : repliSuperCat) {
321 repliSuperCat.setLabel(nameIdx.first) ;
322 const string superLabel = superIndex->getCurrentLabel() ;
323 RooAbsPdf* compPdf = citem.second.simPdf->getPdf(stype.first);
324 if (compPdf) {
325 out->addPdf(*compPdf,superLabel);
326 oocxcoutD(static_cast<RooAbsArg*>(nullptr), InputArguments) << msgPrefix
327 << "assigning pdf " << compPdf->GetName() << "(member of " << citem.second.pdf->GetName()
328 << ") to super label " << superLabel << std::endl ;
329 } else {
330 oocoutW(nullptr, InputArguments) << msgPrefix << "WARNING: No p.d.f. associated with label "
331 << stype.second << " for component RooSimultaneous p.d.f " << citem.second.pdf->GetName()
332 << "which is associated with master index label " << citem.first << std::endl ;
333 }
334 }
335 }
336 }
337 }
338 }
339
340 return out;
341}
342
343
344////////////////////////////////////////////////////////////////////////////////
345/// Copy constructor
346
349 _plotCoefNormSet("!plotCoefNormSet",this,other._plotCoefNormSet),
350 _plotCoefNormRange(other._plotCoefNormRange),
351 _partIntMgr(other._partIntMgr,this),
352 _indexCat("indexCat",this,other._indexCat),
353 _numPdf(other._numPdf)
354{
355 // Copy proxy list
356 for(auto* proxy : static_range_cast<RooRealProxy*>(other._pdfProxyList)) {
357 _pdfProxyList.Add(new RooRealProxy(proxy->GetName(),this,*proxy)) ;
358 }
359}
360
361
362
363////////////////////////////////////////////////////////////////////////////////
364/// Destructor
365
370
371
372
373////////////////////////////////////////////////////////////////////////////////
374/// Return the p.d.f associated with the given index category name
375
377{
379 return proxy ? static_cast<RooAbsPdf*>(proxy->absArg()) : nullptr;
380}
381
382
383
384////////////////////////////////////////////////////////////////////////////////
385/// Associate given PDF with index category state label 'catLabel'.
386/// The name state must be already defined in the index category.
387///
388/// RooSimultaneous can function without having a PDF associated
389/// with every single state. The normalization in such cases is taken
390/// from the number of registered PDFs, but getVal() will fail if
391/// called for an unregistered index state.
392///
393/// PDFs may not overlap (i.e. share any variables) with the index category (function).
394/// \param[in] pdf PDF to be added.
395/// \param[in] catLabel Name of the category state to be associated to the PDF.
396/// \return `true` in case of failure.
397
398bool RooSimultaneous::addPdf(const RooAbsPdf& pdf, const char* catLabel)
399{
400 // PDFs cannot overlap with the index category
401 if (pdf.dependsOn(_indexCat.arg())) {
402 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName() << "): PDF '" << pdf.GetName()
403 << "' overlaps with index category '" << _indexCat.arg().GetName() << "'."<< std::endl ;
404 return true ;
405 }
406
407 // Each index state can only have one PDF associated with it
409 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName() << "): index state '"
410 << catLabel << "' has already an associated PDF." << std::endl ;
411 return true ;
412 }
413
414 const RooSimultaneous* simPdf = dynamic_cast<const RooSimultaneous*>(&pdf) ;
415 if (simPdf) {
416
417 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName()
418 << ") ERROR: you cannot add a RooSimultaneous component to a RooSimultaneous using addPdf()."
419 << " Use the constructor with RooArgList if input p.d.f.s or the map<string,RooAbsPdf&> instead." << std::endl ;
420 return true ;
421
422 } else {
423
424 // Create a proxy named after the associated index state
425 TObject* proxy = new RooRealProxy(catLabel,catLabel,this,const_cast<RooAbsPdf&>(pdf));
427 _numPdf += 1 ;
428 }
429
430 return false ;
431}
432
433////////////////////////////////////////////////////////////////////////////////
434/// Examine the pdf components and check if one of them can be extended or must be extended.
435/// It is enough to have one component that can be extended or must be extended to return the flag in
436/// the total simultaneous pdf.
437
439{
440 bool anyCanExtend = false;
441
443 auto &pdf = static_cast<RooAbsPdf const&>(proxy->arg());
444 if (pdf.mustBeExtended())
445 return MustBeExtended;
446 anyCanExtend |= pdf.canBeExtended();
447 }
449}
450
451////////////////////////////////////////////////////////////////////////////////
452/// Return the current value:
453/// the value of the PDF associated with the current index category state
454
456{
457 // Retrieve the proxy by index name
459 if(!proxy) {
460 return 0;
461 }
462
463 double nEvtTot = 1.0;
464 double nEvtCat = 1.0;
465
466 // Calculate relative weighting factor for sim-pdfs of all extendable components
467 if (canBeExtended()) {
468
469 nEvtTot = 0;
470 nEvtCat = 0;
471
473 auto &pdf2 = static_cast<RooAbsPdf const &>(proxy2->arg());
474 if(!pdf2.canBeExtended()) {
475 // If one of the pdfs can't be expected, reset the normalization
476 // factor to one and break out of the loop.
477 nEvtTot = 1.0;
478 nEvtCat = 1.0;
479 break;
480 }
481 const double nEvt = pdf2.expectedEvents(_normSet);
482 nEvtTot += nEvt;
483 if (proxy == proxy2) {
484 // Matching by proxy by pointer rather than pdfs, because it's
485 // possible to have the same pdf used in different states.
486 nEvtCat += nEvt;
487 }
488 }
489 }
490 double catFrac = nEvtCat / nEvtTot;
491
492 // Return the selected PDF value, normalized by the relative number of
493 // expected events if applicable.
494 return *proxy * catFrac;
495}
496
497////////////////////////////////////////////////////////////////////////////////
498/// Return the number of expected events: If the index is in nset,
499/// then return the sum of the expected events of all components,
500/// otherwise return the number of expected events of the PDF
501/// associated with the current index category state
502
504{
505 if (nset->contains(_indexCat.arg())) {
506
507 double sum(0) ;
508
510 sum += (static_cast<RooAbsPdf*>(proxy->absArg()))->expectedEvents(nset) ;
511 }
512
513 return sum ;
514
515 } else {
516
517 // Retrieve the proxy by index name
519
520 //assert(proxy!=0) ;
521 if (proxy==nullptr) return 0 ;
522
523 // Return the selected PDF value, normalized by the number of index states
524 return (static_cast<RooAbsPdf*>(proxy->absArg()))->expectedEvents(nset);
525 }
526}
527
528
529
530////////////////////////////////////////////////////////////////////////////////
531/// Forward determination of analytical integration capabilities to component p.d.f.s
532/// A unique code is assigned to the combined integration capabilities of all associated
533/// p.d.f.s
534
536 const RooArgSet* normSet, const char* rangeName) const
537{
538 // Declare that we can analytically integrate all requested observables
539 analVars.add(allVars) ;
540
541 // Retrieve (or create) the required partial integral list
542 Int_t code ;
543
544 // Check if this configuration was created before
546 if (cache) {
547 code = _partIntMgr.lastIndex() ;
548 return code+1 ;
549 }
550 cache = new CacheElem ;
551
552 // Create the partial integral set for this request
554 cache->_partIntList.addOwned(std::unique_ptr<RooAbsReal>{proxy->arg().createIntegral(analVars,normSet,nullptr,rangeName)});
555 }
556
557 // Store the partial integral list and return the assigned code ;
559
560 return code+1 ;
561}
562
563
564
565////////////////////////////////////////////////////////////////////////////////
566/// Return analytical integration defined by given code
567
568double RooSimultaneous::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* /*rangeName*/) const
569{
570 // No integration scenario
571 if (code==0) {
572 return getVal(normSet) ;
573 }
574
575 // Partial integration scenarios, rangeName already encoded in 'code'
576 CacheElem* cache = static_cast<CacheElem*>(_partIntMgr.getObjByIndex(code-1)) ;
577
580 return (static_cast<RooAbsReal*>(cache->_partIntList.at(idx)))->getVal(normSet) ;
581}
582
583
584
585
586
587
588////////////////////////////////////////////////////////////////////////////////
589/// Back-end for plotOn() implementation on RooSimultaneous which
590/// needs special handling because a RooSimultaneous PDF cannot
591/// project out its index category via integration. plotOn() will
592/// abort if this is requested without providing a projection dataset.
593
595{
596 // Sanity checks
597 if (plotSanityChecks(frame)) return frame ;
598
599 // Extract projection configuration from command list
600 RooCmdConfig pc("RooSimultaneous::plotOn(" + std::string(GetName()) + ")");
601 pc.defineString("sliceCatState","SliceCat",0,"",true) ;
602 pc.defineDouble("scaleFactor","Normalization",0,1.0) ;
603 pc.defineInt("scaleType","Normalization",0,RooAbsPdf::Relative) ;
604 pc.defineObject("sliceCatList","SliceCat",0,nullptr,true) ;
605 // This dummy is needed for plotOn to recognize the "SliceCatMany" command.
606 // It is not used directly, but the "SliceCat" commands are nested in it.
607 // Removing this dummy definition results in "ERROR: unrecognized command: SliceCatMany".
608 pc.defineObject("dummy1","SliceCatMany",0) ;
609 pc.defineSet("projSet","Project",0) ;
610 pc.defineSet("sliceSet","SliceVars",0) ;
611 pc.defineSet("projDataSet","ProjData",0) ;
612 pc.defineObject("projData","ProjData",1) ;
613 pc.defineMutex("Project","SliceVars") ;
614 pc.allowUndefined() ; // there may be commands we don't handle here
615
616 // Process and check varargs
617 pc.process(cmdList) ;
618 if (!pc.ok(true)) {
619 return frame ;
620 }
621
622 RooAbsData* projData = static_cast<RooAbsData*>(pc.getObject("projData")) ;
623 const RooArgSet* projDataSet = pc.getSet("projDataSet");
624 const RooArgSet* sliceSetTmp = pc.getSet("sliceSet") ;
625 std::unique_ptr<RooArgSet> sliceSet( sliceSetTmp ? (static_cast<RooArgSet*>(sliceSetTmp->Clone())) : nullptr );
626 const RooArgSet* projSet = pc.getSet("projSet") ;
627 double scaleFactor = pc.getDouble("scaleFactor") ;
628 ScaleType stype = (ScaleType) pc.getInt("scaleType") ;
629
630
631 // Look for category slice arguments and add them to the master slice list if found
632 const char* sliceCatState = pc.getString("sliceCatState",nullptr,true) ;
633 const RooLinkedList& sliceCatList = pc.getObjectList("sliceCatList") ;
634 if (sliceCatState) {
635
636 // Make the master slice set if it doesnt exist
637 if (!sliceSet) {
638 sliceSet = std::make_unique<RooArgSet>();
639 }
640
641 // Prepare comma separated label list for parsing
643
644 // Loop over all categories provided by (multiple) Slice() arguments
645 unsigned int tokenIndex = 0;
647 const char* slabel = tokenIndex >= catTokens.size() ? nullptr : catTokens[tokenIndex++].c_str();
648
649 if (slabel) {
650 // Set the slice position to the value indicated by slabel
651 scat->setLabel(slabel) ;
652 // Add the slice category to the master slice set
653 sliceSet->add(*scat,false) ;
654 }
655 }
656 }
657
658 // Check if we have a projection dataset
659 if (!projData) {
660 coutE(InputArguments) << "RooSimultaneous::plotOn(" << GetName() << ") ERROR: must have a projection dataset for index category" << std::endl ;
661 return frame ;
662 }
663
664 // Make list of variables to be projected
666 if (sliceSet) {
667 makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,true) ;
668
669 // Take out the sliced variables
670 for (const auto sliceArg : *sliceSet) {
671 RooAbsArg* arg = projectedVars.find(sliceArg->GetName()) ;
672 if (arg) {
673 projectedVars.remove(*arg) ;
674 } else {
675 coutI(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") slice variable "
676 << sliceArg->GetName() << " was not projected anyway" << std::endl ;
677 }
678 }
679 } else if (projSet) {
680 makeProjectionSet(frame->getPlotVar(),projSet,projectedVars,false) ;
681 } else {
682 makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,true) ;
683 }
684
685 bool projIndex(false) ;
686
687 if (!_indexCat.arg().isDerived()) {
688 // *** Error checking for a fundamental index category ***
689 //cout << "RooSim::plotOn: index is fundamental" << std::endl ;
690
691 // Check that the provided projection dataset contains our index variable
692 if (!projData->get()->find(_indexCat.arg().GetName())) {
693 coutE(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") ERROR: Projection over index category "
694 << "requested, but projection data set doesn't contain index category" << std::endl ;
695 return frame ;
696 }
697
698 if (projectedVars.find(_indexCat.arg().GetName())) {
700 }
701
702 } else {
703 // *** Error checking for a composite index category ***
704
705 // Determine if any servers of the index category are in the projectedVars
707 bool anyServers(false) ;
708 for (const auto server : flattenedCatList()) {
709 if (projectedVars.find(server->GetName())) {
711 projIdxServers.add(*server) ;
712 }
713 }
714
715 // Check that the projection dataset contains all the
716 // index category components we're projecting over
717
718 // Determine if all projected servers of the index category are in the projection dataset
719 bool allServers(true) ;
720 std::string missing;
721 for (const auto server : projIdxServers) {
722 if (!projData->get()->find(server->GetName())) {
724 missing = server->GetName();
725 }
726 }
727
728 if (!allServers) {
729 coutE(Plotting) << "RooSimultaneous::plotOn(" << GetName()
730 << ") ERROR: Projection dataset doesn't contain complete set of index categories to do projection."
731 << "\n\tcategory " << missing << " is missing." << std::endl ;
732 return frame ;
733 }
734
735 if (anyServers) {
736 projIndex = true ;
737 }
738 }
739
740 // Calculate relative weight fractions of components
741 std::unique_ptr<Roo1DTable> wTable( projData->table(_indexCat.arg()) );
742
743 // Clone the index category to be able to cycle through the category states for plotting without
744 // affecting the category state of our instance
746 RooArgSet(*_indexCat).snapshot(idxCloneSet, true);
747 auto idxCatClone = static_cast<RooAbsCategoryLValue*>(idxCloneSet.find(_indexCat->GetName()) );
749
750 // Make list of category columns to exclude from projection data
751 std::unique_ptr<RooArgSet> idxCompSliceSet( idxCatClone->getObservables(frame->getNormVars()) );
752
753 // If we don't project over the index, just do the regular plotOn
754 if (!projIndex) {
755
756 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
757 << " represents a slice in the index category (" << _indexCat.arg().GetName() << ")" << std::endl ;
758
759 // Reduce projData: take out fitCat (component) columns and entries that don't match selected slice
760 // Construct cut string to only select projection data event that match the current slice
761
762 // Make cut string to exclude rows from projection data
763 if (sliceSet) {
765 if (auto* slicedComponent = static_cast<const RooAbsCategory*>(sliceSet->find(*idxComp))) {
766 idxComp->setIndex(slicedComponent->getCurrentIndex(), false);
767 }
768 }
769 }
771
772 // Make temporary projData without RooSim index category components
773 RooArgSet projDataVars(*projData->get()) ;
774 projDataVars.remove(*idxCompSliceSet,true,true) ;
775
776 std::unique_ptr<RooAbsData>
778
779 // Override normalization and projection dataset
781 RooFit::Normalization(scaleFactor * wTable->get(idxCatClone->getCurrentLabel()), RooAbsReal::NumEvent);
783
784 // WVE -- do not adjust normalization for asymmetry plots
786 if (!cmdList.find("Asymmetry")) {
788 }
790
791 // Plot single component
792 RooPlot* retFrame = getPdf(idxCatClone->getCurrentLabel())->plotOn(frame,cmdList2);
793 return retFrame ;
794 }
795
796 // If we project over the index, plot using a temporary RooAddPdf
797 // using the weights from the data as coefficients
798
799 // Build the list of indexCat components that are sliced
800 idxCompSliceSet->remove(projectedVars,true,true) ;
801
802 // Make a new expression that is the weighted sum of requested components
805//RooAbsPdf* pdf ;
806 double sumWeight(0) ;
808
809 idxCatClone->setLabel(proxy->name()) ;
810
811 // Determine if this component is the current slice (if we slice)
812 bool skip(false) ;
813 for (const auto idxSliceCompArg : *idxCompSliceSet) {
814 const auto idxSliceComp = static_cast<RooAbsCategory*>(idxSliceCompArg);
815 RooAbsCategory* idxComp = static_cast<RooAbsCategory*>(idxCloneSet.find(idxSliceComp->GetName())) ;
816 if (idxComp->getCurrentIndex()!=idxSliceComp->getCurrentIndex()) {
817 skip=true ;
818 break ;
819 }
820 }
821 if (skip) continue ;
822
823 // Instantiate a RRV holding this pdfs weight
824 wgtCompList.addOwned(std::make_unique<RooRealVar>(proxy->name(),"coef",wTable->get(proxy->name())));
825 sumWeight += wTable->getFrac(proxy->name()) ;
826
827 // Add the PDF to list list
828 pdfCompList.add(proxy->arg()) ;
829 }
830
832 RooAddPdf plotVar{plotVarName,"weighted sum of RS components",pdfCompList,wgtCompList};
833
834 // Fix appropriate coefficient normalization in plot function
835 if (!_plotCoefNormSet.empty()) {
836 plotVar.fixAddCoefNormalization(_plotCoefNormSet) ;
837 }
838
839 std::unique_ptr<RooAbsData> projDataTmp;
841 if (projData) {
842
843 // Construct cut string to only select projection data event that match the current slice
845
846 // Make temporary projData without RooSim index category components
847 RooArgSet projDataVars(*projData->get()) ;
849 _indexCat.arg().getObservables(frame->getNormVars(), idxCatServers) ;
850
851 projDataVars.remove(idxCatServers,true,true) ;
852
853 projDataTmp = std::unique_ptr<RooAbsData>{projData->reduce(RooFit::SelectVars(projDataVars), RooFit::Cut(cutString.c_str()))};
854
855
856
857 if (projSet) {
858 projSetTmp.add(*projSet) ;
859 projSetTmp.remove(idxCatServers,true,true);
860 }
861 }
862
863
864 if (_indexCat.arg().isDerived() && !idxCompSliceSet->empty()) {
865 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
866 << " represents a slice in index category components " << *idxCompSliceSet << std::endl ;
867
869 _indexCat.arg().getObservables(frame->getNormVars(), idxCompProjSet) ;
870 idxCompProjSet.remove(*idxCompSliceSet,true,true) ;
871 if (!idxCompProjSet.empty()) {
872 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
873 << " averages with data index category components " << idxCompProjSet << std::endl ;
874 }
875 } else {
876 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
877 << " averages with data index category (" << _indexCat.arg().GetName() << ")" << std::endl ;
878 }
879
880
881 // Override normalization and projection dataset
883
884 RooCmdArg tmp1 = RooFit::Normalization(scaleFactor*sumWeight,stype) ;
886 // WVE -- do not adjust normalization for asymmetry plots
887 if (!cmdList.find("Asymmetry")) {
889 }
891
892 RooPlot* frame2 ;
893 if (!projSetTmp.empty()) {
894 // Plot temporary function
897 frame2 = plotVar.plotOn(frame,cmdList2) ;
898 } else {
899 // Plot temporary function
900 frame2 = plotVar.plotOn(frame,cmdList2) ;
901 }
902
903 return frame2 ;
904}
905
906
907////////////////////////////////////////////////////////////////////////////////
908/// Interface function used by test statistics to freeze choice of observables
909/// for interpretation of fraction coefficients. Needed here because a RooSimultaneous
910/// works like a RooAddPdf when plotted
911
917
918
919////////////////////////////////////////////////////////////////////////////////
920/// Interface function used by test statistics to freeze choice of range
921/// for interpretation of fraction coefficients. Needed here because a RooSimultaneous
922/// works like a RooAddPdf when plotted
923
928
929
930
931
932////////////////////////////////////////////////////////////////////////////////
933
935 const RooArgSet* auxProto, bool verbose, bool autoBinned, const char* binnedTag) const
936{
937 const char* idxCatName = _indexCat.arg().GetName() ;
938
939 if (vars.find(idxCatName) && prototype==nullptr
940 && (auxProto==nullptr || auxProto->empty())
941 && (autoBinned || (binnedTag && strlen(binnedTag)))) {
942
943 // Return special generator config that can also do binned generation for selected states
944 return new RooSimSplitGenContext(*this,vars,verbose,autoBinned,binnedTag) ;
945
946 } else {
947
948 // Return regular generator config ;
949 return genContext(vars,prototype,auxProto,verbose) ;
950 }
951}
952
953
954
955////////////////////////////////////////////////////////////////////////////////
956/// Return specialized generator context for simultaneous p.d.f.s
957
959 const RooArgSet* auxProto, bool verbose) const
960{
961 RooArgSet allVars{vars};
962 if(prototype) allVars.add(*prototype->get());
963
966
967 // Not generating index cat: we better error out because it's not clear what
968 // the user expects here. Does the user want to generate according to the
969 // currently-selected pdf? Or does the user want to generate global
970 // observable values according to the union of all category pdfs?
971 // Print an error and tell the user what to do to explicitly.
972 if(catsAmongAllVars.empty()) {
973 coutE(InputArguments) << "RooSimultaneous::generateSimGlobal(" << GetName()
974 << ") asking to generate without the index category!\n"
975 << "It's not clear what to do. you probably want to either:\n"
976 << "\n"
977 << " 1. Generate according to the currently-selected pdf.\n"
978 << " Please do this explicitly with:\n"
979 << " simpdf->getPdf(simpdf->indexCat().getCurrentLabel())->generate(vars, ...)\n"
980 << "\n"
981 << " 1. Generate global observable values according to the union of all component pdfs.\n"
982 << " For this, please use simpdf->generateSimGlobal(vars, ...)\n"
983 << std::endl;
984 return nullptr;
985 }
986
988 if(prototype) {
989 prototype->get()->selectCommon(flattenedCatList(), catsAmongProtoVars);
990
991 if(!catsAmongProtoVars.empty() && catsAmongProtoVars.size() != flattenedCatList().size()) {
992 // Abort if we have only part of the servers
993 coutE(Plotting) << "RooSimultaneous::genContext: ERROR: prototype must include either all "
994 << " components of the RooSimultaneous index category or none " << std::endl;
995 return nullptr;
996 }
997 }
998
999 return new RooSimGenContext(*this,vars,prototype,auxProto,verbose) ;
1000}
1001
1002
1003
1004
1005////////////////////////////////////////////////////////////////////////////////
1006
1008 const RooArgSet* nset,
1009 double scaleFactor,
1011 bool showProgress) const
1012{
1013 if (RooAbsReal::fillDataHist (hist, nset, scaleFactor,
1014 correctForBinVolume, showProgress) == nullptr)
1015 return nullptr;
1016
1017 const double sum = hist->sumEntries();
1018 if (sum != 0) {
1019 for (int i=0 ; i<hist->numEntries() ; i++) {
1020 hist->set(i, hist->weight(i) / sum, 0.);
1021 }
1022 }
1023
1024 return hist;
1025}
1026
1027
1028
1029
1030////////////////////////////////////////////////////////////////////////////////
1031/// Special generator interface for generation of 'global observables' -- for RooStats tools.
1032///
1033/// \note Why one can't just use RooAbsPdf::generate()? That's because when
1034/// using the regular generate() method, a specific component pdf is selected
1035/// for each generated entry according to the index category value. However,
1036/// global observable values are independent of the current index category,
1037/// which can best be illustrated with the case where a global observable
1038/// corresponds to a nuisance parameter that is relevant for multiple channels.
1039/// So the interpretation of what is an entry in the generated dataset is very
1040/// different, hence the separate function.
1041
1043{
1044 // Generating the index category together with the global observables doesn't make any sense.
1047 if(!catsAmongAllVars.empty()) {
1048 coutE(InputArguments) << "RooSimultaneous::generateSimGlobal(" << GetName()
1049 << ") asking to generate global obserables at the same time as the index category!\n"
1050 << "This doesn't make any sense: global observables are generally not related to a specific channel.\n"
1051 << std::endl;
1052 return nullptr;
1053 }
1054
1055 // Make set with clone of variables (placeholder for output)
1057 whatVars.snapshot(globClone);
1058
1059 auto data = std::make_unique<RooDataSet>("gensimglobal","gensimglobal",whatVars);
1060
1061 for (Int_t i=0 ; i<nEvents ; i++) {
1062 for (const auto& nameIdx : indexCat()) {
1063
1064 // Get pdf associated with state from simpdf
1065 RooAbsPdf* pdftmp = getPdf(nameIdx.first);
1066
1068 pdftmp->getObservables(&whatVars, globtmp) ;
1069
1070 // If there are any, generate only global variables defined by the pdf
1071 // associated with this state and transfer values to output placeholder.
1072 if (!globtmp.empty()) {
1073 globClone.assign(*std::unique_ptr<RooDataSet>{pdftmp->generate(globtmp,1)}->get(0)) ;
1074 }
1075 }
1076 data->add(globClone) ;
1077 }
1078
1079 return RooFit::makeOwningPtr(std::move(data));
1080}
1081
1082
1083/// Wraps the components of this RooSimultaneous in RooBinSamplingPdfs.
1084/// \param[in] data The dataset to be used in the eventual fit, used to figure
1085/// out the observables and whether the dataset is binned.
1086/// \param[in] precision Precision argument for all created RooBinSamplingPdfs.
1088
1089 if (precision < 0.) return;
1090
1092
1093 for (auto const &item : this->indexCat()) {
1094
1095 auto const &catName = item.first;
1096 auto &pdf = *this->getPdf(catName);
1097
1098 if (auto newSamplingPdf = RooBinSamplingPdf::create(pdf, data, precision)) {
1099 // Set the "ORIGNAME" attribute the indicate to
1100 // RooAbsArg::redirectServers() which pdf should be replaced by this
1101 // RooBinSamplingPdf in the RooSimultaneous.
1102 newSamplingPdf->setAttribute(
1103 (std::string("ORIGNAME:") + pdf.GetName()).c_str());
1104 newSamplingPdfs.addOwned(std::move(newSamplingPdf));
1105 }
1106 }
1107
1108 this->redirectServers(newSamplingPdfs, false, true);
1109 this->addOwnedComponents(std::move(newSamplingPdfs));
1110}
1111
1112
1113/// Wraps the components of this RooSimultaneous in RooBinSamplingPdfs, with a
1114/// different precision parameter for each component.
1115/// \param[in] data The dataset to be used in the eventual fit, used to figure
1116/// out the observables and whether the dataset is binned.
1117/// \param[in] precisions The map that gives the precision argument for each
1118/// component in the RooSimultaneous. The keys are the pdf names. If
1119/// there is no value for a given component, it will not use the bin
1120/// integration. Otherwise, the value has the same meaning than in
1121/// the IntegrateBins() command argument for RooAbsPdf::fitTo().
1122/// \param[in] useCategoryNames If this flag is set, the category names will be
1123/// used to look up the precision in the precisions map instead of
1124/// the pdf names.
1126 std::map<std::string, double> const& precisions,
1127 bool useCategoryNames /*=false*/) {
1128
1129 constexpr double defaultPrecision = -1.;
1130
1132
1133 for (auto const &item : this->indexCat()) {
1134
1135 auto const &catName = item.first;
1136 auto &pdf = *this->getPdf(catName);
1137 std::string pdfName = pdf.GetName();
1138
1139 auto found = precisions.find(useCategoryNames ? catName : pdfName);
1140 const double precision =
1141 found != precisions.end() ? found->second : defaultPrecision;
1142 if (precision < 0.)
1143 continue;
1144
1145 if (auto newSamplingPdf = RooBinSamplingPdf::create(pdf, data, precision)) {
1146 // Set the "ORIGNAME" attribute the indicate to
1147 // RooAbsArg::redirectServers() which pdf should be replaced by this
1148 // RooBinSamplingPdf in the RooSimultaneous.
1149 newSamplingPdf->setAttribute(
1150 (std::string("ORIGNAME:") + pdf.GetName()).c_str());
1151 newSamplingPdfs.addOwned(std::move(newSamplingPdf));
1152 }
1153 }
1154
1155 this->redirectServers(newSamplingPdfs, false, true);
1156 this->addOwnedComponents(std::move(newSamplingPdfs));
1157}
1158
1159/// Internal utility function to get a list of all category components for this
1160/// RooSimultaneous. The output contains only the index category if it is a
1161/// RooCategory, or the list of all category components if it is a
1162/// RooSuperCategory.
1164{
1165 // Note that the index category of a RooSimultaneous can only be of type
1166 // RooCategory or RooSuperCategory, because these are the only classes that
1167 // inherit from RooAbsCategoryLValue.
1168 if (auto superCat = dynamic_cast<RooSuperCategory const*>(&_indexCat.arg())) {
1169 return superCat->inputCatList();
1170 }
1171
1172 if(!_indexCatSet) {
1173 _indexCatSet = std::make_unique<RooArgSet>(_indexCat.arg());
1174 }
1175 return *_indexCatSet;
1176}
1177
1178namespace {
1179
1180void markObs(RooAbsArg *arg, std::string const &prefix, RooArgSet const &normSet)
1181{
1182 for (RooAbsArg *server : arg->servers()) {
1183 if (server->isFundamental() && normSet.find(*server)) {
1184 markObs(server, prefix, normSet);
1185 server->setAttribute("__obs__");
1186 } else if (!server->isFundamental()) {
1187 markObs(server, prefix, normSet);
1188 }
1189 }
1190}
1191
1192void prefixArgs(RooAbsArg *arg, std::string const &prefix, RooArgSet const &normSet)
1193{
1194 if (!arg->getStringAttribute("__prefix__")) {
1195 arg->SetName((prefix + arg->GetName()).c_str());
1196 arg->setStringAttribute("__prefix__", prefix.c_str());
1197 }
1198 for (RooAbsArg *server : arg->servers()) {
1199 if (server->isFundamental() && normSet.find(*server)) {
1200 prefixArgs(server, prefix, normSet);
1201 } else if (!server->isFundamental()) {
1202 prefixArgs(server, prefix, normSet);
1203 }
1204 }
1205}
1206
1207} // namespace
1208
1209std::unique_ptr<RooAbsArg>
1211{
1212 std::unique_ptr<RooSimultaneous> newSimPdf{static_cast<RooSimultaneous *>(this->Clone())};
1213
1214 const char *rangeName = this->getStringAttribute("RangeName");
1215 bool splitRange = this->getAttribute("SplitRange");
1216
1218 std::vector<std::string> catNames;
1219
1220 for (auto *proxy : static_range_cast<RooRealProxy *>(newSimPdf->_pdfProxyList)) {
1221 catNames.emplace_back(proxy->GetName());
1222 std::string const &catName = catNames.back();
1223 const std::string prefix = "_" + catName + "_";
1224
1225 const std::string origname = proxy->arg().GetName();
1226
1227 auto pdfClone = RooHelpers::cloneTreeWithSameParameters(static_cast<RooAbsPdf const &>(proxy->arg()), &normSet);
1228
1229 markObs(pdfClone.get(), prefix, normSet);
1230
1231 std::unique_ptr<RooArgSet> pdfNormSet{
1232 std::unique_ptr<RooArgSet>(pdfClone->getVariables())->selectByAttrib("__obs__", true)};
1233 std::unique_ptr<RooArgSet> condVarSet{
1234 std::unique_ptr<RooArgSet>(pdfClone->getVariables())->selectByAttrib("__conditional__", true)};
1235
1236 pdfNormSet->remove(*condVarSet, true, true);
1237
1238 if (rangeName) {
1240 }
1241
1243 pdfContext.setLikelihoodMode(ctx.likelihoodMode());
1244 auto *pdfFinal = pdfContext.compile(*pdfClone, *newSimPdf, *pdfNormSet);
1245
1246 // We can only prefix the observables after everything related the
1247 // compiling of the compute graph for the normalization set is done. This
1248 // is because of a subtlety in conditional RooProdPdfs, which stores the
1249 // normalization sets for the individual pdfs in RooArgSets that are
1250 // disconnected from the computation graph, so we have no control over
1251 // them. An alternative would be to use recursive server re-direction,
1252 // but this has more performance overhead.
1253 prefixArgs(pdfFinal, prefix, normSet);
1254
1255 pdfFinal->fixAddCoefNormalization(*pdfNormSet, false);
1256
1257 pdfClone->SetName((std::string("_") + pdfClone->GetName()).c_str());
1258 pdfFinal->addOwnedComponents(std::move(pdfClone));
1259
1260 pdfFinal->setAttribute(("ORIGNAME:" + origname).c_str());
1261 newPdfs.add(*pdfFinal);
1262
1263 // We will remove the old pdf server because we will fill the new ones by
1264 // hand via the creation of new proxies.
1265 newSimPdf->removeServer(const_cast<RooAbsReal &>(proxy->arg()), true);
1266 }
1267
1268 // Replace pdfs with compiled pdfs. Don't use RooAbsArg::redirectServers()
1269 // here, because it doesn't support replacing two servers with the same name
1270 // (it can happen in a RooSimultaneous that two pdfs have the same name).
1271
1272 // First delete old proxies (we have already removed the servers before).
1273 newSimPdf->_pdfProxyList.Delete();
1274
1275 // Recreate the _pdfProxyList with the compiled pdfs
1276 for (std::size_t i = 0; i < newPdfs.size(); ++i) {
1277 const char *label = catNames[i].c_str();
1278 newSimPdf->_pdfProxyList.Add(
1279 new RooRealProxy(label, label, newSimPdf.get(), *static_cast<RooAbsReal *>(newPdfs[i])));
1280 }
1281
1282 ctx.compileServers(*newSimPdf, normSet); // to trigger compiling also the index category
1283
1284 return newSimPdf;
1285}
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define coutI(a)
#define oocoutW(o, a)
#define oocxcoutD(o, a)
#define oocoutI(o, a)
#define coutE(a)
RooTemplateProxy< RooAbsReal > RooRealProxy
Compatibility typedef replacing the old RooRealProxy class.
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char name[80]
Definition TGX11.cxx:110
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
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.
void setStringAttribute(const Text_t *key, const Text_t *value)
Associate string 'value' to this object under key 'key'.
void SetName(const char *name) override
Set the name of the TNamed.
bool addOwnedComponents(const RooAbsCollection &comps)
Take ownership of the contents of 'comps'.
const Text_t * getStringAttribute(const Text_t *key) const
Get string attribute mapped under key 'key'.
bool redirectServers(const RooAbsCollection &newServerList, bool mustReplaceAll=false, bool nameChange=false, bool isRecursionStep=false)
Replace all direct servers of this object with the new servers in newServerList.
bool getAttribute(const Text_t *name) const
Check if a named attribute is set. By default, all attributes are unset.
TObject * Clone(const char *newname=nullptr) const override
Make a clone of an object using the Streamer facility.
Definition RooAbsArg.h:88
Abstract base class for objects that represent a discrete value that can be set from the outside,...
A space to attach TBranches.
bool contains(const char *name) const
Check if collection contains an argument with a specific name.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
virtual bool addOwned(RooAbsArg &var, bool silent=false)
Add an argument and transfer the ownership to the collection.
RooAbsArg * find(const char *name) const
Find object with given name in list.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:56
virtual const RooArgSet * get() const
Definition RooAbsData.h:100
virtual Roo1DTable * table(const RooArgSet &catSet, const char *cuts="", const char *opts="") const
Construct table for product of categories in catSet.
RooFit::OwningPtr< RooAbsData > reduce(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 a reduced copy of this dataset.
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
Abstract base class for generator contexts of RooAbsPdf objects.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
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
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition RooAbsPdf.h:214
@ CanBeExtended
Definition RooAbsPdf.h:208
@ MustBeExtended
Definition RooAbsPdf.h:208
@ CanNotBeExtended
Definition RooAbsPdf.h:208
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.
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 makeProjectionSet(const RooAbsArg *plotVar, const RooArgSet *allVars, RooArgSet &projectedVars, bool silent) const
Utility function for plotOn() that constructs the set of observables to project when plotting ourselv...
Efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:33
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:159
RooArgSet * selectCommon(const RooAbsCollection &refColl) const
Use RooAbsCollection::selecCommon(), but return as RooArgSet.
Definition RooArgSet.h:154
static std::unique_ptr< RooAbsPdf > create(RooAbsPdf &pdf, RooAbsData const &data, double precision)
Creates a wrapping RooBinSamplingPdf if appropriate.
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=nullptr)
Setter function without integration set.
T * getObjByIndex(Int_t index) const
Retrieve payload object by slot index.
Int_t lastIndex() const
Return index of slot used in last get or set operation.
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=nullptr, const TNamed *isetRangeName=nullptr)
Getter function without integration set.
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.
double getDouble(const char *name, double defaultValue=0.0) const
Return double property registered with name 'name'.
bool defineDouble(const char *name, const char *argName, int doubleNum, double defValue=0.0)
Define double property name 'name' mapped to double in slot 'doubleNum' in RooCmdArg with name argNam...
RooArgSet * getSet(const char *name, RooArgSet *set=nullptr) const
Return RooArgSet property registered with name 'name'.
bool defineSet(const char *name, const char *argName, int setNum, const RooArgSet *set=nullptr)
Define TObject property name 'name' mapped to object in slot 'setNum' in RooCmdArg with name argName ...
bool ok(bool verbose) const
Return true of parsing was successful.
bool defineObject(const char *name, const char *argName, int setNum, const TObject *obj=nullptr, bool isArray=false)
Define TObject property name 'name' mapped to object in slot 'setNum' in RooCmdArg with name argName ...
const char * getString(const char *name, const char *defaultValue="", bool convEmptyToNull=false) const
Return string property registered with name 'name'.
bool defineString(const char *name, const char *argName, int stringNum, const char *defValue="", bool appendMode=false)
Define double property name 'name' mapped to double in slot 'stringNum' in RooCmdArg with name argNam...
const RooLinkedList & getObjectList(const char *name) const
Return list of objects registered with name 'name'.
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'.
void removeAll() override
Remove all argument inset using remove(const RooAbsArg&).
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
double weight(std::size_t i) const
Return weight of i-th bin.
void set(std::size_t binNumber, double weight, double wgtErr)
Set bin content of bin that was last loaded with get(std::size_t).
double sumEntries() const override
Sum the weights of all bins.
Container class to hold unbinned data.
Definition RooDataSet.h:32
void compileServers(RooAbsArg &arg, RooArgSet const &normSet)
Collection class for internal use, storing a collection of RooAbsArg pointers in a doubly linked list...
static const TNamed * ptr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
const RooArgSet * getNormVars() const
Definition RooPlot.h:146
RooAbsRealLValue * getPlotVar() const
Definition RooPlot.h:137
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
void selectNormalization(const RooArgSet *depSet=nullptr, bool force=false) override
Interface function used by test statistics to freeze choice of observables for interpretation of frac...
double evaluate() const override
Return the current value: the value of the PDF associated with the current index category state.
Int_t _numPdf
Number of registered PDFs.
void selectNormalizationRange(const char *rangeName=nullptr, bool force=false) override
Interface function used by test statistics to freeze choice of range for interpretation of fraction c...
TList _pdfProxyList
List of PDF proxies (named after applicable category state)
RooObjCacheManager _partIntMgr
! Component normalization manager
~RooSimultaneous() override
Destructor.
RooFit::OwningPtr< RooDataSet > generateSimGlobal(const RooArgSet &whatVars, Int_t nEvents) override
Special generator interface for generation of 'global observables' – for RooStats tools.
RooArgSet const & flattenedCatList() const
Internal utility function to get a list of all category components for this RooSimultaneous.
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
ExtendMode extendMode() const override
Examine the pdf components and check if one of them can be extended or must be extended.
RooCategoryProxy _indexCat
Index category.
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Return analytical integration defined by given code.
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &numVars, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Forward determination of analytical integration capabilities to component p.d.f.s A unique code is as...
double expectedEvents(const RooArgSet *nset) const override
Return the number of expected events: If the index is in nset, then return the sum of the expected ev...
friend class RooSimGenContext
RooAbsGenContext * autoGenContext(const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool verbose=false, bool autoBinned=true, const char *binnedTag="") const override
RooAbsPdf * getPdf(RooStringView catName) const
Return the p.d.f associated with the given index category name.
std::unique_ptr< RooAbsArg > compileForNormSet(RooArgSet const &normSet, RooFit::Detail::CompileContext &ctx) const override
static std::unique_ptr< RooSimultaneous::InitializationOutput > initialize(std::string const &name, RooAbsCategoryLValue &inIndexCat, std::map< std::string, RooAbsPdf * > const &pdfMap)
void wrapPdfsInBinSamplingPdfs(RooAbsData const &data, double precision)
Wraps the components of this RooSimultaneous in RooBinSamplingPdfs.
const TNamed * _plotCoefNormRange
RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool verbose=false) const override
Return specialized generator context for simultaneous p.d.f.s.
std::unique_ptr< RooArgSet > _indexCatSet
! Index category wrapped in a RooArgSet if needed internally
bool addPdf(const RooAbsPdf &pdf, const char *catLabel)
Associate given PDF with index category state label 'catLabel'.
RooSetProxy _plotCoefNormSet
const RooAbsCategoryLValue & indexCat() const
friend class RooSimSplitGenContext
virtual RooDataHist * fillDataHist(RooDataHist *hist, const RooArgSet *nset, double scaleFactor, bool correctForBinVolume=false, bool showProgress=false) const
The RooStringView is a wrapper around a C-style string that can also be constructed from a std::strin...
Joins several RooAbsCategoryLValue objects into a single category.
const char * label() const
Get the label of the current category state. This function only makes sense for category proxies.
const T & arg() const
Return reference to object held in proxy.
TObject * FindObject(const char *name) const override
Find an object in this list using its name.
Definition TList.cxx:708
void Add(TObject *obj) override
Definition TList.h:81
void Delete(Option_t *option="") override
Remove all objects from the list AND delete all heap based objects.
Definition TList.cxx:600
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
Mother of all ROOT objects.
Definition TObject.h:41
virtual const char * GetName() const
Returns name of object.
Definition TObject.cxx:457
virtual Int_t IndexOf(const TObject *obj) const
Return index of object in collection.
Basic string class.
Definition TString.h:138
RooCmdArg SelectVars(const RooArgSet &vars)
RooCmdArg ProjWData(const RooAbsData &projData, bool binData=false)
RooCmdArg Project(const RooArgSet &projSet)
RooCmdArg Normalization(double scaleFactor)
RooCmdArg Cut(const char *cutSpec)
std::vector< std::string > Split(std::string_view str, std::string_view delims, bool skipEmpty=false)
Splits a string at each character in delims.
std::string makeSliceCutString(RooArgSet const &sliceDataSet)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:67
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
std::unique_ptr< T > cloneTreeWithSameParameters(T const &arg, RooArgSet const *observables=nullptr)
Clone RooAbsArg object and reattach to original parameters.
std::string getRangeNameForSimComponent(std::string const &rangeName, bool splitRange, std::string const &catName)
Internal struct used for initialization.
std::vector< RooAbsPdf const * > finalPdfs
std::vector< std::string > finalCatLabels
void addPdf(const RooAbsPdf &pdf, std::string const &catLabel)
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2339