Logo ROOT  
Reference Guide
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
22RooSimultaneous facilitates simultaneous fitting of multiple PDFs
23to subsets of a given dataset.
24The class takes an index category, which is used as a selector
25for PDFs, and a list of PDFs, each associated
26with a state of the index category. RooSimultaneous always returns
27the value of the PDF that is associated with the current value
28of the index category.
29
30Extended likelihood fitting is supported if all components support
31extended likelihood mode. The expected number of events by a RooSimultaneous
32is that of the component p.d.f. selected by the index category.
33
34The index category can be accessed using indexCategory().
35
36###Generating events
37When generating events from a RooSimultaneous, the index category has to be added to
38the dataset. Further, the PDF needs to know the relative probabilities of each category, i.e.,
39how many events are in which category. This can be achieved in two ways:
40- Generating with proto data that have category entries: An event from the same category as
41in the proto data is created for each event in the proto data.
42See RooAbsPdf::generate(const RooArgSet&,const RooDataSet&,Int_t,bool,bool,bool) const.
43- No proto data: A category is chosen randomly.
44\note This requires that the PDFs building the simultaneous are extended. In this way,
45the relative probability of each category can be calculated from the number of events
46in each category.
47**/
48
49#include "RooSimultaneous.h"
51#include "RooPlot.h"
52#include "RooCurve.h"
53#include "RooRealVar.h"
54#include "RooAddPdf.h"
55#include "RooAbsData.h"
56#include "Roo1DTable.h"
57#include "RooSimGenContext.h"
59#include "RooDataSet.h"
60#include "RooCmdConfig.h"
61#include "RooNameReg.h"
62#include "RooGlobalFunc.h"
63#include "RooMsgService.h"
64#include "RooCategory.h"
65#include "RooSuperCategory.h"
66#include "RooDataHist.h"
67#include "RooRandom.h"
68#include "RooArgSet.h"
69#include "RooBinSamplingPdf.h"
70
71#include "ROOT/StringUtils.hxx"
72
73#include <iostream>
74
75using namespace std;
76
78
79
80
81////////////////////////////////////////////////////////////////////////////////
82/// Constructor with index category. PDFs associated with indexCat
83/// states can be added after construction with the addPdf() function.
84///
85/// RooSimultaneous can function without having a PDF associated
86/// with every single state. The normalization in such cases is taken
87/// from the number of registered PDFs, but getVal() will assert if
88/// when called for an unregistered index state.
89
90RooSimultaneous::RooSimultaneous(const char *name, const char *title,
91 RooAbsCategoryLValue& inIndexCat) :
92 RooAbsPdf(name,title),
93 _plotCoefNormSet("!plotCoefNormSet","plotCoefNormSet",this,false,false),
94 _plotCoefNormRange(0),
95 _partIntMgr(this,10),
96 _indexCat("indexCat","Index category",this,inIndexCat),
97 _numPdf(0)
98{
99}
100
101
102
103////////////////////////////////////////////////////////////////////////////////
104/// Constructor from index category and full list of PDFs.
105/// In this constructor form, a PDF must be supplied for each indexCat state
106/// to avoid ambiguities. The PDFs are associated with the states of the
107/// index category as they appear when iterating through the category states
108/// with RooAbsCategory::begin() and RooAbsCategory::end(). This usually means
109/// they are associated by ascending index numbers.
110///
111/// PDFs may not overlap (i.e. share any variables) with the index category (function)
112
113RooSimultaneous::RooSimultaneous(const char *name, const char *title,
114 const RooArgList& inPdfList, RooAbsCategoryLValue& inIndexCat) :
115 RooAbsPdf(name,title),
116 _plotCoefNormSet("!plotCoefNormSet","plotCoefNormSet",this,false,false),
117 _plotCoefNormRange(0),
118 _partIntMgr(this,10),
119 _indexCat("indexCat","Index category",this,inIndexCat),
120 _numPdf(0)
121{
122 if (inPdfList.size() != inIndexCat.size()) {
123 coutE(InputArguments) << "RooSimultaneous::ctor(" << GetName()
124 << " ERROR: Number PDF list entries must match number of index category states, no PDFs added" << endl ;
125 return ;
126 }
127
128 map<string,RooAbsPdf*> pdfMap ;
129 auto indexCatIt = inIndexCat.begin();
130 for (unsigned int i=0; i < inPdfList.size(); ++i) {
131 auto pdf = static_cast<RooAbsPdf*>(&inPdfList[i]);
132 const auto& nameIdx = (*indexCatIt++);
133 pdfMap[nameIdx.first] = pdf;
134 }
135
136 initialize(inIndexCat,pdfMap) ;
137}
138
139
140////////////////////////////////////////////////////////////////////////////////
141
142RooSimultaneous::RooSimultaneous(const char *name, const char *title,
143 map<string,RooAbsPdf*> pdfMap, RooAbsCategoryLValue& inIndexCat) :
144 RooAbsPdf(name,title),
145 _plotCoefNormSet("!plotCoefNormSet","plotCoefNormSet",this,false,false),
146 _plotCoefNormRange(0),
147 _partIntMgr(this,10),
148 _indexCat("indexCat","Index category",this,inIndexCat),
149 _numPdf(0)
150{
151 initialize(inIndexCat,pdfMap) ;
152}
153
154
155
156
157// This class cannot be locally defined in initialize as it cannot be
158// used as a template argument in that case
160 struct CompInfo {
165 } ;
166}
167
168void RooSimultaneous::initialize(RooAbsCategoryLValue& inIndexCat, std::map<std::string,RooAbsPdf*> pdfMap)
169{
170 // First see if there are any RooSimultaneous input components
171 bool simComps(false) ;
172 for (map<string,RooAbsPdf*>::iterator iter=pdfMap.begin() ; iter!=pdfMap.end() ; ++iter) {
173 if (dynamic_cast<RooSimultaneous*>(iter->second)) {
174 simComps = true ;
175 break ;
176 }
177 }
178
179 // If there are no simultaneous component p.d.f. do simple processing through addPdf()
180 if (!simComps) {
181 bool failure = false;
182 for (map<string,RooAbsPdf*>::iterator iter=pdfMap.begin() ; iter!=pdfMap.end() ; ++iter) {
183 failure |= addPdf(*iter->second,iter->first.c_str()) ;
184 }
185
186 if (failure) {
187 throw std::invalid_argument(std::string("At least one of the PDFs of the RooSimultaneous ")
188 + GetName() + " is invalid.");
189 }
190 return ;
191 }
192
193 // Issue info message that we are about to do some rearraning
194 coutI(InputArguments) << "RooSimultaneous::initialize(" << GetName() << ") INFO: one or more input component of simultaneous p.d.f.s are"
195 << " simultaneous p.d.f.s themselves, rewriting composite expressions as one-level simultaneous p.d.f. in terms of"
196 << " final constituents and extended index category" << endl ;
197
198
199 RooArgSet allAuxCats ;
200 map<string,RooSimultaneousAux::CompInfo> compMap ;
201 for (map<string,RooAbsPdf*>::iterator iter=pdfMap.begin() ; iter!=pdfMap.end() ; ++iter) {
203 ci.pdf = iter->second ;
204 RooSimultaneous* simComp = dynamic_cast<RooSimultaneous*>(iter->second) ;
205 if (simComp) {
206 ci.simPdf = simComp ;
207 ci.subIndex = &simComp->indexCat() ;
208 ci.subIndexComps = simComp->indexCat().isFundamental() ? new RooArgSet(simComp->indexCat()) : simComp->indexCat().getVariables() ;
209 allAuxCats.add(*(ci.subIndexComps),true) ;
210 } else {
211 ci.simPdf = 0 ;
212 ci.subIndex = 0 ;
213 ci.subIndexComps = 0 ;
214 }
215 compMap[iter->first] = ci ;
216 }
217
218 // Construct the 'superIndex' from the nominal index category and all auxiliary components
219 RooArgSet allCats(inIndexCat) ;
220 allCats.add(allAuxCats) ;
221 string siname = Form("%s_index",GetName()) ;
222 RooSuperCategory* superIndex = new RooSuperCategory(siname.c_str(),siname.c_str(),allCats) ;
223 bool failure = false;
224
225 // Now process each of original pdf/state map entries
226 for (map<string,RooSimultaneousAux::CompInfo>::iterator citer = compMap.begin() ; citer != compMap.end() ; ++citer) {
227
228 RooArgSet repliCats(allAuxCats) ;
229 if (citer->second.subIndexComps) {
230 repliCats.remove(*citer->second.subIndexComps) ;
231 delete citer->second.subIndexComps ;
232 }
233 inIndexCat.setLabel(citer->first.c_str()) ;
234
235 if (!citer->second.simPdf) {
236
237 // Entry is a plain p.d.f. assign it to every state permutation of the repliCats set
238 RooSuperCategory repliSuperCat("tmp","tmp",repliCats) ;
239
240 // Iterator over all states of repliSuperCat
241 for (const auto& nameIdx : repliSuperCat) {
242 // Set value
243 repliSuperCat.setLabel(nameIdx.first) ;
244 // Retrieve corresponding label of superIndex
245 string superLabel = superIndex->getCurrentLabel() ;
246 failure |= addPdf(*citer->second.pdf,superLabel.c_str()) ;
247 cxcoutD(InputArguments) << "RooSimultaneous::initialize(" << GetName()
248 << ") assigning pdf " << citer->second.pdf->GetName() << " to super label " << superLabel << endl ;
249 }
250 } else {
251
252 // Entry is a simultaneous p.d.f
253
254 if (repliCats.empty()) {
255
256 // Case 1 -- No replication of components of RooSim component are required
257
258 for (const auto& type : *citer->second.subIndex) {
259 const_cast<RooAbsCategoryLValue*>(citer->second.subIndex)->setLabel(type.first.c_str());
260 string superLabel = superIndex->getCurrentLabel() ;
261 RooAbsPdf* compPdf = citer->second.simPdf->getPdf(type.first.c_str());
262 if (compPdf) {
263 failure |= addPdf(*compPdf,superLabel.c_str()) ;
264 cxcoutD(InputArguments) << "RooSimultaneous::initialize(" << GetName()
265 << ") assigning pdf " << compPdf->GetName() << "(member of " << citer->second.pdf->GetName()
266 << ") to super label " << superLabel << endl ;
267 } else {
268 coutW(InputArguments) << "RooSimultaneous::initialize(" << GetName() << ") WARNING: No p.d.f. associated with label "
269 << type.second << " for component RooSimultaneous p.d.f " << citer->second.pdf->GetName()
270 << "which is associated with master index label " << citer->first << endl ;
271 }
272 }
273
274 } else {
275
276 // Case 2 -- Replication of components of RooSim component are required
277
278 // Make replication supercat
279 RooSuperCategory repliSuperCat("tmp","tmp",repliCats) ;
280
281 for (const auto& stype : *citer->second.subIndex) {
282 const_cast<RooAbsCategoryLValue*>(citer->second.subIndex)->setLabel(stype.first.c_str());
283
284 for (const auto& nameIdx : repliSuperCat) {
285 repliSuperCat.setLabel(nameIdx.first) ;
286 const string superLabel = superIndex->getCurrentLabel() ;
287 RooAbsPdf* compPdf = citer->second.simPdf->getPdf(stype.first.c_str());
288 if (compPdf) {
289 failure |= addPdf(*compPdf,superLabel.c_str()) ;
290 cxcoutD(InputArguments) << "RooSimultaneous::initialize(" << GetName()
291 << ") assigning pdf " << compPdf->GetName() << "(member of " << citer->second.pdf->GetName()
292 << ") to super label " << superLabel << endl ;
293 } else {
294 coutW(InputArguments) << "RooSimultaneous::initialize(" << GetName() << ") WARNING: No p.d.f. associated with label "
295 << stype.second << " for component RooSimultaneous p.d.f " << citer->second.pdf->GetName()
296 << "which is associated with master index label " << citer->first << endl ;
297 }
298 }
299 }
300 }
301 }
302 }
303
304 if (failure) {
305 throw std::invalid_argument(std::string("Failed to initialise RooSimultaneous ") + GetName());
306 }
307
308 // Change original master index to super index and take ownership of it
309 _indexCat.setArg(*superIndex) ;
310 addOwnedComponents(*superIndex) ;
311
312}
313
314
315
316////////////////////////////////////////////////////////////////////////////////
317/// Copy constructor
318
320 RooAbsPdf(other,name),
321 _plotCoefNormSet("!plotCoefNormSet",this,other._plotCoefNormSet),
322 _plotCoefNormRange(other._plotCoefNormRange),
323 _partIntMgr(other._partIntMgr,this),
324 _indexCat("indexCat",this,other._indexCat),
325 _numPdf(other._numPdf)
326{
327 // Copy proxy list
328 for(auto* proxy : static_range_cast<RooRealProxy*>(other._pdfProxyList)) {
329 _pdfProxyList.Add(new RooRealProxy(proxy->GetName(),this,*proxy)) ;
330 }
331}
332
333
334
335////////////////////////////////////////////////////////////////////////////////
336/// Destructor
337
339{
341}
342
343
344
345////////////////////////////////////////////////////////////////////////////////
346/// Return the p.d.f associated with the given index category name
347
348RooAbsPdf* RooSimultaneous::getPdf(const char* catName) const
349{
351 return proxy ? ((RooAbsPdf*)proxy->absArg()) : 0 ;
352}
353
354
355
356////////////////////////////////////////////////////////////////////////////////
357/// Associate given PDF with index category state label 'catLabel'.
358/// The name state must be already defined in the index category.
359///
360/// RooSimultaneous can function without having a PDF associated
361/// with every single state. The normalization in such cases is taken
362/// from the number of registered PDFs, but getVal() will fail if
363/// called for an unregistered index state.
364///
365/// PDFs may not overlap (i.e. share any variables) with the index category (function).
366/// \param[in] pdf PDF to be added.
367/// \param[in] catLabel Name of the category state to be associated to the PDF.
368/// \return `true` in case of failure.
369
370bool RooSimultaneous::addPdf(const RooAbsPdf& pdf, const char* catLabel)
371{
372 // PDFs cannot overlap with the index category
373 if (pdf.dependsOn(_indexCat.arg())) {
374 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName() << "): PDF '" << pdf.GetName()
375 << "' overlaps with index category '" << _indexCat.arg().GetName() << "'."<< endl ;
376 return true ;
377 }
378
379 // Each index state can only have one PDF associated with it
380 if (_pdfProxyList.FindObject(catLabel)) {
381 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName() << "): index state '"
382 << catLabel << "' has already an associated PDF." << endl ;
383 return true ;
384 }
385
386 const RooSimultaneous* simPdf = dynamic_cast<const RooSimultaneous*>(&pdf) ;
387 if (simPdf) {
388
389 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName()
390 << ") ERROR: you cannot add a RooSimultaneous component to a RooSimultaneous using addPdf()."
391 << " Use the constructor with RooArgList if input p.d.f.s or the map<string,RooAbsPdf&> instead." << endl ;
392 return true ;
393
394 } else {
395
396 // Create a proxy named after the associated index state
397 TObject* proxy = new RooRealProxy(catLabel,catLabel,this,(RooAbsPdf&)pdf) ;
398 _pdfProxyList.Add(proxy) ;
399 _numPdf += 1 ;
400 }
401
402 return false ;
403}
404
405
406
407
408
409////////////////////////////////////////////////////////////////////////////////
410/// Examine the pdf components and check if one of them can be extended or must be extended
411/// It is enough to have one component that can be exteded or must be extended to return the flag in
412/// the total simultaneous pdf
413
415{
416 bool anyCanExtend(false) ;
417 bool anyMustExtend(false) ;
418
419 for (Int_t i=0 ; i<_numPdf ; i++) {
421 if (proxy) {
422 RooAbsPdf* pdf = (RooAbsPdf*) proxy->absArg() ;
423 //cout << " now processing pdf " << pdf->GetName() << endl;
424 if (pdf->canBeExtended()) {
425 //cout << "RooSim::extendedMode(" << GetName() << ") component " << pdf->GetName() << " can be extended"
426 // << endl;
427 anyCanExtend = true;
428 }
429 if (pdf->mustBeExtended()) {
430 //cout << "RooSim::extendedMode(" << GetName() << ") component " << pdf->GetName() << " MUST be extended" << endl;
431 anyMustExtend = true;
432 }
433 }
434 }
435 if (anyMustExtend) {
436 //cout << "RooSim::extendedMode(" << GetName() << ") returning MustBeExtended" << endl ;
437 return MustBeExtended ;
438 }
439 if (anyCanExtend) {
440 //cout << "RooSim::extendedMode(" << GetName() << ") returning CanBeExtended" << endl ;
441 return CanBeExtended ;
442 }
443 //cout << "RooSim::extendedMode(" << GetName() << ") returning CanNotBeExtended" << endl ;
444 return CanNotBeExtended ;
445}
446
447
448
449
450////////////////////////////////////////////////////////////////////////////////
451/// Return the current value:
452/// the value of the PDF associated with the current index category state
453
455{
456 // Retrieve the proxy by index name
458
459 //assert(proxy!=0) ;
460 if (proxy==0) return 0 ;
461
462 // Calculate relative weighting factor for sim-pdfs of all extendable components
463 double catFrac(1) ;
464 if (canBeExtended()) {
465 double nEvtCat = ((RooAbsPdf*)(proxy->absArg()))->expectedEvents(_normSet) ;
466
467 double nEvtTot(0) ;
468 for(auto * proxy2 : static_range_cast<RooRealProxy*>(_pdfProxyList)) {
469 nEvtTot += ((RooAbsPdf*)(proxy2->absArg()))->expectedEvents(_normSet) ;
470 }
471 catFrac=nEvtCat/nEvtTot ;
472 }
473
474 // Return the selected PDF value, normalized by the number of index states
475 return ((RooAbsPdf*)(proxy->absArg()))->getVal(_normSet)*catFrac ;
476}
477
478
479
480////////////////////////////////////////////////////////////////////////////////
481/// Return the number of expected events: If the index is in nset,
482/// then return the sum of the expected events of all components,
483/// otherwise return the number of expected events of the PDF
484/// associated with the current index category state
485
487{
488 if (nset->contains(_indexCat.arg())) {
489
490 double sum(0) ;
491
492 for(auto * proxy : static_range_cast<RooRealProxy*>(_pdfProxyList)) {
493 sum += ((RooAbsPdf*)(proxy->absArg()))->expectedEvents(nset) ;
494 }
495
496 return sum ;
497
498 } else {
499
500 // Retrieve the proxy by index name
502
503 //assert(proxy!=0) ;
504 if (proxy==0) return 0 ;
505
506 // Return the selected PDF value, normalized by the number of index states
507 return ((RooAbsPdf*)(proxy->absArg()))->expectedEvents(nset);
508 }
509}
510
511
512
513////////////////////////////////////////////////////////////////////////////////
514/// Forward determination of analytical integration capabilities to component p.d.f.s
515/// A unique code is assigned to the combined integration capabilities of all associated
516/// p.d.f.s
517
519 const RooArgSet* normSet, const char* rangeName) const
520{
521 // Declare that we can analytically integrate all requested observables
522 analVars.add(allVars) ;
523
524 // Retrieve (or create) the required partial integral list
525 Int_t code ;
526
527 // Check if this configuration was created before
528 CacheElem* cache = (CacheElem*) _partIntMgr.getObj(normSet,&analVars,0,RooNameReg::ptr(rangeName)) ;
529 if (cache) {
530 code = _partIntMgr.lastIndex() ;
531 return code+1 ;
532 }
533 cache = new CacheElem ;
534
535 // Create the partial integral set for this request
536 for(auto * proxy : static_range_cast<RooRealProxy*>(_pdfProxyList)) {
537 RooAbsReal* pdfInt = proxy->arg().createIntegral(analVars,normSet,0,rangeName) ;
538 cache->_partIntList.addOwned(*pdfInt) ;
539 }
540
541 // Store the partial integral list and return the assigned code ;
542 code = _partIntMgr.setObj(normSet,&analVars,cache,RooNameReg::ptr(rangeName)) ;
543
544 return code+1 ;
545}
546
547
548
549////////////////////////////////////////////////////////////////////////////////
550/// Return analytical integration defined by given code
551
552double RooSimultaneous::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* /*rangeName*/) const
553{
554 // No integration scenario
555 if (code==0) {
556 return getVal(normSet) ;
557 }
558
559 // Partial integration scenarios, rangeName already encoded in 'code'
560 CacheElem* cache = (CacheElem*) _partIntMgr.getObjByIndex(code-1) ;
561
563 Int_t idx = _pdfProxyList.IndexOf(proxy) ;
564 return ((RooAbsReal*)cache->_partIntList.at(idx))->getVal(normSet) ;
565}
566
567
568
569
570
571
572////////////////////////////////////////////////////////////////////////////////
573/// Back-end for plotOn() implementation on RooSimultaneous which
574/// needs special handling because a RooSimultaneous PDF cannot
575/// project out its index category via integration. plotOn() will
576/// abort if this is requested without providing a projection dataset.
577
579{
580 // Sanity checks
581 if (plotSanityChecks(frame)) return frame ;
582
583 // Extract projection configuration from command list
584 RooCmdConfig pc(Form("RooSimultaneous::plotOn(%s)",GetName())) ;
585 pc.defineString("sliceCatState","SliceCat",0,"",true) ;
586 pc.defineDouble("scaleFactor","Normalization",0,1.0) ;
587 pc.defineInt("scaleType","Normalization",0,RooAbsPdf::Relative) ;
588 pc.defineObject("sliceCatList","SliceCat",0,0,true) ;
589 // This dummy is needed for plotOn to recognize the "SliceCatMany" command.
590 // It is not used directly, but the "SliceCat" commands are nested in it.
591 // Removing this dummy definition results in "ERROR: unrecognized command: SliceCatMany".
592 pc.defineObject("dummy1","SliceCatMany",0) ;
593 pc.defineSet("projSet","Project",0) ;
594 pc.defineSet("sliceSet","SliceVars",0) ;
595 pc.defineSet("projDataSet","ProjData",0) ;
596 pc.defineObject("projData","ProjData",1) ;
597 pc.defineMutex("Project","SliceVars") ;
598 pc.allowUndefined() ; // there may be commands we don't handle here
599
600 // Process and check varargs
601 pc.process(cmdList) ;
602 if (!pc.ok(true)) {
603 return frame ;
604 }
605
606 const RooAbsData* projData = (const RooAbsData*) pc.getObject("projData") ;
607 const RooArgSet* projDataSet = pc.getSet("projDataSet");
608 const RooArgSet* sliceSetTmp = pc.getSet("sliceSet") ;
609 std::unique_ptr<RooArgSet> sliceSet( sliceSetTmp ? ((RooArgSet*) sliceSetTmp->Clone()) : nullptr );
610 const RooArgSet* projSet = pc.getSet("projSet") ;
611 double scaleFactor = pc.getDouble("scaleFactor") ;
612 ScaleType stype = (ScaleType) pc.getInt("scaleType") ;
613
614
615 // Look for category slice arguments and add them to the master slice list if found
616 const char* sliceCatState = pc.getString("sliceCatState",0,true) ;
617 const RooLinkedList& sliceCatList = pc.getObjectList("sliceCatList") ;
618 if (sliceCatState) {
619
620 // Make the master slice set if it doesnt exist
621 if (!sliceSet) {
622 sliceSet.reset(new RooArgSet);
623 }
624
625 // Prepare comma separated label list for parsing
626 auto catTokens = ROOT::Split(sliceCatState, ",");
627
628 // Loop over all categories provided by (multiple) Slice() arguments
629 unsigned int tokenIndex = 0;
630 for(auto * scat : static_range_cast<RooCategory*>(sliceCatList)) {
631 const char* slabel = tokenIndex >= catTokens.size() ? nullptr : catTokens[tokenIndex++].c_str();
632
633 if (slabel) {
634 // Set the slice position to the value indicated by slabel
635 scat->setLabel(slabel) ;
636 // Add the slice category to the master slice set
637 sliceSet->add(*scat,false) ;
638 }
639 }
640 }
641
642 // Check if we have a projection dataset
643 if (!projData) {
644 coutE(InputArguments) << "RooSimultaneous::plotOn(" << GetName() << ") ERROR: must have a projection dataset for index category" << endl ;
645 return frame ;
646 }
647
648 // Make list of variables to be projected
649 RooArgSet projectedVars ;
650 if (sliceSet) {
651 makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,true) ;
652
653 // Take out the sliced variables
654 for (const auto sliceArg : *sliceSet) {
655 RooAbsArg* arg = projectedVars.find(sliceArg->GetName()) ;
656 if (arg) {
657 projectedVars.remove(*arg) ;
658 } else {
659 coutI(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") slice variable "
660 << sliceArg->GetName() << " was not projected anyway" << endl ;
661 }
662 }
663 } else if (projSet) {
664 makeProjectionSet(frame->getPlotVar(),projSet,projectedVars,false) ;
665 } else {
666 makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,true) ;
667 }
668
669 bool projIndex(false) ;
670
671 if (!_indexCat.arg().isDerived()) {
672 // *** Error checking for a fundamental index category ***
673 //cout << "RooSim::plotOn: index is fundamental" << endl ;
674
675 // Check that the provided projection dataset contains our index variable
676 if (!projData->get()->find(_indexCat.arg().GetName())) {
677 coutE(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") ERROR: Projection over index category "
678 << "requested, but projection data set doesn't contain index category" << endl ;
679 return frame ;
680 }
681
682 if (projectedVars.find(_indexCat.arg().GetName())) {
683 projIndex=true ;
684 }
685
686 } else {
687 // *** Error checking for a composite index category ***
688
689 // Determine if any servers of the index category are in the projectedVars
690 RooArgSet projIdxServers ;
691 bool anyServers(false) ;
692 for (const auto server : _indexCat->servers()) {
693 if (projectedVars.find(server->GetName())) {
694 anyServers=true ;
695 projIdxServers.add(*server) ;
696 }
697 }
698
699 // Check that the projection dataset contains all the
700 // index category components we're projecting over
701
702 // Determine if all projected servers of the index category are in the projection dataset
703 bool allServers(true) ;
704 std::string missing;
705 for (const auto server : projIdxServers) {
706 if (!projData->get()->find(server->GetName())) {
707 allServers=false ;
708 missing = server->GetName();
709 }
710 }
711
712 if (!allServers) {
713 coutE(Plotting) << "RooSimultaneous::plotOn(" << GetName()
714 << ") ERROR: Projection dataset doesn't contain complete set of index categories to do projection."
715 << "\n\tcategory " << missing << " is missing." << endl ;
716 return frame ;
717 }
718
719 if (anyServers) {
720 projIndex = true ;
721 }
722 }
723
724 // Calculate relative weight fractions of components
725 std::unique_ptr<Roo1DTable> wTable( projData->table(_indexCat.arg()) );
726
727 // Clone the index category to be able to cycle through the category states for plotting without
728 // affecting the category state of our instance
729 std::unique_ptr<RooArgSet> idxCloneSet( RooArgSet(*_indexCat).snapshot(true) );
730 auto idxCatClone = static_cast<RooAbsCategoryLValue*>( idxCloneSet->find(_indexCat->GetName()) );
731 assert(idxCatClone);
732
733 // Make list of category columns to exclude from projection data
734 std::unique_ptr<RooArgSet> idxCompSliceSet( idxCatClone->getObservables(frame->getNormVars()) );
735
736 // If we don't project over the index, just do the regular plotOn
737 if (!projIndex) {
738
739 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
740 << " represents a slice in the index category (" << _indexCat.arg().GetName() << ")" << endl ;
741
742 // Reduce projData: take out fitCat (component) columns and entries that don't match selected slice
743 // Construct cut string to only select projection data event that match the current slice
744
745 // Make cut string to exclude rows from projection data
746 TString cutString ;
747 bool first(true) ;
748 for (const auto arg : *idxCompSliceSet) {
749 auto idxComp = static_cast<RooCategory*>(arg);
750 RooAbsArg* slicedComponent = nullptr;
751 if (sliceSet && (slicedComponent = sliceSet->find(*idxComp)) != nullptr) {
752 auto theCat = static_cast<const RooAbsCategory*>(slicedComponent);
753 idxComp->setIndex(theCat->getCurrentIndex(), false);
754 }
755
756 if (!first) {
757 cutString.Append("&&") ;
758 } else {
759 first=false ;
760 }
761 cutString.Append(Form("%s==%d",idxComp->GetName(),idxComp->getCurrentIndex())) ;
762 }
763
764 // Make temporary projData without RooSim index category components
765 RooArgSet projDataVars(*projData->get()) ;
766 projDataVars.remove(*idxCompSliceSet,true,true) ;
767
768 std::unique_ptr<RooAbsData> projDataTmp( const_cast<RooAbsData*>(projData)->reduce(projDataVars,cutString) );
769
770 // Override normalization and projection dataset
771 RooCmdArg tmp1 = RooFit::Normalization(scaleFactor*wTable->getFrac(idxCatClone->getCurrentLabel()),stype) ;
772 RooCmdArg tmp2 = RooFit::ProjWData(*projDataSet,*projDataTmp) ;
773
774 // WVE -- do not adjust normalization for asymmetry plots
775 RooLinkedList cmdList2(cmdList) ;
776 if (!cmdList.find("Asymmetry")) {
777 cmdList2.Add(&tmp1) ;
778 }
779 cmdList2.Add(&tmp2) ;
780
781 // Plot single component
782 RooPlot* retFrame = getPdf(idxCatClone->getCurrentLabel())->plotOn(frame,cmdList2);
783 return retFrame ;
784 }
785
786 // If we project over the index, plot using a temporary RooAddPdf
787 // using the weights from the data as coefficients
788
789 // Build the list of indexCat components that are sliced
790 idxCompSliceSet->remove(projectedVars,true,true) ;
791
792 // Make a new expression that is the weighted sum of requested components
793 RooArgList pdfCompList ;
794 RooArgList wgtCompList ;
795//RooAbsPdf* pdf ;
796 double sumWeight(0) ;
797 for(auto * proxy : static_range_cast<RooRealProxy*>(_pdfProxyList)) {
798
799 idxCatClone->setLabel(proxy->name()) ;
800
801 // Determine if this component is the current slice (if we slice)
802 bool skip(false) ;
803 for (const auto idxSliceCompArg : *idxCompSliceSet) {
804 const auto idxSliceComp = static_cast<RooAbsCategory*>(idxSliceCompArg);
805 RooAbsCategory* idxComp = (RooAbsCategory*) idxCloneSet->find(idxSliceComp->GetName()) ;
806 if (idxComp->getCurrentIndex()!=idxSliceComp->getCurrentIndex()) {
807 skip=true ;
808 break ;
809 }
810 }
811 if (skip) continue ;
812
813 // Instantiate a RRV holding this pdfs weight fraction
814 RooRealVar *wgtVar = new RooRealVar(proxy->name(),"coef",wTable->getFrac(proxy->name())) ;
815 wgtCompList.addOwned(*wgtVar) ;
816 sumWeight += wTable->getFrac(proxy->name()) ;
817
818 // Add the PDF to list list
819 pdfCompList.add(proxy->arg()) ;
820 }
821
822 TString plotVarName(GetName()) ;
823 RooAddPdf *plotVar = new RooAddPdf(plotVarName,"weighted sum of RS components",pdfCompList,wgtCompList) ;
824
825 // Fix appropriate coefficient normalization in plot function
826 if (_plotCoefNormSet.getSize()>0) {
828 }
829
830 std::unique_ptr<RooAbsData> projDataTmp;
831 RooArgSet projSetTmp ;
832 if (projData) {
833
834 // Construct cut string to only select projection data event that match the current slice
835 TString cutString ;
836 if (idxCompSliceSet->getSize()>0) {
837 bool first(true) ;
838 for (const auto idxSliceCompArg : *idxCompSliceSet) {
839 const auto idxSliceComp = static_cast<RooAbsCategory*>(idxSliceCompArg);
840 if (!first) {
841 cutString.Append("&&") ;
842 } else {
843 first=false ;
844 }
845 cutString.Append(Form("%s==%d",idxSliceComp->GetName(),idxSliceComp->getCurrentIndex())) ;
846 }
847 }
848
849 // Make temporary projData without RooSim index category components
850 RooArgSet projDataVars(*projData->get()) ;
851 RooArgSet* idxCatServers = _indexCat.arg().getObservables(frame->getNormVars()) ;
852
853 projDataVars.remove(*idxCatServers,true,true) ;
854
855 if (idxCompSliceSet->getSize()>0) {
856 projDataTmp.reset( const_cast<RooAbsData*>(projData)->reduce(projDataVars,cutString) );
857 } else {
858 projDataTmp.reset( const_cast<RooAbsData*>(projData)->reduce(projDataVars) );
859 }
860
861
862
863 if (projSet) {
864 projSetTmp.add(*projSet) ;
865 projSetTmp.remove(*idxCatServers,true,true);
866 }
867
868
869 delete idxCatServers ;
870 }
871
872
873 if (_indexCat.arg().isDerived() && idxCompSliceSet->getSize()>0) {
874 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
875 << " represents a slice in index category components " << *idxCompSliceSet << endl ;
876
877 RooArgSet* idxCompProjSet = _indexCat.arg().getObservables(frame->getNormVars()) ;
878 idxCompProjSet->remove(*idxCompSliceSet,true,true) ;
879 if (idxCompProjSet->getSize()>0) {
880 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
881 << " averages with data index category components " << *idxCompProjSet << endl ;
882 }
883 delete idxCompProjSet ;
884 } else {
885 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
886 << " averages with data index category (" << _indexCat.arg().GetName() << ")" << endl ;
887 }
888
889
890 // Override normalization and projection dataset
891 RooLinkedList cmdList2(cmdList) ;
892
893 RooCmdArg tmp1 = RooFit::Normalization(scaleFactor*sumWeight,stype) ;
894 RooCmdArg tmp2 = RooFit::ProjWData(*projDataSet,*projDataTmp) ;
895 // WVE -- do not adjust normalization for asymmetry plots
896 if (!cmdList.find("Asymmetry")) {
897 cmdList2.Add(&tmp1) ;
898 }
899 cmdList2.Add(&tmp2) ;
900
901 RooPlot* frame2 ;
902 if (projSetTmp.getSize()>0) {
903 // Plot temporary function
904 RooCmdArg tmp3 = RooFit::Project(projSetTmp) ;
905 cmdList2.Add(&tmp3) ;
906 frame2 = plotVar->plotOn(frame,cmdList2) ;
907 } else {
908 // Plot temporary function
909 frame2 = plotVar->plotOn(frame,cmdList2) ;
910 }
911
912 // Cleanup
913 delete plotVar ;
914
915 return frame2 ;
916}
917
918
919
920////////////////////////////////////////////////////////////////////////////////
921/// OBSOLETE -- Retained for backward compatibility
922
923RooPlot* RooSimultaneous::plotOn(RooPlot *frame, Option_t* drawOptions, double scaleFactor,
924 ScaleType stype, const RooAbsData* projData, const RooArgSet* projSet,
925 double /*precision*/, bool /*shiftToZero*/, const RooArgSet* /*projDataSet*/,
926 double /*rangeLo*/, double /*rangeHi*/, RooCurve::WingMode /*wmode*/) const
927{
928 // Make command list
929 RooLinkedList cmdList ;
930 cmdList.Add(new RooCmdArg(RooFit::DrawOption(drawOptions))) ;
931 cmdList.Add(new RooCmdArg(RooFit::Normalization(scaleFactor,stype))) ;
932 if (projData) cmdList.Add(new RooCmdArg(RooFit::ProjWData(*projData))) ;
933 if (projSet) cmdList.Add(new RooCmdArg(RooFit::Project(*projSet))) ;
934
935 // Call new method
936 RooPlot* ret = plotOn(frame,cmdList) ;
937
938 // Cleanup
939 cmdList.Delete() ;
940 return ret ;
941}
942
943
944
945////////////////////////////////////////////////////////////////////////////////
946/// Interface function used by test statistics to freeze choice of observables
947/// for interpretation of fraction coefficients. Needed here because a RooSimultaneous
948/// works like a RooAddPdf when plotted
949
950void RooSimultaneous::selectNormalization(const RooArgSet* normSet, bool /*force*/)
951{
953 if (normSet) _plotCoefNormSet.add(*normSet) ;
954}
955
956
957////////////////////////////////////////////////////////////////////////////////
958/// Interface function used by test statistics to freeze choice of range
959/// for interpretation of fraction coefficients. Needed here because a RooSimultaneous
960/// works like a RooAddPdf when plotted
961
962void RooSimultaneous::selectNormalizationRange(const char* normRange2, bool /*force*/)
963{
964 _plotCoefNormRange = RooNameReg::ptr(normRange2) ;
965}
966
967
968
969
970////////////////////////////////////////////////////////////////////////////////
971
973 const RooArgSet* auxProto, bool verbose, bool autoBinned, const char* binnedTag) const
974{
975 const char* idxCatName = _indexCat.arg().GetName() ;
976
977 if (vars.find(idxCatName) && prototype==0
978 && (auxProto==0 || auxProto->empty())
979 && (autoBinned || (binnedTag && strlen(binnedTag)))) {
980
981 // Return special generator config that can also do binned generation for selected states
982 return new RooSimSplitGenContext(*this,vars,verbose,autoBinned,binnedTag) ;
983
984 } else {
985
986 // Return regular generator config ;
987 return genContext(vars,prototype,auxProto,verbose) ;
988 }
989}
990
991
992
993////////////////////////////////////////////////////////////////////////////////
994/// Return specialized generator context for simultaneous p.d.f.s
995
997 const RooArgSet* auxProto, bool verbose) const
998{
999 const char* idxCatName = _indexCat.arg().GetName() ;
1000 const RooArgSet* protoVars = prototype ? prototype->get() : 0 ;
1001
1002 if (vars.find(idxCatName) || (protoVars && protoVars->find(idxCatName))) {
1003
1004 // Generating index category: return special sim-context
1005 return new RooSimGenContext(*this,vars,prototype,auxProto,verbose) ;
1006
1007 } else if (_indexCat.arg().isDerived()) {
1008 // Generating dependents of a derived index category
1009
1010 // Determine if we none,any or all servers
1011 bool anyServer(false), allServers(true) ;
1012 if (prototype) {
1013 for(RooAbsArg * server : _indexCat.arg().servers()) {
1014 if (prototype->get()->find(server->GetName())) {
1015 anyServer=true ;
1016 } else {
1017 allServers=false ;
1018 }
1019 }
1020 } else {
1021 allServers=true ;
1022 }
1023
1024 if (allServers) {
1025 // Use simcontext if we have all servers
1026
1027 return new RooSimGenContext(*this,vars,prototype,auxProto,verbose) ;
1028 } else if (!allServers && anyServer) {
1029 // Abort if we have only part of the servers
1030 coutE(Plotting) << "RooSimultaneous::genContext: ERROR: prototype must include either all "
1031 << " components of the RooSimultaneous index category or none " << endl ;
1032 return 0 ;
1033 }
1034 // Otherwise make single gencontext for current state
1035 }
1036
1037 // Not generating index cat: return context for pdf associated with present index state
1039 if (!proxy) {
1040 coutE(InputArguments) << "RooSimultaneous::genContext(" << GetName()
1041 << ") ERROR: no PDF associated with current state ("
1042 << _indexCat.arg().GetName() << "=" << _indexCat.arg().getCurrentLabel() << ")" << endl ;
1043 return 0 ;
1044 }
1045 return ((RooAbsPdf*)proxy->absArg())->genContext(vars,prototype,auxProto,verbose) ;
1046}
1047
1048
1049
1050
1051////////////////////////////////////////////////////////////////////////////////
1052
1054 const RooArgSet* nset,
1055 double scaleFactor,
1056 bool correctForBinVolume,
1057 bool showProgress) const
1058{
1059 if (RooAbsReal::fillDataHist (hist, nset, scaleFactor,
1060 correctForBinVolume, showProgress) == 0)
1061 return 0;
1062
1063 const double sum = hist->sumEntries();
1064 if (sum != 0) {
1065 for (int i=0 ; i<hist->numEntries() ; i++) {
1066 hist->set(i, hist->weight(i) / sum, 0.);
1067 }
1068 }
1069
1070 return hist;
1071}
1072
1073
1074
1075
1076////////////////////////////////////////////////////////////////////////////////
1077/// Special generator interface for generation of 'global observables' -- for RooStats tools
1078
1080{
1081 // Make set with clone of variables (placeholder for output)
1082 RooArgSet* globClone = (RooArgSet*) whatVars.snapshot() ;
1083
1084 RooDataSet* data = new RooDataSet("gensimglobal","gensimglobal",whatVars) ;
1085
1086 for (Int_t i=0 ; i<nEvents ; i++) {
1087 for (const auto& nameIdx : indexCat()) {
1088
1089 // Get pdf associated with state from simpdf
1090 RooAbsPdf* pdftmp = getPdf(nameIdx.first.c_str());
1091
1092 // Generate only global variables defined by the pdf associated with this state
1093 RooArgSet* globtmp = pdftmp->getObservables(whatVars) ;
1094 RooDataSet* tmp = pdftmp->generate(*globtmp,1) ;
1095
1096 // Transfer values to output placeholder
1097 globClone->assign(*tmp->get(0)) ;
1098
1099 // Cleanup
1100 delete globtmp ;
1101 delete tmp ;
1102 }
1103 data->add(*globClone) ;
1104 }
1105
1106 delete globClone ;
1107 return data ;
1108}
1109
1110
1111/// Wraps the components of this RooSimultaneous in RooBinSamplingPdfs.
1112/// \param[in] data The dataset to be used in the eventual fit, used to figure
1113/// out the observables and whether the dataset is binned.
1114/// \param[in] precision Precision argument for all created RooBinSamplingPdfs.
1116
1117 if (precision < 0.) return;
1118
1119 RooArgSet newSamplingPdfs;
1120
1121 for (auto const &item : this->indexCat()) {
1122
1123 auto const &catName = item.first;
1124 auto &pdf = *this->getPdf(catName.c_str());
1125
1126 if (auto newSamplingPdf = RooBinSamplingPdf::create(pdf, data, precision)) {
1127 // Set the "ORIGNAME" attribute the indicate to
1128 // RooAbsArg::redirectServers() wich pdf should be replaced by this
1129 // RooBinSamplingPdf in the RooSimultaneous.
1130 newSamplingPdf->setAttribute(
1131 (std::string("ORIGNAME:") + pdf.GetName()).c_str());
1132 newSamplingPdfs.addOwned(std::move(newSamplingPdf));
1133 }
1134 }
1135
1136 this->redirectServers(newSamplingPdfs, false, true);
1137 this->addOwnedComponents(std::move(newSamplingPdfs));
1138}
1139
1140
1141/// Wraps the components of this RooSimultaneous in RooBinSamplingPdfs, with a
1142/// different precision parameter for each component.
1143/// \param[in] data The dataset to be used in the eventual fit, used to figure
1144/// out the observables and whether the dataset is binned.
1145/// \param[in] precisions The map that gives the precision argument for each
1146/// component in the RooSimultaneous. The keys are the pdf names. If
1147/// there is no value for a given component, it will not use the bin
1148/// integration. Otherwise, the value has the same meaning than in
1149/// the IntegrateBins() command argument for RooAbsPdf::fitTo().
1150/// \param[in] useCategoryNames If this flag is set, the category names will be
1151/// used to look up the precision in the precisions map instead of
1152/// the pdf names.
1154 std::map<std::string, double> const& precisions,
1155 bool useCategoryNames /*=false*/) {
1156
1157 constexpr double defaultPrecision = -1.;
1158
1159 RooArgSet newSamplingPdfs;
1160
1161 for (auto const &item : this->indexCat()) {
1162
1163 auto const &catName = item.first;
1164 auto &pdf = *this->getPdf(catName.c_str());
1165 std::string pdfName = pdf.GetName();
1166
1167 auto found = precisions.find(useCategoryNames ? catName : pdfName);
1168 const double precision =
1169 found != precisions.end() ? found->second : defaultPrecision;
1170 if (precision < 0.)
1171 continue;
1172
1173 if (auto newSamplingPdf = RooBinSamplingPdf::create(pdf, data, precision)) {
1174 // Set the "ORIGNAME" attribute the indicate to
1175 // RooAbsArg::redirectServers() wich pdf should be replaced by this
1176 // RooBinSamplingPdf in the RooSimultaneous.
1177 newSamplingPdf->setAttribute(
1178 (std::string("ORIGNAME:") + pdf.GetName()).c_str());
1179 newSamplingPdfs.addOwned(std::move(newSamplingPdf));
1180 }
1181 }
1182
1183 this->redirectServers(newSamplingPdfs, false, true);
1184 this->addOwnedComponents(std::move(newSamplingPdfs));
1185}
#define coutI(a)
Definition: RooMsgService.h:34
#define cxcoutD(a)
Definition: RooMsgService.h:85
#define coutW(a)
Definition: RooMsgService.h:36
#define coutE(a)
Definition: RooMsgService.h:37
RooTemplateProxy< RooAbsReal > RooRealProxy
Compatibility typedef replacing the old RooRealProxy class.
Definition: RooRealProxy.h:23
const char Option_t
Definition: RtypesCore.h:66
#define ClassImp(name)
Definition: Rtypes.h:375
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:2456
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition: RooAbsArg.h:71
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.
Definition: RooAbsArg.cxx:805
friend class RooDataSet
Definition: RooAbsArg.h:665
RooArgSet * getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
Definition: RooAbsArg.h:293
bool addOwnedComponents(const RooAbsCollection &comps)
Take ownership of the contents of 'comps'.
Definition: RooAbsArg.cxx:2185
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.
Definition: RooAbsArg.cxx:999
const RefCountList_t & servers() const
List of all servers of this object.
Definition: RooAbsArg.h:198
RooArgSet * getVariables(bool stripDisconnected=true) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
Definition: RooAbsArg.cxx:2057
virtual bool isDerived() const
Does value or shape of this arg depend on any other arg?
Definition: RooAbsArg.h:91
virtual bool isFundamental() const
Is this object a fundamental type that can be added to a dataset? Fundamental-type subclasses overrid...
Definition: RooAbsArg.h:241
RooAbsCategoryLValue is the common abstract base class for objects that represent a discrete value th...
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.
bool empty() const
Int_t getSize() const
Return the number of elements in the collection.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
TObject * Clone(const char *newname=nullptr) const override
Make a clone of an object using the Streamer facility.
void assign(const RooAbsCollection &other) const
Sets the value, cache and constant attribute of any argument in our set that also appears in the othe...
Storage_t::size_type size() const
RooAbsArg * first() const
virtual bool addOwned(RooAbsArg &var, bool silent=false)
Add an argument and transfer the ownership to the collection.
RooAbsArg * find(const char *name) const
Find object with given name in list.
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:62
virtual const RooArgSet * get() const
Definition: RooAbsData.h:106
virtual Roo1DTable * table(const RooArgSet &catSet, const char *cuts="", const char *opts="") const
Construct table for product of categories in catSet.
Definition: RooAbsData.cxx:880
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
Definition: RooAbsData.cxx:374
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
bool mustBeExtended() const
If true PDF must provide extended likelihood term.
Definition: RooAbsPdf.h:274
RooArgSet const * _normSet
Normalization integral (owned by _normMgr)
Definition: RooAbsPdf.h:365
RooDataSet * generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none())
See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,...
Definition: RooAbsPdf.h:61
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition: RooAbsPdf.h:270
@ CanBeExtended
Definition: RooAbsPdf.h:264
@ MustBeExtended
Definition: RooAbsPdf.h:264
@ CanNotBeExtended
Definition: RooAbsPdf.h:264
RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none(), const RooCmdArg &arg9=RooCmdArg::none(), const RooCmdArg &arg10=RooCmdArg::none()) const override
Helper calling plotOn(RooPlot*, RooLinkedList&) const.
Definition: RooAbsPdf.h:127
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:62
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:91
virtual void fixAddCoefNormalization(const RooArgSet &addNormSet=RooArgSet(), bool force=true)
Fix the interpretation of the coefficient of any RooAddPdf component in the expression tree headed by...
RooAbsReal * createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables std::liste...
Definition: RooAbsReal.cxx:522
bool plotSanityChecks(RooPlot *frame) const
Utility function for plotOn(), perform general sanity check on frame to ensure safe plotting operatio...
friend class RooAddPdf
Definition: RooAbsReal.h:551
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...
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:34
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:56
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition: RooArgSet.h:179
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.
RooCategory is an object to represent discrete states.
Definition: RooCategory.h:28
RooCmdArg is a named container for two doubles, two integers two object points and three string point...
Definition: RooCmdArg.h:26
Class RooCmdConfig is a configurable parser for RooCmdArg named arguments.
Definition: RooCmdConfig.h:31
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...
The RooDataHist is a 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.
Definition: RooDataHist.h:104
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.
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:55
const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
RooLinkedList is an collection class for internal use, storing a collection of RooAbsArg pointers in ...
Definition: RooLinkedList.h:38
void Delete(Option_t *o=nullptr) override
Remove all elements in collection and delete all elements NB: Collection does not own elements,...
TObject * find(const char *name) const
Return pointer to object with given name in collection.
virtual void Add(TObject *arg)
Definition: RooLinkedList.h:67
static const TNamed * ptr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
Definition: RooNameReg.cxx:81
A RooPlot is a 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:146
RooAbsRealLValue * getPlotVar() const
Definition: RooPlot.h:137
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:40
RooSimultaneous 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
RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg(), const RooCmdArg &arg9=RooCmdArg(), const RooCmdArg &arg10=RooCmdArg()) const override
RooDataSet * generateSimGlobal(const RooArgSet &whatVars, Int_t nEvents) override
Special generator interface for generation of 'global observables' – for RooStats tools.
~RooSimultaneous() override
Destructor.
ExtendMode extendMode() const override
Examine the pdf components and check if one of them can be extended or must be extended It is enough ...
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
void initialize(RooAbsCategoryLValue &inIndexCat, std::map< std::string, RooAbsPdf * > 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.
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
RooAbsPdf * getPdf(const char *catName) const
Return the p.d.f associated with the given index category name.
The RooSuperCategory can join 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.
bool setArg(T &newRef)
Change object held in proxy into newRef.
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:578
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:470
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition: TList.cxx:357
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:136
TString & Append(const char *cs)
Definition: TString.h:564
RooCmdArg ProjWData(const RooAbsData &projData, bool binData=false)
RooCmdArg DrawOption(const char *opt)
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.
Definition: StringUtils.cxx:23
@ InputArguments
Definition: RooGlobalFunc.h:62
static constexpr double pc
Definition: first.py:1
const RooAbsCategoryLValue * subIndex
static uint64_t sum(uint64_t i)
Definition: Factory.cxx:2345