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