Logo ROOT   6.18/05
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#include "Riostream.h"
51
52#include "TObjString.h"
53#include "RooSimultaneous.h"
55#include "RooPlot.h"
56#include "RooCurve.h"
57#include "RooRealVar.h"
58#include "RooAddPdf.h"
59#include "RooAbsData.h"
60#include "Roo1DTable.h"
61#include "RooSimGenContext.h"
63#include "RooDataSet.h"
64#include "RooCmdConfig.h"
65#include "RooNameReg.h"
66#include "RooGlobalFunc.h"
67#include "RooNameReg.h"
68#include "RooMsgService.h"
69#include "RooCategory.h"
70#include "RooSuperCategory.h"
71#include "RooDataHist.h"
72#include "RooRandom.h"
73#include "RooArgSet.h"
74
75using namespace std ;
76
78;
79
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 in order with the state of the
109/// index category as listed by the index categories type iterator.
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.getSize() != inIndexCat.numTypes()) {
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 // Iterator over PDFs and index cat states and add each pair
130 TIterator* pIter = inPdfList.createIterator() ;
131 TIterator* cIter = inIndexCat.typeIterator() ;
132 RooAbsPdf* pdf ;
133 RooCatType* type(0) ;
134 while ((pdf=(RooAbsPdf*)pIter->Next())) {
135 type = (RooCatType*) cIter->Next() ;
136 pdfMap[string(type->GetName())] = pdf ;
137 }
138 delete pIter ;
139 delete cIter ;
140
141 initialize(inIndexCat,pdfMap) ;
142}
143
144
145////////////////////////////////////////////////////////////////////////////////
146
147RooSimultaneous::RooSimultaneous(const char *name, const char *title,
148 map<string,RooAbsPdf*> pdfMap, RooAbsCategoryLValue& inIndexCat) :
149 RooAbsPdf(name,title),
150 _plotCoefNormSet("!plotCoefNormSet","plotCoefNormSet",this,kFALSE,kFALSE),
151 _plotCoefNormRange(0),
152 _partIntMgr(this,10),
153 _indexCat("indexCat","Index category",this,inIndexCat),
154 _numPdf(0)
155{
156 initialize(inIndexCat,pdfMap) ;
157}
158
159
160
161
162// This class cannot be locally defined in initialize as it cannot be
163// used as a template argument in that case
165 struct CompInfo {
166 RooAbsPdf* pdf ;
167 RooSimultaneous* simPdf ;
168 const RooAbsCategoryLValue* subIndex ;
169 RooArgSet* subIndexComps ;
170 } ;
171}
172
173void RooSimultaneous::initialize(RooAbsCategoryLValue& inIndexCat, std::map<std::string,RooAbsPdf*> pdfMap)
174{
175 // First see if there are any RooSimultaneous input components
176 Bool_t simComps(kFALSE) ;
177 for (map<string,RooAbsPdf*>::iterator iter=pdfMap.begin() ; iter!=pdfMap.end() ; ++iter) {
178 if (dynamic_cast<RooSimultaneous*>(iter->second)) {
179 simComps = kTRUE ;
180 break ;
181 }
182 }
183
184 // If there are no simultaneous component p.d.f. do simple processing through addPdf()
185 if (!simComps) {
186 for (map<string,RooAbsPdf*>::iterator iter=pdfMap.begin() ; iter!=pdfMap.end() ; ++iter) {
187 addPdf(*iter->second,iter->first.c_str()) ;
188 }
189 return ;
190 }
191
192 // Issue info message that we are about to do some rearraning
193 coutI(InputArguments) << "RooSimultaneous::initialize(" << GetName() << ") INFO: one or more input component of simultaneous p.d.f.s are"
194 << " simultaneous p.d.f.s themselves, rewriting composite expressions as one-level simultaneous p.d.f. in terms of"
195 << " final constituents and extended index category" << endl ;
196
197
198 RooArgSet allAuxCats ;
199 map<string,RooSimultaneousAux::CompInfo> compMap ;
200 for (map<string,RooAbsPdf*>::iterator iter=pdfMap.begin() ; iter!=pdfMap.end() ; ++iter) {
201 RooSimultaneousAux::CompInfo ci ;
202 ci.pdf = iter->second ;
203 RooSimultaneous* simComp = dynamic_cast<RooSimultaneous*>(iter->second) ;
204 if (simComp) {
205 ci.simPdf = simComp ;
206 ci.subIndex = &simComp->indexCat() ;
207 ci.subIndexComps = simComp->indexCat().isFundamental() ? new RooArgSet(simComp->indexCat()) : simComp->indexCat().getVariables() ;
208 allAuxCats.add(*(ci.subIndexComps),kTRUE) ;
209 } else {
210 ci.simPdf = 0 ;
211 ci.subIndex = 0 ;
212 ci.subIndexComps = 0 ;
213 }
214 compMap[iter->first] = ci ;
215 }
216
217 // Construct the 'superIndex' from the nominal index category and all auxiliary components
218 RooArgSet allCats(inIndexCat) ;
219 allCats.add(allAuxCats) ;
220 string siname = Form("%s_index",GetName()) ;
221 RooSuperCategory* superIndex = new RooSuperCategory(siname.c_str(),siname.c_str(),allCats) ;
222
223 // Now process each of original pdf/state map entries
224 for (map<string,RooSimultaneousAux::CompInfo>::iterator citer = compMap.begin() ; citer != compMap.end() ; ++citer) {
225
226 RooArgSet repliCats(allAuxCats) ;
227 if (citer->second.subIndexComps) {
228 repliCats.remove(*citer->second.subIndexComps) ;
229 delete citer->second.subIndexComps ;
230 }
231 inIndexCat.setLabel(citer->first.c_str()) ;
232
233
234 if (!citer->second.simPdf) {
235
236 // Entry is a plain p.d.f. assign it to every state permutation of the repliCats set
237 RooSuperCategory repliSuperCat("tmp","tmp",repliCats) ;
238
239 // Iterator over all states of repliSuperCat
240 TIterator* titer = repliSuperCat.typeIterator() ;
242 while ((type=(RooCatType*)titer->Next())) {
243 // Set value
244 repliSuperCat.setLabel(type->GetName()) ;
245 // Retrieve corresponding label of superIndex
246 string superLabel = superIndex->getLabel() ;
247 addPdf(*citer->second.pdf,superLabel.c_str()) ;
248 cxcoutD(InputArguments) << "RooSimultaneous::initialize(" << GetName()
249 << ") assigning pdf " << citer->second.pdf->GetName() << " to super label " << superLabel << endl ;
250 }
251 } else {
252
253 // Entry is a simultaneous p.d.f
254
255 if (repliCats.getSize()==0) {
256
257 // Case 1 -- No replication of components of RooSim component are required
258
259 TIterator* titer = citer->second.subIndex->typeIterator() ;
261 while ((type=(RooCatType*)titer->Next())) {
262 const_cast<RooAbsCategoryLValue*>(citer->second.subIndex)->setLabel(type->GetName()) ;
263 string superLabel = superIndex->getLabel() ;
264 RooAbsPdf* compPdf = citer->second.simPdf->getPdf(type->GetName()) ;
265 if (compPdf) {
266 addPdf(*compPdf,superLabel.c_str()) ;
267 cxcoutD(InputArguments) << "RooSimultaneous::initialize(" << GetName()
268 << ") assigning pdf " << compPdf->GetName() << "(member of " << citer->second.pdf->GetName()
269 << ") to super label " << superLabel << endl ;
270 } else {
271 coutW(InputArguments) << "RooSimultaneous::initialize(" << GetName() << ") WARNING: No p.d.f. associated with label "
272 << type->GetName() << " for component RooSimultaneous p.d.f " << citer->second.pdf->GetName()
273 << "which is associated with master index label " << citer->first << endl ;
274 }
275 }
276 delete titer ;
277
278 } else {
279
280 // Case 2 -- Replication of components of RooSim component are required
281
282 // Make replication supercat
283 RooSuperCategory repliSuperCat("tmp","tmp",repliCats) ;
284 TIterator* triter = repliSuperCat.typeIterator() ;
285
286 TIterator* tsiter = citer->second.subIndex->typeIterator() ;
287 RooCatType* stype, *rtype ;
288 while ((stype=(RooCatType*)tsiter->Next())) {
289 const_cast<RooAbsCategoryLValue*>(citer->second.subIndex)->setLabel(stype->GetName()) ;
290 triter->Reset() ;
291 while ((rtype=(RooCatType*)triter->Next())) {
292 repliSuperCat.setLabel(rtype->GetName()) ;
293 string superLabel = superIndex->getLabel() ;
294 RooAbsPdf* compPdf = citer->second.simPdf->getPdf(stype->GetName()) ;
295 if (compPdf) {
296 addPdf(*compPdf,superLabel.c_str()) ;
297 cxcoutD(InputArguments) << "RooSimultaneous::initialize(" << GetName()
298 << ") assigning pdf " << compPdf->GetName() << "(member of " << citer->second.pdf->GetName()
299 << ") to super label " << superLabel << endl ;
300 } else {
301 coutW(InputArguments) << "RooSimultaneous::initialize(" << GetName() << ") WARNING: No p.d.f. associated with label "
302 << stype->GetName() << " for component RooSimultaneous p.d.f " << citer->second.pdf->GetName()
303 << "which is associated with master index label " << citer->first << endl ;
304 }
305 }
306 }
307
308 delete tsiter ;
309 delete triter ;
310
311 }
312 }
313 }
314
315 // Change original master index to super index and take ownership of it
316 _indexCat.setArg(*superIndex) ;
317 addOwnedComponents(*superIndex) ;
318
319}
320
321
322
323////////////////////////////////////////////////////////////////////////////////
324/// Copy constructor
325
327 RooAbsPdf(other,name),
328 _plotCoefNormSet("!plotCoefNormSet",this,other._plotCoefNormSet),
329 _plotCoefNormRange(other._plotCoefNormRange),
330 _partIntMgr(other._partIntMgr,this),
331 _indexCat("indexCat",this,other._indexCat),
332 _numPdf(other._numPdf)
333{
334 // Copy proxy list
335 TIterator* pIter = other._pdfProxyList.MakeIterator() ;
336 RooRealProxy* proxy ;
337 while ((proxy=(RooRealProxy*)pIter->Next())) {
338 _pdfProxyList.Add(new RooRealProxy(proxy->GetName(),this,*proxy)) ;
339 }
340 delete pIter ;
341}
342
343
344
345////////////////////////////////////////////////////////////////////////////////
346/// Destructor
347
349{
351}
352
353
354
355////////////////////////////////////////////////////////////////////////////////
356/// Return the p.d.f associated with the given index category name
357
358RooAbsPdf* RooSimultaneous::getPdf(const char* catName) const
359{
361 return proxy ? ((RooAbsPdf*)proxy->absArg()) : 0 ;
362}
363
364
365
366////////////////////////////////////////////////////////////////////////////////
367/// Associate given PDF with index category state label 'catLabel'.
368/// The name state must be already defined in the index category.
369///
370/// RooSimultaneous can function without having a PDF associated
371/// with every single state. The normalization in such cases is taken
372/// from the number of registered PDFs, but getVal() will assert if
373/// when called for an unregistered index state.
374///
375/// PDFs may not overlap (i.e. share any variables) with the index category (function).
376
377Bool_t RooSimultaneous::addPdf(const RooAbsPdf& pdf, const char* catLabel)
378{
379 // PDFs cannot overlap with the index category
380 if (pdf.dependsOn(_indexCat.arg())) {
381 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName() << "): ERROR, PDF " << pdf.GetName()
382 << " overlaps with index category " << _indexCat.arg().GetName() << endl ;
383 return kTRUE ;
384 }
385
386 // Each index state can only have one PDF associated with it
387 if (_pdfProxyList.FindObject(catLabel)) {
388 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName() << "): ERROR, index state "
389 << catLabel << " has already an associated PDF" << endl ;
390 return kTRUE ;
391 }
392
393 const RooSimultaneous* simPdf = dynamic_cast<const RooSimultaneous*>(&pdf) ;
394 if (simPdf) {
395
396 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName()
397 << ") ERROR: you cannot add a RooSimultaneous component to a RooSimultaneous using addPdf()."
398 << " Use the constructor with RooArgList if input p.d.f.s or the map<string,RooAbsPdf&> instead." << endl ;
399 return kTRUE ;
400
401 } else {
402
403 // Create a proxy named after the associated index state
404 TObject* proxy = new RooRealProxy(catLabel,catLabel,this,(RooAbsPdf&)pdf) ;
405 _pdfProxyList.Add(proxy) ;
406 _numPdf += 1 ;
407 }
408
409 return kFALSE ;
410}
411
412
413
414
415
416////////////////////////////////////////////////////////////////////////////////
417/// WVE NEEDS FIX
418
420{
421 Bool_t allCanExtend(kTRUE) ;
422 Bool_t anyMustExtend(kFALSE) ;
423
424 for (Int_t i=0 ; i<_numPdf ; i++) {
426 if (proxy) {
427// cout << " now processing pdf " << pdf->GetName() << endl ;
428 RooAbsPdf* pdf = (RooAbsPdf*) proxy->absArg() ;
429 if (!pdf->canBeExtended()) {
430// cout << "RooSim::extendedMode(" << GetName() << ") component " << pdf->GetName() << " cannot be extended" << endl ;
431 allCanExtend=kFALSE ;
432 }
433 if (pdf->mustBeExtended()) {
434 anyMustExtend=kTRUE;
435 }
436 }
437 }
438 if (anyMustExtend) {
439// cout << "RooSim::extendedMode(" << GetName() << ") returning MustBeExtended" << endl ;
440 return MustBeExtended ;
441 }
442 if (allCanExtend) {
443// cout << "RooSim::extendedMode(" << GetName() << ") returning CanBeExtended" << endl ;
444 return CanBeExtended ;
445 }
446// cout << "RooSim::extendedMode(" << GetName() << ") returning CanNotBeExtended" << endl ;
447 return CanNotBeExtended ;
448}
449
450
451
452
453////////////////////////////////////////////////////////////////////////////////
454/// Return the current value:
455/// the value of the PDF associated with the current index category state
456
458{
459 // Retrieve the proxy by index name
461
462 //assert(proxy!=0) ;
463 if (proxy==0) return 0 ;
464
465 // Calculate relative weighting factor for sim-pdfs of all extendable components
466 Double_t catFrac(1) ;
467 if (canBeExtended()) {
468 Double_t nEvtCat = ((RooAbsPdf*)(proxy->absArg()))->expectedEvents(_normSet) ;
469
470 Double_t nEvtTot(0) ;
472 RooRealProxy* proxy2 ;
473 while((proxy2=(RooRealProxy*)iter->Next())) {
474 nEvtTot += ((RooAbsPdf*)(proxy2->absArg()))->expectedEvents(_normSet) ;
475 }
476 delete iter ;
477 catFrac=nEvtCat/nEvtTot ;
478 }
479
480 // Return the selected PDF value, normalized by the number of index states
481 return ((RooAbsPdf*)(proxy->absArg()))->getVal(_normSet)*catFrac ;
482}
483
484
485
486////////////////////////////////////////////////////////////////////////////////
487/// Return the number of expected events: If the index is in nset,
488/// then return the sum of the expected events of all components,
489/// otherwise return the number of expected events of the PDF
490/// associated with the current index category state
491
493{
494 if (nset->contains(_indexCat.arg())) {
495
496 Double_t sum(0) ;
497
499 RooRealProxy* proxy ;
500 while((proxy=(RooRealProxy*)iter->Next())) {
501 sum += ((RooAbsPdf*)(proxy->absArg()))->expectedEvents(nset) ;
502 }
503 delete iter ;
504
505 return sum ;
506
507 } else {
508
509 // Retrieve the proxy by index name
511
512 //assert(proxy!=0) ;
513 if (proxy==0) return 0 ;
514
515 // Return the selected PDF value, normalized by the number of index states
516 return ((RooAbsPdf*)(proxy->absArg()))->expectedEvents(nset);
517 }
518}
519
520
521
522////////////////////////////////////////////////////////////////////////////////
523/// Forward determination of analytical integration capabilities to component p.d.f.s
524/// A unique code is assigned to the combined integration capabilities of all associated
525/// p.d.f.s
526
528 const RooArgSet* normSet, const char* rangeName) const
529{
530 // Declare that we can analytically integrate all requested observables
531 analVars.add(allVars) ;
532
533 // Retrieve (or create) the required partial integral list
534 Int_t code ;
535
536 // Check if this configuration was created before
537 CacheElem* cache = (CacheElem*) _partIntMgr.getObj(normSet,&analVars,0,RooNameReg::ptr(rangeName)) ;
538 if (cache) {
539 code = _partIntMgr.lastIndex() ;
540 return code+1 ;
541 }
542 cache = new CacheElem ;
543
544 // Create the partial integral set for this request
546 RooRealProxy* proxy ;
547 while((proxy=(RooRealProxy*)iter->Next())) {
548 RooAbsReal* pdfInt = proxy->arg().createIntegral(analVars,normSet,0,rangeName) ;
549 cache->_partIntList.addOwned(*pdfInt) ;
550 }
551 delete iter ;
552
553 // Store the partial integral list and return the assigned code ;
554 code = _partIntMgr.setObj(normSet,&analVars,cache,RooNameReg::ptr(rangeName)) ;
555
556 return code+1 ;
557}
558
559
560
561////////////////////////////////////////////////////////////////////////////////
562/// Return analytical integration defined by given code
563
564Double_t RooSimultaneous::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* /*rangeName*/) const
565{
566 // No integration scenario
567 if (code==0) {
568 return getVal(normSet) ;
569 }
570
571 // Partial integration scenarios, rangeName already encoded in 'code'
572 CacheElem* cache = (CacheElem*) _partIntMgr.getObjByIndex(code-1) ;
573
575 Int_t idx = _pdfProxyList.IndexOf(proxy) ;
576 return ((RooAbsReal*)cache->_partIntList.at(idx))->getVal(normSet) ;
577}
578
579
580
581
582
583
584////////////////////////////////////////////////////////////////////////////////
585/// Back-end for plotOn() implementation on RooSimultaneous which
586/// needs special handling because a RooSimultaneous PDF cannot
587/// project out its index category via integration, plotOn() will
588/// abort if this is requested without providing a projection dataset
589
591{
592 // Sanity checks
593 if (plotSanityChecks(frame)) return frame ;
594
595 // Extract projection configuration from command list
596 RooCmdConfig pc(Form("RooSimultaneous::plotOn(%s)",GetName())) ;
597 pc.defineString("sliceCatState","SliceCat",0,"",kTRUE) ;
598 pc.defineDouble("scaleFactor","Normalization",0,1.0) ;
599 pc.defineInt("scaleType","Normalization",0,RooAbsPdf::Relative) ;
600 pc.defineObject("sliceCatList","SliceCat",0,0,kTRUE) ;
601 pc.defineObject("projSet","Project",0) ;
602 pc.defineObject("sliceSet","SliceVars",0) ;
603 pc.defineObject("projDataSet","ProjData",0) ;
604 pc.defineObject("projData","ProjData",1) ;
605 pc.defineMutex("Project","SliceVars") ;
606 pc.allowUndefined() ; // there may be commands we don't handle here
607
608 // Process and check varargs
609 pc.process(cmdList) ;
610 if (!pc.ok(kTRUE)) {
611 return frame ;
612 }
613
614 const RooAbsData* projData = (const RooAbsData*) pc.getObject("projData") ;
615 const RooArgSet* projDataSet = (const RooArgSet*) pc.getObject("projDataSet") ;
616 const RooArgSet* sliceSetTmp = (const RooArgSet*) pc.getObject("sliceSet") ;
617 RooArgSet* sliceSet = sliceSetTmp ? ((RooArgSet*) sliceSetTmp->Clone()) : 0 ;
618 const RooArgSet* projSet = (const RooArgSet*) pc.getObject("projSet") ;
619 Double_t scaleFactor = pc.getDouble("scaleFactor") ;
620 ScaleType stype = (ScaleType) pc.getInt("scaleType") ;
621
622
623 // Look for category slice arguments and add them to the master slice list if found
624 const char* sliceCatState = pc.getString("sliceCatState",0,kTRUE) ;
625 const RooLinkedList& sliceCatList = pc.getObjectList("sliceCatList") ;
626 if (sliceCatState) {
627
628 // Make the master slice set if it doesnt exist
629 if (!sliceSet) {
630 sliceSet = new RooArgSet ;
631 }
632
633 // Prepare comma separated label list for parsing
634 char buf[1024] ;
635 strlcpy(buf,sliceCatState,1024) ;
636 const char* slabel = strtok(buf,",") ;
637
638 // Loop over all categories provided by (multiple) Slice() arguments
639 TIterator* iter = sliceCatList.MakeIterator() ;
640 RooCategory* scat ;
641 while((scat=(RooCategory*)iter->Next())) {
642 if (slabel) {
643 // Set the slice position to the value indicate by slabel
644 scat->setLabel(slabel) ;
645 // Add the slice category to the master slice set
646 sliceSet->add(*scat,kFALSE) ;
647 }
648 slabel = strtok(0,",") ;
649 }
650 delete iter ;
651 }
652
653 // Check if we have a projection dataset
654 if (!projData) {
655 coutE(InputArguments) << "RooSimultaneous::plotOn(" << GetName() << ") ERROR: must have a projection dataset for index category" << endl ;
656 return frame ;
657 }
658
659 // Make list of variables to be projected
660 RooArgSet projectedVars ;
661 if (sliceSet) {
662 //cout << "frame->getNormVars() = " ; frame->getNormVars()->Print("1") ;
663
664 makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,kTRUE) ;
665
666 // Take out the sliced variables
667 TIterator* iter = sliceSet->createIterator() ;
668 RooAbsArg* sliceArg ;
669 while((sliceArg=(RooAbsArg*)iter->Next())) {
670 RooAbsArg* arg = projectedVars.find(sliceArg->GetName()) ;
671 if (arg) {
672 projectedVars.remove(*arg) ;
673 } else {
674 coutI(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") slice variable "
675 << sliceArg->GetName() << " was not projected anyway" << endl ;
676 }
677 }
678 delete iter ;
679 } else if (projSet) {
680 makeProjectionSet(frame->getPlotVar(),projSet,projectedVars,kFALSE) ;
681 } else {
682 makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,kTRUE) ;
683 }
684
685 Bool_t projIndex(kFALSE) ;
686
687 if (!_indexCat.arg().isDerived()) {
688 // *** Error checking for a fundamental index category ***
689 //cout << "RooSim::plotOn: index is fundamental" << endl ;
690
691 // Check that the provided projection dataset contains our index variable
692 if (!projData->get()->find(_indexCat.arg().GetName())) {
693 coutE(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") ERROR: Projection over index category "
694 << "requested, but projection data set doesn't contain index category" << endl ;
695 return frame ;
696 }
697
698 if (projectedVars.find(_indexCat.arg().GetName())) {
699 projIndex=kTRUE ;
700 }
701
702 } else {
703 // *** Error checking for a composite index category ***
704
705 // Determine if any servers of the index category are in the projectedVars
707 RooAbsArg* server ;
708 RooArgSet projIdxServers ;
709 Bool_t anyServers(kFALSE) ;
710 while((server=(RooAbsArg*)sIter->Next())) {
711 if (projectedVars.find(server->GetName())) {
712 anyServers=kTRUE ;
713 projIdxServers.add(*server) ;
714 }
715 }
716 delete sIter ;
717
718 // Check that the projection dataset contains all the
719 // index category components we're projecting over
720
721 // Determine if all projected servers of the index category are in the projection dataset
722 sIter = projIdxServers.createIterator() ;
723 Bool_t allServers(kTRUE) ;
724 while((server=(RooAbsArg*)sIter->Next())) {
725 if (!projData->get()->find(server->GetName())) {
726 allServers=kFALSE ;
727 }
728 }
729 delete sIter ;
730
731 if (!allServers) {
732 coutE(Plotting) << "RooSimultaneous::plotOn(" << GetName()
733 << ") ERROR: Projection dataset doesn't contain complete set of index category dependents" << endl ;
734 return frame ;
735 }
736
737 if (anyServers) {
738 projIndex = kTRUE ;
739 }
740 }
741
742 // Calculate relative weight fractions of components
743 Roo1DTable* wTable = projData->table(_indexCat.arg()) ;
744
745 // If we don't project over the index, just do the regular plotOn
746 if (!projIndex) {
747
748 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
749 << " represents a slice in the index category (" << _indexCat.arg().GetName() << ")" << endl ;
750
751 // Reduce projData: take out fitCat (component) columns and entries that don't match selected slice
752 // Construct cut string to only select projection data event that match the current slice
753
754 const RooAbsData* projDataTmp(projData) ;
755 if (projData) {
756 // Make list of categories columns to exclude from projection data
757 RooArgSet* indexCatComps = _indexCat.arg().getObservables(frame->getNormVars());
758
759 // Make cut string to exclude rows from projection data
760 TString cutString ;
761 TIterator* compIter = indexCatComps->createIterator() ;
762 RooAbsCategory* idxComp ;
764 while((idxComp=(RooAbsCategory*)compIter->Next())) {
765 if (!first) {
766 cutString.Append("&&") ;
767 } else {
768 first=kFALSE ;
769 }
770 cutString.Append(Form("%s==%d",idxComp->GetName(),idxComp->getIndex())) ;
771 }
772 delete compIter ;
773
774 // Make temporary projData without RooSim index category components
775 RooArgSet projDataVars(*projData->get()) ;
776 projDataVars.remove(*indexCatComps,kTRUE,kTRUE) ;
777
778 projDataTmp = ((RooAbsData*)projData)->reduce(projDataVars,cutString) ;
779 delete indexCatComps ;
780 }
781
782 // Multiply scale factor with fraction of events in current state of index
783
784// RooPlot* retFrame = getPdf(_indexCat.arg().getLabel())->plotOn(frame,drawOptions,
785// scaleFactor*wTable->getFrac(_indexCat.arg().getLabel()),
786// stype,projDataTmp,projSet) ;
787
788 // Override normalization and projection dataset
789 RooCmdArg tmp1 = RooFit::Normalization(scaleFactor*wTable->getFrac(_indexCat.arg().getLabel()),stype) ;
790 RooCmdArg tmp2 = RooFit::ProjWData(*projDataSet,*projDataTmp) ;
791
792 // WVE -- do not adjust normalization for asymmetry plots
793 RooLinkedList cmdList2(cmdList) ;
794 if (!cmdList.find("Asymmetry")) {
795 cmdList2.Add(&tmp1) ;
796 }
797 cmdList2.Add(&tmp2) ;
798
799 // Plot single component
800 RooPlot* retFrame = getPdf(_indexCat.arg().getLabel())->plotOn(frame,cmdList2) ;
801
802 // Delete temporary dataset
803 if (projDataTmp) {
804 delete projDataTmp ;
805 }
806
807 delete wTable ;
808 delete sliceSet ;
809 return retFrame ;
810 }
811
812 // If we project over the index, plot using a temporary RooAddPdf
813 // using the weights from the data as coefficients
814
815 // Make a deep clone of our index category
816 RooArgSet* idxCloneSet = (RooArgSet*) RooArgSet(_indexCat.arg()).snapshot(kTRUE) ;
817 RooAbsCategoryLValue* idxCatClone = (RooAbsCategoryLValue*) idxCloneSet->find(_indexCat.arg().GetName()) ;
818
819 // Build the list of indexCat components that are sliced
820 RooArgSet* idxCompSliceSet = _indexCat.arg().getObservables(frame->getNormVars()) ;
821 idxCompSliceSet->remove(projectedVars,kTRUE,kTRUE) ;
822 TIterator* idxCompSliceIter = idxCompSliceSet->createIterator() ;
823
824 // Make a new expression that is the weighted sum of requested components
825 RooArgList pdfCompList ;
826 RooArgList wgtCompList ;
827//RooAbsPdf* pdf ;
828 RooRealProxy* proxy ;
830 Double_t sumWeight(0) ;
831 while((proxy=(RooRealProxy*)pIter->Next())) {
832
833 idxCatClone->setLabel(proxy->name()) ;
834
835 // Determine if this component is the current slice (if we slice)
836 Bool_t skip(kFALSE) ;
837 idxCompSliceIter->Reset() ;
838 RooAbsCategory* idxSliceComp ;
839 while((idxSliceComp=(RooAbsCategory*)idxCompSliceIter->Next())) {
840 RooAbsCategory* idxComp = (RooAbsCategory*) idxCloneSet->find(idxSliceComp->GetName()) ;
841 if (idxComp->getIndex()!=idxSliceComp->getIndex()) {
842 skip=kTRUE ;
843 break ;
844 }
845 }
846 if (skip) continue ;
847
848 // Instantiate a RRV holding this pdfs weight fraction
849 RooRealVar *wgtVar = new RooRealVar(proxy->name(),"coef",wTable->getFrac(proxy->name())) ;
850 wgtCompList.addOwned(*wgtVar) ;
851 sumWeight += wTable->getFrac(proxy->name()) ;
852
853 // Add the PDF to list list
854 pdfCompList.add(proxy->arg()) ;
855 }
856
857 TString plotVarName(GetName()) ;
858 RooAddPdf *plotVar = new RooAddPdf(plotVarName,"weighted sum of RS components",pdfCompList,wgtCompList) ;
859
860 // Fix appropriate coefficient normalization in plot function
861 if (_plotCoefNormSet.getSize()>0) {
863 }
864
865 RooAbsData* projDataTmp(0) ;
866 RooArgSet projSetTmp ;
867 if (projData) {
868
869 // Construct cut string to only select projection data event that match the current slice
870 TString cutString ;
871 if (idxCompSliceSet->getSize()>0) {
872 idxCompSliceIter->Reset() ;
873 RooAbsCategory* idxSliceComp ;
875 while((idxSliceComp=(RooAbsCategory*)idxCompSliceIter->Next())) {
876 if (!first) {
877 cutString.Append("&&") ;
878 } else {
879 first=kFALSE ;
880 }
881 cutString.Append(Form("%s==%d",idxSliceComp->GetName(),idxSliceComp->getIndex())) ;
882 }
883 }
884
885 // Make temporary projData without RooSim index category components
886 RooArgSet projDataVars(*projData->get()) ;
887 RooArgSet* idxCatServers = _indexCat.arg().getObservables(frame->getNormVars()) ;
888
889 projDataVars.remove(*idxCatServers,kTRUE,kTRUE) ;
890
891 if (idxCompSliceSet->getSize()>0) {
892 projDataTmp = ((RooAbsData*)projData)->reduce(projDataVars,cutString) ;
893 } else {
894 projDataTmp = ((RooAbsData*)projData)->reduce(projDataVars) ;
895 }
896
897
898
899 if (projSet) {
900 projSetTmp.add(*projSet) ;
901 projSetTmp.remove(*idxCatServers,kTRUE,kTRUE);
902 }
903
904
905 delete idxCatServers ;
906 }
907
908
909 if (_indexCat.arg().isDerived() && idxCompSliceSet->getSize()>0) {
910 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
911 << " represents a slice in index category components " << *idxCompSliceSet << endl ;
912
913 RooArgSet* idxCompProjSet = _indexCat.arg().getObservables(frame->getNormVars()) ;
914 idxCompProjSet->remove(*idxCompSliceSet,kTRUE,kTRUE) ;
915 if (idxCompProjSet->getSize()>0) {
916 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
917 << " averages with data index category components " << *idxCompProjSet << endl ;
918 }
919 delete idxCompProjSet ;
920 } else {
921 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
922 << " averages with data index category (" << _indexCat.arg().GetName() << ")" << endl ;
923 }
924
925
926 // Override normalization and projection dataset
927 RooLinkedList cmdList2(cmdList) ;
928
929 RooCmdArg tmp1 = RooFit::Normalization(scaleFactor*sumWeight,stype) ;
930 RooCmdArg tmp2 = RooFit::ProjWData(*projDataSet,*projDataTmp) ;
931 // WVE -- do not adjust normalization for asymmetry plots
932 if (!cmdList.find("Asymmetry")) {
933 cmdList2.Add(&tmp1) ;
934 }
935 cmdList2.Add(&tmp2) ;
936
937 RooPlot* frame2 ;
938 if (projSetTmp.getSize()>0) {
939 // Plot temporary function
940 RooCmdArg tmp3 = RooFit::Project(projSetTmp) ;
941 cmdList2.Add(&tmp3) ;
942 frame2 = plotVar->plotOn(frame,cmdList2) ;
943 } else {
944 // Plot temporary function
945 frame2 = plotVar->plotOn(frame,cmdList2) ;
946 }
947
948 // Cleanup
949 delete sliceSet ;
950 delete pIter ;
951 delete wTable ;
952 delete idxCloneSet ;
953 delete idxCompSliceIter ;
954 delete idxCompSliceSet ;
955 delete plotVar ;
956
957 if (projDataTmp) delete projDataTmp ;
958
959 return frame2 ;
960}
961
962
963
964////////////////////////////////////////////////////////////////////////////////
965/// OBSOLETE -- Retained for backward compatibility
966
967RooPlot* RooSimultaneous::plotOn(RooPlot *frame, Option_t* drawOptions, Double_t scaleFactor,
968 ScaleType stype, const RooAbsData* projData, const RooArgSet* projSet,
969 Double_t /*precision*/, Bool_t /*shiftToZero*/, const RooArgSet* /*projDataSet*/,
970 Double_t /*rangeLo*/, Double_t /*rangeHi*/, RooCurve::WingMode /*wmode*/) const
971{
972 // Make command list
973 RooLinkedList cmdList ;
974 cmdList.Add(new RooCmdArg(RooFit::DrawOption(drawOptions))) ;
975 cmdList.Add(new RooCmdArg(RooFit::Normalization(scaleFactor,stype))) ;
976 if (projData) cmdList.Add(new RooCmdArg(RooFit::ProjWData(*projData))) ;
977 if (projSet) cmdList.Add(new RooCmdArg(RooFit::Project(*projSet))) ;
978
979 // Call new method
980 RooPlot* ret = plotOn(frame,cmdList) ;
981
982 // Cleanup
983 cmdList.Delete() ;
984 return ret ;
985}
986
987
988
989////////////////////////////////////////////////////////////////////////////////
990/// Interface function used by test statistics to freeze choice of observables
991/// for interpretation of fraction coefficients. Needed here because a RooSimultaneous
992/// works like a RooAddPdf when plotted
993
995{
997 if (normSet) _plotCoefNormSet.add(*normSet) ;
998}
999
1000
1001////////////////////////////////////////////////////////////////////////////////
1002/// Interface function used by test statistics to freeze choice of range
1003/// for interpretation of fraction coefficients. Needed here because a RooSimultaneous
1004/// works like a RooAddPdf when plotted
1005
1006void RooSimultaneous::selectNormalizationRange(const char* normRange2, Bool_t /*force*/)
1007{
1008 _plotCoefNormRange = RooNameReg::ptr(normRange2) ;
1009}
1010
1011
1012
1013
1014////////////////////////////////////////////////////////////////////////////////
1015
1017 const RooArgSet* auxProto, Bool_t verbose, Bool_t autoBinned, const char* binnedTag) const
1018{
1019 const char* idxCatName = _indexCat.arg().GetName() ;
1020
1021 if (vars.find(idxCatName) && prototype==0 && (auxProto==0 || auxProto->getSize()==0) && (autoBinned || (binnedTag && strlen(binnedTag)))) {
1022
1023 // Return special generator config that can also do binned generation for selected states
1024 return new RooSimSplitGenContext(*this,vars,verbose,autoBinned,binnedTag) ;
1025
1026 } else {
1027
1028 // Return regular generator config ;
1029 return genContext(vars,prototype,auxProto,verbose) ;
1030 }
1031}
1032
1033
1034
1035////////////////////////////////////////////////////////////////////////////////
1036/// Return specialized generator contenxt for simultaneous p.d.f.s
1037
1039 const RooArgSet* auxProto, Bool_t verbose) const
1040{
1041 const char* idxCatName = _indexCat.arg().GetName() ;
1042 const RooArgSet* protoVars = prototype ? prototype->get() : 0 ;
1043
1044 if (vars.find(idxCatName) || (protoVars && protoVars->find(idxCatName))) {
1045
1046 // Generating index category: return special sim-context
1047 return new RooSimGenContext(*this,vars,prototype,auxProto,verbose) ;
1048
1049 } else if (_indexCat.arg().isDerived()) {
1050 // Generating dependents of a derived index category
1051
1052 // Determine if we none,any or all servers
1053 Bool_t anyServer(kFALSE), allServers(kTRUE) ;
1054 if (prototype) {
1055 TIterator* sIter = _indexCat.arg().serverIterator() ;
1056 RooAbsArg* server ;
1057 while((server=(RooAbsArg*)sIter->Next())) {
1058 if (prototype->get()->find(server->GetName())) {
1059 anyServer=kTRUE ;
1060 } else {
1061 allServers=kFALSE ;
1062 }
1063 }
1064 delete sIter ;
1065 } else {
1066 allServers=kTRUE ;
1067 }
1068
1069 if (allServers) {
1070 // Use simcontext if we have all servers
1071
1072 return new RooSimGenContext(*this,vars,prototype,auxProto,verbose) ;
1073 } else if (!allServers && anyServer) {
1074 // Abort if we have only part of the servers
1075 coutE(Plotting) << "RooSimultaneous::genContext: ERROR: prototype must include either all "
1076 << " components of the RooSimultaneous index category or none " << endl ;
1077 return 0 ;
1078 }
1079 // Otherwise make single gencontext for current state
1080 }
1081
1082 // Not generating index cat: return context for pdf associated with present index state
1084 if (!proxy) {
1085 coutE(InputArguments) << "RooSimultaneous::genContext(" << GetName()
1086 << ") ERROR: no PDF associated with current state ("
1087 << _indexCat.arg().GetName() << "=" << _indexCat.arg().getLabel() << ")" << endl ;
1088 return 0 ;
1089 }
1090 return ((RooAbsPdf*)proxy->absArg())->genContext(vars,prototype,auxProto,verbose) ;
1091}
1092
1093
1094
1095
1096////////////////////////////////////////////////////////////////////////////////
1097
1099 const RooArgSet* nset,
1100 Double_t scaleFactor,
1101 Bool_t correctForBinVolume,
1102 Bool_t showProgress) const
1103{
1104 if (RooAbsReal::fillDataHist (hist, nset, scaleFactor,
1105 correctForBinVolume, showProgress) == 0)
1106 return 0;
1107
1108 Double_t sum = 0;
1109 for (int i=0 ; i<hist->numEntries() ; i++) {
1110 hist->get(i) ;
1111 sum += hist->weight();
1112 }
1113 if (sum != 0) {
1114 for (int i=0 ; i<hist->numEntries() ; i++) {
1115 hist->get(i) ;
1116 hist->set (hist->weight() / sum);
1117 }
1118 }
1119
1120 return hist;
1121}
1122
1123
1124
1125
1126////////////////////////////////////////////////////////////////////////////////
1127/// Special generator interface for generation of 'global observables' -- for RooStats tools
1128
1130{
1131 // Make set with clone of variables (placeholder for output)
1132 RooArgSet* globClone = (RooArgSet*) whatVars.snapshot() ;
1133
1134 RooDataSet* data = new RooDataSet("gensimglobal","gensimglobal",whatVars) ;
1135
1136 // Construct iterator over index types
1137 TIterator* iter = indexCat().typeIterator() ;
1138
1139 for (Int_t i=0 ; i<nEvents ; i++) {
1140 iter->Reset() ;
1141 RooCatType* tt ;
1142 while((tt=(RooCatType*) iter->Next())) {
1143
1144 // Get pdf associated with state from simpdf
1145 RooAbsPdf* pdftmp = getPdf(tt->GetName()) ;
1146
1147 // Generate only global variables defined by the pdf associated with this state
1148 RooArgSet* globtmp = pdftmp->getObservables(whatVars) ;
1149 RooDataSet* tmp = pdftmp->generate(*globtmp,1) ;
1150
1151 // Transfer values to output placeholder
1152 *globClone = *tmp->get(0) ;
1153
1154 // Cleanup
1155 delete globtmp ;
1156 delete tmp ;
1157 }
1158 data->add(*globClone) ;
1159 }
1160
1161
1162 delete iter ;
1163 delete globClone ;
1164 return data ;
1165}
1166
1167
1168
1169
1170
1171
1172
1173
1174
#define coutI(a)
Definition: RooMsgService.h:31
#define cxcoutD(a)
Definition: RooMsgService.h:79
#define coutW(a)
Definition: RooMsgService.h:33
#define coutE(a)
Definition: RooMsgService.h:34
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
const char Option_t
Definition: RtypesCore.h:62
#define ClassImp(name)
Definition: Rtypes.h:365
char name[80]
Definition: TGX11.cxx:109
int type
Definition: TGX11.cxx:120
char * Form(const char *fmt,...)
Roo1DTable implements a one-dimensional table.
Definition: Roo1DTable.h:24
Double_t getFrac(const char *label, Bool_t silent=kFALSE) const
Return the fraction of entries in the table contained in the slot named 'label'.
Definition: Roo1DTable.cxx:301
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:70
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Return the observables of this pdf given a set of observables.
Definition: RooAbsArg.h:240
friend class RooDataSet
Definition: RooAbsArg.h:583
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:729
friend class RooArgSet
Definition: RooAbsArg.h:516
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:206
RooArgSet * getVariables(Bool_t stripDisconnected=kTRUE) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
Definition: RooAbsArg.cxx:2034
virtual Bool_t isDerived() const
Does value or shape of this arg depend on any other arg?
Definition: RooAbsArg.h:90
Bool_t addOwnedComponents(const RooArgSet &comps)
Take ownership of the contents of 'comps'.
Definition: RooAbsArg.cxx:2232
TIterator * serverIterator() const R__SUGGEST_ALTERNATIVE("Use servers() and begin()
RooAbsCategoryLValue is the common abstract base class for objects that represent a discrete value th...
virtual Bool_t setLabel(const char *label, Bool_t printError=kTRUE)=0
RooAbsCategory is the common abstract base class for objects that represent a discrete value with a f...
TIterator * typeIterator() const
Return iterator over all defined states.
virtual Int_t getIndex() const
Return index number of current state.
virtual const char * getLabel() const
Return label string of current state.
Int_t numTypes(const char *=0) const
Int_t getSize() const
Bool_t contains(const RooAbsArg &var) const
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.
TIterator * createIterator(Bool_t dir=kIterForward) const R__SUGGEST_ALTERNATIVE("begin()
TIterator-style iteration over contained elements.
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:37
virtual const RooArgSet * get() const
Definition: RooAbsData.h:80
virtual Roo1DTable * table(const RooArgSet &catSet, const char *cuts="", const char *opts="") const
Construct table for product of categories in catSet.
Definition: RooAbsData.cxx:754
RooAbsData * reduce(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())
Create a reduced copy of this dataset.
Definition: RooAbsData.cxx:381
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
friend class CacheElem
Definition: RooAbsPdf.h:328
Bool_t canBeExtended() const
Definition: RooAbsPdf.h:230
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:56
@ CanBeExtended
Definition: RooAbsPdf.h:223
@ MustBeExtended
Definition: RooAbsPdf.h:223
@ CanNotBeExtended
Definition: RooAbsPdf.h:223
Bool_t mustBeExtended() const
Definition: RooAbsPdf.h:234
RooArgSet * _normSet
Normalization integral (owned by _normMgr)
Definition: RooAbsPdf.h:316
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:119
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
friend class RooAddPdf
Definition: RooAbsReal.h:491
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:81
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:531
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.
friend class RooRealProxy
Definition: RooAbsReal.h:422
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:88
virtual const char * name() const
Definition: RooArgProxy.h:42
RooAbsArg * absArg() const
Definition: RooArgProxy.h:37
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:134
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
Int_t lastIndex() const
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=0)
RooCatType is an auxilary class for RooAbsCategory and defines a a single category state.
Definition: RooCatType.h:22
virtual const Text_t * GetName() const
Returns name of object.
Definition: RooCatType.h:44
virtual Bool_t setArg(RooAbsCategory &newRef)
Change object held in proxy into newRef.
const RooAbsCategory & arg() const
const char * label() const
RooCategory represents a fundamental (non-derived) discrete value object.
Definition: RooCategory.h:24
virtual Bool_t setLabel(const char *label, Bool_t printError=kTRUE)
Set value by specifying the name of the desired state If printError is set, a message will be printed...
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.
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:98
void set(Double_t weight, Double_t wgtErr=-1)
Increment the weight of the bin enclosing the coordinates given by 'row' by the specified amount.
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:31
virtual const RooArgSet * get(Int_t index) const
Return RooArgSet with coordinates of event 'index'.
virtual void add(const RooArgSet &row, Double_t weight=1.0, Double_t weightError=0)
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:36
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:63
static const TNamed * ptr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
Definition: RooNameReg.cxx:104
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:41
const RooArgSet * getNormVars() const
Definition: RooPlot.h:141
RooAbsRealLValue * getPlotVar() const
Definition: RooPlot.h:132
RooRealProxy is the concrete proxy for RooAbsReal objects A RooRealProxy is the general mechanism to ...
Definition: RooRealProxy.h:23
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
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 contenxt 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.
RooSuperCategory can join several RooAbsCategoryLValue objects into a single category.
virtual Bool_t setLabel(const char *label, Bool_t printError=kTRUE)
Set the value of the super category by specifying the state name by setting the state names of the co...
Iterator abstract base class.
Definition: TIterator.h:30
virtual void Reset()=0
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:575
virtual TIterator * MakeIterator(Bool_t dir=kIterForward) const
Return a list iterator.
Definition: TList.cxx:719
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition: TList.cxx:467
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)
@ InputArguments
Definition: RooGlobalFunc.h:58
RooCmdArg Normalization(Double_t scaleFactor)
static constexpr double pc
Definition: first.py:1
auto * tt
Definition: textangle.C:16
static long int sum(long int i)
Definition: Factory.cxx:2258