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;
81 auto indexCatIt = inIndexCat.begin();
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
90} // namespace
91
93
94void RooSimultaneous::InitializationOutput::addPdf(const RooAbsPdf &pdf, std::string const &catLabel)
95{
96 finalPdfs.push_back(&pdf);
97 finalCatLabels.emplace_back(catLabel);
98}
99
100using std::string, std::endl;
101
103
104
105
106////////////////////////////////////////////////////////////////////////////////
107/// Constructor with index category. PDFs associated with indexCat
108/// states can be added after construction with the addPdf() function.
109///
110/// RooSimultaneous can function without having a PDF associated
111/// with every single state. The normalization in such cases is taken
112/// from the number of registered PDFs, but getVal() will assert if
113/// when called for an unregistered index state.
114
115RooSimultaneous::RooSimultaneous(const char *name, const char *title,
116 RooAbsCategoryLValue& inIndexCat) :
117 RooSimultaneous{name, title, std::map<std::string, RooAbsPdf*>{}, inIndexCat}
118{
119}
120
121
122////////////////////////////////////////////////////////////////////////////////
123/// Constructor from index category and full list of PDFs.
124/// In this constructor form, a PDF must be supplied for each indexCat state
125/// to avoid ambiguities. The PDFs are associated with the states of the
126/// index category as they appear when iterating through the category states
127/// with RooAbsCategory::begin() and RooAbsCategory::end(). This usually means
128/// they are associated by ascending index numbers.
129///
130/// PDFs may not overlap (i.e. share any variables) with the index category (function)
131
132RooSimultaneous::RooSimultaneous(const char *name, const char *title,
133 const RooArgList& inPdfList, RooAbsCategoryLValue& inIndexCat) :
134 RooSimultaneous{name, title, createPdfMap(inPdfList, inIndexCat), inIndexCat}
135{
136 if (inPdfList.size() != inIndexCat.size()) {
137 std::stringstream errMsg;
138 errMsg << "RooSimultaneous::ctor(" << GetName()
139 << " ERROR: Number PDF list entries must match number of index category states, no PDFs added";
140 coutE(InputArguments) << errMsg.str() << std::endl;
141 throw std::invalid_argument(errMsg.str());
142 }
143}
144
145
146////////////////////////////////////////////////////////////////////////////////
147
148RooSimultaneous::RooSimultaneous(const char *name, const char *title, std::map<string, RooAbsPdf *> pdfMap,
149 RooAbsCategoryLValue &inIndexCat)
150 : RooSimultaneous(name, title, std::move(*initialize(name ? name : "", inIndexCat, pdfMap)))
151{
152}
153
154/// For internal use in RooFit.
155RooSimultaneous::RooSimultaneous(const char *name, const char *title,
157 RooAbsCategoryLValue &inIndexCat)
158 : RooSimultaneous(name, title, RooFit::Detail::flatMapToStdMap(pdfMap), inIndexCat)
159{
160}
161
163 : RooAbsPdf(name, title),
164 _plotCoefNormSet("!plotCoefNormSet", "plotCoefNormSet", this, false, false),
165 _partIntMgr(this, 10),
166 _indexCat("indexCat", "Index category", this, *initInfo.indexCat)
167{
168 for (std::size_t i = 0; i < initInfo.finalPdfs.size(); ++i) {
169 addPdf(*initInfo.finalPdfs[i], initInfo.finalCatLabels[i].c_str());
170 }
171
172 // Take ownership of eventual super category
173 if (initInfo.superIndex) {
174 addOwnedComponents(std::move(initInfo.superIndex));
175 }
176}
177
178/// \cond ROOFIT_INTERNAL
179
180// This class cannot be locally defined in initialize as it cannot be
181// used as a template argument in that case
182namespace RooSimultaneousAux {
183 struct CompInfo {
184 RooAbsPdf* pdf ;
185 RooSimultaneous* simPdf ;
186 const RooAbsCategoryLValue* subIndex ;
187 std::unique_ptr<RooArgSet> subIndexComps;
188 } ;
189}
190
191/// \endcond
192
193std::unique_ptr<RooSimultaneous::InitializationOutput>
195 std::map<std::string, RooAbsPdf *> const& pdfMap)
196
197{
198 auto out = std::make_unique<RooSimultaneous::InitializationOutput>();
199 out->indexCat = &inIndexCat;
200
201 // First see if there are any RooSimultaneous input components
202 bool simComps(false) ;
203 for (auto const& item : pdfMap) {
204 if (dynamic_cast<RooSimultaneous*>(item.second)) {
205 simComps = true ;
206 break ;
207 }
208 }
209
210 // If there are no simultaneous component p.d.f. do simple processing through addPdf()
211 if (!simComps) {
212 for (auto const& item : pdfMap) {
213 out->addPdf(*item.second,item.first);
214 }
215 return out;
216 }
217
218 std::string msgPrefix = "RooSimultaneous::initialize(" + name + ") ";
219
220 // Issue info message that we are about to do some rearranging
221 oocoutI(nullptr, InputArguments) << msgPrefix << "INFO: one or more input component of simultaneous p.d.f.s are"
222 << " simultaneous p.d.f.s themselves, rewriting composite expressions as one-level simultaneous p.d.f. in terms of"
223 << " final constituents and extended index category" << std::endl;
224
225
226 RooArgSet allAuxCats ;
227 std::map<string,RooSimultaneousAux::CompInfo> compMap ;
228 for (auto const& item : pdfMap) {
229 RooSimultaneousAux::CompInfo ci ;
230 ci.pdf = item.second ;
231 RooSimultaneous* simComp = dynamic_cast<RooSimultaneous*>(item.second) ;
232 if (simComp) {
233 ci.simPdf = simComp ;
234 ci.subIndex = &simComp->indexCat() ;
235 ci.subIndexComps = simComp->indexCat().isFundamental()
236 ? std::make_unique<RooArgSet>(simComp->indexCat())
237 : std::unique_ptr<RooArgSet>(simComp->indexCat().getVariables());
238 allAuxCats.add(*ci.subIndexComps,true) ;
239 } else {
240 ci.simPdf = nullptr;
241 ci.subIndex = nullptr;
242 }
243 compMap[item.first] = std::move(ci);
244 }
245
246 // Construct the 'superIndex' from the nominal index category and all auxiliary components
247 RooArgSet allCats(inIndexCat) ;
248 allCats.add(allAuxCats) ;
249 std::string siname = name + "_index";
250 out->superIndex = std::make_unique<RooSuperCategory>(siname.c_str(),siname.c_str(),allCats) ;
251 auto *superIndex = out->superIndex.get();
252 out->indexCat = superIndex;
253
254 // Now process each of original pdf/state map entries
255 for (auto const& citem : compMap) {
256
257 RooArgSet repliCats(allAuxCats) ;
258 if (citem.second.subIndexComps) {
259 repliCats.remove(*citem.second.subIndexComps) ;
260 }
261 inIndexCat.setLabel(citem.first.c_str()) ;
262
263 if (!citem.second.simPdf) {
264
265 // Entry is a plain p.d.f. assign it to every state permutation of the repliCats set
266 RooSuperCategory repliSuperCat("tmp","tmp",repliCats) ;
267
268 // Iterator over all states of repliSuperCat
269 for (const auto& nameIdx : repliSuperCat) {
270 // Set value
271 repliSuperCat.setLabel(nameIdx.first) ;
272 // Retrieve corresponding label of superIndex
273 string superLabel = superIndex->getCurrentLabel() ;
274 out->addPdf(*citem.second.pdf,superLabel);
275 oocxcoutD(static_cast<RooAbsArg*>(nullptr), InputArguments) << msgPrefix
276 << "assigning pdf " << citem.second.pdf->GetName() << " to super label " << superLabel << endl ;
277 }
278 } else {
279
280 // Entry is a simultaneous p.d.f
281
282 if (repliCats.empty()) {
283
284 // Case 1 -- No replication of components of RooSim component are required
285
286 for (const auto& type : *citem.second.subIndex) {
287 const_cast<RooAbsCategoryLValue*>(citem.second.subIndex)->setLabel(type.first.c_str());
288 string superLabel = superIndex->getCurrentLabel() ;
289 RooAbsPdf* compPdf = citem.second.simPdf->getPdf(type.first);
290 if (compPdf) {
291 out->addPdf(*compPdf,superLabel);
292 oocxcoutD(static_cast<RooAbsArg*>(nullptr), InputArguments) << msgPrefix
293 << "assigning pdf " << compPdf->GetName() << "(member of " << citem.second.pdf->GetName()
294 << ") to super label " << superLabel << endl ;
295 } else {
296 oocoutW(nullptr, InputArguments) << msgPrefix << "WARNING: No p.d.f. associated with label "
297 << type.second << " for component RooSimultaneous p.d.f " << citem.second.pdf->GetName()
298 << "which is associated with master index label " << citem.first << endl ;
299 }
300 }
301
302 } else {
303
304 // Case 2 -- Replication of components of RooSim component are required
305
306 // Make replication supercat
307 RooSuperCategory repliSuperCat("tmp","tmp",repliCats) ;
308
309 for (const auto& stype : *citem.second.subIndex) {
310 const_cast<RooAbsCategoryLValue*>(citem.second.subIndex)->setLabel(stype.first.c_str());
311
312 for (const auto& nameIdx : repliSuperCat) {
313 repliSuperCat.setLabel(nameIdx.first) ;
314 const string superLabel = superIndex->getCurrentLabel() ;
315 RooAbsPdf* compPdf = citem.second.simPdf->getPdf(stype.first);
316 if (compPdf) {
317 out->addPdf(*compPdf,superLabel);
318 oocxcoutD(static_cast<RooAbsArg*>(nullptr), InputArguments) << msgPrefix
319 << "assigning pdf " << compPdf->GetName() << "(member of " << citem.second.pdf->GetName()
320 << ") to super label " << superLabel << endl ;
321 } else {
322 oocoutW(nullptr, InputArguments) << msgPrefix << "WARNING: No p.d.f. associated with label "
323 << stype.second << " for component RooSimultaneous p.d.f " << citem.second.pdf->GetName()
324 << "which is associated with master index label " << citem.first << endl ;
325 }
326 }
327 }
328 }
329 }
330 }
331
332 return out;
333}
334
335
336////////////////////////////////////////////////////////////////////////////////
337/// Copy constructor
338
340 RooAbsPdf(other,name),
341 _plotCoefNormSet("!plotCoefNormSet",this,other._plotCoefNormSet),
342 _plotCoefNormRange(other._plotCoefNormRange),
343 _partIntMgr(other._partIntMgr,this),
344 _indexCat("indexCat",this,other._indexCat),
345 _numPdf(other._numPdf)
346{
347 // Copy proxy list
348 for(auto* proxy : static_range_cast<RooRealProxy*>(other._pdfProxyList)) {
349 _pdfProxyList.Add(new RooRealProxy(proxy->GetName(),this,*proxy)) ;
350 }
351}
352
353
354
355////////////////////////////////////////////////////////////////////////////////
356/// Destructor
357
359{
361}
362
363
364
365////////////////////////////////////////////////////////////////////////////////
366/// Return the p.d.f associated with the given index category name
367
369{
370 RooRealProxy* proxy = static_cast<RooRealProxy*>(_pdfProxyList.FindObject(catName.c_str()));
371 return proxy ? static_cast<RooAbsPdf*>(proxy->absArg()) : nullptr;
372}
373
374
375
376////////////////////////////////////////////////////////////////////////////////
377/// Associate given PDF with index category state label 'catLabel'.
378/// The name state must be already defined in the index category.
379///
380/// RooSimultaneous can function without having a PDF associated
381/// with every single state. The normalization in such cases is taken
382/// from the number of registered PDFs, but getVal() will fail if
383/// called for an unregistered index state.
384///
385/// PDFs may not overlap (i.e. share any variables) with the index category (function).
386/// \param[in] pdf PDF to be added.
387/// \param[in] catLabel Name of the category state to be associated to the PDF.
388/// \return `true` in case of failure.
389
390bool RooSimultaneous::addPdf(const RooAbsPdf& pdf, const char* catLabel)
391{
392 // PDFs cannot overlap with the index category
393 if (pdf.dependsOn(_indexCat.arg())) {
394 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName() << "): PDF '" << pdf.GetName()
395 << "' overlaps with index category '" << _indexCat.arg().GetName() << "'."<< endl ;
396 return true ;
397 }
398
399 // Each index state can only have one PDF associated with it
400 if (_pdfProxyList.FindObject(catLabel)) {
401 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName() << "): index state '"
402 << catLabel << "' has already an associated PDF." << endl ;
403 return true ;
404 }
405
406 const RooSimultaneous* simPdf = dynamic_cast<const RooSimultaneous*>(&pdf) ;
407 if (simPdf) {
408
409 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName()
410 << ") ERROR: you cannot add a RooSimultaneous component to a RooSimultaneous using addPdf()."
411 << " Use the constructor with RooArgList if input p.d.f.s or the map<string,RooAbsPdf&> instead." << endl ;
412 return true ;
413
414 } else {
415
416 // Create a proxy named after the associated index state
417 TObject* proxy = new RooRealProxy(catLabel,catLabel,this,const_cast<RooAbsPdf&>(pdf));
418 _pdfProxyList.Add(proxy) ;
419 _numPdf += 1 ;
420 }
421
422 return false ;
423}
424
425
426
427
428
429////////////////////////////////////////////////////////////////////////////////
430/// Examine the pdf components and check if one of them can be extended or must be extended.
431/// It is enough to have one component that can be extended or must be extended to return the flag in
432/// the total simultaneous pdf.
433
435{
436 bool anyCanExtend(false) ;
437 bool anyMustExtend(false) ;
438
439 for (Int_t i=0 ; i<_numPdf ; i++) {
440 RooRealProxy* proxy = static_cast<RooRealProxy*>(_pdfProxyList.At(i));
441 if (proxy) {
442 RooAbsPdf* pdf = static_cast<RooAbsPdf*>(proxy->absArg()) ;
443 //cout << " now processing pdf " << pdf->GetName() << endl;
444 if (pdf->canBeExtended()) {
445 //cout << "RooSim::extendedMode(" << GetName() << ") component " << pdf->GetName() << " can be extended"
446 // << endl;
447 anyCanExtend = true;
448 }
449 if (pdf->mustBeExtended()) {
450 //cout << "RooSim::extendedMode(" << GetName() << ") component " << pdf->GetName() << " MUST be extended" << endl;
451 anyMustExtend = true;
452 }
453 }
454 }
455 if (anyMustExtend) {
456 //cout << "RooSim::extendedMode(" << GetName() << ") returning MustBeExtended" << endl ;
457 return MustBeExtended ;
458 }
459 if (anyCanExtend) {
460 //cout << "RooSim::extendedMode(" << GetName() << ") returning CanBeExtended" << endl ;
461 return CanBeExtended ;
462 }
463 //cout << "RooSim::extendedMode(" << GetName() << ") returning CanNotBeExtended" << endl ;
464 return CanNotBeExtended ;
465}
466
467
468
469
470////////////////////////////////////////////////////////////////////////////////
471/// Return the current value:
472/// the value of the PDF associated with the current index category state
473
475{
476 // Retrieve the proxy by index name
478
479 //assert(proxy!=0) ;
480 if (proxy==nullptr) return 0 ;
481
482 // Calculate relative weighting factor for sim-pdfs of all extendable components
483 double catFrac(1) ;
484 if (canBeExtended()) {
485 double nEvtCat = (static_cast<RooAbsPdf*>(proxy->absArg()))->expectedEvents(_normSet) ;
486
487 double nEvtTot(0) ;
488 for(auto * proxy2 : static_range_cast<RooRealProxy*>(_pdfProxyList)) {
489 nEvtTot += (static_cast<RooAbsPdf*>(proxy2->absArg()))->expectedEvents(_normSet) ;
490 }
491 catFrac=nEvtCat/nEvtTot ;
492 }
493
494 // Return the selected PDF value, normalized by the number of index states
495 return (static_cast<RooAbsPdf*>(proxy->absArg()))->getVal(_normSet)*catFrac ;
496}
497
498
499
500////////////////////////////////////////////////////////////////////////////////
501/// Return the number of expected events: If the index is in nset,
502/// then return the sum of the expected events of all components,
503/// otherwise return the number of expected events of the PDF
504/// associated with the current index category state
505
507{
508 if (nset->contains(_indexCat.arg())) {
509
510 double sum(0) ;
511
512 for(auto * proxy : static_range_cast<RooRealProxy*>(_pdfProxyList)) {
513 sum += (static_cast<RooAbsPdf*>(proxy->absArg()))->expectedEvents(nset) ;
514 }
515
516 return sum ;
517
518 } else {
519
520 // Retrieve the proxy by index name
522
523 //assert(proxy!=0) ;
524 if (proxy==nullptr) return 0 ;
525
526 // Return the selected PDF value, normalized by the number of index states
527 return (static_cast<RooAbsPdf*>(proxy->absArg()))->expectedEvents(nset);
528 }
529}
530
531
532
533////////////////////////////////////////////////////////////////////////////////
534/// Forward determination of analytical integration capabilities to component p.d.f.s
535/// A unique code is assigned to the combined integration capabilities of all associated
536/// p.d.f.s
537
539 const RooArgSet* normSet, const char* rangeName) const
540{
541 // Declare that we can analytically integrate all requested observables
542 analVars.add(allVars) ;
543
544 // Retrieve (or create) the required partial integral list
545 Int_t code ;
546
547 // Check if this configuration was created before
548 CacheElem* cache = static_cast<CacheElem*>(_partIntMgr.getObj(normSet,&analVars,nullptr,RooNameReg::ptr(rangeName))) ;
549 if (cache) {
550 code = _partIntMgr.lastIndex() ;
551 return code+1 ;
552 }
553 cache = new CacheElem ;
554
555 // Create the partial integral set for this request
556 for(auto * proxy : static_range_cast<RooRealProxy*>(_pdfProxyList)) {
557 cache->_partIntList.addOwned(std::unique_ptr<RooAbsReal>{proxy->arg().createIntegral(analVars,normSet,nullptr,rangeName)});
558 }
559
560 // Store the partial integral list and return the assigned code ;
561 code = _partIntMgr.setObj(normSet,&analVars,cache,RooNameReg::ptr(rangeName)) ;
562
563 return code+1 ;
564}
565
566
567
568////////////////////////////////////////////////////////////////////////////////
569/// Return analytical integration defined by given code
570
571double RooSimultaneous::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* /*rangeName*/) const
572{
573 // No integration scenario
574 if (code==0) {
575 return getVal(normSet) ;
576 }
577
578 // Partial integration scenarios, rangeName already encoded in 'code'
579 CacheElem* cache = static_cast<CacheElem*>(_partIntMgr.getObjByIndex(code-1)) ;
580
582 Int_t idx = _pdfProxyList.IndexOf(proxy) ;
583 return (static_cast<RooAbsReal*>(cache->_partIntList.at(idx)))->getVal(normSet) ;
584}
585
586
587
588
589
590
591////////////////////////////////////////////////////////////////////////////////
592/// Back-end for plotOn() implementation on RooSimultaneous which
593/// needs special handling because a RooSimultaneous PDF cannot
594/// project out its index category via integration. plotOn() will
595/// abort if this is requested without providing a projection dataset.
596
598{
599 // Sanity checks
600 if (plotSanityChecks(frame)) return frame ;
601
602 // Extract projection configuration from command list
603 RooCmdConfig pc("RooSimultaneous::plotOn(" + std::string(GetName()) + ")");
604 pc.defineString("sliceCatState","SliceCat",0,"",true) ;
605 pc.defineDouble("scaleFactor","Normalization",0,1.0) ;
606 pc.defineInt("scaleType","Normalization",0,RooAbsPdf::Relative) ;
607 pc.defineObject("sliceCatList","SliceCat",0,nullptr,true) ;
608 // This dummy is needed for plotOn to recognize the "SliceCatMany" command.
609 // It is not used directly, but the "SliceCat" commands are nested in it.
610 // Removing this dummy definition results in "ERROR: unrecognized command: SliceCatMany".
611 pc.defineObject("dummy1","SliceCatMany",0) ;
612 pc.defineSet("projSet","Project",0) ;
613 pc.defineSet("sliceSet","SliceVars",0) ;
614 pc.defineSet("projDataSet","ProjData",0) ;
615 pc.defineObject("projData","ProjData",1) ;
616 pc.defineMutex("Project","SliceVars") ;
617 pc.allowUndefined() ; // there may be commands we don't handle here
618
619 // Process and check varargs
620 pc.process(cmdList) ;
621 if (!pc.ok(true)) {
622 return frame ;
623 }
624
625 const RooAbsData* projData = static_cast<const RooAbsData*>(pc.getObject("projData")) ;
626 const RooArgSet* projDataSet = pc.getSet("projDataSet");
627 const RooArgSet* sliceSetTmp = pc.getSet("sliceSet") ;
628 std::unique_ptr<RooArgSet> sliceSet( sliceSetTmp ? (static_cast<RooArgSet*>(sliceSetTmp->Clone())) : nullptr );
629 const RooArgSet* projSet = pc.getSet("projSet") ;
630 double scaleFactor = pc.getDouble("scaleFactor") ;
631 ScaleType stype = (ScaleType) pc.getInt("scaleType") ;
632
633
634 // Look for category slice arguments and add them to the master slice list if found
635 const char* sliceCatState = pc.getString("sliceCatState",nullptr,true) ;
636 const RooLinkedList& sliceCatList = pc.getObjectList("sliceCatList") ;
637 if (sliceCatState) {
638
639 // Make the master slice set if it doesnt exist
640 if (!sliceSet) {
641 sliceSet = std::make_unique<RooArgSet>();
642 }
643
644 // Prepare comma separated label list for parsing
645 auto catTokens = ROOT::Split(sliceCatState, ",");
646
647 // Loop over all categories provided by (multiple) Slice() arguments
648 unsigned int tokenIndex = 0;
649 for(auto * scat : static_range_cast<RooCategory*>(sliceCatList)) {
650 const char* slabel = tokenIndex >= catTokens.size() ? nullptr : catTokens[tokenIndex++].c_str();
651
652 if (slabel) {
653 // Set the slice position to the value indicated by slabel
654 scat->setLabel(slabel) ;
655 // Add the slice category to the master slice set
656 sliceSet->add(*scat,false) ;
657 }
658 }
659 }
660
661 // Check if we have a projection dataset
662 if (!projData) {
663 coutE(InputArguments) << "RooSimultaneous::plotOn(" << GetName() << ") ERROR: must have a projection dataset for index category" << endl ;
664 return frame ;
665 }
666
667 // Make list of variables to be projected
668 RooArgSet projectedVars ;
669 if (sliceSet) {
670 makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,true) ;
671
672 // Take out the sliced variables
673 for (const auto sliceArg : *sliceSet) {
674 RooAbsArg* arg = projectedVars.find(sliceArg->GetName()) ;
675 if (arg) {
676 projectedVars.remove(*arg) ;
677 } else {
678 coutI(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") slice variable "
679 << sliceArg->GetName() << " was not projected anyway" << endl ;
680 }
681 }
682 } else if (projSet) {
683 makeProjectionSet(frame->getPlotVar(),projSet,projectedVars,false) ;
684 } else {
685 makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,true) ;
686 }
687
688 bool projIndex(false) ;
689
690 if (!_indexCat.arg().isDerived()) {
691 // *** Error checking for a fundamental index category ***
692 //cout << "RooSim::plotOn: index is fundamental" << endl ;
693
694 // Check that the provided projection dataset contains our index variable
695 if (!projData->get()->find(_indexCat.arg().GetName())) {
696 coutE(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") ERROR: Projection over index category "
697 << "requested, but projection data set doesn't contain index category" << endl ;
698 return frame ;
699 }
700
701 if (projectedVars.find(_indexCat.arg().GetName())) {
702 projIndex=true ;
703 }
704
705 } else {
706 // *** Error checking for a composite index category ***
707
708 // Determine if any servers of the index category are in the projectedVars
709 RooArgSet projIdxServers ;
710 bool anyServers(false) ;
711 for (const auto server : flattenedCatList()) {
712 if (projectedVars.find(server->GetName())) {
713 anyServers=true ;
714 projIdxServers.add(*server) ;
715 }
716 }
717
718 // Check that the projection dataset contains all the
719 // index category components we're projecting over
720
721 // Determine if all projected servers of the index category are in the projection dataset
722 bool allServers(true) ;
723 std::string missing;
724 for (const auto server : projIdxServers) {
725 if (!projData->get()->find(server->GetName())) {
726 allServers=false ;
727 missing = server->GetName();
728 }
729 }
730
731 if (!allServers) {
732 coutE(Plotting) << "RooSimultaneous::plotOn(" << GetName()
733 << ") ERROR: Projection dataset doesn't contain complete set of index categories to do projection."
734 << "\n\tcategory " << missing << " is missing." << endl ;
735 return frame ;
736 }
737
738 if (anyServers) {
739 projIndex = true ;
740 }
741 }
742
743 // Calculate relative weight fractions of components
744 std::unique_ptr<Roo1DTable> wTable( projData->table(_indexCat.arg()) );
745
746 // Clone the index category to be able to cycle through the category states for plotting without
747 // affecting the category state of our instance
748 RooArgSet idxCloneSet;
749 RooArgSet(*_indexCat).snapshot(idxCloneSet, true);
750 auto idxCatClone = static_cast<RooAbsCategoryLValue*>(idxCloneSet.find(_indexCat->GetName()) );
751 assert(idxCatClone);
752
753 // Make list of category columns to exclude from projection data
754 std::unique_ptr<RooArgSet> idxCompSliceSet( idxCatClone->getObservables(frame->getNormVars()) );
755
756 // If we don't project over the index, just do the regular plotOn
757 if (!projIndex) {
758
759 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
760 << " represents a slice in the index category (" << _indexCat.arg().GetName() << ")" << endl ;
761
762 // Reduce projData: take out fitCat (component) columns and entries that don't match selected slice
763 // Construct cut string to only select projection data event that match the current slice
764
765 // Make cut string to exclude rows from projection data
766 TString cutString ;
767 bool first(true) ;
768 for (const auto arg : *idxCompSliceSet) {
769 auto idxComp = static_cast<RooCategory*>(arg);
770 RooAbsArg* slicedComponent = nullptr;
771 if (sliceSet && (slicedComponent = sliceSet->find(*idxComp)) != nullptr) {
772 auto theCat = static_cast<const RooAbsCategory*>(slicedComponent);
773 idxComp->setIndex(theCat->getCurrentIndex(), false);
774 }
775
776 if (!first) {
777 cutString.Append("&&") ;
778 } else {
779 first=false ;
780 }
781 cutString.Append(Form("%s==%d",idxComp->GetName(),idxComp->getCurrentIndex())) ;
782 }
783
784 // Make temporary projData without RooSim index category components
785 RooArgSet projDataVars(*projData->get()) ;
786 projDataVars.remove(*idxCompSliceSet,true,true) ;
787
788 std::unique_ptr<RooAbsData> projDataTmp( const_cast<RooAbsData*>(projData)->reduce(projDataVars,cutString) );
789
790 // Override normalization and projection dataset
791 RooCmdArg tmp1 = RooFit::Normalization(scaleFactor*wTable->getFrac(idxCatClone->getCurrentLabel()),stype) ;
792 RooCmdArg tmp2 = RooFit::ProjWData(*projDataSet,*projDataTmp) ;
793
794 // WVE -- do not adjust normalization for asymmetry plots
795 RooLinkedList cmdList2(cmdList) ;
796 if (!cmdList.find("Asymmetry")) {
797 cmdList2.Add(&tmp1) ;
798 }
799 cmdList2.Add(&tmp2) ;
800
801 // Plot single component
802 RooPlot* retFrame = getPdf(idxCatClone->getCurrentLabel())->plotOn(frame,cmdList2);
803 return retFrame ;
804 }
805
806 // If we project over the index, plot using a temporary RooAddPdf
807 // using the weights from the data as coefficients
808
809 // Build the list of indexCat components that are sliced
810 idxCompSliceSet->remove(projectedVars,true,true) ;
811
812 // Make a new expression that is the weighted sum of requested components
813 RooArgList pdfCompList ;
814 RooArgList wgtCompList ;
815//RooAbsPdf* pdf ;
816 double sumWeight(0) ;
817 for(auto * proxy : static_range_cast<RooRealProxy*>(_pdfProxyList)) {
818
819 idxCatClone->setLabel(proxy->name()) ;
820
821 // Determine if this component is the current slice (if we slice)
822 bool skip(false) ;
823 for (const auto idxSliceCompArg : *idxCompSliceSet) {
824 const auto idxSliceComp = static_cast<RooAbsCategory*>(idxSliceCompArg);
825 RooAbsCategory* idxComp = static_cast<RooAbsCategory*>(idxCloneSet.find(idxSliceComp->GetName())) ;
826 if (idxComp->getCurrentIndex()!=idxSliceComp->getCurrentIndex()) {
827 skip=true ;
828 break ;
829 }
830 }
831 if (skip) continue ;
832
833 // Instantiate a RRV holding this pdfs weight fraction
834 wgtCompList.addOwned(std::make_unique<RooRealVar>(proxy->name(),"coef",wTable->getFrac(proxy->name())));
835 sumWeight += wTable->getFrac(proxy->name()) ;
836
837 // Add the PDF to list list
838 pdfCompList.add(proxy->arg()) ;
839 }
840
841 TString plotVarName(GetName()) ;
842 RooAddPdf plotVar{plotVarName,"weighted sum of RS components",pdfCompList,wgtCompList};
843
844 // Fix appropriate coefficient normalization in plot function
845 if (!_plotCoefNormSet.empty()) {
846 plotVar.fixAddCoefNormalization(_plotCoefNormSet) ;
847 }
848
849 std::unique_ptr<RooAbsData> projDataTmp;
850 RooArgSet projSetTmp ;
851 if (projData) {
852
853 // Construct cut string to only select projection data event that match the current slice
854 TString cutString ;
855 if (!idxCompSliceSet->empty()) {
856 bool first(true) ;
857 for (const auto idxSliceCompArg : *idxCompSliceSet) {
858 const auto idxSliceComp = static_cast<RooAbsCategory*>(idxSliceCompArg);
859 if (!first) {
860 cutString.Append("&&") ;
861 } else {
862 first=false ;
863 }
864 cutString.Append(Form("%s==%d",idxSliceComp->GetName(),idxSliceComp->getCurrentIndex())) ;
865 }
866 }
867
868 // Make temporary projData without RooSim index category components
869 RooArgSet projDataVars(*projData->get()) ;
870 RooArgSet idxCatServers;
871 _indexCat.arg().getObservables(frame->getNormVars(), idxCatServers) ;
872
873 projDataVars.remove(idxCatServers,true,true) ;
874
875 if (!idxCompSliceSet->empty()) {
876 projDataTmp = std::unique_ptr<RooAbsData>{const_cast<RooAbsData*>(projData)->reduce(projDataVars,cutString)};
877 } else {
878 projDataTmp = std::unique_ptr<RooAbsData>{const_cast<RooAbsData*>(projData)->reduce(projDataVars)};
879 }
880
881
882
883 if (projSet) {
884 projSetTmp.add(*projSet) ;
885 projSetTmp.remove(idxCatServers,true,true);
886 }
887 }
888
889
890 if (_indexCat.arg().isDerived() && !idxCompSliceSet->empty()) {
891 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
892 << " represents a slice in index category components " << *idxCompSliceSet << endl ;
893
894 RooArgSet idxCompProjSet;
895 _indexCat.arg().getObservables(frame->getNormVars(), idxCompProjSet) ;
896 idxCompProjSet.remove(*idxCompSliceSet,true,true) ;
897 if (!idxCompProjSet.empty()) {
898 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
899 << " averages with data index category components " << idxCompProjSet << endl ;
900 }
901 } else {
902 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
903 << " averages with data index category (" << _indexCat.arg().GetName() << ")" << endl ;
904 }
905
906
907 // Override normalization and projection dataset
908 RooLinkedList cmdList2(cmdList) ;
909
910 RooCmdArg tmp1 = RooFit::Normalization(scaleFactor*sumWeight,stype) ;
911 RooCmdArg tmp2 = RooFit::ProjWData(*projDataSet,*projDataTmp) ;
912 // WVE -- do not adjust normalization for asymmetry plots
913 if (!cmdList.find("Asymmetry")) {
914 cmdList2.Add(&tmp1) ;
915 }
916 cmdList2.Add(&tmp2) ;
917
918 RooPlot* frame2 ;
919 if (!projSetTmp.empty()) {
920 // Plot temporary function
921 RooCmdArg tmp3 = RooFit::Project(projSetTmp) ;
922 cmdList2.Add(&tmp3) ;
923 frame2 = plotVar.plotOn(frame,cmdList2) ;
924 } else {
925 // Plot temporary function
926 frame2 = plotVar.plotOn(frame,cmdList2) ;
927 }
928
929 return frame2 ;
930}
931
932
933////////////////////////////////////////////////////////////////////////////////
934/// Interface function used by test statistics to freeze choice of observables
935/// for interpretation of fraction coefficients. Needed here because a RooSimultaneous
936/// works like a RooAddPdf when plotted
937
938void RooSimultaneous::selectNormalization(const RooArgSet* normSet, bool /*force*/)
939{
941 if (normSet) _plotCoefNormSet.add(*normSet) ;
942}
943
944
945////////////////////////////////////////////////////////////////////////////////
946/// Interface function used by test statistics to freeze choice of range
947/// for interpretation of fraction coefficients. Needed here because a RooSimultaneous
948/// works like a RooAddPdf when plotted
949
950void RooSimultaneous::selectNormalizationRange(const char* normRange2, bool /*force*/)
951{
952 _plotCoefNormRange = RooNameReg::ptr(normRange2) ;
953}
954
955
956
957
958////////////////////////////////////////////////////////////////////////////////
959
961 const RooArgSet* auxProto, bool verbose, bool autoBinned, const char* binnedTag) const
962{
963 const char* idxCatName = _indexCat.arg().GetName() ;
964
965 if (vars.find(idxCatName) && prototype==nullptr
966 && (auxProto==nullptr || auxProto->empty())
967 && (autoBinned || (binnedTag && strlen(binnedTag)))) {
968
969 // Return special generator config that can also do binned generation for selected states
970 return new RooSimSplitGenContext(*this,vars,verbose,autoBinned,binnedTag) ;
971
972 } else {
973
974 // Return regular generator config ;
975 return genContext(vars,prototype,auxProto,verbose) ;
976 }
977}
978
979
980
981////////////////////////////////////////////////////////////////////////////////
982/// Return specialized generator context for simultaneous p.d.f.s
983
985 const RooArgSet* auxProto, bool verbose) const
986{
987 RooArgSet allVars{vars};
988 if(prototype) allVars.add(*prototype->get());
989
990 RooArgSet catsAmongAllVars;
991 allVars.selectCommon(flattenedCatList(), catsAmongAllVars);
992
993 // Not generating index cat: return context for pdf associated with present index state
994 if(catsAmongAllVars.empty()) {
995 auto* proxy = static_cast<RooRealProxy*>(_pdfProxyList.FindObject(_indexCat->getCurrentLabel()));
996 if (!proxy) {
997 coutE(InputArguments) << "RooSimultaneous::genContext(" << GetName()
998 << ") ERROR: no PDF associated with current state ("
999 << _indexCat.arg().GetName() << "=" << _indexCat.arg().getCurrentLabel() << ")" << endl ;
1000 return nullptr;
1001 }
1002 return static_cast<RooAbsPdf*>(proxy->absArg())->genContext(vars,prototype,auxProto,verbose) ;
1003 }
1004
1005 RooArgSet catsAmongProtoVars;
1006 if(prototype) {
1007 prototype->get()->selectCommon(flattenedCatList(), catsAmongProtoVars);
1008
1009 if(!catsAmongProtoVars.empty() && catsAmongProtoVars.size() != flattenedCatList().size()) {
1010 // Abort if we have only part of the servers
1011 coutE(Plotting) << "RooSimultaneous::genContext: ERROR: prototype must include either all "
1012 << " components of the RooSimultaneous index category or none " << std::endl;
1013 return nullptr;
1014 }
1015 }
1016
1017 return new RooSimGenContext(*this,vars,prototype,auxProto,verbose) ;
1018}
1019
1020
1021
1022
1023////////////////////////////////////////////////////////////////////////////////
1024
1026 const RooArgSet* nset,
1027 double scaleFactor,
1028 bool correctForBinVolume,
1029 bool showProgress) const
1030{
1031 if (RooAbsReal::fillDataHist (hist, nset, scaleFactor,
1032 correctForBinVolume, showProgress) == nullptr)
1033 return nullptr;
1034
1035 const double sum = hist->sumEntries();
1036 if (sum != 0) {
1037 for (int i=0 ; i<hist->numEntries() ; i++) {
1038 hist->set(i, hist->weight(i) / sum, 0.);
1039 }
1040 }
1041
1042 return hist;
1043}
1044
1045
1046
1047
1048////////////////////////////////////////////////////////////////////////////////
1049/// Special generator interface for generation of 'global observables' -- for RooStats tools
1050
1052{
1053 // Make set with clone of variables (placeholder for output)
1054 RooArgSet globClone;
1055 whatVars.snapshot(globClone);
1056
1057 auto data = std::make_unique<RooDataSet>("gensimglobal","gensimglobal",whatVars);
1058
1059 for (Int_t i=0 ; i<nEvents ; i++) {
1060 for (const auto& nameIdx : indexCat()) {
1061
1062 // Get pdf associated with state from simpdf
1063 RooAbsPdf* pdftmp = getPdf(nameIdx.first);
1064
1065 RooArgSet globtmp;
1066 pdftmp->getObservables(&whatVars, globtmp) ;
1067
1068 // If there are any, generate only global variables defined by the pdf
1069 // associated with this state and transfer values to output placeholder.
1070 if (!globtmp.empty()) {
1071 globClone.assign(*std::unique_ptr<RooDataSet>{pdftmp->generate(globtmp,1)}->get(0)) ;
1072 }
1073 }
1074 data->add(globClone) ;
1075 }
1076
1077 return RooFit::makeOwningPtr(std::move(data));
1078}
1079
1080
1081/// Wraps the components of this RooSimultaneous in RooBinSamplingPdfs.
1082/// \param[in] data The dataset to be used in the eventual fit, used to figure
1083/// out the observables and whether the dataset is binned.
1084/// \param[in] precision Precision argument for all created RooBinSamplingPdfs.
1086
1087 if (precision < 0.) return;
1088
1089 RooArgSet newSamplingPdfs;
1090
1091 for (auto const &item : this->indexCat()) {
1092
1093 auto const &catName = item.first;
1094 auto &pdf = *this->getPdf(catName);
1095
1096 if (auto newSamplingPdf = RooBinSamplingPdf::create(pdf, data, precision)) {
1097 // Set the "ORIGNAME" attribute the indicate to
1098 // RooAbsArg::redirectServers() which pdf should be replaced by this
1099 // RooBinSamplingPdf in the RooSimultaneous.
1100 newSamplingPdf->setAttribute(
1101 (std::string("ORIGNAME:") + pdf.GetName()).c_str());
1102 newSamplingPdfs.addOwned(std::move(newSamplingPdf));
1103 }
1104 }
1105
1106 this->redirectServers(newSamplingPdfs, false, true);
1107 this->addOwnedComponents(std::move(newSamplingPdfs));
1108}
1109
1110
1111/// Wraps the components of this RooSimultaneous in RooBinSamplingPdfs, with a
1112/// different precision parameter for each component.
1113/// \param[in] data The dataset to be used in the eventual fit, used to figure
1114/// out the observables and whether the dataset is binned.
1115/// \param[in] precisions The map that gives the precision argument for each
1116/// component in the RooSimultaneous. The keys are the pdf names. If
1117/// there is no value for a given component, it will not use the bin
1118/// integration. Otherwise, the value has the same meaning than in
1119/// the IntegrateBins() command argument for RooAbsPdf::fitTo().
1120/// \param[in] useCategoryNames If this flag is set, the category names will be
1121/// used to look up the precision in the precisions map instead of
1122/// the pdf names.
1124 std::map<std::string, double> const& precisions,
1125 bool useCategoryNames /*=false*/) {
1126
1127 constexpr double defaultPrecision = -1.;
1128
1129 RooArgSet newSamplingPdfs;
1130
1131 for (auto const &item : this->indexCat()) {
1132
1133 auto const &catName = item.first;
1134 auto &pdf = *this->getPdf(catName);
1135 std::string pdfName = pdf.GetName();
1136
1137 auto found = precisions.find(useCategoryNames ? catName : pdfName);
1138 const double precision =
1139 found != precisions.end() ? found->second : defaultPrecision;
1140 if (precision < 0.)
1141 continue;
1142
1143 if (auto newSamplingPdf = RooBinSamplingPdf::create(pdf, data, precision)) {
1144 // Set the "ORIGNAME" attribute the indicate to
1145 // RooAbsArg::redirectServers() which pdf should be replaced by this
1146 // RooBinSamplingPdf in the RooSimultaneous.
1147 newSamplingPdf->setAttribute(
1148 (std::string("ORIGNAME:") + pdf.GetName()).c_str());
1149 newSamplingPdfs.addOwned(std::move(newSamplingPdf));
1150 }
1151 }
1152
1153 this->redirectServers(newSamplingPdfs, false, true);
1154 this->addOwnedComponents(std::move(newSamplingPdfs));
1155}
1156
1157/// Internal utility function to get a list of all category components for this
1158/// RooSimultaneous. The output contains only the index category if it is a
1159/// RooCategory, or the list of all category components if it is a
1160/// RooSuperCategory.
1162{
1163 // Note that the index category of a RooSimultaneous can only be of type
1164 // RooCategory or RooSuperCategory, because these are the only classes that
1165 // inherit from RooAbsCategoryLValue.
1166 if (auto superCat = dynamic_cast<RooSuperCategory const*>(&_indexCat.arg())) {
1167 return superCat->inputCatList();
1168 }
1169
1170 if(!_indexCatSet) {
1171 _indexCatSet = std::make_unique<RooArgSet>(_indexCat.arg());
1172 }
1173 return *_indexCatSet;
1174}
1175
1176namespace {
1177
1178void markObs(RooAbsArg *arg, std::string const &prefix, RooArgSet const &normSet)
1179{
1180 for (RooAbsArg *server : arg->servers()) {
1181 if (server->isFundamental() && normSet.find(*server)) {
1182 markObs(server, prefix, normSet);
1183 server->setAttribute("__obs__");
1184 } else if (!server->isFundamental()) {
1185 markObs(server, prefix, normSet);
1186 }
1187 }
1188}
1189
1190void prefixArgs(RooAbsArg *arg, std::string const &prefix, RooArgSet const &normSet)
1191{
1192 if (!arg->getStringAttribute("__prefix__")) {
1193 arg->SetName((prefix + arg->GetName()).c_str());
1194 arg->setStringAttribute("__prefix__", prefix.c_str());
1195 }
1196 for (RooAbsArg *server : arg->servers()) {
1197 if (server->isFundamental() && normSet.find(*server)) {
1198 prefixArgs(server, prefix, normSet);
1199 } else if (!server->isFundamental()) {
1200 prefixArgs(server, prefix, normSet);
1201 }
1202 }
1203}
1204
1205} // namespace
1206
1207std::unique_ptr<RooAbsArg>
1209{
1210 std::unique_ptr<RooSimultaneous> newSimPdf{static_cast<RooSimultaneous *>(this->Clone())};
1211
1212 const char *rangeName = this->getStringAttribute("RangeName");
1213 bool splitRange = this->getAttribute("SplitRange");
1214
1215 RooArgSet newPdfs;
1216 std::vector<std::string> catNames;
1217
1218 for (auto *proxy : static_range_cast<RooRealProxy *>(newSimPdf->_pdfProxyList)) {
1219 catNames.emplace_back(proxy->GetName());
1220 std::string const &catName = catNames.back();
1221 const std::string prefix = "_" + catName + "_";
1222
1223 const std::string origname = proxy->arg().GetName();
1224
1225 auto pdfClone = RooHelpers::cloneTreeWithSameParameters(static_cast<RooAbsPdf const &>(proxy->arg()), &normSet);
1226
1227 markObs(pdfClone.get(), prefix, normSet);
1228
1229 std::unique_ptr<RooArgSet> pdfNormSet(
1230 static_cast<RooArgSet *>(std::unique_ptr<RooArgSet>(pdfClone->getVariables())->selectByAttrib("__obs__", true)));
1231
1232 if (rangeName) {
1233 pdfClone->setNormRange(RooHelpers::getRangeNameForSimComponent(rangeName, splitRange, catName).c_str());
1234 }
1235
1236 RooFit::Detail::CompileContext pdfContext{*pdfNormSet};
1237 pdfContext.setLikelihoodMode(ctx.likelihoodMode());
1238 auto *pdfFinal = pdfContext.compile(*pdfClone, *newSimPdf, *pdfNormSet);
1239
1240 // We can only prefix the observables after everything related the
1241 // compiling of the compute graph for the normalization set is done. This
1242 // is because of a subtlety in conditional RooProdPdfs, which stores the
1243 // normalization sets for the individual pdfs in RooArgSets that are
1244 // disconnected from the computation graph, so we have no control over
1245 // them. An alternative would be to use recursive server re-direction,
1246 // but this has more performance overhead.
1247 prefixArgs(pdfFinal, prefix, normSet);
1248
1249 pdfFinal->fixAddCoefNormalization(*pdfNormSet, false);
1250
1251 pdfClone->SetName((std::string("_") + pdfClone->GetName()).c_str());
1252 pdfFinal->addOwnedComponents(std::move(pdfClone));
1253
1254 pdfFinal->setAttribute(("ORIGNAME:" + origname).c_str());
1255 newPdfs.add(*pdfFinal);
1256
1257 // We will remove the old pdf server because we will fill the new ones by
1258 // hand via the creation of new proxies.
1259 newSimPdf->removeServer(const_cast<RooAbsReal &>(proxy->arg()), true);
1260 }
1261
1262 // Replace pdfs with compiled pdfs. Don't use RooAbsArg::redirectServers()
1263 // here, because it doesn't support replacing two servers with the same name
1264 // (it can happen in a RooSimultaneous that two pdfs have the same name).
1265
1266 // First delete old proxies (we have already removed the servers before).
1267 newSimPdf->_pdfProxyList.Delete();
1268
1269 // Recreate the _pdfProxyList with the compiled pdfs
1270 for (std::size_t i = 0; i < newPdfs.size(); ++i) {
1271 const char *label = catNames[i].c_str();
1272 newSimPdf->_pdfProxyList.Add(
1273 new RooRealProxy(label, label, newSimPdf.get(), *static_cast<RooAbsReal *>(newPdfs[i])));
1274 }
1275
1276 ctx.compileServers(*newSimPdf, normSet); // to trigger compiling also the index category
1277
1278 return newSimPdf;
1279}
bool initialize()
#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.
#define ClassImp(name)
Definition Rtypes.h:377
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
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
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'.
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
void 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.
virtual bool isDerived() const
Does value or shape of this arg depend on any other arg?
Definition RooAbsArg.h:97
bool getAttribute(const Text_t *name) const
Check if a named attribute is set. By default, all attributes are unset.
RooFit::OwningPtr< RooArgSet > getVariables(bool stripDisconnected=true) const
Return RooArgSet with all variables (tree leaf nodes of expression tree)
TObject * Clone(const char *newname=nullptr) const override
Make a clone of an object using the Streamer facility.
Definition RooAbsArg.h:89
virtual bool isFundamental() const
Is this object a fundamental type that can be added to a dataset? Fundamental-type subclasses overrid...
Definition RooAbsArg.h:249
Abstract base class for objects that represent a discrete value that can be set from the outside,...
virtual bool setLabel(const char *label, bool printError=true)=0
Change category state by specifying a state name.
A space to attach TBranches.
virtual value_type getCurrentIndex() const
Return index number of current state.
virtual const char * getCurrentLabel() const
Return label string of current state.
std::map< std::string, value_type >::const_iterator begin() const
Iterator for category state names. Points to pairs of index and name.
std::size_t size() const
Number of states defined.
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
bool contains(const RooAbsArg &var) const
Check if collection contains an argument with the same name as var.
const char * GetName() const override
Returns name of object.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
TObject * Clone(const char *newname=nullptr) const override
Make a clone of an object using the Streamer facility.
void assign(const RooAbsCollection &other) const
Sets the value, cache and constant attribute of any argument in our set that also appears in the othe...
Storage_t::size_type size() const
RooAbsArg * first() const
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:57
virtual const RooArgSet * get() const
Definition RooAbsData.h:101
virtual Roo1DTable * table(const RooArgSet &catSet, const char *cuts="", const char *opts="") const
Construct table for product of categories in catSet.
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:40
bool mustBeExtended() const
If true PDF must provide extended likelihood term.
Definition RooAbsPdf.h:223
RooArgSet const * _normSet
Normalization integral (owned by _normMgr)
Definition RooAbsPdf.h:321
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:124
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition RooAbsPdf.h:219
@ CanBeExtended
Definition RooAbsPdf.h:213
@ MustBeExtended
Definition RooAbsPdf.h:213
@ CanNotBeExtended
Definition RooAbsPdf.h:213
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:57
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
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:103
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:55
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:191
RooArgSet * selectCommon(const RooAbsCollection &refColl) const
Use RooAbsCollection::selecCommon(), but return as RooArgSet.
Definition RooArgSet.h:186
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.
Object to represent discrete states.
Definition RooCategory.h:28
Named container for two doubles, two integers two object points and three string pointers that can be...
Definition RooCmdArg.h:26
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:57
const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
void compileServers(RooAbsArg &arg, RooArgSet const &normSet)
Collection class for internal use, storing a collection of RooAbsArg pointers in a doubly linked list...
TObject * find(const char *name) const
Return pointer to object with given name in collection.
virtual void Add(TObject *arg)
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:45
void remove(const char *name=nullptr, bool deleteToo=true)
Remove object with given name, or last object added if no name is given.
Definition RooPlot.cxx:868
const RooArgSet * getNormVars() const
Definition RooPlot.h:152
RooAbsRealLValue * getPlotVar() const
Definition RooPlot.h:143
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...
const char * c_str() const
Joins several RooAbsCategoryLValue objects into a single category.
bool setLabel(const char *label, bool printError=true) override
Set the value of the super category by specifying the state name.
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:576
void Add(TObject *obj) override
Definition TList.h:83
void Delete(Option_t *option="") override
Remove all objects from the list AND delete all heap based objects.
Definition TList.cxx:468
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:355
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Mother of all ROOT objects.
Definition TObject.h:41
virtual Int_t IndexOf(const TObject *obj) const
Return index of object in collection.
Basic string class.
Definition TString.h:139
TString & Append(const char *cs)
Definition TString.h:572
RooCmdArg ProjWData(const RooAbsData &projData, bool binData=false)
RooCmdArg Project(const RooArgSet &projSet)
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.
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
OwningPtr< T > makeOwningPtr(std::unique_ptr< T > &&ptr)
Internal helper to turn a std::unique_ptr<T> into an OwningPtr.
Definition Config.h:40
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
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:2345