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