#ifndef RooStats_FeldmanCousins
#include "RooStats/FeldmanCousins.h"
#endif
#ifndef RooStats_RooStatsUtils
#include "RooStats/RooStatsUtils.h"
#endif
#ifndef RooStats_PointSetInterval
#include "RooStats/PointSetInterval.h"
#endif
#include "RooStats/ModelConfig.h"
#include "RooStats/SamplingDistribution.h"
#include "RooStats/ProfileLikelihoodTestStat.h"
#include "RooStats/NeymanConstruction.h"
#include "RooStats/RooStatsUtils.h"
#include "RooDataSet.h"
#include "RooDataHist.h"
#include "RooGlobalFunc.h"
#include "RooMsgService.h"
#include "TFile.h"
#include "TTree.h"
ClassImp(RooStats::FeldmanCousins) ;
using namespace RooFit;
using namespace RooStats;
FeldmanCousins::FeldmanCousins(RooAbsData& data, ModelConfig& model) :
fSize(0.05),
fModel(model),
fData(data),
fTestStatSampler(0),
fPointsToTest(0),
fConfBelt(0),
fAdaptiveSampling(false),
fAdditionalNToysFactor(1.),
fNbins(10),
fFluctuateData(true),
fDoProfileConstruction(true),
fSaveBeltToFile(false),
fCreateBelt(false)
{
}
FeldmanCousins::~FeldmanCousins() {
if(fPointsToTest) delete fPointsToTest;
if(fTestStatSampler) delete fTestStatSampler;
}
void FeldmanCousins::SetModel(const ModelConfig & model) {
fModel = model;
}
TestStatSampler* FeldmanCousins::GetTestStatSampler() const{
if(!fTestStatSampler)
this->CreateTestStatSampler();
return fTestStatSampler;
}
void FeldmanCousins::CreateTestStatSampler() const{
ProfileLikelihoodTestStat* testStatistic = new ProfileLikelihoodTestStat(*fModel.GetPdf());
fTestStatSampler = new ToyMCSampler(*testStatistic,int(fAdditionalNToysFactor*50./fSize)) ;
fTestStatSampler->SetParametersForTestStat(*fModel.GetParametersOfInterest() );
if(fModel.GetObservables())
fTestStatSampler->SetObservables(*fModel.GetObservables());
fTestStatSampler->SetPdf(*fModel.GetPdf());
if(!fAdaptiveSampling){
ooccoutP(&fModel,Generation) << "FeldmanCousins: ntoys per point = " << int (fAdditionalNToysFactor*50./fSize )<< endl;
} else{
ooccoutP(&fModel,Generation) << "FeldmanCousins: ntoys per point: adaptive" << endl;
}
if(fFluctuateData){
ooccoutP(&fModel,Generation) << "FeldmanCousins: nEvents per toy will fluctuate about expectation" << endl;
} else{
ooccoutP(&fModel,Generation) << "FeldmanCousins: nEvents per toy will not fluctuate, will always be " << fData.numEntries() << endl;
fTestStatSampler->SetNEventsPerToy(fData.numEntries());
}
}
void FeldmanCousins::CreateParameterPoints() const{
RooAbsPdf* pdf = fModel.GetPdf();
if (!pdf ){
ooccoutE(&fModel,Generation) << "FeldmanCousins: ModelConfig has no PDF" << endl;
return;
}
RooArgSet* parameters = new RooArgSet(*fModel.GetParametersOfInterest());
if(fModel.GetNuisanceParameters())
parameters->add(*fModel.GetNuisanceParameters());
if( fModel.GetNuisanceParameters() && ! fModel.GetParametersOfInterest()->equals(*parameters) && fDoProfileConstruction) {
ooccoutP(&fModel,Generation) << "FeldmanCousins: Model has nuisance parameters, will do profile construction" << endl;
TIter it2 = fModel.GetParametersOfInterest()->createIterator();
RooRealVar *myarg2;
while ((myarg2 = dynamic_cast<RooRealVar*>(it2.Next()))) {
myarg2->setBins(fNbins);
}
RooDataHist* parameterScan = new RooDataHist("parameterScan", "", *fModel.GetParametersOfInterest());
ooccoutP(&fModel,Generation) << "FeldmanCousins: # points to test = " << parameterScan->numEntries() << endl;
RooFit::MsgLevel previous = RooMsgService::instance().globalKillBelow();
RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL) ;
RooAbsReal* nll = pdf->createNLL(fData,RooFit::CloneData(false));
RooAbsReal* profile = nll->createProfile(*fModel.GetParametersOfInterest());
RooDataSet* profileConstructionPoints = new RooDataSet("profileConstruction",
"profileConstruction",
*parameters);
for(Int_t i=0; i<parameterScan->numEntries(); ++i){
*parameters = *parameterScan->get(i);
profile->getVal();
profileConstructionPoints->add(*parameters);
}
RooMsgService::instance().setGlobalKillBelow(previous) ;
delete profile;
delete nll;
delete parameterScan;
fPointsToTest = profileConstructionPoints;
} else{
ooccoutP(&fModel,Generation) << "FeldmanCousins: Model has no nuisance parameters" << endl;
TIter it = parameters->createIterator();
RooRealVar *myarg;
while ((myarg = dynamic_cast<RooRealVar*>(it.Next()))) {
myarg->setBins(fNbins);
}
RooDataHist* parameterScan = new RooDataHist("parameterScan", "", *parameters);
ooccoutP(&fModel,Generation) << "FeldmanCousins: # points to test = " << parameterScan->numEntries() << endl;
fPointsToTest = parameterScan;
}
delete parameters;
}
PointSetInterval* FeldmanCousins::GetInterval() const {
fModel.GuessObsAndNuisance(fData);
if(!fTestStatSampler)
this->CreateTestStatSampler();
fTestStatSampler->SetObservables(*fModel.GetObservables());
if(!fFluctuateData)
fTestStatSampler->SetNEventsPerToy(fData.numEntries());
this->CreateParameterPoints();
NeymanConstruction nc(fData,fModel);
nc.SetTestStatSampler(*fTestStatSampler);
nc.SetTestSize(fSize);
nc.SetParameterPointsToTest( *fPointsToTest );
nc.SetLeftSideTailFraction(0.);
nc.SetData(fData);
nc.UseAdaptiveSampling(fAdaptiveSampling);
nc.AdditionalNToysFactor(fAdditionalNToysFactor);
nc.SaveBeltToFile(fSaveBeltToFile);
nc.CreateConfBelt(fCreateBelt);
fConfBelt = nc.GetConfidenceBelt();
return nc.GetInterval();
}