// RooAbsTestStatistic is the abstract base class for all test
// statistics. Test statistics that evaluate the PDF at each data
// point should inherit from the RooAbsOptTestStatistic class which
// implements several generic optimizations that can be done for such
// quantities.
//
// This test statistic base class organizes calculation of test
// statistic values for RooSimultaneous PDF as a combination of test
// statistic values for the PDF components of the simultaneous PDF and
// organizes multi-processor parallel calculation of test statistic
// values. For the latter, the test statistic value is calculated in
// partitions in parallel executing processes and a posteriori
// combined in the main thread.
// END_HTML
#include "RooFit.h"
#include "Riostream.h"
#include "RooAbsTestStatistic.h"
#include "RooAbsPdf.h"
#include "RooSimultaneous.h"
#include "RooAbsData.h"
#include "RooArgSet.h"
#include "RooRealVar.h"
#include "RooNLLVar.h"
#include "RooRealMPFE.h"
#include "RooErrorHandler.h"
#include "RooMsgService.h"
ClassImp(RooAbsTestStatistic)
;
RooAbsTestStatistic::RooAbsTestStatistic()
{
_init = kFALSE ;
_gofArray = 0 ;
_mpfeArray = 0 ;
}
RooAbsTestStatistic::RooAbsTestStatistic(const char *name, const char *title, RooAbsReal& real, RooAbsData& data,
const RooArgSet& projDeps, const char* rangeName, const char* addCoefRangeName,
Int_t nCPU, Bool_t interleave, Bool_t verbose, Bool_t splitCutRange) :
RooAbsReal(name,title),
_paramSet("paramSet","Set of parameters",this),
_func(&real),
_data(&data),
_projDeps((RooArgSet*)projDeps.Clone()),
_rangeName(rangeName?rangeName:""),
_addCoefRangeName(addCoefRangeName?addCoefRangeName:""),
_splitRange(splitCutRange),
_simCount(1),
_verbose(verbose),
_nGof(0),
_gofArray(0),
_nCPU(nCPU),
_mpfeArray(0),
_mpinterl(interleave)
{
RooArgSet* params = real.getParameters(&data) ;
_paramSet.add(*params) ;
delete params ;
if (_nCPU>1) {
_gofOpMode = MPMaster ;
} else {
Bool_t simMode = dynamic_cast<RooSimultaneous*>(&real)?kTRUE:kFALSE ;
if (simMode) {
_gofOpMode = SimMaster ;
} else {
_gofOpMode = Slave ;
}
}
_setNum = 0 ;
_numSets = 1 ;
_init = kFALSE ;
_nEvents = data.numEntries() ;
}
RooAbsTestStatistic::RooAbsTestStatistic(const RooAbsTestStatistic& other, const char* name) :
RooAbsReal(other,name),
_paramSet("paramSet","Set of parameters",this),
_func(other._func),
_data(other._data),
_projDeps((RooArgSet*)other._projDeps->Clone()),
_rangeName(other._rangeName),
_addCoefRangeName(other._addCoefRangeName),
_splitRange(other._splitRange),
_simCount(1),
_verbose(other._verbose),
_nGof(0),
_gofArray(0),
_nCPU(other._nCPU),
_mpfeArray(0),
_mpinterl(other._mpinterl)
{
_paramSet.add(other._paramSet) ;
if (_nCPU>1) {
_gofOpMode = MPMaster ;
} else {
Bool_t simMode = dynamic_cast<RooSimultaneous*>(_func)?kTRUE:kFALSE ;
if (simMode) {
_gofOpMode = SimMaster ;
} else {
_gofOpMode = Slave ;
}
}
_setNum = 0 ;
_numSets = 1 ;
_init = kFALSE ;
_nEvents = _data->numEntries() ;
}
RooAbsTestStatistic::~RooAbsTestStatistic()
{
if (_gofOpMode==MPMaster && _init) {
Int_t i ;
for (i=0 ; i<_nCPU ; i++) {
delete _mpfeArray[i] ;
}
delete[] _mpfeArray ;
}
if (_gofOpMode==SimMaster && _init) {
Int_t i ;
for (i=0 ; i<_nGof ; i++) {
delete _gofArray[i] ;
}
delete[] _gofArray ;
}
delete _projDeps ;
}
Double_t RooAbsTestStatistic::evaluate() const
{
if (!_init) {
const_cast<RooAbsTestStatistic*>(this)->initialize() ;
}
if (_gofOpMode==SimMaster) {
Double_t ret = combinedValue((RooAbsReal**)_gofArray,_nGof) ;
if (numSets()==1) {
ret /= globalNormalization() ;
}
return ret ;
} else if (_gofOpMode==MPMaster) {
Int_t i ;
for (i=0 ; i<_nCPU ; i++) {
_mpfeArray[i]->calculate() ;
}
Double_t ret = combinedValue((RooAbsReal**)_mpfeArray,_nCPU)/globalNormalization() ;
return ret ;
} else {
Int_t nFirst, nLast, nStep ;
if (_mpinterl) {
nFirst = _setNum ;
nLast = _nEvents ;
nStep = _numSets ;
} else {
nFirst = _nEvents * _setNum / _numSets ;
nLast = _nEvents * (_setNum+1) / _numSets ;
nStep = 1 ;
}
Double_t ret = evaluatePartition(nFirst,nLast,nStep) ;
if (numSets()==1) {
ret /= globalNormalization() ;
}
return ret ;
}
}
Bool_t RooAbsTestStatistic::initialize()
{
if (_init) return kFALSE ;
if (_gofOpMode==MPMaster) {
initMPMode(_func,_data,_projDeps,_rangeName.size()?_rangeName.c_str():0,_addCoefRangeName.size()?_addCoefRangeName.c_str():0) ;
} else if (_gofOpMode==SimMaster) {
initSimMode((RooSimultaneous*)_func,_data,_projDeps,_rangeName.size()?_rangeName.c_str():0,_addCoefRangeName.size()?_addCoefRangeName.c_str():0) ;
}
_init = kTRUE ;
return kFALSE ;
}
Bool_t RooAbsTestStatistic::redirectServersHook(const RooAbsCollection& newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t)
{
if (_gofOpMode==SimMaster && _gofArray) {
Int_t i ;
for (i=0 ; i<_nGof ; i++) {
if (_gofArray[i]) {
_gofArray[i]->recursiveRedirectServers(newServerList,mustReplaceAll,nameChange) ;
}
}
} else if (_gofOpMode==MPMaster && _mpfeArray) {
Int_t i ;
for (i=0 ; i<_nCPU ; i++) {
if (_mpfeArray[i]) {
_mpfeArray[i]->recursiveRedirectServers(newServerList,mustReplaceAll,nameChange) ;
cout << "redirecting servers on " << _mpfeArray[i]->GetName() << endl ;
}
}
}
return kFALSE ;
}
void RooAbsTestStatistic::printCompactTreeHook(ostream& os, const char* indent)
{
if (_gofOpMode==SimMaster) {
Int_t i ;
os << indent << "RooAbsTestStatistic begin GOF contents" << endl ;
for (i=0 ; i<_nGof ; i++) {
if (_gofArray[i]) {
TString indent2(indent) ;
indent2 += Form("[%d] ",i) ;
_gofArray[i]->printCompactTreeHook(os,indent2) ;
}
}
os << indent << "RooAbsTestStatistic end GOF contents" << endl ;
} else if (_gofOpMode==MPMaster) {
}
}
void RooAbsTestStatistic::constOptimizeTestStatistic(ConstOpCode opcode)
{
Int_t i ;
initialize() ;
if (_gofOpMode==SimMaster) {
for (i=0 ; i<_nGof ; i++) {
if (_gofArray[i]) _gofArray[i]->constOptimizeTestStatistic(opcode) ;
}
} else if (_gofOpMode==MPMaster) {
for (i=0 ; i<_nCPU ; i++) {
_mpfeArray[i]->constOptimizeTestStatistic(opcode) ;
}
}
}
void RooAbsTestStatistic::setMPSet(Int_t inSetNum, Int_t inNumSets)
{
_setNum = inSetNum ; _numSets = inNumSets ;
if (_gofOpMode==SimMaster) {
initialize() ;
Int_t i ;
for (i=0 ; i<_nGof ; i++) {
if (_gofArray[i]) _gofArray[i]->setMPSet(inSetNum,inNumSets) ;
}
}
}
void RooAbsTestStatistic::initMPMode(RooAbsReal* real, RooAbsData* data, const RooArgSet* projDeps, const char* rangeName, const char* addCoefRangeName)
{
Int_t i ;
_mpfeArray = new pRooRealMPFE[_nCPU] ;
RooAbsTestStatistic* gof = create(GetName(),GetTitle(),*real,*data,*projDeps,rangeName,addCoefRangeName,1,_mpinterl,_verbose,_splitRange) ;
gof->recursiveRedirectServers(_paramSet) ;
for (i=0 ; i<_nCPU ; i++) {
gof->setMPSet(i,_nCPU) ;
gof->SetName(Form("%s_GOF%d",GetName(),i)) ;
gof->SetTitle(Form("%s_GOF%d",GetTitle(),i)) ;
Bool_t doInline = (i==_nCPU-1) ;
if (!doInline) coutI(Eval) << "RooAbsTestStatistic::initMPMode: starting remote server process #" << i << endl ;
_mpfeArray[i] = new RooRealMPFE(Form("%s_%x_MPFE%d",GetName(),this,i),Form("%s_%x_MPFE%d",GetTitle(),this,i),*gof,doInline) ;
_mpfeArray[i]->initialize() ;
}
return ;
}
void RooAbsTestStatistic::initSimMode(RooSimultaneous* simpdf, RooAbsData* data,
const RooArgSet* projDeps, const char* rangeName, const char* addCoefRangeName)
{
RooAbsCategoryLValue& simCat = (RooAbsCategoryLValue&) simpdf->indexCat() ;
TString simCatName(simCat.GetName()) ;
TList* dsetList = const_cast<RooAbsData*>(data)->split(simCat) ;
if (!dsetList) {
coutE(Fitting) << "RooAbsTestStatistic::initSimMode(" << GetName() << ") unable to split dataset, abort" << endl ;
RooErrorHandler::softAbort() ;
}
Int_t n(0) ;
_nGof = 0 ;
RooCatType* type ;
TIterator* catIter = simCat.typeIterator() ;
while((type=(RooCatType*)catIter->Next())){
RooAbsPdf* pdf = simpdf->getPdf(type->GetName()) ;
RooAbsData* dset = (RooAbsData*) dsetList->FindObject(type->GetName()) ;
if (pdf && dset && (dset->numEntries(kTRUE)!=0. || processEmptyDataSets() )) {
_nGof++ ;
}
}
_gofArray = new pRooAbsTestStatistic[_nGof] ;
catIter->Reset() ;
while((type=(RooCatType*)catIter->Next())){
RooAbsPdf* pdf = simpdf->getPdf(type->GetName()) ;
RooAbsData* dset = (RooAbsData*) dsetList->FindObject(type->GetName()) ;
if (pdf && dset && (dset->numEntries(kTRUE)!=0. || processEmptyDataSets())) {
coutI(Fitting) << "RooAbsTestStatistic::initSimMode: creating slave calculator #" << n << " for state " << type->GetName()
<< " (" << dset->numEntries() << " dataset entries)" << endl ;
if (_splitRange && rangeName) {
_gofArray[n] = create(type->GetName(),type->GetName(),*pdf,*dset,*projDeps,
Form("%s_%s",rangeName,type->GetName()),addCoefRangeName,_nCPU*(_mpinterl?-1:1),_mpinterl,_verbose,_splitRange) ;
} else {
_gofArray[n] = create(type->GetName(),type->GetName(),*pdf,*dset,*projDeps,
rangeName,addCoefRangeName,_nCPU,_mpinterl,_verbose,_splitRange) ;
}
_gofArray[n]->setSimCount(_nGof) ;
_gofArray[n]->recursiveRedirectServers(_paramSet) ;
n++ ;
} else {
if ((!dset || (dset->numEntries(kTRUE)==0. && !processEmptyDataSets()) ) && pdf) {
if (_verbose) {
coutI(Fitting) << "RooAbsTestStatistic::initSimMode: state " << type->GetName()
<< " has no data entries, no slave calculator created" << endl ;
}
}
}
}
dsetList->Delete() ;
delete dsetList ;
delete catIter ;
}
Last change: Fri Oct 31 17:53:53 2008
Last generated: 2008-10-31 17:53
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.