#include "RooDataHist.h"
#include "RooDataSet.h"
#include "RooGlobalFunc.h"
#include "RooNLLVar.h"
#include "RooRealVar.h"
#include "RooAbsData.h"
#include "RooWorkspace.h"
#include "TH1.h"
#include "RooStats/HybridCalculator.h"
ClassImp(RooStats::HybridCalculator)
using namespace RooStats;
HybridCalculator::HybridCalculator(const char *name) :
TNamed(name,name),
fSbModel(0),
fBModel(0),
fObservables(0),
fNuisanceParameters(0),
fPriorPdf(0),
fData(0),
fUsePriorPdf(false)
{
SetTestStatistic(1);
SetNumberOfToys(1000);
}
HybridCalculator::HybridCalculator( RooAbsPdf& sbModel,
RooAbsPdf& bModel,
RooArgList& observables,
const RooArgSet* nuisance_parameters,
RooAbsPdf* priorPdf ,
bool GenerateBinned,
int testStatistics,
int numToys) :
fSbModel(&sbModel),
fBModel(&bModel),
fNuisanceParameters(nuisance_parameters),
fPriorPdf(priorPdf),
fData(0),
fGenerateBinned(GenerateBinned),
fUsePriorPdf(false)
{
fObservables = new RooArgList(observables);
SetTestStatistic(testStatistics);
SetNumberOfToys(numToys);
if (priorPdf) UseNuisance(true);
}
HybridCalculator::HybridCalculator( RooAbsData & data,
RooAbsPdf& sbModel,
RooAbsPdf& bModel,
const RooArgSet* nuisance_parameters,
RooAbsPdf* priorPdf,
bool GenerateBinned,
int testStatistics,
int numToys) :
fSbModel(&sbModel),
fBModel(&bModel),
fObservables(0),
fNuisanceParameters(nuisance_parameters),
fPriorPdf(priorPdf),
fData(&data),
fGenerateBinned(GenerateBinned),
fUsePriorPdf(false)
{
SetTestStatistic(testStatistics);
SetNumberOfToys(numToys);
if (priorPdf) UseNuisance(true);
}
HybridCalculator::HybridCalculator( RooAbsData& data,
const ModelConfig& sbModel,
const ModelConfig& bModel,
bool GenerateBinned,
int testStatistics,
int numToys) :
fSbModel(sbModel.GetPdf()),
fBModel(bModel.GetPdf()),
fObservables(0),
fNuisanceParameters((sbModel.GetNuisanceParameters()) ? sbModel.GetNuisanceParameters() : bModel.GetNuisanceParameters()),
fPriorPdf((sbModel.GetPriorPdf()) ? sbModel.GetPriorPdf() : bModel.GetPriorPdf()),
fData(&data),
fGenerateBinned(GenerateBinned),
fUsePriorPdf(false)
{
if (fPriorPdf) UseNuisance(true);
SetTestStatistic(testStatistics);
SetNumberOfToys(numToys);
}
HybridCalculator::~HybridCalculator()
{
if (fObservables) delete fObservables;
}
void HybridCalculator::SetNullModel(const ModelConfig& model)
{
fBModel = model.GetPdf();
if (!fPriorPdf) fPriorPdf = model.GetPriorPdf();
if (!fNuisanceParameters) fNuisanceParameters = model.GetNuisanceParameters();
}
void HybridCalculator::SetAlternateModel(const ModelConfig& model)
{
fSbModel = model.GetPdf();
fPriorPdf = model.GetPriorPdf();
fNuisanceParameters = model.GetNuisanceParameters();
}
void HybridCalculator::SetTestStatistic(int index)
{
fTestStatisticsIdx = index;
}
HybridResult* HybridCalculator::Calculate(TH1& data, unsigned int nToys, bool usePriors) const
{
TString dataHistName = GetName(); dataHistName += "_roodatahist";
RooDataHist dataHist(dataHistName,"Data distribution as RooDataHist converted from TH1",*fObservables,&data);
HybridResult* result = Calculate(dataHist,nToys,usePriors);
return result;
}
HybridResult* HybridCalculator::Calculate(RooAbsData& data, unsigned int nToys, bool usePriors) const
{
double testStatData = 0;
if ( fTestStatisticsIdx==2 ) {
double nEvents = data.sumEntries();
testStatData = nEvents;
} else if ( fTestStatisticsIdx==3 ) {
RooNLLVar sb_nll("sb_nll","sb_nll",*fSbModel,data,RooFit::Extended());
fSbModel->fitTo(data);
double sb_nll_val = sb_nll.getVal();
RooNLLVar b_nll("b_nll","b_nll",*fBModel,data,RooFit::Extended());
fBModel->fitTo(data);
double b_nll_val = b_nll.getVal();
double m2lnQ = 2*(sb_nll_val-b_nll_val);
testStatData = m2lnQ;
} else if ( fTestStatisticsIdx==1 ) {
RooNLLVar sb_nll("sb_nll","sb_nll",*fSbModel,data,RooFit::Extended());
RooNLLVar b_nll("b_nll","b_nll",*fBModel,data,RooFit::Extended());
double m2lnQ = 2*(sb_nll.getVal()-b_nll.getVal());
testStatData = m2lnQ;
}
HybridResult* result = Calculate(nToys,usePriors);
result->SetDataTestStatistics(testStatData);
return result;
}
HybridResult* HybridCalculator::Calculate(unsigned int nToys, bool usePriors) const
{
std::vector<double> bVals;
bVals.reserve(nToys);
std::vector<double> sbVals;
sbVals.reserve(nToys);
RunToys(bVals,sbVals,nToys,usePriors);
HybridResult* result;
TString name = "HybridResult_" + TString(GetName() );
if ( fTestStatisticsIdx==2 )
result = new HybridResult(name,sbVals,bVals,false);
else
result = new HybridResult(name,sbVals,bVals);
return result;
}
void HybridCalculator::RunToys(std::vector<double>& bVals, std::vector<double>& sbVals, unsigned int nToys, bool usePriors) const
{
std::cout << "HybridCalculator: run " << nToys << " toy-MC experiments\n";
std::cout << "with test statistics index: " << fTestStatisticsIdx << "\n";
if (usePriors) std::cout << "marginalize nuisance parameters \n";
assert(nToys > 0);
assert(fBModel);
assert(fSbModel);
if (usePriors) {
assert(fPriorPdf);
assert(fNuisanceParameters);
}
std::vector<double> parameterValues;
int nParameters = (fNuisanceParameters) ? fNuisanceParameters->getSize() : 0;
RooArgList parametersList("parametersList");
if (usePriors && nParameters>0) {
parametersList.add(*fNuisanceParameters);
parameterValues.resize(nParameters);
for (int iParameter=0; iParameter<nParameters; iParameter++) {
RooRealVar* oneParam = (RooRealVar*) parametersList.at(iParameter);
parameterValues[iParameter] = oneParam->getVal();
}
}
for (unsigned int iToy=0; iToy<nToys; iToy++) {
if ( iToy%500==0 ) {
std::cout << "....... toy number " << iToy << " / " << nToys << std::endl;
}
if (usePriors && nParameters>0) {
RooDataSet* tmpValues = (RooDataSet*) fPriorPdf->generate(*fNuisanceParameters,1);
for (int iParameter=0; iParameter<nParameters; iParameter++) {
RooRealVar* oneParam = (RooRealVar*) parametersList.at(iParameter);
oneParam->setVal(tmpValues->get()->getRealValue(oneParam->GetName()));
}
delete tmpValues;
}
RooAbsData* bData;
if (fGenerateBinned)
bData = static_cast<RooAbsData*> (fBModel->generateBinned(*fObservables,RooFit::Extended()));
else
bData = static_cast<RooAbsData*> (fBModel->generate(*fObservables,RooFit::Extended()));
bool bIsEmpty = false;
if (bData==NULL) {
bIsEmpty = true;
RooDataSet* bDataDummy=new RooDataSet("bDataDummy","empty dataset",*fObservables);
bData = static_cast<RooAbsData*>(new RooDataHist ("bDataEmpty","",*fObservables,*bDataDummy));
delete bDataDummy;
}
RooAbsData* sbData;
if (fGenerateBinned)
sbData = static_cast<RooAbsData*> (fSbModel->generateBinned(*fObservables,RooFit::Extended()));
else
sbData = static_cast<RooAbsData*> (fSbModel->generate(*fObservables,RooFit::Extended()));
bool sbIsEmpty = false;
if (sbData==NULL) {
sbIsEmpty = true;
RooDataSet* sbDataDummy=new RooDataSet("sbDataDummy","empty dataset",*fObservables);
sbData = static_cast<RooAbsData*>(new RooDataHist ("sbDataEmpty","",*fObservables,*sbDataDummy));
delete sbDataDummy;
}
if (usePriors && nParameters>0) {
for (int iParameter=0; iParameter<nParameters; iParameter++) {
RooRealVar* oneParam = (RooRealVar*) parametersList.at(iParameter);
oneParam->setVal(parameterValues[iParameter]);
}
}
if ( fTestStatisticsIdx==2 ) {
double nEvents = 0;
if ( !sbIsEmpty ) nEvents = sbData->numEntries();
sbVals.push_back(nEvents);
} else if ( fTestStatisticsIdx==3 ) {
RooNLLVar sb_nll("sb_nll","sb_nll",*fSbModel,*sbData,RooFit::Extended());
fSbModel->fitTo(*sbData);
double sb_nll_val = sb_nll.getVal();
RooNLLVar b_nll("b_nll","b_nll",*fBModel,*sbData,RooFit::Extended());
fBModel->fitTo(*sbData);
double b_nll_val = b_nll.getVal();
double m2lnQ = 2*(sb_nll_val-b_nll_val);
sbVals.push_back(m2lnQ);
} else if ( fTestStatisticsIdx==1 ) {
RooNLLVar sb_nll("sb_nll","sb_nll",*fSbModel,*sbData,RooFit::Extended());
RooNLLVar b_nll("b_nll","b_nll",*fBModel,*sbData,RooFit::Extended());
double m2lnQ = 2*(sb_nll.getVal()-b_nll.getVal());
sbVals.push_back(m2lnQ);
}
if ( fTestStatisticsIdx==2 ) {
double nEvents = 0;
if ( !bIsEmpty ) nEvents = bData->numEntries();
bVals.push_back(nEvents);
} else if ( fTestStatisticsIdx==3 ) {
RooNLLVar sb_nll("sb_nll","sb_nll",*fSbModel,*bData,RooFit::Extended());
fSbModel->fitTo(*bData);
double sb_nll_val = sb_nll.getVal();
RooNLLVar b_nll("b_nll","b_nll",*fBModel,*bData,RooFit::Extended());
fBModel->fitTo(*bData);
double b_nll_val = b_nll.getVal();
double m2lnQ = 2*(sb_nll_val-b_nll_val);
bVals.push_back(m2lnQ);
} else if ( fTestStatisticsIdx==1 ) {
RooNLLVar sb_nll("sb_nll","sb_nll",*fSbModel,*bData,RooFit::Extended());
RooNLLVar b_nll("b_nll","b_nll",*fBModel,*bData,RooFit::Extended());
double m2lnQ = 2*(sb_nll.getVal()-b_nll.getVal());
bVals.push_back(m2lnQ);
}
delete sbData;
delete bData;
}
if (usePriors && nParameters>0) {
for (int iParameter=0; iParameter<nParameters; iParameter++) {
RooRealVar* oneParam = (RooRealVar*) parametersList.at(iParameter);
oneParam->setVal(parameterValues[iParameter]);
}
}
return;
}
void HybridCalculator::PrintMore(const char* options) const
{
if (fSbModel) {
std::cout << "Signal plus background model:\n";
fSbModel->Print(options);
}
if (fBModel) {
std::cout << "\nBackground model:\n";
fBModel->Print(options);
}
if (fObservables) {
std::cout << "\nObservables:\n";
fObservables->Print(options);
}
if (fNuisanceParameters) {
std::cout << "\nParameters being integrated:\n";
fNuisanceParameters->Print(options);
}
if (fPriorPdf) {
std::cout << "\nPrior PDF model for integration:\n";
fPriorPdf->Print(options);
}
return;
}
HybridResult* HybridCalculator::GetHypoTest() const {
if (!DoCheckInputs()) return 0;
RooAbsData * treeData = dynamic_cast<RooAbsData *> (fData);
if (!treeData) {
std::cerr << "Error in HybridCalculator::GetHypoTest - invalid data type - return NULL" << std::endl;
return 0;
}
bool usePrior = (fUsePriorPdf && fPriorPdf );
return Calculate( *treeData, fNToys, usePrior);
}
bool HybridCalculator::DoCheckInputs() const {
if (!fData) {
std::cerr << "Error in HybridCalculator - data have not been set" << std::endl;
return false;
}
if (!fObservables && fData->get() ) fObservables = new RooArgList( *fData->get() );
if (!fObservables) {
std::cerr << "Error in HybridCalculator - no observables" << std::endl;
return false;
}
if (!fSbModel) {
std::cerr << "Error in HybridCalculator - S+B pdf has not been set " << std::endl;
return false;
}
if (!fBModel) {
std::cerr << "Error in HybridCalculator - B pdf has not been set" << std::endl;
return false;
}
if (fUsePriorPdf && !fNuisanceParameters) {
std::cerr << "Error in HybridCalculator - nuisance parameters have not been set " << std::endl;
return false;
}
if (fUsePriorPdf && !fPriorPdf) {
std::cerr << "Error in HybridCalculator - prior pdf has not been set " << std::endl;
return false;
}
return true;
}