RooSimultaneous.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitCore                                                       *
00004  * @(#)root/roofitcore:$Id: RooSimultaneous.cxx 36230 2010-10-09 20:21:02Z wouter $
00005  * Authors:                                                                  *
00006  *   WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu       *
00007  *   DK, David Kirkby,    UC Irvine,         dkirkby@uci.edu                 *
00008  *                                                                           *
00009  * Copyright (c) 2000-2005, Regents of the University of California          *
00010  *                          and Stanford University. All rights reserved.    *
00011  *                                                                           *
00012  * Redistribution and use in source and binary forms,                        *
00013  * with or without modification, are permitted according to the terms        *
00014  * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             *
00015  *****************************************************************************/
00016 
00017 //////////////////////////////////////////////////////////////////////////////
00018 //
00019 // BEGIN_HTML
00020 // RooSimultaneous facilitates simultaneous fitting of multiple PDFs
00021 // to subsets of a given dataset.
00022 // <p>
00023 // The class takes an index category, which is interpreted as
00024 // the data subset indicator, and a list of PDFs, each associated
00025 // with a state of the index category. RooSimultaneous always returns
00026 // the value of the PDF that is associated with the current value
00027 // of the index category
00028 // <p>
00029 // Extended likelihood fitting is supported if all components support
00030 // extended likelihood mode. The expected number of events by a RooSimultaneous
00031 // is that of the component p.d.f. selected by the index category
00032 // END_HTML
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   // Constructor with index category. PDFs associated with indexCat
00077   // states can be added after construction with the addPdf() function.
00078   // 
00079   // RooSimultaneous can function without having a PDF associated
00080   // with every single state. The normalization in such cases is taken
00081   // from the number of registered PDFs, but getVal() will assert if
00082   // when called for an unregistered index state.
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   // Constructor from index category and full list of PDFs. 
00098   // In this constructor form, a PDF must be supplied for each indexCat state
00099   // to avoid ambiguities. The PDFS are associated in order with the state of the
00100   // index category as listed by the index categories type iterator.
00101   //
00102   // PDFs may not overlap (i.e. share any variables) with the index category (function)
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   // Iterator over PDFs and index cat states and add each pair
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 // This class cannot be locally defined in initialize as it cannot be
00144 // used as a template argument in that case
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   // First see if there are any RooSimultaneous input components
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   // If there are no simultaneous component p.d.f. do simple processing through addPdf()
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   // Issue info message that we are about to do some rearraning
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   // Construct the 'superIndex' from the nominal index category and all auxiliary components
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   // Now process each of original pdf/state map entries
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       // Entry is a plain p.d.f. assign it to every state permutation of the repliCats set
00218       RooSuperCategory repliSuperCat("tmp","tmp",repliCats) ;
00219 
00220       // Iterator over all states of repliSuperCat
00221       TIterator* titer = repliSuperCat.typeIterator() ;
00222       RooCatType* type ;
00223       while ((type=(RooCatType*)titer->Next())) {
00224         // Set value 
00225         repliSuperCat.setLabel(type->GetName()) ;
00226         // Retrieve corresponding label of superIndex 
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       // Entry is a simultaneous p.d.f
00235 
00236       if (repliCats.getSize()==0) {
00237 
00238         // Case 1 -- No replication of components of RooSim component are required
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         // Case 2 -- Replication of components of RooSim component are required
00262 
00263         // Make replication supercat
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   // Change original master index to super index and take ownership of it
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   // Copy constructor
00314 
00315   // Copy proxy list 
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   // Destructor
00330 
00331   _pdfProxyList.Delete() ;
00332 }
00333 
00334 
00335 
00336 //_____________________________________________________________________________
00337 RooAbsPdf* RooSimultaneous::getPdf(const char* catName) const 
00338 {
00339   // Return the p.d.f associated with the given index category name
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   // Associate given PDF with index category state label 'catLabel'.
00351   // The names state must be already defined in the index category
00352   //
00353   // RooSimultaneous can function without having a PDF associated
00354   // with every single state. The normalization in such cases is taken
00355   // from the number of registered PDFs, but getVal() will assert if
00356   // when called for an unregistered index state.
00357   //
00358   // PDFs may not overlap (i.e. share any variables) with the index category (function)
00359 
00360   // PDFs cannot overlap with the index category
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   // Each index state can only have one PDF associated with it
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     // Create a proxy named after the associated index state
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   // WVE NEEDS FIX
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 //       cout << " now processing pdf " << pdf->GetName() << endl ;
00408       RooAbsPdf* pdf = (RooAbsPdf*) proxy->absArg() ;
00409       if (!pdf->canBeExtended()) {
00410 //      cout << "RooSim::extendedMode(" << GetName() << ") component " << pdf->GetName() << " cannot be extended" << endl ;
00411         allCanExtend=kFALSE ;
00412       }
00413       if (pdf->mustBeExtended()) {
00414         anyMustExtend=kTRUE;
00415       }
00416     }
00417   }
00418   if (anyMustExtend) {
00419 //     cout << "RooSim::extendedMode(" << GetName() << ") returning MustBeExtended" << endl ;
00420     return MustBeExtended ;
00421   }
00422   if (allCanExtend) {
00423 //     cout << "RooSim::extendedMode(" << GetName() << ") returning CanBeExtended" << endl ;
00424     return CanBeExtended ;
00425   }
00426 //   cout << "RooSim::extendedMode(" << GetName() << ") returning CanNotBeExtended" << endl ;
00427   return CanNotBeExtended ; 
00428 }
00429 
00430 
00431 
00432 
00433 //_____________________________________________________________________________
00434 Double_t RooSimultaneous::evaluate() const
00435 {  
00436   // Return the current value: 
00437   // the value of the PDF associated with the current index category state
00438 
00439   // Retrieve the proxy by index name
00440   RooRealProxy* proxy = (RooRealProxy*) _pdfProxyList.FindObject(_indexCat.label()) ;
00441   
00442   //assert(proxy!=0) ;
00443   if (proxy==0) return 0 ;
00444 
00445   // Calculate relative weighting factor for sim-pdfs of all extendable components
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   // Return the selected PDF value, normalized by the number of index states  
00461   return ((RooAbsPdf*)(proxy->absArg()))->getVal(_normSet)*catFrac ; 
00462 }
00463 
00464 
00465 
00466 //_____________________________________________________________________________
00467 Double_t RooSimultaneous::expectedEvents(const RooArgSet* nset) const 
00468 {
00469   // Return the number of expected events: If the index is in nset,
00470   // then return the sum of the expected events of all components,
00471   // otherwise return the number of expected events of the PDF
00472   // associated with the current index category state
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     // Retrieve the proxy by index name
00490     RooRealProxy* proxy = (RooRealProxy*) _pdfProxyList.FindObject(_indexCat.label()) ;
00491     
00492     //assert(proxy!=0) ;
00493     if (proxy==0) return 0 ;
00494 
00495     // Return the selected PDF value, normalized by the number of index states
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   // Forward determination of analytical integration capabilities to component p.d.f.s
00507   // A unique code is assigned to the combined integration capabilities of all associated
00508   // p.d.f.s
00509   
00510   // Declare that we can analytically integrate all requested observables
00511   analVars.add(allVars) ;
00512 
00513   // Retrieve (or create) the required partial integral list
00514   Int_t code ;
00515 
00516   // Check if this configuration was created before
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   // Create the partial integral set for this request
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   // Store the partial integral list and return the assigned code ;
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* /*rangeName*/) const 
00543 {
00544   // Return analytical integration defined by given code
00545 
00546   // No integration scenario
00547   if (code==0) {
00548     return getVal(normSet) ;
00549   }
00550 
00551   // Partial integration scenarios, rangeName already encoded in 'code'
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   // Back-end for plotOn() implementation on RooSimultaneous which
00568   // needs special handling because a RooSimultaneous PDF cannot
00569   // project out its index category via integration, plotOn() will
00570   // abort if this is requested without providing a projection dataset
00571 
00572   // Sanity checks
00573   if (plotSanityChecks(frame)) return frame ;
00574   
00575   // Extract projection configuration from command list
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() ; // there may be commands we don't handle here
00587 
00588   // Process and check varargs 
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   // Look for category slice arguments and add them to the master slice list if found
00604   const char* sliceCatState = pc.getString("sliceCatState",0,kTRUE) ;
00605   const RooLinkedList& sliceCatList = pc.getObjectList("sliceCatList") ;
00606   if (sliceCatState) {
00607 
00608     // Make the master slice set if it doesnt exist
00609     if (!sliceSet) {
00610       sliceSet = new RooArgSet ;
00611     }
00612 
00613     // Prepare comma separated label list for parsing
00614     char buf[1024] ;
00615     strlcpy(buf,sliceCatState,1024) ;
00616     const char* slabel = strtok(buf,",") ;
00617 
00618     // Loop over all categories provided by (multiple) Slice() arguments
00619     TIterator* iter = sliceCatList.MakeIterator() ;
00620     RooCategory* scat ;
00621     while((scat=(RooCategory*)iter->Next())) {
00622       if (slabel) {
00623         // Set the slice position to the value indicate by slabel
00624         scat->setLabel(slabel) ;
00625         // Add the slice category to the master slice set
00626         sliceSet->add(*scat,kFALSE) ;
00627       }
00628       slabel = strtok(0,",") ;
00629     }
00630     delete iter ;
00631   }
00632 
00633   // Check if we have a projection dataset
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   // Make list of variables to be projected
00640   RooArgSet projectedVars ;
00641   if (sliceSet) {
00642     //cout << "frame->getNormVars() = " ; frame->getNormVars()->Print("1") ;
00643 
00644     makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,kTRUE) ;
00645     
00646     // Take out the sliced variables
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     // *** Error checking for a fundamental index category ***
00669     //cout << "RooSim::plotOn: index is fundamental" << endl ;
00670       
00671     // Check that the provided projection dataset contains our index variable
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     // *** Error checking for a composite index category ***
00684 
00685     // Determine if any servers of the index category are in the projectedVars
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     // Check that the projection dataset contains all the 
00699     // index category components we're projecting over
00700 
00701     // Determine if all projected servers of the index category are in the projection dataset
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   // Calculate relative weight fractions of components
00723   Roo1DTable* wTable = projData->table(_indexCat.arg()) ;
00724 
00725   // If we don't project over the index, just do the regular plotOn
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     // Reduce projData: take out fitCat (component) columns and entries that don't match selected slice
00732     // Construct cut string to only select projection data event that match the current slice
00733 
00734     const RooAbsData* projDataTmp(projData) ;
00735     if (projData) {
00736       // Make list of categories columns to exclude from projection data      
00737       RooArgSet* indexCatComps = _indexCat.arg().getObservables(frame->getNormVars());
00738 
00739       // Make cut string to exclude rows from projection data
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       // Make temporary projData without RooSim index category components
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     // Multiply scale factor with fraction of events in current state of index
00763 
00764 //     RooPlot* retFrame =  getPdf(_indexCat.arg().getLabel())->plotOn(frame,drawOptions,
00765 //                                         scaleFactor*wTable->getFrac(_indexCat.arg().getLabel()),
00766 //                                         stype,projDataTmp,projSet) ;
00767 
00768     // Override normalization and projection dataset
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     // WVE -- do not adjust normalization for asymmetry plots
00774     if (!cmdList.find("Asymmetry")) {
00775       cmdList2.Add(&tmp1) ;
00776     }
00777     cmdList2.Add(&tmp2) ;
00778 
00779     // Plot single component
00780     RooPlot* retFrame =  getPdf(_indexCat.arg().getLabel())->plotOn(frame,cmdList2) ;
00781 
00782     // Delete temporary dataset
00783     if (projDataTmp) {
00784       delete projDataTmp ;
00785     }
00786 
00787     delete wTable ;
00788     delete sliceSet ;
00789     return retFrame ;
00790   }
00791 
00792   // If we project over the index, plot using a temporary RooAddPdf
00793   // using the weights from the data as coefficients
00794 
00795   // Make a deep clone of our index category
00796   RooArgSet* idxCloneSet = (RooArgSet*) RooArgSet(_indexCat.arg()).snapshot(kTRUE) ;
00797   RooAbsCategoryLValue* idxCatClone = (RooAbsCategoryLValue*) idxCloneSet->find(_indexCat.arg().GetName()) ;
00798 
00799   // Build the list of indexCat components that are sliced
00800   RooArgSet* idxCompSliceSet = _indexCat.arg().getObservables(frame->getNormVars()) ;
00801   idxCompSliceSet->remove(projectedVars,kTRUE,kTRUE) ;
00802   TIterator* idxCompSliceIter = idxCompSliceSet->createIterator() ;
00803 
00804   // Make a new expression that is the weighted sum of requested components
00805   RooArgList pdfCompList ;
00806   RooArgList wgtCompList ;
00807 //RooAbsPdf* pdf ;
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     // Determine if this component is the current slice (if we slice)
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     // Instantiate a RRV holding this pdfs weight fraction
00829     RooRealVar *wgtVar = new RooRealVar(proxy->name(),"coef",wTable->getFrac(proxy->name())) ;
00830     wgtCompList.addOwned(*wgtVar) ;
00831     sumWeight += wTable->getFrac(proxy->name()) ;
00832 
00833     // Add the PDF to list list
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   // Fix appropriate coefficient normalization in plot function
00841   if (_plotCoefNormSet.getSize()>0) {
00842     plotVar->fixAddCoefNormalization(_plotCoefNormSet) ;
00843   }
00844 
00845   RooAbsData* projDataTmp(0) ;
00846   RooArgSet projSetTmp ;
00847   if (projData) {
00848     
00849     // Construct cut string to only select projection data event that match the current slice
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     // Make temporary projData without RooSim index category components
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   // Override normalization and projection dataset
00907   RooLinkedList cmdList2(cmdList) ;
00908 
00909   RooCmdArg tmp1 = RooFit::Normalization(scaleFactor*sumWeight,stype) ;
00910   RooCmdArg tmp2 = RooFit::ProjWData(*projDataSet,*projDataTmp) ;
00911   // WVE -- do not adjust normalization for asymmetry plots
00912   if (!cmdList.find("Asymmetry")) {
00913     cmdList2.Add(&tmp1) ;
00914   }
00915   cmdList2.Add(&tmp2) ;
00916 
00917   RooPlot* frame2 ;
00918   if (projSetTmp.getSize()>0) {
00919     // Plot temporary function  
00920     RooCmdArg tmp3 = RooFit::Project(projSetTmp) ;
00921     cmdList2.Add(&tmp3) ;
00922     frame2 = plotVar->plotOn(frame,cmdList2) ;
00923   } else {
00924     // Plot temporary function  
00925     frame2 = plotVar->plotOn(frame,cmdList2) ;
00926   }
00927 
00928   // Cleanup
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 /*precision*/, Bool_t /*shiftToZero*/, const RooArgSet* /*projDataSet*/,
00948                                  Double_t /*rangeLo*/, Double_t /*rangeHi*/, RooCurve::WingMode /*wmode*/) const
00949 {
00950   // OBSOLETE -- Retained for backward compatibility
00951 
00952   // Make command list
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   // Call new method
00960   RooPlot* ret = plotOn(frame,cmdList) ;
00961 
00962   // Cleanup
00963   cmdList.Delete() ;
00964   return ret ;  
00965 }
00966 
00967 
00968 
00969 //_____________________________________________________________________________
00970 void RooSimultaneous::selectNormalization(const RooArgSet* normSet, Bool_t /*force*/) 
00971 {
00972   // Interface function used by test statistics to freeze choice of observables
00973   // for interpretation of fraction coefficients. Needed here because a RooSimultaneous
00974   // works like a RooAddPdf when plotted
00975   
00976   _plotCoefNormSet.removeAll() ;
00977   if (normSet) _plotCoefNormSet.add(*normSet) ;
00978 }
00979 
00980 
00981 //_____________________________________________________________________________
00982 void RooSimultaneous::selectNormalizationRange(const char* normRange2, Bool_t /*force*/) 
00983 {
00984   // Interface function used by test statistics to freeze choice of range
00985   // for interpretation of fraction coefficients. Needed here because a RooSimultaneous
00986   // works like a RooAddPdf when plotted
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   // Return specialized generator contenxt for simultaneous p.d.f.s
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     // Generating index category: return special sim-context
01005     return new RooSimGenContext(*this,vars,prototype,auxProto,verbose) ;
01006   } else if (_indexCat.arg().isDerived()) {
01007     // Generating dependents of a derived index category
01008 
01009     // Determine if we none,any or all servers
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       // Use simcontext if we have all servers
01028       return new RooSimGenContext(*this,vars,prototype,auxProto,verbose) ;
01029     } else if (!allServers && anyServer) {
01030       // Abort if we have only part of the servers
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     // Otherwise make single gencontext for current state
01036   } 
01037 
01038   // Not generating index cat: return context for pdf associated with present index state
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 

Generated on Tue Jul 5 15:07:32 2011 for ROOT_528-00b_version by  doxygen 1.5.1