#include <BCRooInterface.h>


Public Member Functions | |
| BCRooInterface () | |
| BCRooInterface (const char *name) | |
| void | DefineParameters () |
| void | Initialize (const char *rootFile, const char *wsName="batWS", const char *dataName="data", const char *modelName="model", const char *priorName="priorPOI", const char *priorNuisanceName="priorNuisance", const char *paramsName="parameters", const char *listPOIName="POI") |
| void | Initialize (RooAbsData &data, RooAbsPdf &model, RooAbsPdf &prior, const RooArgSet *params, const RooArgSet &listPOI) |
| double | LogAPrioriProbability (std::vector< double > parameters) |
| double | LogLikelihood (std::vector< double > parameters) |
| ~BCRooInterface () | |
Private Attributes | |
| RooAbsData * | fData |
| RooAbsPdf * | fModel |
| RooAbsReal * | fNll |
| RooArgSet * | fObservables |
| RooArgList * | fParams |
| RooAbsPdf * | fPrior |
Definition at line 17 of file BCRooInterface.h.
| BCRooInterface::BCRooInterface | ( | ) |
Definition at line 118 of file BCRooInterface.cxx.
: BCModel() { // default constructor }
| BCRooInterface::BCRooInterface | ( | const char * | name | ) |
Definition at line 124 of file BCRooInterface.cxx.
: BCModel(name) { // another constructor }
| BCRooInterface::~BCRooInterface | ( | ) |
Definition at line 130 of file BCRooInterface.cxx.
{ // default destructor
}
| void BCRooInterface::DefineParameters | ( | ) |
Definition at line 136 of file BCRooInterface.cxx.
{ // define for BAT the list of parameters, range and plot binning
int default_nbins = 500;
int nParams = fParams->getSize();
for (int iParam=0; iParam<nParams; iParam++) {
RooRealVar* ipar = (RooRealVar*) fParams->at(iParam);
this->AddParameter(ipar->GetName(),ipar->getMin(),ipar->getMax());
this->SetNbins(ipar->GetName(),default_nbins);
std::cout << "added parameter: " << ipar->GetName() << " defined in range [ " << ipar->getMin() << " - " << ipar->getMax() << " ]\n";
}
}
| void BCRooInterface::Initialize | ( | const char * | rootFile, | |
| const char * | wsName = "batWS", |
|||
| const char * | dataName = "data", |
|||
| const char * | modelName = "model", |
|||
| const char * | priorName = "priorPOI", |
|||
| const char * | priorNuisanceName = "priorNuisance", |
|||
| const char * | paramsName = "parameters", |
|||
| const char * | listPOIName = "POI" | |||
| ) |
Definition at line 51 of file BCRooInterface.cxx.
{
// retrieve the RooFit inputs from the ROOT file
/*
// hard coded names in the workspace
char* rootFile = "bat_workspace.root";
char* wsName= "batWS";
char* dataName= "data";
char* modelName= "model";
char* priorName= "priorPOI";
char* priorNuisanceName= "priorNuisance";
char* paramsName= "parameters";
char* listPOIName= "POI";
*/
std::cout << "Opening " << rootFile << std::endl;
TFile* file = new TFile(rootFile);
std::cout << "content :\n";
file->ls();
RooWorkspace* bat_ws = (RooWorkspace*) file->Get(wsName);
bat_ws->Print("v");
fData = (RooAbsData*) bat_ws->data(dataName);
fModel = (RooAbsPdf*) bat_ws->function(modelName);
// make the product of both priors to get the full prior probability function
RooAbsPdf* priorPOI = (RooAbsPdf*) bat_ws->function(priorName);
RooAbsPdf* priorNuisance = (RooAbsPdf*) bat_ws->pdf(priorNuisanceName);
if (priorNuisance!=0 && priorPOI!=0) {
fPrior = new RooProdPdf("fPrior","complete prior",*priorPOI,*priorNuisance);
}
else {
if ( priorNuisance!=0 )
fPrior=priorNuisance;
else if ( priorPOI!=0 )
fPrior = priorPOI;
else
std::cout << "No prior PDF: the program will crash\n";
}
std::cout << "Imported parameters:\n";
fParams = new RooArgList(*(bat_ws->set(listPOIName)));
RooArgSet* paramsTmp = (RooArgSet*) bat_ws->set(paramsName);
if (paramsTmp!=0) fParams->add(*paramsTmp);
fParams->Print("v");
// create the log-likelihood function
//fNll = new RooNLLVar("fNll","",*fModel,*fData,true/*extended*/);
RooArgSet* constrainedParams = fModel->getParameters(*fData);
fNll = fModel->createNLL(*fData, RooFit::Constrain(*constrainedParams) );
file->Close();
DefineParameters();
}
| void BCRooInterface::Initialize | ( | RooAbsData & | data, | |
| RooAbsPdf & | model, | |||
| RooAbsPdf & | prior, | |||
| const RooArgSet * | params, | |||
| const RooArgSet & | listPOI | |||
| ) |
Definition at line 14 of file BCRooInterface.cxx.
{
fData = &data;
fModel = &model;
// make the product of both priors to get the full prior probability function
RooAbsPdf* prior_total = &prior_trans;
if (prior_total!=0 ) {
fPrior = prior_total;
}
else {
std::cout << "No prior PDF: the program will crash\n";
}
std::cout << "Imported parameters:\n";
fParams = new RooArgList(listPOI);
const RooArgSet* paramsTmp = params;
if (paramsTmp!=0)
fParams->add(*paramsTmp);
fParams->Print("v");
// create the log-likelihood function
// fNll = new RooNLLVar("fNll","",*fModel,*fData,true/*extended*/);
RooArgSet* constrainedParams = fModel->getParameters(*fData);
fNll = fModel->createNLL(*fData, RooFit::Constrain(*constrainedParams) );
DefineParameters();
}
| double BCRooInterface::LogAPrioriProbability | ( | std::vector< double > | parameters | ) | [virtual] |
Returns natural logarithm of the prior probability. Method needs to be overloaded by the user.
| parameters | A set of parameter values |
Reimplemented from BCModel.
Definition at line 167 of file BCRooInterface.cxx.
{ // this method returs the logarithm of the prior probability for the parameters: p(parameters).
// retrieve the values of the parameters to be tested
int nParams = fParams->getSize();
for (int iParam=0; iParam<nParams; iParam++) {
RooRealVar* ipar = (RooRealVar*) fParams->at(iParam);
ipar->setVal(parameters.at(iParam));
}
// compute the log of the prior function
RooArgSet* tmpArgSet = new RooArgSet(*fParams);
double prob = fPrior->getVal(tmpArgSet);
delete tmpArgSet;
if (prob<1e-300)
prob = 1e-300;
return log(prob);
}
| double BCRooInterface::LogLikelihood | ( | std::vector< double > | parameter | ) | [virtual] |
Calculates natural logarithm of the likelihood. Method needs to be overloaded by the user.
| parameters | A set of parameter values |
Reimplemented from BCModel.
Definition at line 151 of file BCRooInterface.cxx.
{ // this methods returns the logarithm of the conditional probability: p(data|parameters)
// retrieve the values of the parameters to be tested
int nParams = fParams->getSize();
for (int iParam=0; iParam<nParams; iParam++) {
RooRealVar* ipar = (RooRealVar*) fParams->at(iParam);
ipar->setVal(parameters.at(iParam));
}
// compute the log of the likelihood function
double logprob = -fNll->getVal();
return logprob;
}
RooAbsData* BCRooInterface::fData [private] |
Definition at line 52 of file BCRooInterface.h.
RooAbsPdf* BCRooInterface::fModel [private] |
Definition at line 53 of file BCRooInterface.h.
RooAbsReal* BCRooInterface::fNll [private] |
Definition at line 54 of file BCRooInterface.h.
RooArgSet* BCRooInterface::fObservables [private] |
Definition at line 55 of file BCRooInterface.h.
RooArgList* BCRooInterface::fParams [private] |
Definition at line 56 of file BCRooInterface.h.
RooAbsPdf* BCRooInterface::fPrior [private] |
Definition at line 57 of file BCRooInterface.h.
1.7.1