RooHistPdf.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitCore                                                       *
00004  * @(#)root/roofit:$Id: RooHistPdf.cxx 37219 2010-12-03 12:51:36Z 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 // RooHistPdf implements a probablity density function sampled from a 
00021 // multidimensional histogram. The histogram distribution is explicitly
00022 // normalized by RooHistPdf and can have an arbitrary number of real or 
00023 // discrete dimensions.
00024 // END_HTML
00025 //
00026 
00027 #include "RooFit.h"
00028 #include "Riostream.h"
00029 
00030 #include "RooHistPdf.h"
00031 #include "RooDataHist.h"
00032 #include "RooMsgService.h"
00033 #include "RooRealVar.h"
00034 #include "RooCategory.h"
00035 #include "RooWorkspace.h"
00036 
00037 
00038 
00039 ClassImp(RooHistPdf)
00040 ;
00041 
00042 
00043 
00044 //_____________________________________________________________________________
00045 RooHistPdf::RooHistPdf() : _dataHist(0), _totVolume(0), _unitNorm(kFALSE)
00046 {
00047   // Default constructor
00048   // coverity[UNINIT_CTOR]
00049   _histObsIter = _histObsList.createIterator() ;
00050   _pdfObsIter = _pdfObsList.createIterator() ;
00051 }
00052 
00053 
00054 //_____________________________________________________________________________
00055 RooHistPdf::RooHistPdf(const char *name, const char *title, const RooArgSet& vars, 
00056                        const RooDataHist& dhist, Int_t intOrder) :
00057   RooAbsPdf(name,title), 
00058   _pdfObsList("pdfObs","List of p.d.f. observables",this),
00059   _dataHist((RooDataHist*)&dhist), 
00060   _codeReg(10),
00061   _intOrder(intOrder),
00062   _cdfBoundaries(kFALSE),
00063   _totVolume(0),
00064   _unitNorm(kFALSE)
00065 {
00066   // Constructor from a RooDataHist. RooDataHist dimensions
00067   // can be either real or discrete. See RooDataHist::RooDataHist for details on the binning.
00068   // RooHistPdf neither owns or clone 'dhist' and the user must ensure the input histogram exists
00069   // for the entire life span of this PDF.
00070 
00071   _histObsList.addClone(vars) ;
00072   _pdfObsList.add(vars) ;
00073 
00074   // Verify that vars and dhist.get() have identical contents
00075   const RooArgSet* dvars = dhist.get() ;
00076   if (vars.getSize()!=dvars->getSize()) {
00077     coutE(InputArguments) << "RooHistPdf::ctor(" << GetName() 
00078                           << ") ERROR variable list and RooDataHist must contain the same variables." << endl ;
00079     assert(0) ;
00080   }
00081   TIterator* iter = vars.createIterator() ;
00082   RooAbsArg* arg ;
00083   while((arg=(RooAbsArg*)iter->Next())) {
00084     if (!dvars->find(arg->GetName())) {
00085       coutE(InputArguments) << "RooHistPdf::ctor(" << GetName() 
00086                             << ") ERROR variable list and RooDataHist must contain the same variables." << endl ;
00087       assert(0) ;
00088     }
00089   }
00090   delete iter ;
00091 
00092   _histObsIter = _histObsList.createIterator() ;
00093   _pdfObsIter = _pdfObsList.createIterator() ;
00094 }
00095 
00096 
00097 
00098 
00099 //_____________________________________________________________________________
00100 RooHistPdf::RooHistPdf(const char *name, const char *title, const RooArgList& pdfObs, 
00101                        const RooArgList& histObs, const RooDataHist& dhist, Int_t intOrder) :
00102   RooAbsPdf(name,title), 
00103   _pdfObsList("pdfObs","List of p.d.f. observables",this),
00104   _dataHist((RooDataHist*)&dhist), 
00105   _codeReg(10),
00106   _intOrder(intOrder),
00107   _cdfBoundaries(kFALSE),
00108   _totVolume(0),
00109   _unitNorm(kFALSE)
00110 {
00111   // Constructor from a RooDataHist. The first list of observables are the p.d.f.
00112   // observables, which may any RooAbsReal (function or variable). The second list
00113   // are the corresponding observables in the RooDataHist which must be of type
00114   // RooRealVar or RooCategory This constructor thus allows to apply a coordinate transformation
00115   // on the histogram data to be applied.
00116 
00117   _histObsList.addClone(histObs) ;
00118   _pdfObsList.add(pdfObs) ;
00119 
00120   // Verify that vars and dhist.get() have identical contents
00121   const RooArgSet* dvars = dhist.get() ;
00122   if (histObs.getSize()!=dvars->getSize()) {
00123     coutE(InputArguments) << "RooHistPdf::ctor(" << GetName() 
00124                           << ") ERROR histogram variable list and RooDataHist must contain the same variables." << endl ;
00125     throw(string("RooHistPdf::ctor() ERROR: histogram variable list and RooDataHist must contain the same variables")) ;
00126   }
00127   TIterator* iter = histObs.createIterator() ;
00128   RooAbsArg* arg ;
00129   while((arg=(RooAbsArg*)iter->Next())) {
00130     if (!dvars->find(arg->GetName())) {
00131       coutE(InputArguments) << "RooHistPdf::ctor(" << GetName() 
00132                             << ") ERROR variable list and RooDataHist must contain the same variables." << endl ;
00133       throw(string("RooHistPdf::ctor() ERROR: histogram variable list and RooDataHist must contain the same variables")) ;
00134     }
00135     if (!arg->isFundamental()) {
00136       coutE(InputArguments) << "RooHistPdf::ctor(" << GetName() 
00137                             << ") ERROR all elements of histogram observables set must be of type RooRealVar or RooCategory." << endl ;
00138       throw(string("RooHistPdf::ctor() ERROR all elements of histogram observables set must be of type RooRealVar or RooCategory.")) ;
00139     }
00140   }
00141   delete iter ;
00142 
00143   _histObsIter = _histObsList.createIterator() ;
00144   _pdfObsIter = _pdfObsList.createIterator() ;
00145 }
00146 
00147 
00148 
00149 //_____________________________________________________________________________
00150 RooHistPdf::RooHistPdf(const RooHistPdf& other, const char* name) :
00151   RooAbsPdf(other,name), 
00152   _pdfObsList("pdfObs",this,other._pdfObsList),
00153   _dataHist(other._dataHist),
00154   _codeReg(other._codeReg),
00155   _intOrder(other._intOrder),
00156   _cdfBoundaries(other._cdfBoundaries),
00157   _totVolume(other._totVolume),
00158   _unitNorm(other._unitNorm)
00159 {
00160   // Copy constructor
00161 
00162   _histObsList.addClone(other._histObsList) ;
00163 
00164   _histObsIter = _histObsList.createIterator() ;
00165   _pdfObsIter = _pdfObsList.createIterator() ;
00166 }
00167 
00168 
00169 
00170 
00171 //_____________________________________________________________________________
00172 RooHistPdf::~RooHistPdf()
00173 {
00174   // Destructor
00175 
00176   delete _histObsIter ;
00177   delete _pdfObsIter ;
00178 }
00179 
00180 
00181 
00182 
00183 
00184 //_____________________________________________________________________________
00185 Double_t RooHistPdf::evaluate() const
00186 {
00187   // Return the current value: The value of the bin enclosing the current coordinates
00188   // of the observables, normalized by the histograms contents. Interpolation
00189   // is applied if the RooHistPdf is configured to do that
00190   
00191   // Transfer values from   
00192   if (_pdfObsList.getSize()>0) {
00193     _histObsIter->Reset() ;
00194     _pdfObsIter->Reset() ;
00195     RooAbsArg* harg, *parg ;
00196     while((harg=(RooAbsArg*)_histObsIter->Next())) {
00197       parg = (RooAbsArg*)_pdfObsIter->Next() ;
00198       if (harg != parg) {
00199         parg->syncCache() ;
00200         harg->copyCache(parg,kTRUE) ;
00201       }
00202     }
00203   }
00204 
00205   Double_t ret =  _dataHist->weight(_histObsList,_intOrder,_unitNorm?kFALSE:kTRUE,_cdfBoundaries) ;  
00206   if (ret<0) {
00207     ret=0 ;
00208   }
00209   
00210 //   cout << "RooHistPdf::evaluate(" << GetName() << ") ret = " << ret << " at " << ((RooAbsReal*)_histObsList.first())->getVal() << endl ;
00211   return ret ;
00212 }
00213 
00214 
00215 //_____________________________________________________________________________
00216 Double_t RooHistPdf::totVolume() const
00217 {
00218   // Return the total volume spanned by the observables of the RooHistPdf
00219 
00220   // Return previously calculated value, if any
00221   if (_totVolume>0) {
00222     return _totVolume ;
00223   }
00224   _totVolume = 1. ;
00225   TIterator* iter = _histObsList.createIterator() ;
00226   RooAbsArg* arg ;
00227   while((arg=(RooAbsArg*)iter->Next())) {
00228     RooRealVar* real = dynamic_cast<RooRealVar*>(arg) ;
00229     if (real) {
00230       _totVolume *= (real->getMax()-real->getMin()) ;
00231     } else {
00232       RooCategory* cat = dynamic_cast<RooCategory*>(arg) ;
00233       if (cat) {
00234         _totVolume *= cat->numTypes() ;
00235       }
00236     }
00237   }
00238   delete iter ;
00239   return _totVolume ;
00240 }
00241 
00242 
00243 
00244 //_____________________________________________________________________________
00245 Int_t RooHistPdf::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const 
00246 {
00247   // Determine integration scenario. If no interpolation is used,
00248   // RooHistPdf can perform all integrals over its dependents
00249   // analytically via partial or complete summation of the input
00250   // histogram. If interpolation is used on the integral over
00251   // all histogram observables is supported
00252 
00253   // Only analytical integrals over the full range are defined
00254   if (rangeName!=0) {
00255     return 0 ;
00256   }
00257 
00258   // First make list of pdf observables to histogram observables
00259   RooArgList hobsl(_histObsList),pobsl(_pdfObsList) ;
00260   RooArgSet allVarsHist ;
00261   TIterator* iter = allVars.createIterator() ;
00262   RooAbsArg* pdfobs ;
00263   while((pdfobs=(RooAbsArg*)iter->Next())) {
00264     Int_t idx = pobsl.index(pdfobs) ;
00265     if (idx>=0) {
00266       RooAbsArg* hobs = hobsl.at(idx) ;
00267       if (hobs) {
00268         allVarsHist.add(*hobs) ;
00269       }
00270     }
00271   }
00272   delete iter ;
00273 
00274   // Simplest scenario, integrate over all dependents
00275   RooAbsCollection *allVarsCommon = allVarsHist.selectCommon(_histObsList) ;  
00276   Bool_t intAllObs = (allVarsCommon->getSize()==_histObsList.getSize()) ;
00277   delete allVarsCommon ;
00278   if (intAllObs) {
00279     analVars.add(allVars) ;
00280     return 1000 ;
00281   }
00282 
00283   // Disable partial analytical integrals if interpolation is used
00284 //   if (_intOrder>0) {
00285 //     return 0 ;
00286 //   }
00287 
00288   // Find subset of _histObsList that integration is requested over
00289   RooArgSet* allVarsSel = (RooArgSet*) allVarsHist.selectCommon(_histObsList) ;
00290   if (allVarsSel->getSize()==0) {
00291     delete allVarsSel ;
00292     return 0 ;
00293   }
00294 
00295   // Partial integration scenarios.
00296   // Build unique code from bit mask of integrated variables in depList
00297   Int_t code(0),n(0) ;
00298   iter = _histObsList.createIterator() ;
00299   RooAbsArg* arg ;
00300   while((arg=(RooAbsArg*)iter->Next())) {
00301     if (allVars.find(arg->GetName())) {
00302       code |= (1<<n) ;
00303       analVars.add(*pobsl.at(n)) ;
00304     }
00305     n++ ;
00306   }
00307   delete iter ;
00308 
00309   return code ;
00310 
00311 }
00312 
00313 
00314 
00315 //_____________________________________________________________________________
00316 Double_t RooHistPdf::analyticalIntegral(Int_t code, const char* /*rangeName*/) const 
00317 {
00318   // Return integral identified by 'code'. The actual integration
00319   // is deferred to RooDataHist::sum() which implements partial
00320   // or complete summation over the histograms contents
00321 
00322   // WVE needs adaptation for rangeName feature
00323   // Simplest scenario, integration over all dependents
00324   if (code==1000) {    
00325     return _dataHist->sum(kFALSE) ;
00326   }
00327 
00328   // Partial integration scenario, retrieve set of variables, calculate partial sum
00329   RooArgSet intSet ;
00330   TIterator* iter = _histObsList.createIterator() ;
00331   RooAbsArg* arg ;
00332   Int_t n(0) ;
00333   while((arg=(RooAbsArg*)iter->Next())) {
00334     if (code & (1<<n)) {
00335       intSet.add(*arg) ;
00336     }
00337     n++ ;
00338   }
00339   delete iter ;
00340 
00341   // WVE must sync hist slice list values to pdf slice list
00342   // Transfer values from   
00343   if (_pdfObsList.getSize()>0) {
00344     _histObsIter->Reset() ;
00345     _pdfObsIter->Reset() ;
00346     RooAbsArg* harg, *parg ;
00347     while((harg=(RooAbsArg*)_histObsIter->Next())) {
00348       parg = (RooAbsArg*)_pdfObsIter->Next() ;
00349       if (harg != parg) {
00350         parg->syncCache() ;
00351         harg->copyCache(parg,kTRUE) ;
00352       }
00353     }
00354   }  
00355 
00356 
00357   Double_t ret =  _dataHist->sum(intSet,_histObsList,kTRUE) ;
00358 //   cout << "RooHistPdf::ai(" << GetName() << ") code = " << code << " ret = " << ret << endl ;
00359 //   cout << "intSet = " << intSet << endl ;
00360 //   cout << "slice position = " << endl ;
00361 //   _histObsList.Print("v") ;
00362   return ret ;
00363 }
00364 
00365 
00366 
00367 //_____________________________________________________________________________
00368 list<Double_t>* RooHistPdf::plotSamplingHint(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
00369 {
00370   // Return sampling hint for making curves of (projections) of this function
00371   // as the recursive division strategy of RooCurve cannot deal efficiently
00372   // with the vertical lines that occur in a non-interpolated histogram
00373 
00374   // No hints are required when interpolation is used
00375   if (_intOrder>0) {
00376     return 0 ;
00377   }
00378 
00379   // Check that observable is in dataset, if not no hint is generated
00380   RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dataHist->get()->find(obs.GetName())) ;
00381   if (!lvarg) {
00382     return 0 ;
00383   }
00384 
00385   // Retrieve position of all bin boundaries
00386   const RooAbsBinning* binning = lvarg->getBinningPtr(0) ;
00387   Double_t* boundaries = binning->array() ;
00388 
00389   list<Double_t>* hint = new list<Double_t> ;
00390 
00391   // Widen range slighty
00392   xlo = xlo - 0.01*(xhi-xlo) ;
00393   xhi = xhi + 0.01*(xhi-xlo) ;
00394 
00395   Double_t delta = (xhi-xlo)*1e-8 ;
00396  
00397   // Construct array with pairs of points positioned epsilon to the left and
00398   // right of the bin boundaries
00399   for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
00400     if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
00401       hint->push_back(boundaries[i]-delta) ;
00402       hint->push_back(boundaries[i]+delta) ;
00403     }
00404   }
00405 
00406   return hint ;
00407 }
00408 
00409 
00410 
00411 //_____________________________________________________________________________
00412 Int_t RooHistPdf::getMaxVal(const RooArgSet& vars) const 
00413 {
00414   // Only handle case of maximum in all variables
00415   RooAbsCollection* common = _pdfObsList.selectCommon(vars) ;
00416   if (common->getSize()==_pdfObsList.getSize()) {
00417     delete common ;
00418     return 1;
00419   }
00420   delete common ;
00421   return 0 ;
00422 }
00423 
00424 
00425 //_____________________________________________________________________________
00426 Double_t RooHistPdf::maxVal(Int_t code) const 
00427 {
00428   assert(code==1) ;
00429 
00430   Double_t max(-1) ;
00431   for (Int_t i=0 ; i<_dataHist->numEntries() ; i++) {
00432     _dataHist->get(i) ;
00433     Double_t wgt = _dataHist->weight() ;
00434     if (wgt>max) max=wgt ;
00435   }
00436 
00437   return max*1.05 ;
00438 }
00439 
00440 
00441 //_____________________________________________________________________________
00442 Bool_t RooHistPdf::importWorkspaceHook(RooWorkspace& ws) 
00443 {  
00444   // Check if our datahist is already in the workspace
00445   std::list<RooAbsData*> allData = ws.allData() ;
00446   std::list<RooAbsData*>::const_iterator iter ;
00447   for (iter = allData.begin() ; iter != allData.end() ; ++iter) {
00448     // If your dataset is already in this workspace nothing needs to be done
00449     if (*iter == _dataHist) {
00450       return kFALSE ;
00451     }
00452   }
00453 
00454   // Check if dataset with given name already exists
00455   RooAbsData* wsdata = ws.data(_dataHist->GetName()) ;
00456   if (wsdata) {
00457     if (wsdata->InheritsFrom(RooDataHist::Class())) {
00458       // Exists and is of correct type -- adjust internal pointer
00459       _dataHist = (RooDataHist*) wsdata ;
00460       return kFALSE ;
00461     } else {
00462       // Exists and is NOT of correct type -- abort
00463       
00464     }
00465   }
00466 
00467   // We need to import our datahist into the workspace
00468   Bool_t flag = ws.import(*_dataHist) ;
00469   if (flag) {
00470     coutE(ObjectHandling) << "RooHistPdf::importWorkspaceHook(" << GetName() 
00471                           << ") error importing RooDataHist into workspace: dataset of different type with same name already exists." << endl ;
00472     return kTRUE ;
00473   }
00474 
00475   // Redirect our internal pointer to the copy in the workspace
00476   _dataHist = (RooDataHist*) ws.data(_dataHist->GetName()) ;
00477   return kFALSE ;
00478 }

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