// RooHistFunc implements a real-valued function sampled from a
// multidimensional histogram. The histogram can have an arbitrary number of real or
// discrete dimensions and may have negative values
// END_HTML
#include "RooFit.h"
#include "Riostream.h"
#include "RooHistFunc.h"
#include "RooDataHist.h"
#include "RooMsgService.h"
#include "RooRealVar.h"
#include "RooCategory.h"
#include "RooWorkspace.h"
#include "TError.h"
using namespace std;
ClassImp(RooHistFunc)
;
RooHistFunc::RooHistFunc() :
_dataHist(0),
_intOrder(0),
_cdfBoundaries(kFALSE),
_totVolume(0),
_unitNorm(kFALSE)
{
}
RooHistFunc::RooHistFunc(const char *name, const char *title, const RooArgSet& vars,
const RooDataHist& dhist, Int_t intOrder) :
RooAbsReal(name,title),
_depList("depList","List of dependents",this),
_dataHist((RooDataHist*)&dhist),
_codeReg(10),
_intOrder(intOrder),
_cdfBoundaries(kFALSE),
_totVolume(0),
_unitNorm(kFALSE)
{
_depList.add(vars) ;
const RooArgSet* dvars = dhist.get() ;
if (vars.getSize()!=dvars->getSize()) {
coutE(InputArguments) << "RooHistFunc::ctor(" << GetName()
<< ") ERROR variable list and RooDataHist must contain the same variables." << endl ;
assert(0) ;
}
TIterator* iter = vars.createIterator() ;
RooAbsArg* arg ;
while((arg=(RooAbsArg*)iter->Next())) {
if (!dvars->find(arg->GetName())) {
coutE(InputArguments) << "RooHistFunc::ctor(" << GetName()
<< ") ERROR variable list and RooDataHist must contain the same variables." << endl ;
assert(0) ;
}
}
delete iter ;
}
RooHistFunc::RooHistFunc(const RooHistFunc& other, const char* name) :
RooAbsReal(other,name),
_depList("depList",this,other._depList),
_dataHist(other._dataHist),
_codeReg(other._codeReg),
_intOrder(other._intOrder),
_cdfBoundaries(other._cdfBoundaries),
_totVolume(other._totVolume),
_unitNorm(other._unitNorm)
{
}
Double_t RooHistFunc::evaluate() const
{
Double_t ret = _dataHist->weight(_depList,_intOrder,kFALSE,_cdfBoundaries) ;
return ret ;
}
Int_t RooHistFunc::getMaxVal(const RooArgSet& vars) const
{
RooAbsCollection* common = _depList.selectCommon(vars) ;
if (common->getSize()==_depList.getSize()) {
delete common ;
return 1;
}
delete common ;
return 0 ;
}
Double_t RooHistFunc::maxVal(Int_t code) const
{
R__ASSERT(code==1) ;
Double_t max(-1) ;
for (Int_t i=0 ; i<_dataHist->numEntries() ; i++) {
_dataHist->get(i) ;
Double_t wgt = _dataHist->weight() ;
if (wgt>max) max=wgt ;
}
return max*1.05 ;
}
Double_t RooHistFunc::totVolume() const
{
if (_totVolume>0) {
return _totVolume ;
}
_totVolume = 1. ;
TIterator* iter = _depList.createIterator() ;
RooAbsArg* arg ;
while((arg=(RooAbsArg*)iter->Next())) {
RooRealVar* real = dynamic_cast<RooRealVar*>(arg) ;
if (real) {
_totVolume *= (real->getMax()-real->getMin()) ;
} else {
RooCategory* cat = dynamic_cast<RooCategory*>(arg) ;
if (cat) {
_totVolume *= cat->numTypes() ;
}
}
}
delete iter ;
return _totVolume ;
}
Int_t RooHistFunc::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
{
if (rangeName!=0) {
return 0 ;
}
RooAbsCollection *allVarsCommon = allVars.selectCommon(_depList) ;
Bool_t intAllObs = (allVarsCommon->getSize()==_depList.getSize()) ;
delete allVarsCommon ;
if (intAllObs && matchArgs(allVars,analVars,_depList)) {
return 1000 ;
}
if (_intOrder>0) {
return 0 ;
}
RooArgSet* allVarsSel = (RooArgSet*) allVars.selectCommon(_depList) ;
if (allVarsSel->getSize()==0) {
delete allVarsSel ;
return 0 ;
}
Int_t code(0),n(0) ;
TIterator* iter = _depList.createIterator() ;
RooAbsArg* arg ;
while((arg=(RooAbsArg*)iter->Next())) {
if (allVars.find(arg->GetName())) code |= (1<<n) ;
n++ ;
}
delete iter ;
analVars.add(*allVarsSel) ;
return code ;
}
Double_t RooHistFunc::analyticalIntegral(Int_t code, const char* ) const
{
if (code==1000) {
return _dataHist->sum(kTRUE) ;
}
RooArgSet intSet ;
TIterator* iter = _depList.createIterator() ;
RooAbsArg* arg ;
Int_t n(0) ;
while((arg=(RooAbsArg*)iter->Next())) {
if (code & (1<<n)) {
intSet.add(*arg) ;
}
n++ ;
}
delete iter ;
Double_t ret = _dataHist->sum(intSet,_depList,kTRUE) ;
return ret ;
}
list<Double_t>* RooHistFunc::plotSamplingHint(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
{
if (_intOrder>1) {
return 0 ;
}
RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dataHist->get()->find(obs.GetName())) ;
if (!lvarg) {
return 0 ;
}
const RooAbsBinning* binning = lvarg->getBinningPtr(0) ;
Double_t* boundaries = binning->array() ;
list<Double_t>* hint = new list<Double_t> ;
xlo = xlo - 0.01*(xhi-xlo) ;
xhi = xhi + 0.01*(xhi-xlo) ;
Double_t delta = (xhi-xlo)*1e-8 ;
for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
hint->push_back(boundaries[i]-delta) ;
hint->push_back(boundaries[i]+delta) ;
}
}
return hint ;
}
std::list<Double_t>* RooHistFunc::binBoundaries(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
{
if (_intOrder>1) {
return 0 ;
}
RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dataHist->get()->find(obs.GetName())) ;
if (!lvarg) {
return 0 ;
}
const RooAbsBinning* binning = lvarg->getBinningPtr(0) ;
Double_t* boundaries = binning->array() ;
list<Double_t>* hint = new list<Double_t> ;
for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
hint->push_back(boundaries[i]) ;
}
}
return hint ;
}
Bool_t RooHistFunc::importWorkspaceHook(RooWorkspace& ws)
{
std::list<RooAbsData*> allData = ws.allEmbeddedData() ;
std::list<RooAbsData*>::const_iterator iter ;
for (iter = allData.begin() ; iter != allData.end() ; ++iter) {
if (*iter == _dataHist) {
return kFALSE ;
}
}
RooAbsData* wsdata = ws.embeddedData(_dataHist->GetName()) ;
if (wsdata) {
if (wsdata->InheritsFrom(RooDataHist::Class())) {
if (areIdentical((RooDataHist&)*wsdata,*_dataHist)) {
_dataHist = (RooDataHist*) wsdata ;
} else {
TString uniqueName = Form("%s_%s",_dataHist->GetName(),GetName()) ;
Bool_t flag = ws.import(*_dataHist,RooFit::Rename(uniqueName.Data()),RooFit::Embedded()) ;
if (flag) {
coutE(ObjectHandling) << " RooHistPdf::importWorkspaceHook(" << GetName() << ") unable to import clone of underlying RooDataHist with unique name " << uniqueName << ", abort" << endl ;
return kTRUE ;
}
_dataHist = (RooDataHist*) ws.embeddedData(uniqueName.Data()) ;
}
} else {
TString uniqueName = Form("%s_%s",_dataHist->GetName(),GetName()) ;
Bool_t flag = ws.import(*_dataHist,RooFit::Rename(uniqueName.Data()),RooFit::Embedded()) ;
if (flag) {
coutE(ObjectHandling) << " RooHistPdf::importWorkspaceHook(" << GetName() << ") unable to import clone of underlying RooDataHist with unique name " << uniqueName << ", abort" << endl ;
return kTRUE ;
}
_dataHist = (RooDataHist*) ws.embeddedData(uniqueName.Data()) ;
}
return kFALSE ;
}
ws.import(*_dataHist,RooFit::Embedded()) ;
_dataHist = (RooDataHist*) ws.embeddedData(_dataHist->GetName()) ;
return kFALSE ;
}
Bool_t RooHistFunc::areIdentical(const RooDataHist& dh1, const RooDataHist& dh2)
{
if (fabs(dh1.sumEntries()-dh2.sumEntries())>1e-8) return kFALSE ;
if (dh1.numEntries() != dh2.numEntries()) return kFALSE ;
for (int i=0 ; i < dh1.numEntries() ; i++) {
dh1.get(i) ;
dh2.get(i) ;
if (fabs(dh1.weight()-dh2.weight())>1e-8) return kFALSE ;
}
return kTRUE ;
}
void RooHistFunc::Streamer(TBuffer &R__b)
{
if (R__b.IsReading()) {
R__b.ReadClassBuffer(RooHistFunc::Class(),this);
_proxyList.Clear() ;
registerProxy(_depList) ;
} else {
R__b.WriteClassBuffer(RooHistFunc::Class(),this);
}
}