#include "RooFit.h"
#include "Riostream.h"
#include <string.h>
#include "RooAbsOptTestStatistic.h"
#include "RooMsgService.h"
#include "RooAbsPdf.h"
#include "RooAbsData.h"
#include "RooArgSet.h"
#include "RooRealVar.h"
#include "RooErrorHandler.h"
#include "RooGlobalFunc.h"
ClassImp(RooAbsOptTestStatistic)
;
RooAbsOptTestStatistic:: RooAbsOptTestStatistic()
{
_normSet = 0 ;
_pdfCloneSet = 0 ;
_dataClone = 0 ;
_pdfClone = 0 ;
_projDeps = 0 ;
}
RooAbsOptTestStatistic::RooAbsOptTestStatistic(const char *name, const char *title, RooAbsPdf& pdf, RooAbsData& data,
const RooArgSet& projDeps, const char* rangeName, const char* addCoefRangeName,
Int_t nCPU, Bool_t verbose, Bool_t splitCutRange) :
RooAbsTestStatistic(name,title,pdf,data,projDeps,rangeName, addCoefRangeName, nCPU, verbose, splitCutRange),
_projDeps(0)
{
if (operMode()!=Slave) {
_normSet = 0 ;
return ;
}
RooArgSet obs(*data.get()) ;
obs.remove(projDeps,kTRUE,kTRUE) ;
if (pdf.recursiveCheckObservables(&obs)) {
coutE(InputArguments) << "RooAbsOptTestStatistic: ERROR in PDF dependents, abort" << endl ;
RooErrorHandler::softAbort() ;
return ;
}
RooArgSet* pdfDepSet = pdf.getObservables(&data) ;
const RooArgSet* dataDepSet = data.get() ;
TIterator* iter = pdfDepSet->createIterator() ;
RooAbsArg* arg ;
while((arg=(RooAbsArg*)iter->Next())) {
RooRealVar* pdfReal = dynamic_cast<RooRealVar*>(arg) ;
if (!pdfReal) continue ;
RooRealVar* datReal = dynamic_cast<RooRealVar*>(dataDepSet->find(pdfReal->GetName())) ;
if (!datReal) continue ;
if (pdfReal->getMin()<(datReal->getMin()-1e-6)) {
coutE(InputArguments) << "RooAbsOptTestStatistic: ERROR minimum of PDF variable " << arg->GetName()
<< "(" << pdfReal->getMin() << ") is smaller than that of "
<< arg->GetName() << " in the dataset (" << datReal->getMin() << ")" << endl ;
RooErrorHandler::softAbort() ;
return ;
}
if (pdfReal->getMax()>(datReal->getMax()+1e-6)) {
coutE(InputArguments) << "RooAbsOptTestStatistic: ERROR maximum of PDF variable " << arg->GetName()
<< " is smaller than that of " << arg->GetName() << " in the dataset" << endl ;
RooErrorHandler::softAbort() ;
return ;
}
}
if (rangeName && strlen(rangeName)) {
_dataClone = ((RooAbsData&)data).reduce(RooFit::SelectVars(*pdfDepSet),RooFit::CutRange(rangeName)) ;
} else {
_dataClone = ((RooAbsData&)data).reduce(RooFit::SelectVars(*pdfDepSet)) ;
}
if (rangeName && strlen(rangeName)) {
cxcoutI(Fitting) << "RooAbsOptTestStatistic::ctor(" << GetName() << ") constructing likelihood for sub-range named " << rangeName << endl ;
TIterator* iter2 = _dataClone->get()->createIterator() ;
while((arg=(RooAbsArg*)iter2->Next())) {
RooRealVar* pdfReal = dynamic_cast<RooRealVar*>(arg) ;
if (pdfReal) {
if (!(addCoefRangeName && strlen(addCoefRangeName))) {
pdfReal->setRange(Form("NormalizationRangeFor%s",rangeName),pdfReal->getMin(),pdfReal->getMax()) ;
}
pdfReal->setRange(pdfReal->getMin(rangeName),pdfReal->getMax(rangeName)) ;
}
}
if (!_splitRange) {
iter->Reset() ;
while((arg=(RooAbsArg*)iter->Next())) {
RooRealVar* pdfReal = dynamic_cast<RooRealVar*>(arg) ;
if (pdfReal) {
pdfReal->setRange("fit",pdfReal->getMin(rangeName),pdfReal->getMax(rangeName)) ;
}
}
}
delete iter2 ;
}
delete iter ;
delete pdfDepSet ;
setEventCount(_dataClone->numEntries()) ;
RooArgSet tmp("PdfBranchNodeList") ;
pdf.branchNodeServerList(&tmp) ;
_pdfCloneSet = (RooArgSet*) tmp.snapshot(kFALSE) ;
_pdfClone = (RooAbsPdf*) _pdfCloneSet->find(pdf.GetName()) ;
if (rangeName) {
_pdfClone->fixAddCoefNormalization(*_dataClone->get()) ;
if (addCoefRangeName) {
cxcoutI(Fitting) << "RooAbsOptTestStatistic::ctor(" << GetName() << ") fixing interpretation of coefficients of any RooAddPdf component to range " << addCoefRangeName << endl ;
_pdfClone->fixAddCoefRange(addCoefRangeName,kFALSE) ;
} else {
cxcoutI(Fitting) << "RooAbsOptTestStatistic::ctor(" << GetName() << ") fixing interpretation of coefficients of any RooAddPdf to full domain of observables " << endl ;
_pdfClone->fixAddCoefRange(Form("NormalizationRangeFor%s",rangeName),kFALSE) ;
}
}
_pdfClone->attachDataSet(*_dataClone) ;
_normSet = (RooArgSet*) data.get()->snapshot(kFALSE) ;
if (projDeps.getSize()>0) {
_projDeps = (RooArgSet*) projDeps.snapshot(kFALSE) ;
RooArgSet* tobedel = (RooArgSet*) _normSet->selectCommon(*_projDeps) ;
_normSet->remove(*_projDeps,kTRUE,kTRUE) ;
TIterator* ii = tobedel->createIterator() ;
RooAbsArg* aa ;
while((aa=(RooAbsArg*)ii->Next())) {
delete aa ;
}
delete ii ;
delete tobedel ;
RooArgSet *projDataDeps = (RooArgSet*) _dataClone->get()->selectCommon(*_projDeps) ;
projDataDeps->setAttribAll("projectedDependent") ;
delete projDataDeps ;
}
RooArgSet* params = _pdfClone->getParameters(_dataClone) ;
_paramSet.add(*params) ;
delete params ;
if (_projDeps) {
RooArgSet *projDataDeps = (RooArgSet*) _dataClone->get()->selectCommon(*_projDeps) ;
projDataDeps->setAttribAll("projectedDependent") ;
delete projDataDeps ;
}
coutI(Optimization) << "RooAbsOptTestStatistic::ctor(" << GetName() << ") optimizing internal clone of p.d.f for likelihood evaluation."
<< "Lazy evaluation and associated change tracking will disabled for all nodes that depend on observables" << endl ;
optimizeCaching() ;
}
RooAbsOptTestStatistic::RooAbsOptTestStatistic(const RooAbsOptTestStatistic& other, const char* name) :
RooAbsTestStatistic(other,name)
{
if (operMode()!=Slave) {
_projDeps = 0 ;
_normSet = other._normSet ? ((RooArgSet*) other._normSet->snapshot()) : 0 ;
return ;
}
_pdfCloneSet = (RooArgSet*) other._pdfCloneSet->snapshot(kFALSE) ;
_pdfClone = (RooAbsPdf*) _pdfCloneSet->find(other._pdfClone->GetName()) ;
TIterator* iter = _pdfCloneSet->createIterator() ;
RooAbsArg* branch ;
while((branch=(RooAbsArg*)iter->Next())) {
branch->setOperMode(other._pdfCloneSet->find(branch->GetName())->operMode()) ;
}
delete iter ;
_dataClone = (RooAbsData*) other._dataClone->cacheClone(_pdfCloneSet) ;
_pdfClone->attachDataSet(*_dataClone) ;
_normSet = (RooArgSet*) other._normSet->snapshot() ;
if (other._projDeps) {
_projDeps = (RooArgSet*) other._projDeps->snapshot() ;
} else {
_projDeps = 0 ;
}
}
RooAbsOptTestStatistic::~RooAbsOptTestStatistic()
{
if (operMode()==Slave) {
delete _pdfCloneSet ;
delete _dataClone ;
delete _projDeps ;
}
delete _normSet ;
}
Double_t RooAbsOptTestStatistic::combinedValue(RooAbsReal** array, Int_t n) const
{
Double_t sum(0) ;
Int_t i ;
for (i=0 ; i<n ; i++) {
Double_t tmp = array[i]->getVal() ;
if (tmp==0) return 0 ;
sum += tmp ;
}
return sum ;
}
Bool_t RooAbsOptTestStatistic::redirectServersHook(const RooAbsCollection& newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t isRecursive)
{
RooAbsTestStatistic::redirectServersHook(newServerList,mustReplaceAll,nameChange,isRecursive) ;
if (operMode()!=Slave) return kFALSE ;
Bool_t ret = _pdfClone->recursiveRedirectServers(newServerList,kFALSE,nameChange) ;
return ret ;
}
void RooAbsOptTestStatistic::printCompactTreeHook(ostream& os, const char* indent)
{
RooAbsTestStatistic::printCompactTreeHook(os,indent) ;
if (operMode()!=Slave) return ;
TString indent2(indent) ;
indent2 += ">>" ;
_pdfClone->printCompactTree(os,indent2) ;
}
void RooAbsOptTestStatistic::constOptimizeTestStatistic(ConstOpCode opcode)
{
RooAbsTestStatistic::constOptimizeTestStatistic(opcode);
if (operMode()!=Slave) return ;
switch(opcode) {
case Activate:
cxcoutI(Optimization) << "RooAbsOptTestStatistic::constOptimize(" << GetName() << ") optimizing evaluation of test statistic by finding all nodes in p.d.f that depend exclusively"
<< " on observables and constant parameters and precalculating their values" << endl ;
optimizeConstantTerms(kTRUE) ;
break ;
case DeActivate:
cxcoutI(Optimization) << "RooAbsOptTestStatistic::constOptimize(" << GetName() << ") deactivating optimization of constant terms in test statistic" << endl ;
optimizeConstantTerms(kFALSE) ;
break ;
case ConfigChange:
cxcoutI(Optimization) << "RooAbsOptTestStatistic::constOptimize(" << GetName() << ") one ore more parameter were changed from constant to floating or vice versa, "
<< "re-evaluating constant term optimization" << endl ;
optimizeConstantTerms(kFALSE) ;
optimizeConstantTerms(kTRUE) ;
break ;
case ValueChange:
cxcoutI(Optimization) << "RooAbsOptTestStatistic::constOptimize(" << GetName() << ") the value of one ore more constant parameter were changed re-evaluating constant term optimization" << endl ;
optimizeConstantTerms(kFALSE) ;
optimizeConstantTerms(kTRUE) ;
break ;
}
}
void RooAbsOptTestStatistic::optimizeCaching()
{
_pdfClone->getVal(_normSet) ;
_pdfClone->optimizeCacheMode(*_dataClone->get()) ;
_dataClone->setDirtyProp(kFALSE) ;
_dataClone->optimizeReadingWithCaching(*_pdfClone, RooArgSet()) ;
}
void RooAbsOptTestStatistic::optimizeConstantTerms(Bool_t activate)
{
if(activate) {
_pdfClone->getVal(_normSet) ;
RooArgSet cacheableNodes ;
_pdfClone->findConstantNodes(*_dataClone->get(),cacheableNodes) ;
_dataClone->cacheArgs(cacheableNodes,_normSet) ;
TIterator* cIter = cacheableNodes.createIterator() ;
RooAbsArg *cacheArg ;
while((cacheArg=(RooAbsArg*)cIter->Next())){
cacheArg->setOperMode(RooAbsArg::AClean) ;
}
delete cIter ;
_dataClone->optimizeReadingWithCaching(*_pdfClone, cacheableNodes) ;
} else {
_dataClone->resetCache() ;
_dataClone->setArgStatus(*_dataClone->get(),kTRUE) ;
optimizeCaching() ;
}
}
Last update: Thu Jan 17 08:43:37 2008
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.