00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035 #include "RooFit.h"
00036 #include "Riostream.h"
00037
00038 #include "TObjString.h"
00039 #include "RooSimultaneous.h"
00040 #include "RooAbsCategoryLValue.h"
00041 #include "RooPlot.h"
00042 #include "RooCurve.h"
00043 #include "RooRealVar.h"
00044 #include "RooAddPdf.h"
00045 #include "RooAbsData.h"
00046 #include "Roo1DTable.h"
00047 #include "RooSimGenContext.h"
00048 #include "RooDataSet.h"
00049 #include "RooCmdConfig.h"
00050 #include "RooNameReg.h"
00051 #include "RooGlobalFunc.h"
00052 #include "RooNameReg.h"
00053 #include "RooMsgService.h"
00054 #include "RooCategory.h"
00055 #include "RooSuperCategory.h"
00056 #include "RooArgSet.h"
00057
00058 using namespace std ;
00059
00060 ClassImp(RooSimultaneous)
00061 ;
00062
00063
00064
00065
00066
00067 RooSimultaneous::RooSimultaneous(const char *name, const char *title,
00068 RooAbsCategoryLValue& inIndexCat) :
00069 RooAbsPdf(name,title),
00070 _plotCoefNormSet("!plotCoefNormSet","plotCoefNormSet",this,kFALSE,kFALSE),
00071 _plotCoefNormRange(0),
00072 _partIntMgr(this,10),
00073 _indexCat("indexCat","Index category",this,inIndexCat),
00074 _numPdf(0)
00075 {
00076
00077
00078
00079
00080
00081
00082
00083 }
00084
00085
00086
00087
00088 RooSimultaneous::RooSimultaneous(const char *name, const char *title,
00089 const RooArgList& inPdfList, RooAbsCategoryLValue& inIndexCat) :
00090 RooAbsPdf(name,title),
00091 _plotCoefNormSet("!plotCoefNormSet","plotCoefNormSet",this,kFALSE,kFALSE),
00092 _plotCoefNormRange(0),
00093 _partIntMgr(this,10),
00094 _indexCat("indexCat","Index category",this,inIndexCat),
00095 _numPdf(0)
00096 {
00097
00098
00099
00100
00101
00102
00103
00104 if (inPdfList.getSize() != inIndexCat.numTypes()) {
00105 coutE(InputArguments) << "RooSimultaneous::ctor(" << GetName()
00106 << " ERROR: Number PDF list entries must match number of index category states, no PDFs added" << endl ;
00107 return ;
00108 }
00109
00110 map<string,RooAbsPdf*> pdfMap ;
00111
00112 TIterator* pIter = inPdfList.createIterator() ;
00113 TIterator* cIter = inIndexCat.typeIterator() ;
00114 RooAbsPdf* pdf ;
00115 RooCatType* type(0) ;
00116 while ((pdf=(RooAbsPdf*)pIter->Next())) {
00117 type = (RooCatType*) cIter->Next() ;
00118 pdfMap[string(type->GetName())] = pdf ;
00119 }
00120 delete pIter ;
00121 delete cIter ;
00122
00123 initialize(inIndexCat,pdfMap) ;
00124 }
00125
00126
00127
00128 RooSimultaneous::RooSimultaneous(const char *name, const char *title,
00129 map<string,RooAbsPdf*> pdfMap, RooAbsCategoryLValue& inIndexCat) :
00130 RooAbsPdf(name,title),
00131 _plotCoefNormSet("!plotCoefNormSet","plotCoefNormSet",this,kFALSE,kFALSE),
00132 _plotCoefNormRange(0),
00133 _partIntMgr(this,10),
00134 _indexCat("indexCat","Index category",this,inIndexCat),
00135 _numPdf(0)
00136 {
00137 initialize(inIndexCat,pdfMap) ;
00138 }
00139
00140
00141
00142
00143
00144
00145 namespace RooSimultaneousAux {
00146 struct CompInfo {
00147 RooAbsPdf* pdf ;
00148 RooSimultaneous* simPdf ;
00149 const RooAbsCategoryLValue* subIndex ;
00150 RooArgSet* subIndexComps ;
00151 } ;
00152 }
00153
00154 void RooSimultaneous::initialize(RooAbsCategoryLValue& inIndexCat, std::map<std::string,RooAbsPdf*> pdfMap)
00155 {
00156
00157 Bool_t simComps(kFALSE) ;
00158 for (map<string,RooAbsPdf*>::iterator iter=pdfMap.begin() ; iter!=pdfMap.end() ; iter++) {
00159 if (dynamic_cast<RooSimultaneous*>(iter->second)) {
00160 simComps = kTRUE ;
00161 break ;
00162 }
00163 }
00164
00165
00166 if (!simComps) {
00167 for (map<string,RooAbsPdf*>::iterator iter=pdfMap.begin() ; iter!=pdfMap.end() ; iter++) {
00168 addPdf(*iter->second,iter->first.c_str()) ;
00169 }
00170 return ;
00171 }
00172
00173
00174 coutI(InputArguments) << "RooSimultaneous::initialize(" << GetName() << ") INFO: one or more input component of simultaneous p.d.f.s are"
00175 << " simultaneous p.d.f.s themselves, rewriting composite expressions as one-level simultaneous p.d.f. in terms of"
00176 << " final constituents and extended index category" << endl ;
00177
00178
00179 RooArgSet allAuxCats ;
00180 map<string,RooSimultaneousAux::CompInfo> compMap ;
00181 for (map<string,RooAbsPdf*>::iterator iter=pdfMap.begin() ; iter!=pdfMap.end() ; iter++) {
00182 RooSimultaneousAux::CompInfo ci ;
00183 ci.pdf = iter->second ;
00184 RooSimultaneous* simComp = dynamic_cast<RooSimultaneous*>(iter->second) ;
00185 if (simComp) {
00186 ci.simPdf = simComp ;
00187 ci.subIndex = &simComp->indexCat() ;
00188 ci.subIndexComps = simComp->indexCat().isFundamental() ? new RooArgSet(simComp->indexCat()) : simComp->indexCat().getVariables() ;
00189 allAuxCats.add(*(ci.subIndexComps),kTRUE) ;
00190 } else {
00191 ci.simPdf = 0 ;
00192 ci.subIndex = 0 ;
00193 ci.subIndexComps = 0 ;
00194 }
00195 compMap[iter->first] = ci ;
00196 }
00197
00198
00199 RooArgSet allCats(inIndexCat) ;
00200 allCats.add(allAuxCats) ;
00201 string siname = Form("%s_index",GetName()) ;
00202 RooSuperCategory* superIndex = new RooSuperCategory(siname.c_str(),siname.c_str(),allCats) ;
00203
00204
00205 for (map<string,RooSimultaneousAux::CompInfo>::iterator citer = compMap.begin() ; citer != compMap.end() ; citer++) {
00206
00207 RooArgSet repliCats(allAuxCats) ;
00208 if (citer->second.subIndexComps) {
00209 repliCats.remove(*citer->second.subIndexComps) ;
00210 delete citer->second.subIndexComps ;
00211 }
00212 inIndexCat.setLabel(citer->first.c_str()) ;
00213
00214
00215 if (!citer->second.simPdf) {
00216
00217
00218 RooSuperCategory repliSuperCat("tmp","tmp",repliCats) ;
00219
00220
00221 TIterator* titer = repliSuperCat.typeIterator() ;
00222 RooCatType* type ;
00223 while ((type=(RooCatType*)titer->Next())) {
00224
00225 repliSuperCat.setLabel(type->GetName()) ;
00226
00227 string superLabel = superIndex->getLabel() ;
00228 addPdf(*citer->second.pdf,superLabel.c_str()) ;
00229 cxcoutD(InputArguments) << "RooSimultaneous::initialize(" << GetName()
00230 << ") assigning pdf " << citer->second.pdf->GetName() << " to super label " << superLabel << endl ;
00231 }
00232 } else {
00233
00234
00235
00236 if (repliCats.getSize()==0) {
00237
00238
00239
00240 TIterator* titer = citer->second.subIndex->typeIterator() ;
00241 RooCatType* type ;
00242 while ((type=(RooCatType*)titer->Next())) {
00243 const_cast<RooAbsCategoryLValue*>(citer->second.subIndex)->setLabel(type->GetName()) ;
00244 string superLabel = superIndex->getLabel() ;
00245 RooAbsPdf* compPdf = citer->second.simPdf->getPdf(type->GetName()) ;
00246 if (compPdf) {
00247 addPdf(*compPdf,superLabel.c_str()) ;
00248 cxcoutD(InputArguments) << "RooSimultaneous::initialize(" << GetName()
00249 << ") assigning pdf " << compPdf->GetName() << "(member of " << citer->second.pdf->GetName()
00250 << ") to super label " << superLabel << endl ;
00251 } else {
00252 coutW(InputArguments) << "RooSimultaneous::initialize(" << GetName() << ") WARNING: No p.d.f. associated with label "
00253 << type->GetName() << " for component RooSimultaneous p.d.f " << citer->second.pdf->GetName()
00254 << "which is associated with master index label " << citer->first << endl ;
00255 }
00256 }
00257 delete titer ;
00258
00259 } else {
00260
00261
00262
00263
00264 RooSuperCategory repliSuperCat("tmp","tmp",repliCats) ;
00265 TIterator* triter = repliSuperCat.typeIterator() ;
00266
00267 TIterator* tsiter = citer->second.subIndex->typeIterator() ;
00268 RooCatType* stype, *rtype ;
00269 while ((stype=(RooCatType*)tsiter->Next())) {
00270 const_cast<RooAbsCategoryLValue*>(citer->second.subIndex)->setLabel(stype->GetName()) ;
00271 triter->Reset() ;
00272 while ((rtype=(RooCatType*)triter->Next())) {
00273 repliSuperCat.setLabel(rtype->GetName()) ;
00274 string superLabel = superIndex->getLabel() ;
00275 RooAbsPdf* compPdf = citer->second.simPdf->getPdf(stype->GetName()) ;
00276 if (compPdf) {
00277 addPdf(*compPdf,superLabel.c_str()) ;
00278 cxcoutD(InputArguments) << "RooSimultaneous::initialize(" << GetName()
00279 << ") assigning pdf " << compPdf->GetName() << "(member of " << citer->second.pdf->GetName()
00280 << ") to super label " << superLabel << endl ;
00281 } else {
00282 coutW(InputArguments) << "RooSimultaneous::initialize(" << GetName() << ") WARNING: No p.d.f. associated with label "
00283 << stype->GetName() << " for component RooSimultaneous p.d.f " << citer->second.pdf->GetName()
00284 << "which is associated with master index label " << citer->first << endl ;
00285 }
00286 }
00287 }
00288
00289 delete tsiter ;
00290 delete triter ;
00291
00292 }
00293 }
00294 }
00295
00296
00297 _indexCat.setArg(*superIndex) ;
00298 addOwnedComponents(*superIndex) ;
00299
00300 }
00301
00302
00303
00304
00305 RooSimultaneous::RooSimultaneous(const RooSimultaneous& other, const char* name) :
00306 RooAbsPdf(other,name),
00307 _plotCoefNormSet("!plotCoefNormSet",this,other._plotCoefNormSet),
00308 _plotCoefNormRange(other._plotCoefNormRange),
00309 _partIntMgr(other._partIntMgr,this),
00310 _indexCat("indexCat",this,other._indexCat),
00311 _numPdf(other._numPdf)
00312 {
00313
00314
00315
00316 TIterator* pIter = other._pdfProxyList.MakeIterator() ;
00317 RooRealProxy* proxy ;
00318 while ((proxy=(RooRealProxy*)pIter->Next())) {
00319 _pdfProxyList.Add(new RooRealProxy(proxy->GetName(),this,*proxy)) ;
00320 }
00321 delete pIter ;
00322 }
00323
00324
00325
00326
00327 RooSimultaneous::~RooSimultaneous()
00328 {
00329
00330
00331 _pdfProxyList.Delete() ;
00332 }
00333
00334
00335
00336
00337 RooAbsPdf* RooSimultaneous::getPdf(const char* catName) const
00338 {
00339
00340
00341 RooRealProxy* proxy = (RooRealProxy*) _pdfProxyList.FindObject(catName) ;
00342 return proxy ? ((RooAbsPdf*)proxy->absArg()) : 0 ;
00343 }
00344
00345
00346
00347
00348 Bool_t RooSimultaneous::addPdf(const RooAbsPdf& pdf, const char* catLabel)
00349 {
00350
00351
00352
00353
00354
00355
00356
00357
00358
00359
00360
00361 if (pdf.dependsOn(_indexCat.arg())) {
00362 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName() << "): ERROR, PDF " << pdf.GetName()
00363 << " overlaps with index category " << _indexCat.arg().GetName() << endl ;
00364 return kTRUE ;
00365 }
00366
00367
00368 if (_pdfProxyList.FindObject(catLabel)) {
00369 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName() << "): ERROR, index state "
00370 << catLabel << " has already an associated PDF" << endl ;
00371 return kTRUE ;
00372 }
00373
00374 const RooSimultaneous* simPdf = dynamic_cast<const RooSimultaneous*>(&pdf) ;
00375 if (simPdf) {
00376
00377 coutE(InputArguments) << "RooSimultaneous::addPdf(" << GetName()
00378 << ") ERROR: you cannot add a RooSimultaneous component to a RooSimultaneous using addPdf()."
00379 << " Use the constructor with RooArgList if input p.d.f.s or the map<string,RooAbsPdf&> instead." << endl ;
00380 return kTRUE ;
00381
00382 } else {
00383
00384
00385 TObject* proxy = new RooRealProxy(catLabel,catLabel,this,(RooAbsPdf&)pdf) ;
00386 _pdfProxyList.Add(proxy) ;
00387 _numPdf += 1 ;
00388 }
00389
00390 return kFALSE ;
00391 }
00392
00393
00394
00395
00396
00397
00398 RooAbsPdf::ExtendMode RooSimultaneous::extendMode() const
00399 {
00400
00401 Bool_t allCanExtend(kTRUE) ;
00402 Bool_t anyMustExtend(kFALSE) ;
00403
00404 for (Int_t i=0 ; i<_numPdf ; i++) {
00405 RooRealProxy* proxy = (RooRealProxy*) _pdfProxyList.FindObject(_indexCat.label()) ;
00406 if (proxy) {
00407
00408 RooAbsPdf* pdf = (RooAbsPdf*) proxy->absArg() ;
00409 if (!pdf->canBeExtended()) {
00410
00411 allCanExtend=kFALSE ;
00412 }
00413 if (pdf->mustBeExtended()) {
00414 anyMustExtend=kTRUE;
00415 }
00416 }
00417 }
00418 if (anyMustExtend) {
00419
00420 return MustBeExtended ;
00421 }
00422 if (allCanExtend) {
00423
00424 return CanBeExtended ;
00425 }
00426
00427 return CanNotBeExtended ;
00428 }
00429
00430
00431
00432
00433
00434 Double_t RooSimultaneous::evaluate() const
00435 {
00436
00437
00438
00439
00440 RooRealProxy* proxy = (RooRealProxy*) _pdfProxyList.FindObject(_indexCat.label()) ;
00441
00442
00443 if (proxy==0) return 0 ;
00444
00445
00446 Double_t catFrac(1) ;
00447 if (canBeExtended()) {
00448 Double_t nEvtCat = ((RooAbsPdf*)(proxy->absArg()))->expectedEvents(_normSet) ;
00449
00450 Double_t nEvtTot(0) ;
00451 TIterator* iter = _pdfProxyList.MakeIterator() ;
00452 RooRealProxy* proxy2 ;
00453 while((proxy2=(RooRealProxy*)iter->Next())) {
00454 nEvtTot += ((RooAbsPdf*)(proxy2->absArg()))->expectedEvents(_normSet) ;
00455 }
00456 delete iter ;
00457 catFrac=nEvtCat/nEvtTot ;
00458 }
00459
00460
00461 return ((RooAbsPdf*)(proxy->absArg()))->getVal(_normSet)*catFrac ;
00462 }
00463
00464
00465
00466
00467 Double_t RooSimultaneous::expectedEvents(const RooArgSet* nset) const
00468 {
00469
00470
00471
00472
00473
00474 if (nset->contains(_indexCat.arg())) {
00475
00476 Double_t sum(0) ;
00477
00478 TIterator* iter = _pdfProxyList.MakeIterator() ;
00479 RooRealProxy* proxy ;
00480 while((proxy=(RooRealProxy*)iter->Next())) {
00481 sum += ((RooAbsPdf*)(proxy->absArg()))->expectedEvents(nset) ;
00482 }
00483 delete iter ;
00484
00485 return sum ;
00486
00487 } else {
00488
00489
00490 RooRealProxy* proxy = (RooRealProxy*) _pdfProxyList.FindObject(_indexCat.label()) ;
00491
00492
00493 if (proxy==0) return 0 ;
00494
00495
00496 return ((RooAbsPdf*)(proxy->absArg()))->expectedEvents(nset);
00497 }
00498 }
00499
00500
00501
00502
00503 Int_t RooSimultaneous::getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars,
00504 const RooArgSet* normSet, const char* rangeName) const
00505 {
00506
00507
00508
00509
00510
00511 analVars.add(allVars) ;
00512
00513
00514 Int_t code ;
00515
00516
00517 CacheElem* cache = (CacheElem*) _partIntMgr.getObj(normSet,&analVars,0,RooNameReg::ptr(rangeName)) ;
00518 if (cache) {
00519 code = _partIntMgr.lastIndex() ;
00520 return code+1 ;
00521 }
00522 cache = new CacheElem ;
00523
00524
00525 TIterator* iter = _pdfProxyList.MakeIterator() ;
00526 RooRealProxy* proxy ;
00527 while((proxy=(RooRealProxy*)iter->Next())) {
00528 RooAbsReal* pdfInt = proxy->arg().createIntegral(analVars,normSet,0,rangeName) ;
00529 cache->_partIntList.addOwned(*pdfInt) ;
00530 }
00531 delete iter ;
00532
00533
00534 code = _partIntMgr.setObj(normSet,&analVars,cache,RooNameReg::ptr(rangeName)) ;
00535
00536 return code+1 ;
00537 }
00538
00539
00540
00541
00542 Double_t RooSimultaneous::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* ) const
00543 {
00544
00545
00546
00547 if (code==0) {
00548 return getVal(normSet) ;
00549 }
00550
00551
00552 CacheElem* cache = (CacheElem*) _partIntMgr.getObjByIndex(code-1) ;
00553
00554 RooRealProxy* proxy = (RooRealProxy*) _pdfProxyList.FindObject(_indexCat.label()) ;
00555 Int_t idx = _pdfProxyList.IndexOf(proxy) ;
00556 return ((RooAbsReal*)cache->_partIntList.at(idx))->getVal(normSet) ;
00557 }
00558
00559
00560
00561
00562
00563
00564
00565 RooPlot* RooSimultaneous::plotOn(RooPlot *frame, RooLinkedList& cmdList) const
00566 {
00567
00568
00569
00570
00571
00572
00573 if (plotSanityChecks(frame)) return frame ;
00574
00575
00576 RooCmdConfig pc(Form("RooSimultaneous::plotOn(%s)",GetName())) ;
00577 pc.defineString("sliceCatState","SliceCat",0,"",kTRUE) ;
00578 pc.defineDouble("scaleFactor","Normalization",0,1.0) ;
00579 pc.defineInt("scaleType","Normalization",0,RooAbsPdf::Relative) ;
00580 pc.defineObject("sliceCatList","SliceCat",0,0,kTRUE) ;
00581 pc.defineObject("projSet","Project",0) ;
00582 pc.defineObject("sliceSet","SliceVars",0) ;
00583 pc.defineObject("projDataSet","ProjData",0) ;
00584 pc.defineObject("projData","ProjData",1) ;
00585 pc.defineMutex("Project","SliceVars") ;
00586 pc.allowUndefined() ;
00587
00588
00589 pc.process(cmdList) ;
00590 if (!pc.ok(kTRUE)) {
00591 return frame ;
00592 }
00593
00594 const RooAbsData* projData = (const RooAbsData*) pc.getObject("projData") ;
00595 const RooArgSet* projDataSet = (const RooArgSet*) pc.getObject("projDataSet") ;
00596 const RooArgSet* sliceSetTmp = (const RooArgSet*) pc.getObject("sliceSet") ;
00597 RooArgSet* sliceSet = sliceSetTmp ? ((RooArgSet*) sliceSetTmp->Clone()) : 0 ;
00598 const RooArgSet* projSet = (const RooArgSet*) pc.getObject("projSet") ;
00599 Double_t scaleFactor = pc.getDouble("scaleFactor") ;
00600 ScaleType stype = (ScaleType) pc.getInt("scaleType") ;
00601
00602
00603
00604 const char* sliceCatState = pc.getString("sliceCatState",0,kTRUE) ;
00605 const RooLinkedList& sliceCatList = pc.getObjectList("sliceCatList") ;
00606 if (sliceCatState) {
00607
00608
00609 if (!sliceSet) {
00610 sliceSet = new RooArgSet ;
00611 }
00612
00613
00614 char buf[1024] ;
00615 strlcpy(buf,sliceCatState,1024) ;
00616 const char* slabel = strtok(buf,",") ;
00617
00618
00619 TIterator* iter = sliceCatList.MakeIterator() ;
00620 RooCategory* scat ;
00621 while((scat=(RooCategory*)iter->Next())) {
00622 if (slabel) {
00623
00624 scat->setLabel(slabel) ;
00625
00626 sliceSet->add(*scat,kFALSE) ;
00627 }
00628 slabel = strtok(0,",") ;
00629 }
00630 delete iter ;
00631 }
00632
00633
00634 if (!projData) {
00635 coutE(InputArguments) << "RooSimultaneous::plotOn(" << GetName() << ") ERROR: must have a projection dataset for index category" << endl ;
00636 return frame ;
00637 }
00638
00639
00640 RooArgSet projectedVars ;
00641 if (sliceSet) {
00642
00643
00644 makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,kTRUE) ;
00645
00646
00647 TIterator* iter = sliceSet->createIterator() ;
00648 RooAbsArg* sliceArg ;
00649 while((sliceArg=(RooAbsArg*)iter->Next())) {
00650 RooAbsArg* arg = projectedVars.find(sliceArg->GetName()) ;
00651 if (arg) {
00652 projectedVars.remove(*arg) ;
00653 } else {
00654 coutI(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") slice variable "
00655 << sliceArg->GetName() << " was not projected anyway" << endl ;
00656 }
00657 }
00658 delete iter ;
00659 } else if (projSet) {
00660 makeProjectionSet(frame->getPlotVar(),projSet,projectedVars,kFALSE) ;
00661 } else {
00662 makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,kTRUE) ;
00663 }
00664
00665 Bool_t projIndex(kFALSE) ;
00666
00667 if (!_indexCat.arg().isDerived()) {
00668
00669
00670
00671
00672 if (!projData->get()->find(_indexCat.arg().GetName())) {
00673 coutE(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") ERROR: Projection over index category "
00674 << "requested, but projection data set doesn't contain index category" << endl ;
00675 return frame ;
00676 }
00677
00678 if (projectedVars.find(_indexCat.arg().GetName())) {
00679 projIndex=kTRUE ;
00680 }
00681
00682 } else {
00683
00684
00685
00686 TIterator* sIter = _indexCat.arg().serverIterator() ;
00687 RooAbsArg* server ;
00688 RooArgSet projIdxServers ;
00689 Bool_t anyServers(kFALSE) ;
00690 while((server=(RooAbsArg*)sIter->Next())) {
00691 if (projectedVars.find(server->GetName())) {
00692 anyServers=kTRUE ;
00693 projIdxServers.add(*server) ;
00694 }
00695 }
00696 delete sIter ;
00697
00698
00699
00700
00701
00702 sIter = projIdxServers.createIterator() ;
00703 Bool_t allServers(kTRUE) ;
00704 while((server=(RooAbsArg*)sIter->Next())) {
00705 if (!projData->get()->find(server->GetName())) {
00706 allServers=kFALSE ;
00707 }
00708 }
00709 delete sIter ;
00710
00711 if (!allServers) {
00712 coutE(Plotting) << "RooSimultaneous::plotOn(" << GetName()
00713 << ") ERROR: Projection dataset doesn't contain complete set of index category dependents" << endl ;
00714 return frame ;
00715 }
00716
00717 if (anyServers) {
00718 projIndex = kTRUE ;
00719 }
00720 }
00721
00722
00723 Roo1DTable* wTable = projData->table(_indexCat.arg()) ;
00724
00725
00726 if (!projIndex) {
00727
00728 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
00729 << " represents a slice in the index category (" << _indexCat.arg().GetName() << ")" << endl ;
00730
00731
00732
00733
00734 const RooAbsData* projDataTmp(projData) ;
00735 if (projData) {
00736
00737 RooArgSet* indexCatComps = _indexCat.arg().getObservables(frame->getNormVars());
00738
00739
00740 TString cutString ;
00741 TIterator* compIter = indexCatComps->createIterator() ;
00742 RooAbsCategory* idxComp ;
00743 Bool_t first(kTRUE) ;
00744 while((idxComp=(RooAbsCategory*)compIter->Next())) {
00745 if (!first) {
00746 cutString.Append("&&") ;
00747 } else {
00748 first=kFALSE ;
00749 }
00750 cutString.Append(Form("%s==%d",idxComp->GetName(),idxComp->getIndex())) ;
00751 }
00752 delete compIter ;
00753
00754
00755 RooArgSet projDataVars(*projData->get()) ;
00756 projDataVars.remove(*indexCatComps,kTRUE,kTRUE) ;
00757
00758 projDataTmp = ((RooAbsData*)projData)->reduce(projDataVars,cutString) ;
00759 delete indexCatComps ;
00760 }
00761
00762
00763
00764
00765
00766
00767
00768
00769 RooLinkedList cmdList2(cmdList) ;
00770 RooCmdArg tmp1 = RooFit::Normalization(scaleFactor*wTable->getFrac(_indexCat.arg().getLabel()),stype) ;
00771 RooCmdArg tmp2 = RooFit::ProjWData(*projDataSet,*projDataTmp) ;
00772
00773
00774 if (!cmdList.find("Asymmetry")) {
00775 cmdList2.Add(&tmp1) ;
00776 }
00777 cmdList2.Add(&tmp2) ;
00778
00779
00780 RooPlot* retFrame = getPdf(_indexCat.arg().getLabel())->plotOn(frame,cmdList2) ;
00781
00782
00783 if (projDataTmp) {
00784 delete projDataTmp ;
00785 }
00786
00787 delete wTable ;
00788 delete sliceSet ;
00789 return retFrame ;
00790 }
00791
00792
00793
00794
00795
00796 RooArgSet* idxCloneSet = (RooArgSet*) RooArgSet(_indexCat.arg()).snapshot(kTRUE) ;
00797 RooAbsCategoryLValue* idxCatClone = (RooAbsCategoryLValue*) idxCloneSet->find(_indexCat.arg().GetName()) ;
00798
00799
00800 RooArgSet* idxCompSliceSet = _indexCat.arg().getObservables(frame->getNormVars()) ;
00801 idxCompSliceSet->remove(projectedVars,kTRUE,kTRUE) ;
00802 TIterator* idxCompSliceIter = idxCompSliceSet->createIterator() ;
00803
00804
00805 RooArgList pdfCompList ;
00806 RooArgList wgtCompList ;
00807
00808 RooRealProxy* proxy ;
00809 TIterator* pIter = _pdfProxyList.MakeIterator() ;
00810 Double_t sumWeight(0) ;
00811 while((proxy=(RooRealProxy*)pIter->Next())) {
00812
00813 idxCatClone->setLabel(proxy->name()) ;
00814
00815
00816 Bool_t skip(kFALSE) ;
00817 idxCompSliceIter->Reset() ;
00818 RooAbsCategory* idxSliceComp ;
00819 while((idxSliceComp=(RooAbsCategory*)idxCompSliceIter->Next())) {
00820 RooAbsCategory* idxComp = (RooAbsCategory*) idxCloneSet->find(idxSliceComp->GetName()) ;
00821 if (idxComp->getIndex()!=idxSliceComp->getIndex()) {
00822 skip=kTRUE ;
00823 break ;
00824 }
00825 }
00826 if (skip) continue ;
00827
00828
00829 RooRealVar *wgtVar = new RooRealVar(proxy->name(),"coef",wTable->getFrac(proxy->name())) ;
00830 wgtCompList.addOwned(*wgtVar) ;
00831 sumWeight += wTable->getFrac(proxy->name()) ;
00832
00833
00834 pdfCompList.add(proxy->arg()) ;
00835 }
00836
00837 TString plotVarName(GetName()) ;
00838 RooAddPdf *plotVar = new RooAddPdf(plotVarName,"weighted sum of RS components",pdfCompList,wgtCompList) ;
00839
00840
00841 if (_plotCoefNormSet.getSize()>0) {
00842 plotVar->fixAddCoefNormalization(_plotCoefNormSet) ;
00843 }
00844
00845 RooAbsData* projDataTmp(0) ;
00846 RooArgSet projSetTmp ;
00847 if (projData) {
00848
00849
00850 TString cutString ;
00851 if (idxCompSliceSet->getSize()>0) {
00852 idxCompSliceIter->Reset() ;
00853 RooAbsCategory* idxSliceComp ;
00854 Bool_t first(kTRUE) ;
00855 while((idxSliceComp=(RooAbsCategory*)idxCompSliceIter->Next())) {
00856 if (!first) {
00857 cutString.Append("&&") ;
00858 } else {
00859 first=kFALSE ;
00860 }
00861 cutString.Append(Form("%s==%d",idxSliceComp->GetName(),idxSliceComp->getIndex())) ;
00862 }
00863 }
00864
00865
00866 RooArgSet projDataVars(*projData->get()) ;
00867 RooArgSet* idxCatServers = _indexCat.arg().getObservables(frame->getNormVars()) ;
00868
00869 projDataVars.remove(*idxCatServers,kTRUE,kTRUE) ;
00870
00871 if (idxCompSliceSet->getSize()>0) {
00872 projDataTmp = ((RooAbsData*)projData)->reduce(projDataVars,cutString) ;
00873 } else {
00874 projDataTmp = ((RooAbsData*)projData)->reduce(projDataVars) ;
00875 }
00876
00877
00878
00879 if (projSet) {
00880 projSetTmp.add(*projSet) ;
00881 projSetTmp.remove(*idxCatServers,kTRUE,kTRUE);
00882 }
00883
00884
00885 delete idxCatServers ;
00886 }
00887
00888
00889 if (_indexCat.arg().isDerived() && idxCompSliceSet->getSize()>0) {
00890 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
00891 << " represents a slice in index category components " << *idxCompSliceSet << endl ;
00892
00893 RooArgSet* idxCompProjSet = _indexCat.arg().getObservables(frame->getNormVars()) ;
00894 idxCompProjSet->remove(*idxCompSliceSet,kTRUE,kTRUE) ;
00895 if (idxCompProjSet->getSize()>0) {
00896 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
00897 << " averages with data index category components " << *idxCompProjSet << endl ;
00898 }
00899 delete idxCompProjSet ;
00900 } else {
00901 coutI(Plotting) << "RooSimultaneous::plotOn(" << GetName() << ") plot on " << frame->getPlotVar()->GetName()
00902 << " averages with data index category (" << _indexCat.arg().GetName() << ")" << endl ;
00903 }
00904
00905
00906
00907 RooLinkedList cmdList2(cmdList) ;
00908
00909 RooCmdArg tmp1 = RooFit::Normalization(scaleFactor*sumWeight,stype) ;
00910 RooCmdArg tmp2 = RooFit::ProjWData(*projDataSet,*projDataTmp) ;
00911
00912 if (!cmdList.find("Asymmetry")) {
00913 cmdList2.Add(&tmp1) ;
00914 }
00915 cmdList2.Add(&tmp2) ;
00916
00917 RooPlot* frame2 ;
00918 if (projSetTmp.getSize()>0) {
00919
00920 RooCmdArg tmp3 = RooFit::Project(projSetTmp) ;
00921 cmdList2.Add(&tmp3) ;
00922 frame2 = plotVar->plotOn(frame,cmdList2) ;
00923 } else {
00924
00925 frame2 = plotVar->plotOn(frame,cmdList2) ;
00926 }
00927
00928
00929 delete sliceSet ;
00930 delete pIter ;
00931 delete wTable ;
00932 delete idxCloneSet ;
00933 delete idxCompSliceIter ;
00934 delete idxCompSliceSet ;
00935 delete plotVar ;
00936
00937 if (projDataTmp) delete projDataTmp ;
00938
00939 return frame2 ;
00940 }
00941
00942
00943
00944
00945 RooPlot* RooSimultaneous::plotOn(RooPlot *frame, Option_t* drawOptions, Double_t scaleFactor,
00946 ScaleType stype, const RooAbsData* projData, const RooArgSet* projSet,
00947 Double_t , Bool_t , const RooArgSet* ,
00948 Double_t , Double_t , RooCurve::WingMode ) const
00949 {
00950
00951
00952
00953 RooLinkedList cmdList ;
00954 cmdList.Add(new RooCmdArg(RooFit::DrawOption(drawOptions))) ;
00955 cmdList.Add(new RooCmdArg(RooFit::Normalization(scaleFactor,stype))) ;
00956 if (projData) cmdList.Add(new RooCmdArg(RooFit::ProjWData(*projData))) ;
00957 if (projSet) cmdList.Add(new RooCmdArg(RooFit::Project(*projSet))) ;
00958
00959
00960 RooPlot* ret = plotOn(frame,cmdList) ;
00961
00962
00963 cmdList.Delete() ;
00964 return ret ;
00965 }
00966
00967
00968
00969
00970 void RooSimultaneous::selectNormalization(const RooArgSet* normSet, Bool_t )
00971 {
00972
00973
00974
00975
00976 _plotCoefNormSet.removeAll() ;
00977 if (normSet) _plotCoefNormSet.add(*normSet) ;
00978 }
00979
00980
00981
00982 void RooSimultaneous::selectNormalizationRange(const char* normRange2, Bool_t )
00983 {
00984
00985
00986
00987
00988 _plotCoefNormRange = RooNameReg::ptr(normRange2) ;
00989 }
00990
00991
00992
00993
00994
00995 RooAbsGenContext* RooSimultaneous::genContext(const RooArgSet &vars, const RooDataSet *prototype,
00996 const RooArgSet* auxProto, Bool_t verbose) const
00997 {
00998
00999
01000 const char* idxCatName = _indexCat.arg().GetName() ;
01001 const RooArgSet* protoVars = prototype ? prototype->get() : 0 ;
01002
01003 if (vars.find(idxCatName) || (protoVars && protoVars->find(idxCatName))) {
01004
01005 return new RooSimGenContext(*this,vars,prototype,auxProto,verbose) ;
01006 } else if (_indexCat.arg().isDerived()) {
01007
01008
01009
01010 Bool_t anyServer(kFALSE), allServers(kTRUE) ;
01011 if (prototype) {
01012 TIterator* sIter = _indexCat.arg().serverIterator() ;
01013 RooAbsArg* server ;
01014 while((server=(RooAbsArg*)sIter->Next())) {
01015 if (prototype->get()->find(server->GetName())) {
01016 anyServer=kTRUE ;
01017 } else {
01018 allServers=kFALSE ;
01019 }
01020 }
01021 delete sIter ;
01022 } else {
01023 allServers=kTRUE ;
01024 }
01025
01026 if (allServers) {
01027
01028 return new RooSimGenContext(*this,vars,prototype,auxProto,verbose) ;
01029 } else if (!allServers && anyServer) {
01030
01031 coutE(Plotting) << "RooSimultaneous::genContext: ERROR: prototype must include either all "
01032 << " components of the RooSimultaneous index category or none " << endl ;
01033 return 0 ;
01034 }
01035
01036 }
01037
01038
01039 RooRealProxy* proxy = (RooRealProxy*) _pdfProxyList.FindObject(_indexCat.arg().getLabel()) ;
01040 if (!proxy) {
01041 coutE(InputArguments) << "RooSimultaneous::genContext(" << GetName()
01042 << ") ERROR: no PDF associated with current state ("
01043 << _indexCat.arg().GetName() << "=" << _indexCat.arg().getLabel() << ")" << endl ;
01044 return 0 ;
01045 }
01046 return ((RooAbsPdf*)proxy->absArg())->genContext(vars,prototype,auxProto,verbose) ;
01047 }
01048