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_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 "RooHelpers.h"
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,kFALSE,kFALSE),
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,kFALSE,kFALSE),
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,kFALSE,kFALSE),
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 {
161 RooAbsPdf* pdf ;
162 RooSimultaneous* simPdf ;
163 const RooAbsCategoryLValue* subIndex ;
164 RooArgSet* subIndexComps ;
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_t simComps(kFALSE) ;
172 for (map<string,RooAbsPdf*>::iterator iter=pdfMap.begin() ; iter!=pdfMap.end() ; ++iter) {
173 if (dynamic_cast<RooSimultaneous*>(iter->second)) {
174 simComps = kTRUE ;
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) {
202 RooSimultaneousAux::CompInfo ci ;
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),kTRUE) ;
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.getSize()==0) {
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 TIterator* pIter = other._pdfProxyList.MakeIterator() ;
329 RooRealProxy* proxy ;
330 while ((proxy=(RooRealProxy*)pIter->Next())) {
331 _pdfProxyList.Add(new RooRealProxy(proxy->GetName(),this,*proxy)) ;
332 }
333 delete pIter ;
334}
335
336
337
338////////////////////////////////////////////////////////////////////////////////
339/// Destructor
340
342{
344}
345
346
347
348////////////////////////////////////////////////////////////////////////////////
349/// Return the p.d.f associated with the given index category name
350
351RooAbsPdf* RooSimultaneous::getPdf(const char* catName) const
352{
354 return proxy ? ((RooAbsPdf*)proxy->absArg()) : 0 ;
355}
356
357
358
359////////////////////////////////////////////////////////////////////////////////
360/// Associate given PDF with index category state label 'catLabel'.
361/// The name state must be already defined in the index category.
362///
363/// RooSimultaneous can function without having a PDF associated
364/// with every single state. The normalization in such cases is taken
365/// from the number of registered PDFs, but getVal() will fail if
366/// called for an unregistered index state.
367///
368/// PDFs may not overlap (i.e. share any variables) with the index category (function).
369/// \param[in] pdf PDF to be added.
370/// \param[in] catLabel Name of the category state to be associated to the PDF.
371/// \return `true` in case of failure.
372
373Bool_t RooSimultaneous::addPdf(const RooAbsPdf& pdf, const char* catLabel)
374{
375 // PDFs cannot overlap with the index category
376 if (pdf.dependsOn(_indexCat.arg())) {
377 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName() << "): PDF '" << pdf.GetName()
378 << "' overlaps with index category '" << _indexCat.arg().GetName() << "'."<< endl ;
379 return kTRUE ;
380 }
381
382 // Each index state can only have one PDF associated with it
383 if (_pdfProxyList.FindObject(catLabel)) {
384 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName() << "): index state '"
385 << catLabel << "' has already an associated PDF." << endl ;
386 return kTRUE ;
387 }
388
389 const RooSimultaneous* simPdf = dynamic_cast<const RooSimultaneous*>(&pdf) ;
390 if (simPdf) {
391
392 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName()
393 << ") ERROR: you cannot add a RooSimultaneous component to a RooSimultaneous using addPdf()."
394 << " Use the constructor with RooArgList if input p.d.f.s or the map<string,RooAbsPdf&> instead." << endl ;
395 return kTRUE ;
396
397 } else {
398
399 // Create a proxy named after the associated index state
400 TObject* proxy = new RooRealProxy(catLabel,catLabel,this,(RooAbsPdf&)pdf) ;
401 _pdfProxyList.Add(proxy) ;
402 _numPdf += 1 ;
403 }
404
405 return kFALSE ;
406}
407
408
409
410
411
412////////////////////////////////////////////////////////////////////////////////
413/// WVE NEEDS FIX
414
416{
417 Bool_t allCanExtend(kTRUE) ;
418 Bool_t anyMustExtend(kFALSE) ;
419
420 for (Int_t i=0 ; i<_numPdf ; i++) {
422 if (proxy) {
423// cout << " now processing pdf " << pdf->GetName() << endl ;
424 RooAbsPdf* pdf = (RooAbsPdf*) proxy->absArg() ;
425 if (!pdf->canBeExtended()) {
426// cout << "RooSim::extendedMode(" << GetName() << ") component " << pdf->GetName() << " cannot be extended" << endl ;
427 allCanExtend=kFALSE ;
428 }
429 if (pdf->mustBeExtended()) {
430 anyMustExtend=kTRUE;
431 }
432 }
433 }
434 if (anyMustExtend) {
435// cout << "RooSim::extendedMode(" << GetName() << ") returning MustBeExtended" << endl ;
436 return MustBeExtended ;
437 }
438 if (allCanExtend) {
439// cout << "RooSim::extendedMode(" << GetName() << ") returning CanBeExtended" << endl ;
440 return CanBeExtended ;
441 }
442// cout << "RooSim::extendedMode(" << GetName() << ") returning CanNotBeExtended" << endl ;
443 return CanNotBeExtended ;
444}
445
446
447
448
449////////////////////////////////////////////////////////////////////////////////
450/// Return the current value:
451/// the value of the PDF associated with the current index category state
452
454{
455 // Retrieve the proxy by index name
457
458 //assert(proxy!=0) ;
459 if (proxy==0) return 0 ;
460
461 // Calculate relative weighting factor for sim-pdfs of all extendable components
462 Double_t catFrac(1) ;
463 if (canBeExtended()) {
464 Double_t nEvtCat = ((RooAbsPdf*)(proxy->absArg()))->expectedEvents(_normSet) ;
465
466 Double_t nEvtTot(0) ;
468 RooRealProxy* proxy2 ;
469 while((proxy2=(RooRealProxy*)iter->Next())) {
470 nEvtTot += ((RooAbsPdf*)(proxy2->absArg()))->expectedEvents(_normSet) ;
471 }
472 delete iter ;
473 catFrac=nEvtCat/nEvtTot ;
474 }
475
476 // Return the selected PDF value, normalized by the number of index states
477 return ((RooAbsPdf*)(proxy->absArg()))->getVal(_normSet)*catFrac ;
478}
479
480
481
482////////////////////////////////////////////////////////////////////////////////
483/// Return the number of expected events: If the index is in nset,
484/// then return the sum of the expected events of all components,
485/// otherwise return the number of expected events of the PDF
486/// associated with the current index category state
487
489{
490 if (nset->contains(_indexCat.arg())) {
491
492 Double_t sum(0) ;
493
495 RooRealProxy* proxy ;
496 while((proxy=(RooRealProxy*)iter->Next())) {
497 sum += ((RooAbsPdf*)(proxy->absArg()))->expectedEvents(nset) ;
498 }
499 delete iter ;
500
501 return sum ;
502
503 } else {
504
505 // Retrieve the proxy by index name
507
508 //assert(proxy!=0) ;
509 if (proxy==0) return 0 ;
510
511 // Return the selected PDF value, normalized by the number of index states
512 return ((RooAbsPdf*)(proxy->absArg()))->expectedEvents(nset);
513 }
514}
515
516
517
518////////////////////////////////////////////////////////////////////////////////
519/// Forward determination of analytical integration capabilities to component p.d.f.s
520/// A unique code is assigned to the combined integration capabilities of all associated
521/// p.d.f.s
522
524 const RooArgSet* normSet, const char* rangeName) const
525{
526 // Declare that we can analytically integrate all requested observables
527 analVars.add(allVars) ;
528
529 // Retrieve (or create) the required partial integral list
530 Int_t code ;
531
532 // Check if this configuration was created before
533 CacheElem* cache = (CacheElem*) _partIntMgr.getObj(normSet,&analVars,0,RooNameReg::ptr(rangeName)) ;
534 if (cache) {
535 code = _partIntMgr.lastIndex() ;
536 return code+1 ;
537 }
538 cache = new CacheElem ;
539
540 // Create the partial integral set for this request
542 RooRealProxy* proxy ;
543 while((proxy=(RooRealProxy*)iter->Next())) {
544 RooAbsReal* pdfInt = proxy->arg().createIntegral(analVars,normSet,0,rangeName) ;
545 cache->_partIntList.addOwned(*pdfInt) ;
546 }
547 delete iter ;
548
549 // Store the partial integral list and return the assigned code ;
550 code = _partIntMgr.setObj(normSet,&analVars,cache,RooNameReg::ptr(rangeName)) ;
551
552 return code+1 ;
553}
554
555
556
557////////////////////////////////////////////////////////////////////////////////
558/// Return analytical integration defined by given code
559
560Double_t RooSimultaneous::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* /*rangeName*/) const
561{
562 // No integration scenario
563 if (code==0) {
564 return getVal(normSet) ;
565 }
566
567 // Partial integration scenarios, rangeName already encoded in 'code'
568 CacheElem* cache = (CacheElem*) _partIntMgr.getObjByIndex(code-1) ;
569
571 Int_t idx = _pdfProxyList.IndexOf(proxy) ;
572 return ((RooAbsReal*)cache->_partIntList.at(idx))->getVal(normSet) ;
573}
574
575
576
577
578
579
580////////////////////////////////////////////////////////////////////////////////
581/// Back-end for plotOn() implementation on RooSimultaneous which
582/// needs special handling because a RooSimultaneous PDF cannot
583/// project out its index category via integration, plotOn() will
584/// abort if this is requested without providing a projection dataset
585
587{
588 // Sanity checks
589 if (plotSanityChecks(frame)) return frame ;
590
591 // Extract projection configuration from command list
592 RooCmdConfig pc(Form("RooSimultaneous::plotOn(%s)",GetName())) ;
593 pc.defineString("sliceCatState","SliceCat",0,"",kTRUE) ;
594 pc.defineDouble("scaleFactor","Normalization",0,1.0) ;
595 pc.defineInt("scaleType","Normalization",0,RooAbsPdf::Relative) ;
596 pc.defineObject("sliceCatList","SliceCat",0,0,kTRUE) ;
597 pc.defineObject("projSet","Project",0) ;
598 pc.defineObject("sliceSet","SliceVars",0) ;
599 pc.defineObject("projDataSet","ProjData",0) ;
600 pc.defineObject("projData","ProjData",1) ;
601 pc.defineMutex("Project","SliceVars") ;
602 pc.allowUndefined() ; // there may be commands we don't handle here
603
604 // Process and check varargs
605 pc.process(cmdList) ;
606 if (!pc.ok(kTRUE)) {
607 return frame ;
608 }
609
610 const RooAbsData* projData = (const RooAbsData*) pc.getObject("projData") ;
611 const RooArgSet* projDataSet = (const RooArgSet*) pc.getObject("projDataSet") ;
612 const RooArgSet* sliceSetTmp = (const RooArgSet*) pc.getObject("sliceSet") ;
613 std::unique_ptr<RooArgSet> sliceSet( sliceSetTmp ? ((RooArgSet*) sliceSetTmp->Clone()) : nullptr );
614 const RooArgSet* projSet = (const RooArgSet*) pc.getObject("projSet") ;
615 Double_t scaleFactor = pc.getDouble("scaleFactor") ;
616 ScaleType stype = (ScaleType) pc.getInt("scaleType") ;
617
618
619 // Look for category slice arguments and add them to the master slice list if found
620 const char* sliceCatState = pc.getString("sliceCatState",0,kTRUE) ;
621 const RooLinkedList& sliceCatList = pc.getObjectList("sliceCatList") ;
622 if (sliceCatState) {
623
624 // Make the master slice set if it doesnt exist
625 if (!sliceSet) {
626 sliceSet.reset(new RooArgSet);
627 }
628
629 // Prepare comma separated label list for parsing
630 auto catTokens = RooHelpers::tokenise(sliceCatState, ",");
631
632 // Loop over all categories provided by (multiple) Slice() arguments
633 TIterator* iter = sliceCatList.MakeIterator() ;
634 RooCategory* scat ;
635 unsigned int tokenIndex = 0;
636 while((scat=(RooCategory*)iter->Next())) {
637 const char* slabel = tokenIndex >= catTokens.size() ? nullptr : catTokens[tokenIndex++].c_str();
638
639 if (slabel) {
640 // Set the slice position to the value indicate by slabel
641 scat->setLabel(slabel) ;
642 // Add the slice category to the master slice set
643 sliceSet->add(*scat,kFALSE) ;
644 }
645 }
646 delete iter ;
647 }
648
649 // Check if we have a projection dataset
650 if (!projData) {
651 coutE(InputArguments) << "RooSimultaneous::plotOn(" << GetName() << ") ERROR: must have a projection dataset for index category" << endl ;
652 return frame ;
653 }
654
655 // Make list of variables to be projected
656 RooArgSet projectedVars ;
657 if (sliceSet) {
658 makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,kTRUE) ;
659
660 // Take out the sliced variables
661 for (const auto sliceArg : *sliceSet) {
662 RooAbsArg* arg = projectedVars.find(sliceArg->GetName()) ;
663 if (arg) {
664 projectedVars.remove(*arg) ;
665 } else {
666 coutI(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") slice variable "
667 << sliceArg->GetName() << " was not projected anyway" << endl ;
668 }
669 }
670 } else if (projSet) {
671 makeProjectionSet(frame->getPlotVar(),projSet,projectedVars,kFALSE) ;
672 } else {
673 makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,kTRUE) ;
674 }
675
676 Bool_t projIndex(kFALSE) ;
677
678 if (!_indexCat.arg().isDerived()) {
679 // *** Error checking for a fundamental index category ***
680 //cout << "RooSim::plotOn: index is fundamental" << endl ;
681
682 // Check that the provided projection dataset contains our index variable
683 if (!projData->get()->find(_indexCat.arg().GetName())) {
684 coutE(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") ERROR: Projection over index category "
685 << "requested, but projection data set doesn't contain index category" << endl ;
686 return frame ;
687 }
688
689 if (projectedVars.find(_indexCat.arg().GetName())) {
690 projIndex=kTRUE ;
691 }
692
693 } else {
694 // *** Error checking for a composite index category ***
695
696 // Determine if any servers of the index category are in the projectedVars
697 RooArgSet projIdxServers ;
698 Bool_t anyServers(kFALSE) ;
699 for (const auto server : _indexCat->servers()) {
700 if (projectedVars.find(server->GetName())) {
701 anyServers=kTRUE ;
702 projIdxServers.add(*server) ;
703 }
704 }
705
706 // Check that the projection dataset contains all the
707 // index category components we're projecting over
708
709 // Determine if all projected servers of the index category are in the projection dataset
710 Bool_t allServers(kTRUE) ;
711 for (const auto server : projIdxServers) {
712 if (!projData->get()->find(server->GetName())) {
713 allServers=kFALSE ;
714 }
715 }
716
717 if (!allServers) {
718 coutE(Plotting) << "RooSimultaneous::plotOn(" << GetName()
719 << ") ERROR: Projection dataset doesn't contain complete set of index category dependents" << endl ;
720 return frame ;
721 }
722
723 if (anyServers) {
724 projIndex = kTRUE ;
725 }
726 }
727
728 // Calculate relative weight fractions of components
729 std::unique_ptr<Roo1DTable> wTable( projData->table(_indexCat.arg()) );
730
731 // Clone the index category to be able to cycle through the category states for plotting without
732 // affecting the category state of our instance
733 std::unique_ptr<RooArgSet> idxCloneSet( RooArgSet(*_indexCat).snapshot(true) );
734 auto idxCatClone = static_cast<RooAbsCategoryLValue*>( idxCloneSet->find(_indexCat->GetName()) );
735 assert(idxCatClone);
736
737 // Make list of category columns to exclude from projection data
738 std::unique_ptr<RooArgSet> idxCompSliceSet( idxCatClone->getObservables(frame->getNormVars()) );
739
740 // If we don't project over the index, just do the regular plotOn
741 if (!projIndex) {
742
743 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
744 << " represents a slice in the index category (" << _indexCat.arg().GetName() << ")" << endl ;
745
746 // Reduce projData: take out fitCat (component) columns and entries that don't match selected slice
747 // Construct cut string to only select projection data event that match the current slice
748
749 // Make cut string to exclude rows from projection data
750 TString cutString ;
752 for (const auto arg : *idxCompSliceSet) {
753 auto idxComp = static_cast<RooCategory*>(arg);
754 RooAbsArg* slicedComponent = nullptr;
755 if (sliceSet && (slicedComponent = sliceSet->find(*idxComp)) != nullptr) {
756 auto theCat = static_cast<const RooAbsCategory*>(slicedComponent);
757 idxComp->setIndex(theCat->getCurrentIndex(), false);
758 }
759
760 if (!first) {
761 cutString.Append("&&") ;
762 } else {
763 first=kFALSE ;
764 }
765 cutString.Append(Form("%s==%d",idxComp->GetName(),idxComp->getCurrentIndex())) ;
766 }
767
768 // Make temporary projData without RooSim index category components
769 RooArgSet projDataVars(*projData->get()) ;
770 projDataVars.remove(*idxCompSliceSet,kTRUE,kTRUE) ;
771
772 std::unique_ptr<RooAbsData> projDataTmp( const_cast<RooAbsData*>(projData)->reduce(projDataVars,cutString) );
773
774 // Override normalization and projection dataset
775 RooCmdArg tmp1 = RooFit::Normalization(scaleFactor*wTable->getFrac(idxCatClone->getCurrentLabel()),stype) ;
776 RooCmdArg tmp2 = RooFit::ProjWData(*projDataSet,*projDataTmp) ;
777
778 // WVE -- do not adjust normalization for asymmetry plots
779 RooLinkedList cmdList2(cmdList) ;
780 if (!cmdList.find("Asymmetry")) {
781 cmdList2.Add(&tmp1) ;
782 }
783 cmdList2.Add(&tmp2) ;
784
785 // Plot single component
786 RooPlot* retFrame = getPdf(idxCatClone->getCurrentLabel())->plotOn(frame,cmdList2);
787 return retFrame ;
788 }
789
790 // If we project over the index, plot using a temporary RooAddPdf
791 // using the weights from the data as coefficients
792
793 // Build the list of indexCat components that are sliced
794 idxCompSliceSet->remove(projectedVars,kTRUE,kTRUE) ;
795
796 // Make a new expression that is the weighted sum of requested components
797 RooArgList pdfCompList ;
798 RooArgList wgtCompList ;
799//RooAbsPdf* pdf ;
800 RooRealProxy* proxy ;
802 Double_t sumWeight(0) ;
803 while((proxy=(RooRealProxy*)pIter->Next())) {
804
805 idxCatClone->setLabel(proxy->name()) ;
806
807 // Determine if this component is the current slice (if we slice)
808 Bool_t skip(kFALSE) ;
809 for (const auto idxSliceCompArg : *idxCompSliceSet) {
810 const auto idxSliceComp = static_cast<RooAbsCategory*>(idxSliceCompArg);
811 RooAbsCategory* idxComp = (RooAbsCategory*) idxCloneSet->find(idxSliceComp->GetName()) ;
812 if (idxComp->getCurrentIndex()!=idxSliceComp->getCurrentIndex()) {
813 skip=kTRUE ;
814 break ;
815 }
816 }
817 if (skip) continue ;
818
819 // Instantiate a RRV holding this pdfs weight fraction
820 RooRealVar *wgtVar = new RooRealVar(proxy->name(),"coef",wTable->getFrac(proxy->name())) ;
821 wgtCompList.addOwned(*wgtVar) ;
822 sumWeight += wTable->getFrac(proxy->name()) ;
823
824 // Add the PDF to list list
825 pdfCompList.add(proxy->arg()) ;
826 }
827
828 TString plotVarName(GetName()) ;
829 RooAddPdf *plotVar = new RooAddPdf(plotVarName,"weighted sum of RS components",pdfCompList,wgtCompList) ;
830
831 // Fix appropriate coefficient normalization in plot function
832 if (_plotCoefNormSet.getSize()>0) {
834 }
835
836 std::unique_ptr<RooAbsData> projDataTmp;
837 RooArgSet projSetTmp ;
838 if (projData) {
839
840 // Construct cut string to only select projection data event that match the current slice
841 TString cutString ;
842 if (idxCompSliceSet->getSize()>0) {
844 for (const auto idxSliceCompArg : *idxCompSliceSet) {
845 const auto idxSliceComp = static_cast<RooAbsCategory*>(idxSliceCompArg);
846 if (!first) {
847 cutString.Append("&&") ;
848 } else {
849 first=kFALSE ;
850 }
851 cutString.Append(Form("%s==%d",idxSliceComp->GetName(),idxSliceComp->getCurrentIndex())) ;
852 }
853 }
854
855 // Make temporary projData without RooSim index category components
856 RooArgSet projDataVars(*projData->get()) ;
857 RooArgSet* idxCatServers = _indexCat.arg().getObservables(frame->getNormVars()) ;
858
859 projDataVars.remove(*idxCatServers,kTRUE,kTRUE) ;
860
861 if (idxCompSliceSet->getSize()>0) {
862 projDataTmp.reset( const_cast<RooAbsData*>(projData)->reduce(projDataVars,cutString) );
863 } else {
864 projDataTmp.reset( const_cast<RooAbsData*>(projData)->reduce(projDataVars) );
865 }
866
867
868
869 if (projSet) {
870 projSetTmp.add(*projSet) ;
871 projSetTmp.remove(*idxCatServers,kTRUE,kTRUE);
872 }
873
874
875 delete idxCatServers ;
876 }
877
878
879 if (_indexCat.arg().isDerived() && idxCompSliceSet->getSize()>0) {
880 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
881 << " represents a slice in index category components " << *idxCompSliceSet << endl ;
882
883 RooArgSet* idxCompProjSet = _indexCat.arg().getObservables(frame->getNormVars()) ;
884 idxCompProjSet->remove(*idxCompSliceSet,kTRUE,kTRUE) ;
885 if (idxCompProjSet->getSize()>0) {
886 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
887 << " averages with data index category components " << *idxCompProjSet << endl ;
888 }
889 delete idxCompProjSet ;
890 } else {
891 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
892 << " averages with data index category (" << _indexCat.arg().GetName() << ")" << endl ;
893 }
894
895
896 // Override normalization and projection dataset
897 RooLinkedList cmdList2(cmdList) ;
898
899 RooCmdArg tmp1 = RooFit::Normalization(scaleFactor*sumWeight,stype) ;
900 RooCmdArg tmp2 = RooFit::ProjWData(*projDataSet,*projDataTmp) ;
901 // WVE -- do not adjust normalization for asymmetry plots
902 if (!cmdList.find("Asymmetry")) {
903 cmdList2.Add(&tmp1) ;
904 }
905 cmdList2.Add(&tmp2) ;
906
907 RooPlot* frame2 ;
908 if (projSetTmp.getSize()>0) {
909 // Plot temporary function
910 RooCmdArg tmp3 = RooFit::Project(projSetTmp) ;
911 cmdList2.Add(&tmp3) ;
912 frame2 = plotVar->plotOn(frame,cmdList2) ;
913 } else {
914 // Plot temporary function
915 frame2 = plotVar->plotOn(frame,cmdList2) ;
916 }
917
918 // Cleanup
919 delete plotVar ;
920
921 return frame2 ;
922}
923
924
925
926////////////////////////////////////////////////////////////////////////////////
927/// OBSOLETE -- Retained for backward compatibility
928
929RooPlot* RooSimultaneous::plotOn(RooPlot *frame, Option_t* drawOptions, Double_t scaleFactor,
930 ScaleType stype, const RooAbsData* projData, const RooArgSet* projSet,
931 Double_t /*precision*/, Bool_t /*shiftToZero*/, const RooArgSet* /*projDataSet*/,
932 Double_t /*rangeLo*/, Double_t /*rangeHi*/, RooCurve::WingMode /*wmode*/) const
933{
934 // Make command list
935 RooLinkedList cmdList ;
936 cmdList.Add(new RooCmdArg(RooFit::DrawOption(drawOptions))) ;
937 cmdList.Add(new RooCmdArg(RooFit::Normalization(scaleFactor,stype))) ;
938 if (projData) cmdList.Add(new RooCmdArg(RooFit::ProjWData(*projData))) ;
939 if (projSet) cmdList.Add(new RooCmdArg(RooFit::Project(*projSet))) ;
940
941 // Call new method
942 RooPlot* ret = plotOn(frame,cmdList) ;
943
944 // Cleanup
945 cmdList.Delete() ;
946 return ret ;
947}
948
949
950
951////////////////////////////////////////////////////////////////////////////////
952/// Interface function used by test statistics to freeze choice of observables
953/// for interpretation of fraction coefficients. Needed here because a RooSimultaneous
954/// works like a RooAddPdf when plotted
955
957{
959 if (normSet) _plotCoefNormSet.add(*normSet) ;
960}
961
962
963////////////////////////////////////////////////////////////////////////////////
964/// Interface function used by test statistics to freeze choice of range
965/// for interpretation of fraction coefficients. Needed here because a RooSimultaneous
966/// works like a RooAddPdf when plotted
967
968void RooSimultaneous::selectNormalizationRange(const char* normRange2, Bool_t /*force*/)
969{
970 _plotCoefNormRange = RooNameReg::ptr(normRange2) ;
971}
972
973
974
975
976////////////////////////////////////////////////////////////////////////////////
977
979 const RooArgSet* auxProto, Bool_t verbose, Bool_t autoBinned, const char* binnedTag) const
980{
981 const char* idxCatName = _indexCat.arg().GetName() ;
982
983 if (vars.find(idxCatName) && prototype==0
984 && (auxProto==0 || auxProto->getSize()==0)
985 && (autoBinned || (binnedTag && strlen(binnedTag)))) {
986
987 // Return special generator config that can also do binned generation for selected states
988 return new RooSimSplitGenContext(*this,vars,verbose,autoBinned,binnedTag) ;
989
990 } else {
991
992 // Return regular generator config ;
993 return genContext(vars,prototype,auxProto,verbose) ;
994 }
995}
996
997
998
999////////////////////////////////////////////////////////////////////////////////
1000/// Return specialized generator context for simultaneous p.d.f.s
1001
1003 const RooArgSet* auxProto, Bool_t verbose) const
1004{
1005 const char* idxCatName = _indexCat.arg().GetName() ;
1006 const RooArgSet* protoVars = prototype ? prototype->get() : 0 ;
1007
1008 if (vars.find(idxCatName) || (protoVars && protoVars->find(idxCatName))) {
1009
1010 // Generating index category: return special sim-context
1011 return new RooSimGenContext(*this,vars,prototype,auxProto,verbose) ;
1012
1013 } else if (_indexCat.arg().isDerived()) {
1014 // Generating dependents of a derived index category
1015
1016 // Determine if we none,any or all servers
1017 Bool_t anyServer(kFALSE), allServers(kTRUE) ;
1018 if (prototype) {
1019 TIterator* sIter = _indexCat.arg().serverIterator() ;
1020 RooAbsArg* server ;
1021 while((server=(RooAbsArg*)sIter->Next())) {
1022 if (prototype->get()->find(server->GetName())) {
1023 anyServer=kTRUE ;
1024 } else {
1025 allServers=kFALSE ;
1026 }
1027 }
1028 delete sIter ;
1029 } else {
1030 allServers=kTRUE ;
1031 }
1032
1033 if (allServers) {
1034 // Use simcontext if we have all servers
1035
1036 return new RooSimGenContext(*this,vars,prototype,auxProto,verbose) ;
1037 } else if (!allServers && anyServer) {
1038 // Abort if we have only part of the servers
1039 coutE(Plotting) << "RooSimultaneous::genContext: ERROR: prototype must include either all "
1040 << " components of the RooSimultaneous index category or none " << endl ;
1041 return 0 ;
1042 }
1043 // Otherwise make single gencontext for current state
1044 }
1045
1046 // Not generating index cat: return context for pdf associated with present index state
1048 if (!proxy) {
1049 coutE(InputArguments) << "RooSimultaneous::genContext(" << GetName()
1050 << ") ERROR: no PDF associated with current state ("
1051 << _indexCat.arg().GetName() << "=" << _indexCat.arg().getCurrentLabel() << ")" << endl ;
1052 return 0 ;
1053 }
1054 return ((RooAbsPdf*)proxy->absArg())->genContext(vars,prototype,auxProto,verbose) ;
1055}
1056
1057
1058
1059
1060////////////////////////////////////////////////////////////////////////////////
1061
1063 const RooArgSet* nset,
1064 Double_t scaleFactor,
1065 Bool_t correctForBinVolume,
1066 Bool_t showProgress) const
1067{
1068 if (RooAbsReal::fillDataHist (hist, nset, scaleFactor,
1069 correctForBinVolume, showProgress) == 0)
1070 return 0;
1071
1072 Double_t sum = 0;
1073 for (int i=0 ; i<hist->numEntries() ; i++) {
1074 hist->get(i) ;
1075 sum += hist->weight();
1076 }
1077 if (sum != 0) {
1078 for (int i=0 ; i<hist->numEntries() ; i++) {
1079 hist->get(i) ;
1080 hist->set (hist->weight() / sum);
1081 }
1082 }
1083
1084 return hist;
1085}
1086
1087
1088
1089
1090////////////////////////////////////////////////////////////////////////////////
1091/// Special generator interface for generation of 'global observables' -- for RooStats tools
1092
1094{
1095 // Make set with clone of variables (placeholder for output)
1096 RooArgSet* globClone = (RooArgSet*) whatVars.snapshot() ;
1097
1098 RooDataSet* data = new RooDataSet("gensimglobal","gensimglobal",whatVars) ;
1099
1100 for (Int_t i=0 ; i<nEvents ; i++) {
1101 for (const auto& nameIdx : indexCat()) {
1102
1103 // Get pdf associated with state from simpdf
1104 RooAbsPdf* pdftmp = getPdf(nameIdx.first.c_str());
1105
1106 // Generate only global variables defined by the pdf associated with this state
1107 RooArgSet* globtmp = pdftmp->getObservables(whatVars) ;
1108 RooDataSet* tmp = pdftmp->generate(*globtmp,1) ;
1109
1110 // Transfer values to output placeholder
1111 *globClone = *tmp->get(0) ;
1112
1113 // Cleanup
1114 delete globtmp ;
1115 delete tmp ;
1116 }
1117 data->add(*globClone) ;
1118 }
1119
1120 delete globClone ;
1121 return data ;
1122}
1123
1124
1125
1126
1127
1128
1129
1130
1131
#define coutI(a)
Definition: RooMsgService.h:30
#define cxcoutD(a)
Definition: RooMsgService.h:81
#define coutW(a)
Definition: RooMsgService.h:32
#define coutE(a)
Definition: RooMsgService.h:33
RooTemplateProxy< RooAbsReal > RooRealProxy
Compatibility typedef replacing the old RooRealProxy class.
Definition: RooRealProxy.h:23
const Bool_t kFALSE
Definition: RtypesCore.h:90
const Bool_t kTRUE
Definition: RtypesCore.h:89
const char Option_t
Definition: RtypesCore.h:64
#define ClassImp(name)
Definition: Rtypes.h:361
char name[80]
Definition: TGX11.cxx:109
int type
Definition: TGX11.cxx:120
char * Form(const char *fmt,...)
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:73
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Return the observables of this pdf given a set of observables.
Definition: RooAbsArg.h:276
friend class RooDataSet
Definition: RooAbsArg.h:620
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.
Definition: RooAbsArg.cxx:730
friend class RooArgSet
Definition: RooAbsArg.h:572
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:243
RooArgSet * getVariables(Bool_t stripDisconnected=kTRUE) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
Definition: RooAbsArg.cxx:1909
virtual Bool_t isDerived() const
Does value or shape of this arg depend on any other arg?
Definition: RooAbsArg.h:93
const RefCountList_t & servers() const
List of all servers of this object.
Definition: RooAbsArg.h:199
Bool_t addOwnedComponents(const RooArgSet &comps)
Take ownership of the contents of 'comps'.
Definition: RooAbsArg.cxx:2107
TIterator * serverIterator() const
Definition: RooAbsArg.h:141
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 Bool_t addOwned(RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
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.
Storage_t::size_type size() 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:44
virtual const RooArgSet * get() const
Definition: RooAbsData.h:87
virtual Roo1DTable * table(const RooArgSet &catSet, const char *cuts="", const char *opts="") const
Construct table for product of categories in catSet.
Definition: RooAbsData.cxx:748
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
friend class CacheElem
Definition: RooAbsPdf.h:333
Bool_t canBeExtended() const
Definition: RooAbsPdf.h:236
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:55
@ CanBeExtended
Definition: RooAbsPdf.h:229
@ MustBeExtended
Definition: RooAbsPdf.h:229
@ CanNotBeExtended
Definition: RooAbsPdf.h:229
Bool_t mustBeExtended() const
Definition: RooAbsPdf.h:240
RooArgSet * _normSet
Normalization integral (owned by _normMgr)
Definition: RooAbsPdf.h:321
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:118
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:60
friend class RooAddPdf
Definition: RooAbsReal.h:524
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:90
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 ...
Definition: RooAbsReal.cxx:560
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:29
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition: RooArgList.h:74
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition: RooArgSet.h:126
virtual Bool_t add(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling add() for each element in the source coll...
Definition: RooArgSet.h:88
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:23
virtual Bool_t setLabel(const char *label, bool printError=true) override
Set value by specifying the name of the desired state.
RooCmdArg is a named container for two doubles, two integers two object points and three string point...
Definition: RooCmdArg.h:28
Class RooCmdConfig is a configurable parser for RooCmdArg named arguments.
Definition: RooCmdConfig.h:27
The RooDataHist is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
virtual Double_t weight() const
Definition: RooDataHist.h:106
void set(Double_t weight, Double_t wgtErr=-1)
Set the weight and weight error of the bin enclosing the current (i.e.
virtual Int_t numEntries() const
Return the number of bins.
virtual const RooArgSet * get() const
Definition: RooDataHist.h:79
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:33
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 ...
Definition: RooLinkedList.h:35
TIterator * MakeIterator(Bool_t forward=kTRUE) const
Create a TIterator for this list.
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)
Definition: RooLinkedList.h:62
static const TNamed * ptr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
Definition: RooNameReg.cxx:93
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:148
RooAbsRealLValue * getPlotVar() const
Definition: RooPlot.h:139
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:950
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Overloaded RooArgSet::add() method inserts 'var' into set and registers 'var' as server to owner with...
virtual void removeAll()
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
RooCategoryProxy _indexCat
virtual ExtendMode extendMode() const
WVE NEEDS FIX.
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
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
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
Helper calling plotOn(RooPlot*, RooLinkedList&) const.
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.
Iterator abstract base class.
Definition: TIterator.h:30
virtual TObject * Next()=0
virtual void Add(TObject *obj)
Definition: TList.h:87
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition: TList.cxx:577
virtual TIterator * MakeIterator(Bool_t dir=kIterForward) const
Return a list iterator.
Definition: TList.cxx:721
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition: TList.cxx:469
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Mother of all ROOT objects.
Definition: TObject.h:37
virtual Int_t IndexOf(const TObject *obj) const
Return index of object in collection.
Basic string class.
Definition: TString.h:131
TString & Append(const char *cs)
Definition: TString.h:559
RooCmdArg ProjWData(const RooAbsData &projData, Bool_t binData=kFALSE)
RooCmdArg DrawOption(const char *opt)
RooCmdArg Project(const RooArgSet &projSet)
RooCmdArg Normalization(Double_t scaleFactor)
@ InputArguments
Definition: RooGlobalFunc.h:68
std::vector< std::string > tokenise(const std::string &str, const std::string &delims, bool returnEmptyToken=true)
Tokenise the string by splitting at the characters in delims.
Definition: RooHelpers.cxx:62
static constexpr double pc
Definition: first.py:1
static long int sum(long int i)
Definition: Factory.cxx:2275