BCRooInterface Class Reference

#include <BCRooInterface.h>

Inherits BCModel.

Collaboration diagram for BCRooInterface:
Collaboration graph
[legend]

List of all members.

Public Member Functions

 BCRooInterface (const char *name)
 BCRooInterface ()
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")
double LogAPrioriProbability (std::vector< double > parameters)
double LogLikelihood (std::vector< double > parameters)
 ~BCRooInterface ()

Private Attributes

RooAbsData * fData
RooAbsPdf * fModel
RooNLLVar * fNll
RooArgSet * fObservables
RooArgList * fParams
RooAbsPdf * fPrior

Detailed Description

Definition at line 17 of file BCRooInterface.h.


Constructor & Destructor Documentation

BCRooInterface::BCRooInterface (  ) 

Definition at line 74 of file BCRooInterface.cxx.

00074                                : BCModel()
00075 {  // default constructor
00076 
00077 }

BCRooInterface::BCRooInterface ( const char *  name  ) 

Definition at line 80 of file BCRooInterface.cxx.

00080                                                : BCModel(name)
00081 {  // another constructor
00082 
00083 }

BCRooInterface::~BCRooInterface (  ) 

Definition at line 86 of file BCRooInterface.cxx.

00087 {  // default destructor
00088 
00089 }


Member Function Documentation

void BCRooInterface::DefineParameters (  ) 

Definition at line 92 of file BCRooInterface.cxx.

00093 {  // define for BAT the list of parameters, range and plot binning
00094 
00095   int default_nbins = 500;
00096    
00097   int nParams = fParams->getSize();
00098   for (int iParam=0; iParam<nParams; iParam++) {
00099     RooRealVar* ipar = (RooRealVar*) fParams->at(iParam);
00100     this->AddParameter(ipar->GetName(),ipar->getMin(),ipar->getMax());
00101     this->SetNbins(ipar->GetName(),default_nbins);
00102     std::cout << "added parameter: " << ipar->GetName() << " defined in range [ " << ipar->getMin() << " - " << ipar->getMax() << " ]\n";
00103   }
00104 }

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 14 of file BCRooInterface.cxx.

00022 {
00023   // retrieve the RooFit inputs from the ROOT file
00024 
00025   /*
00026   // hard coded names in the workspace
00027   char* rootFile = "bat_workspace.root";
00028   char* wsName= "batWS";
00029   char* dataName= "data";
00030   char* modelName= "model";
00031   char* priorName= "priorPOI";
00032   char* priorNuisanceName= "priorNuisance";
00033   char* paramsName= "parameters";
00034   char* listPOIName= "POI";
00035   */
00036 
00037   std::cout << "Opening " << rootFile << std::endl;
00038   TFile* file = new TFile(rootFile);
00039   std::cout << "content :\n";
00040   file->ls();
00041   
00042   RooWorkspace* bat_ws = (RooWorkspace*) file->Get(wsName);
00043   bat_ws->Print("v");
00044   
00045   fData = (RooAbsData*) bat_ws->data(dataName);
00046   fModel = (RooAbsPdf*) bat_ws->function(modelName);
00047   
00048   // make the product of both priors to get the full prior probability function
00049   RooAbsPdf* priorPOI = (RooAbsPdf*) bat_ws->function(priorName);
00050   RooAbsPdf* priorNuisance = (RooAbsPdf*) bat_ws->pdf(priorNuisanceName);
00051   if (priorNuisance!=0 && priorPOI!=0) {
00052     fPrior = new RooProdPdf("fPrior","complete prior",*priorPOI,*priorNuisance);
00053   } else {
00054     if ( priorNuisance!=0 ) fPrior=priorNuisance;
00055     else if ( priorPOI!=0 ) fPrior = priorPOI;
00056     else std::cout << "No prior PDF: the program will crash\n";
00057   }
00058   
00059   std::cout << "Imported parameters:\n";
00060   fParams  = new RooArgList(*(bat_ws->set(listPOIName)));
00061   RooArgSet* paramsTmp = (RooArgSet*) bat_ws->set(paramsName);
00062   if (paramsTmp!=0) fParams->add(*paramsTmp);
00063   fParams->Print("v");
00064   
00065   // create the log-likelihood function
00066   fNll = new RooNLLVar("fNll","",*fModel,*fData,true/*extended*/);
00067   
00068   file->Close();
00069   
00070   DefineParameters();
00071 }

double BCRooInterface::LogAPrioriProbability ( std::vector< double >  parameters  )  [virtual]

Returns natural logarithm of the prior probability. Method needs to be overloaded by the user.

Parameters:
parameters A set of parameter values
Returns:
The prior probability p(parameters)
See also:
GetPrior(std::vector <double> parameters)

Reimplemented from BCModel.

Definition at line 123 of file BCRooInterface.cxx.

00124 {  // this method returs the logarithm of the prior probability for the parameters: p(parameters).
00125 
00126    // retrieve the values of the parameters to be tested
00127   int nParams = fParams->getSize();
00128   for (int iParam=0; iParam<nParams; iParam++) {
00129     RooRealVar* ipar = (RooRealVar*) fParams->at(iParam);
00130     ipar->setVal(parameters.at(iParam));
00131   }
00132 
00133   // compute the log of the prior function
00134   RooArgSet* tmpArgSet = new RooArgSet(*fParams);
00135   double prob = fPrior->getVal(tmpArgSet);
00136   delete tmpArgSet;
00137   if (prob<1e-300) prob = 1e-300;
00138   return log(prob);
00139 }

double BCRooInterface::LogLikelihood ( std::vector< double >  parameter  )  [virtual]

Calculates natural logarithm of the likelihood. Method needs to be overloaded by the user.

Parameters:
parameters A set of parameter values
Returns:
Natural logarithm of the likelihood

Reimplemented from BCModel.

Definition at line 107 of file BCRooInterface.cxx.

00108 {  // this methods returns the logarithm of the conditional probability: p(data|parameters)
00109 
00110    // retrieve the values of the parameters to be tested
00111   int nParams = fParams->getSize();
00112   for (int iParam=0; iParam<nParams; iParam++) {
00113     RooRealVar* ipar = (RooRealVar*) fParams->at(iParam);
00114     ipar->setVal(parameters.at(iParam));
00115   }
00116 
00117   // compute the log of the likelihood function
00118   double logprob = -fNll->getVal();
00119   return logprob;
00120 }


Member Data Documentation

RooAbsData* BCRooInterface::fData [private]

Definition at line 46 of file BCRooInterface.h.

RooAbsPdf* BCRooInterface::fModel [private]

Definition at line 47 of file BCRooInterface.h.

RooNLLVar* BCRooInterface::fNll [private]

Definition at line 48 of file BCRooInterface.h.

RooArgSet* BCRooInterface::fObservables [private]

Definition at line 49 of file BCRooInterface.h.

RooArgList* BCRooInterface::fParams [private]

Definition at line 50 of file BCRooInterface.h.

RooAbsPdf* BCRooInterface::fPrior [private]

Definition at line 51 of file BCRooInterface.h.


The documentation for this class was generated from the following files:

Generated on Tue Oct 6 09:48:23 2009 for Bayesian Analysis Toolkit by  doxygen 1.6.1