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