BayesianAnalysisToolkit  0.9.3
Public Member Functions | Private Member Functions | Private Attributes | List of all members
BCRooInterface Class Reference

#include <BCRooInterface.h>

Inheritance diagram for BCRooInterface:
Inheritance graph
[legend]
Collaboration diagram for BCRooInterface:
Collaboration graph
[legend]

Public Member Functions

 BCRooInterface (const char *name="", bool fillChain=false)
 
 ~BCRooInterface ()
 
void DefineParameters ()
 
double LogAPrioriProbability (const std::vector< double > &parameters)
 
double LogLikelihood (const std::vector< double > &parameters)
 
void Initialize (RooAbsData &data, RooAbsPdf &model, RooAbsPdf &prior, const RooArgSet *params, const RooArgSet &listPOI)
 
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 SetNumBins (const char *parname, int nbins)
 
void SetNumBins (int nbins)
 
void SetupRooStatsMarkovChain ()
 
void MCMCIterationInterface ()
 
RooStats::MarkovChain * GetRooStatsMarkovChain ()
 
RooArgSet * GetArgSetForMarkovChain ()
 
- Public Member Functions inherited from BCModel
 BCModel (const char *name="model")
 
 BCModel (const BCModel &bcmodel)
 
virtual ~BCModel ()
 
BCModeloperator= (const BCModel &bcmodel)
 
const std::string & GetName () const
 
double GetModelAPrioriProbability () const
 
double GetModelAPosterioriProbability () const
 
BCDataSetGetDataSet () const
 
BCDataPointGetDataPointLowerBoundaries () const
 
BCDataPointGetDataPointUpperBoundaries () const
 
double GetDataPointLowerBoundary (unsigned int index) const
 
double GetDataPointUpperBoundary (unsigned int index) const
 
bool GetFlagBoundaries () const
 
unsigned GetNDataPoints () const
 
BCDataPointGetDataPoint (unsigned int index) const
 
void SetName (const char *name)
 
void SetModelAPrioriProbability (double probability)
 
void SetModelAPosterioriProbability (double probability)
 
virtual int AddParameter (BCParameter *parameter)
 
void SetDataSet (BCDataSet *dataset)
 
void SetSingleDataPoint (BCDataPoint *datapoint)
 
void SetSingleDataPoint (BCDataSet *dataset, unsigned int index)
 
void SetDataBoundaries (unsigned int index, double lowerboundary, double upperboundary, bool fixed=false)
 
void SetDataPointLowerBoundaries (BCDataPoint *datasetlowerboundaries)
 
void SetDataPointUpperBoundaries (BCDataPoint *datasetupperboundaries)
 
void SetDataPointLowerBoundary (int index, double lowerboundary)
 
void SetDataPointUpperBoundary (int index, double upperboundary)
 
int SetPrior (int index, TF1 *f)
 
int SetPrior (const char *name, TF1 *f)
 
int SetPriorDelta (int index, double value)
 
int SetPriorDelta (const char *name, double value)
 
int SetPriorGauss (int index, double mean, double sigma)
 
int SetPriorGauss (const char *name, double mean, double sigma)
 
int SetPriorGauss (int index, double mean, double sigmadown, double sigmaup)
 
int SetPriorGauss (const char *name, double mean, double sigmadown, double sigmaup)
 
int SetPrior (int index, TH1 *h, bool flag=false)
 
int SetPrior (const char *name, TH1 *h, bool flag=false)
 
int SetPriorConstant (int index)
 
int SetPriorConstant (const char *name)
 
int SetPriorConstantAll ()
 
void Copy (const BCModel &bcmodel)
 
double APrioriProbability (const std::vector< double > &parameters)
 
virtual double Likelihood (const std::vector< double > &params)
 
double ProbabilityNN (const std::vector< double > &params)
 
double LogProbabilityNN (const std::vector< double > &parameters)
 
double Probability (const std::vector< double > &parameter)
 
double LogProbability (const std::vector< double > &parameter)
 
virtual double SamplingFunction (const std::vector< double > &parameters)
 
double Eval (const std::vector< double > &parameters)
 
virtual double LogEval (const std::vector< double > &parameters)
 
virtual void CorrelateDataPointValues (std::vector< double > &x)
 
double GetPvalueFromChi2 (const std::vector< double > &par, int sigma_index)
 
double GetPvalueFromKolmogorov (const std::vector< double > &par, int index)
 
double GetPvalueFromChi2NDoF (std::vector< double > par, int sigma_index)
 
BCH1DCalculatePValue (std::vector< double > par, bool flag_histogram=false)
 
double GetPValue ()
 
double GetPValueNDoF ()
 
double GetChi2NDoF ()
 
std::vector< double > GetChi2Runs (int dataIndex, int sigmaIndex)
 
void SetGoFNIterationsMax (int n)
 
void SetGoFNIterationsRun (int n)
 
void SetGoFNChains (int n)
 
double HessianMatrixElement (const BCParameter *parameter1, const BCParameter *parameter2, std::vector< double > point)
 
void PrintSummary ()
 
void PrintResults (const char *file)
 
void PrintShortFitSummary (int chi2flag=0)
 
void PrintHessianMatrix (std::vector< double > parameters)
 
virtual double CDF (const std::vector< double > &, int, bool)
 
- Public Member Functions inherited from BCIntegrate
 BCIntegrate ()
 
 BCIntegrate (const BCIntegrate &bcintegrate)
 
virtual ~BCIntegrate ()
 
BCIntegrateoperator= (const BCIntegrate &bcintegrate)
 
int ReadMarginalizedFromFile (const char *file)
 
BCH1DGetMarginalized (const BCParameter *parameter)
 
BCH1DGetMarginalized (const char *name)
 
BCH1DGetMarginalized (unsigned index)
 
BCH2DGetMarginalized (const BCParameter *parameter1, const BCParameter *parameter2)
 
BCH2DGetMarginalized (const char *name1, const char *name2)
 
BCH2DGetMarginalized (unsigned index1, unsigned index2)
 
int PrintAllMarginalized1D (const char *filebase)
 
int PrintAllMarginalized2D (const char *filebase)
 
int PrintAllMarginalized (const char *file, std::string options1d="BTsiB3CS1D0pdf0Lmeanmode", std::string options2d="BTfB3CS1meangmode", unsigned int hdiv=1, unsigned int ndiv=1)
 
double GetIntegral () const
 
BCIntegrate::BCOptimizationMethod GetOptimizationMethod () const
 
BCIntegrate::BCIntegrationMethod GetIntegrationMethod () const
 
BCIntegrate::BCMarginalizationMethod GetMarginalizationMethod () const
 
BCIntegrate::BCSASchedule GetSASchedule () const
 
void GetRandomVectorUnitHypercube (std::vector< double > &x) const
 
void GetRandomVectorInParameterSpace (std::vector< double > &x) const
 
double GetRandomPoint (std::vector< double > &x)
 
int GetNIterationsMin () const
 
int GetNIterationsMax () const
 
int GetNIterationsPrecisionCheck () const
 
int GetNIterationsOutput () const
 
int GetNIterations () const
 
double GetRelativePrecision () const
 
double GetAbsolutePrecision () const
 
BCCubaMethod GetCubaIntegrationMethod () const
 
const BCCubaOptions::VegasGetCubaVegasOptions () const
 
const BCCubaOptions::SuaveGetCubaSuaveOptions () const
 
const BCCubaOptions::DivonneGetCubaDivonneOptions () const
 
const BCCubaOptions::CuhreGetCubaCuhreOptions () const
 
BCH1DGetSlice (const BCParameter *parameter, const std::vector< double > parameters=std::vector< double >(0), int bins=0, bool flag_norm=true)
 
BCH1DGetSlice (const char *name, const std::vector< double > parameters=std::vector< double >(0), int nbins=0, bool flag_norm=true)
 
BCH2DGetSlice (const BCParameter *parameter1, const BCParameter *parameter2, const std::vector< double > parameters=std::vector< double >(0), int bins=0, bool flag_norm=true)
 
BCH2DGetSlice (const char *name1, const char *name2, const std::vector< double > parameters=std::vector< double >(0), int nbins=0, bool flag_norm=true)
 
BCH2DGetSlice (unsigned index1, unsigned index2, const std::vector< double > parameters=std::vector< double >(0), int nbins=0, bool flag_norm=true)
 
double GetError () const
 
TMinuit * GetMinuit ()
 
int GetMinuitErrorFlag () const
 
double GetSAT0 () const
 
double GetSATmin () const
 
double GetBestFitParameter (unsigned index) const
 
double GetBestFitParameterError (unsigned index) const
 
double GetLogMaximum ()
 
const std::vector< double > & GetBestFitParameters () const
 
const std::vector< double > & GetBestFitParameterErrors () const
 
void SetMinuitArlist (double *arglist)
 
void SetFlagIgnorePrevOptimization (bool flag)
 
void SetOptimizationMethod (BCIntegrate::BCOptimizationMethod method)
 
void SetIntegrationMethod (BCIntegrate::BCIntegrationMethod method)
 
void SetMarginalizationMethod (BCIntegrate::BCMarginalizationMethod method)
 
void SetSASchedule (BCIntegrate::BCSASchedule schedule)
 
void SetNIterationsMin (int niterations)
 
void SetNIterationsMax (int niterations)
 
void SetNIterationsPrecisionCheck (int niterations)
 
void SetNIterationsOutput (int niterations)
 
void SetRelativePrecision (double relprecision)
 
void SetAbsolutePrecision (double absprecision)
 
void SetCubaIntegrationMethod (BCCubaMethod type)
 
void SetCubaOptions (const BCCubaOptions::Vegas &options)
 
void SetCubaOptions (const BCCubaOptions::Suave &options)
 
void SetCubaOptions (const BCCubaOptions::Divonne &options)
 
void SetCubaOptions (const BCCubaOptions::Cuhre &options)
 
void SetSAT0 (double T0)
 
void SetSATmin (double Tmin)
 
void SetFlagWriteSAToFile (bool flag)
 
TTree * GetSATree ()
 
void InitializeSATree ()
 
double Normalize ()
 
double Integrate (BCIntegrationMethod intmethod)
 
double Integrate ()
 
double Integrate (BCIntegrationMethod type, tRandomizer randomizer, tEvaluator evaluator, tIntegralUpdater updater, std::vector< double > &sums)
 
double EvaluatorMC (std::vector< double > &sums, const std::vector< double > &point, bool &accepted)
 
int MarginalizeAll ()
 
int MarginalizeAll (BCMarginalizationMethod margmethod)
 
virtual void MarginalizePreprocess ()
 
virtual void MarginalizePostprocess ()
 
void SAInitialize ()
 
std::vector< double > FindMode (std::vector< double > start=std::vector< double >())
 
std::vector< double > FindMode (BCIntegrate::BCOptimizationMethod optmethod, std::vector< double > start=std::vector< double >())
 
double SATemperature (double t)
 
double SATemperatureBoltzmann (double t)
 
double SATemperatureCauchy (double t)
 
virtual double SATemperatureCustom (double t)
 
std::vector< double > GetProposalPointSA (const std::vector< double > &x, int t)
 
std::vector< double > GetProposalPointSABoltzmann (const std::vector< double > &x, int t)
 
std::vector< double > GetProposalPointSACauchy (const std::vector< double > &x, int t)
 
virtual std::vector< double > GetProposalPointSACustom (const std::vector< double > &x, int t)
 
std::vector< double > SAHelperGetRandomPointOnHypersphere ()
 
double SAHelperGetRadialCauchy ()
 
double SAHelperSinusToNIntegral (int dim, double theta)
 
virtual void ResetResults ()
 
std::string DumpIntegrationMethod (BCIntegrationMethod type)
 
std::string DumpCurrentIntegrationMethod ()
 
std::string DumpUsedIntegrationMethod ()
 
std::string DumpMarginalizationMethod (BCMarginalizationMethod type)
 
std::string DumpCurrentMarginalizationMethod ()
 
std::string DumpUsedMarginalizationMethod ()
 
std::string DumpOptimizationMethod (BCOptimizationMethod type)
 
std::string DumpCurrentOptimizationMethod ()
 
std::string DumpUsedOptimizationMethod ()
 
std::string DumpCubaIntegrationMethod (BCCubaMethod type)
 
std::string DumpCubaIntegrationMethod ()
 
void SetBestFitParameters (const std::vector< double > &x)
 
void SetBestFitParameters (const std::vector< double > &x, const double &new_value, double &old_value)
 
unsigned GetNIntegrationVariables ()
 
double CalculateIntegrationVolume ()
 
bool CheckMarginalizationAvailability (BCMarginalizationMethod type)
 
bool CheckMarginalizationIndices (TH1 *hist, const std::vector< unsigned > &index)
 
- Public Member Functions inherited from BCEngineMCMC
void WriteMarkovChain (bool flag)
 
 BCEngineMCMC ()
 
 BCEngineMCMC (const BCEngineMCMC &enginemcmc)
 
virtual ~BCEngineMCMC ()
 
BCEngineMCMCoperator= (const BCEngineMCMC &engineMCMC)
 
unsigned MCMCGetNChains () const
 
unsigned MCMCGetNLag () const
 
const std::vector< unsigned > & MCMCGetNIterations () const
 
unsigned MCMCGetCurrentIteration () const
 
unsigned MCMCGetCurrentChain () const
 
unsigned MCMCGetNIterationsConvergenceGlobal () const
 
bool MCMCGetFlagConvergenceGlobal () const
 
unsigned MCMCGetNIterationsMax () const
 
unsigned MCMCGetNIterationsRun () const
 
unsigned MCMCGetNIterationsPreRunMin () const
 
unsigned MCMCGetNIterationsUpdate () const
 
unsigned MCMCGetNIterationsUpdateMax () const
 
std::vector< int > MCMCGetNTrialsTrue () const
 
int MCMCGetNTrials () const
 
const std::vector< double > & MCMCGetprobMean () const
 
const std::vector< double > & MCMCGetVariance () const
 
const std::vector< double > & MCMCGetTrialFunctionScaleFactor () const
 
std::vector< double > MCMCGetTrialFunctionScaleFactor (unsigned ichain) const
 
double MCMCGetTrialFunctionScaleFactor (unsigned ichain, unsigned ipar)
 
const std::vector< double > & MCMCGetx () const
 
std::vector< double > MCMCGetx (unsigned ichain)
 
double MCMCGetx (unsigned ichain, unsigned ipar) const
 
const std::vector< double > & MCMCGetLogProbx () const
 
double MCMCGetLogProbx (unsigned ichain)
 
int MCMCGetPhase () const
 
const std::vector< double > & MCMCGetMaximumPoints () const
 
std::vector< double > MCMCGetMaximumPoint (unsigned i) const
 
const std::vector< double > & MCMCGetMaximumLogProb () const
 
int MCMCGetFlagInitialPosition () const
 
double MCMCGetRValueCriterion () const
 
double MCMCGetRValueParametersCriterion () const
 
double MCMCGetRValue () const
 
double MCMCGetRValueParameters (unsigned i)
 
bool MCMCGetRValueStrict () const
 
bool MCMCGetFlagRun () const
 
TTree * MCMCGetMarkovChainTree (unsigned i)
 
BCH1DMCMCGetH1Marginalized (unsigned i)
 
BCH2DMCMCGetH2Marginalized (unsigned i, unsigned j)
 
BCParameterGetParameter (int index) const
 
BCParameterGetParameter (const char *name) const
 
unsigned int GetNParameters () const
 
unsigned int GetNFixedParameters ()
 
unsigned int GetNFreeParameters ()
 
const std::vector< double > & GetBestFitParametersMarginalized () const
 
void MCMCSetTrialFunctionScaleFactor (std::vector< double > scale)
 
void MCMCSetNChains (unsigned n)
 
void MCMCSetNLag (unsigned n)
 
void MCMCSetNIterationsMax (unsigned n)
 
void MCMCSetNIterationsRun (unsigned n)
 
void MCMCSetNIterationsPreRunMin (unsigned n)
 
void MCMCSetNIterationsUpdate (unsigned n)
 
void MCMCSetNIterationsUpdateMax (unsigned n)
 
void MCMCSetMinimumEfficiency (double efficiency)
 
void MCMCSetMaximumEfficiency (double efficiency)
 
void MCMCSetRandomSeed (unsigned seed)
 
void MCMCSetWriteChainToFile (bool flag)
 
void MCMCSetWritePreRunToFile (bool flag)
 
void MCMCSetInitialPositions (const std::vector< double > &x0s)
 
void MCMCSetInitialPositions (std::vector< std::vector< double > > x0s)
 
void MCMCSetFlagInitialPosition (int flag)
 
void MCMCSetFlagOrderParameters (bool flag)
 
void MCMCSetFlagFillHistograms (bool flag)
 
void MCMCSetFlagPreRun (bool flag)
 
void MCMCSetRValueCriterion (double r)
 
void MCMCSetRValueParametersCriterion (double r)
 
void MCMCSetRValueStrict (bool strict=true)
 
void MCMCSetMarkovChainTrees (const std::vector< TTree * > &trees)
 
void MCMCInitializeMarkovChainTrees ()
 
int SetMarginalized (unsigned index, TH1D *h)
 
int SetMarginalized (unsigned index1, unsigned index2, TH2D *h)
 
void MCMCSetValuesDefault ()
 
void MCMCSetValuesQuick ()
 
void MCMCSetValuesDetail ()
 
void MCMCSetPrecision (BCEngineMCMC::Precision precision)
 
void SetNbins (unsigned int nbins)
 
void Copy (const BCEngineMCMC &enginemcmc)
 
virtual int AddParameter (const char *name, double min, double max, const char *latexname="")
 
virtual void MCMCTrialFunction (unsigned ichain, std::vector< double > &x)
 
virtual double MCMCTrialFunctionSingle (unsigned ichain, unsigned ipar)
 
bool MCMCGetProposalPointMetropolis (unsigned chain, std::vector< double > &x)
 
bool MCMCGetProposalPointMetropolis (unsigned chain, unsigned parameter, std::vector< double > &x)
 
bool MCMCGetNewPointMetropolis (unsigned chain=0)
 
bool MCMCGetNewPointMetropolis (unsigned chain, unsigned parameter)
 
void MCMCInChainCheckMaximum ()
 
void MCMCInChainUpdateStatistics ()
 
void MCMCInChainFillHistograms ()
 
void MCMCInChainTestConvergenceAllChains ()
 
void MCMCInChainWriteChains ()
 
int MCMCMetropolis ()
 
int MCMCMetropolisPreRun ()
 
void MCMCResetRunStatistics ()
 
void MCMCInitializeMarkovChains ()
 
int MCMCInitialize ()
 
void ClearParameters (bool hard=false)
 
virtual void MCMCUserIterationInterface ()
 
virtual void MCMCCurrentPointInterface (std::vector< double > &, int, bool)
 

Private Member Functions

void AddToCurrentChainElement (double xij, int chainNum, int poiNum)
 
bool EqualsLastChainElement (int chainNum)
 
double GetWeightForChain (int chainNum)
 

Private Attributes

RooAbsData * fData
 
RooAbsPdf * fModel
 
RooAbsReal * fNll
 
RooArgSet * fObservables
 
RooArgList * fParams
 
RooArgList * fParamsPOI
 
RooAbsPdf * fPrior
 
int _default_nbins
 
RooRealVar * priorhelpvar
 
bool _addeddummyprior
 
bool _fillChain
 
bool fFirstComparison
 
RooStats::MarkovChain * _roostatsMarkovChain
 
RooArgSet _parametersForMarkovChainPrevious
 
RooArgSet _parametersForMarkovChainCurrent
 
std::vector< std::vector
< double > > 
fPreviousStep
 
std::vector< std::vector
< double > > 
fCurrentStep
 
std::vector< double > fVecWeights
 
std::list< std::pair< const
char *, int > > 
_nbins_list
 

Additional Inherited Members

- Public Types inherited from BCIntegrate
enum  BCOptimizationMethod {
  kOptEmpty, kOptSimAnn, kOptMetropolis, kOptMinuit,
  kOptDefault, NOptMethods
}
 
enum  BCIntegrationMethod {
  kIntEmpty, kIntMonteCarlo, kIntCuba, kIntGrid,
  kIntDefault, NIntMethods
}
 
enum  BCMarginalizationMethod {
  kMargEmpty, kMargMetropolis, kMargMonteCarlo, kMargGrid,
  kMargDefault, NMargMethods
}
 
enum  BCSASchedule { kSACauchy, kSABoltzmann, kSACustom, NSAMethods }
 
enum  BCCubaMethod {
  kCubaVegas, kCubaSuave, kCubaDivonne, kCubaCuhre,
  NCubaMethods
}
 
typedef void(BCIntegrate::* tRandomizer )(std::vector< double > &) const
 
typedef double(BCIntegrate::* tEvaluator )(std::vector< double > &, const std::vector< double > &, bool &)
 
typedef void(* tIntegralUpdater )(const std::vector< double > &, const int &, double &, double &)
 
- Static Public Member Functions inherited from BCIntegrate
static void IntegralUpdaterMC (const std::vector< double > &sums, const int &nIterations, double &integral, double &absprecision)
 
static int CubaIntegrand (const int *ndim, const double xx[], const int *ncomp, double ff[], void *userdata)
 
static void FCNLikelihood (int &npar, double *grad, double &fval, double *par, int flag)
 
- Protected Member Functions inherited from BCIntegrate
unsigned IntegrationOutputFrequency () const
 
void LogOutputAtStartOfIntegration (BCIntegrationMethod type, BCCubaMethod cubatype)
 
void LogOutputAtIntegrationStatusUpdate (BCIntegrationMethod type, double integral, double absprecision, int nIterations)
 
void LogOutputAtEndOfIntegration (double integral, double absprecision, double relprecision, int nIterations)
 
void Copy (const BCIntegrate &bcintegrate)
 
- Protected Attributes inherited from BCModel
std::string fName
 
double fModelAPriori
 
double fModelAPosteriori
 
BCDataSetfDataSet
 
BCDataPointfDataPointLowerBoundaries
 
BCDataPointfDataPointUpperBoundaries
 
std::vector< bool > fDataFixedValues
 
double fPValue
 
double fChi2NDoF
 
double fPValueNDoF
 
bool flag_discrete
 
int fGoFNIterationsMax
 
int fGoFNIterationsRun
 
int fGoFNChains
 
std::vector< TNamed * > fPriorContainer
 
bool fPriorConstantAll
 
std::vector< bool > fPriorContainerConstant
 
std::vector< bool > fPriorContainerInterpolate
 

Detailed Description

Definition at line 23 of file BCRooInterface.h.

Constructor & Destructor Documentation

BCRooInterface::BCRooInterface ( const char *  name = "",
bool  fillChain = false 
)

Definition at line 144 of file BCRooInterface.cxx.

144  :
145  BCModel(name),
146  fData(NULL),
147  fModel(NULL),
148  fNll(NULL),
149  fObservables(NULL),
150  fParams(NULL),
151  fParamsPOI(NULL),
152  fPrior(NULL),
153  _default_nbins(150),
154  priorhelpvar(NULL),
155  _addeddummyprior(false),
156  _fillChain(fillChain),
157  fFirstComparison(false),
159 {
160  // todo this interface not ready for grid marginalization yet
162 }
BCRooInterface::~BCRooInterface ( )

Definition at line 165 of file BCRooInterface.cxx.

166 { // default destructor
167  if(_fillChain) {
168  delete _roostatsMarkovChain;
169  }
170 }

Member Function Documentation

void BCRooInterface::AddToCurrentChainElement ( double  xij,
int  chainNum,
int  poiNum 
)
private

Definition at line 355 of file BCRooInterface.cxx.

356 {
357  fCurrentStep[chainNum][poiNum] = xij;
358 }
void BCRooInterface::DefineParameters ( )

Definition at line 173 of file BCRooInterface.cxx.

174 { // define for BAT the list of parameters, range and plot binning
175 
176  int nParams = fParams->getSize();
177  for (int iParam=0; iParam<nParams; iParam++) {
178  RooRealVar* ipar = (RooRealVar*) fParams->at(iParam);
179  BCParameter * bcpar = new BCParameter(ipar->GetName(),ipar->getMin(),ipar->getMax());
180  bcpar->SetNbins(_default_nbins);
181  this->AddParameter(bcpar);
182  std::cout << "added parameter: " << bcpar->GetName() << " defined in range [ " << bcpar->GetLowerLimit() << " - "
183  << bcpar->GetUpperLimit() << " ] with number of bins: " << bcpar->GetNbins() << " \n";
184  }
185 
186  for(std::list< std::pair<const char*,int> >::iterator listiter = _nbins_list.begin(); listiter != _nbins_list.end(); listiter++) {
187  GetParameter((*listiter).first)->SetNbins((*listiter).second);
188  std::cout << "adjusted parameter: " << (*listiter).first << " to number of bins: " << (*listiter).second << " \n";
189  }
190 }
bool BCRooInterface::EqualsLastChainElement ( int  chainNum)
private

Definition at line 360 of file BCRooInterface.cxx.

361 {
362  bool equals = true;
363  std::vector<double>::iterator itPrevious = fPreviousStep[chainNum].begin();
364 
365  if(fFirstComparison == true) {
366  fFirstComparison = false;
368  return true;
369  }
370 
371 
372  for (std::vector<double>::iterator itCurrent = fCurrentStep[chainNum].begin(); itCurrent != fCurrentStep[chainNum].end(); ++itCurrent) {
373  if(*itCurrent != *itPrevious) {
374  equals = false;
375  fPreviousStep[chainNum] = fCurrentStep[chainNum];
376  break;
377  }
378  ++itPrevious;
379  }
380 
381  if(equals == true) {
382  fVecWeights[chainNum] += 1.0;
383  }
384 
385  return equals;
386 
387 }
RooArgSet* BCRooInterface::GetArgSetForMarkovChain ( )
inline

Definition at line 63 of file BCRooInterface.h.

RooStats::MarkovChain* BCRooInterface::GetRooStatsMarkovChain ( )
inline

Definition at line 62 of file BCRooInterface.h.

62 { return _roostatsMarkovChain;}
double BCRooInterface::GetWeightForChain ( int  chainNum)
private

Definition at line 389 of file BCRooInterface.cxx.

390 {
391  double retval = fVecWeights[chainNum];
392  fVecWeights[chainNum]= 1.0 ;
393  return retval;
394 }
void BCRooInterface::Initialize ( RooAbsData &  data,
RooAbsPdf &  model,
RooAbsPdf &  prior,
const RooArgSet *  params,
const RooArgSet &  listPOI 
)

Definition at line 20 of file BCRooInterface.cxx.

25 {
26 
27  fData = &data;
28  fModel = &model;
29 
30  // make the product of both priors to get the full prior probability function
31  RooAbsPdf* prior_total = &prior_trans;
32  if (prior_total!=0 ) {
33  fPrior = prior_total;
34  }
35  else {
36  std::cout << "No prior PDF: without taking action the program would crash\n";
37  std::cout << "No prior PDF: adding dummy uniform pdf on the interval [0..1]\n";
38  priorhelpvar = new RooRealVar("priorhelpvar","",0.0, 1.0 );
39  _addeddummyprior = true;
40  RooUniform* priorhelpdist = new RooUniform("priorhelpdist","", *priorhelpvar);
41  fPrior = priorhelpdist;
42  }
43 
44  std::cout << "Imported parameters:\n";
45  fParams = new RooArgList(listPOI);
46  const RooArgSet* paramsTmp = params;
47  if (paramsTmp!=0)
48  fParams->add(*paramsTmp);
49  fParams->Print("v");
50 
51  fParamsPOI = new RooArgList(listPOI);
52  // create the log-likelihood function
53  // fNll = new RooNLLVar("fNll","",*fModel,*fData,true/*extended*/);
54 
55  RooArgSet* constrainedParams = fModel->getParameters(*fData);
56  fNll = fModel->createNLL(*fData, RooFit::Constrain(*constrainedParams) );
57 
59 
60  if(_fillChain) {
62  }
63 }
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 67 of file BCRooInterface.cxx.

75 {
76  // retrieve the RooFit inputs from the ROOT file
77 
78  /*
79  // hard coded names in the workspace
80  char* rootFile = "bat_workspace.root";
81  char* wsName= "batWS";
82  char* dataName= "data";
83  char* modelName= "model";
84  char* priorName= "priorPOI";
85  char* priorNuisanceName= "priorNuisance";
86  char* paramsName= "parameters";
87  char* listPOIName= "POI";
88  */
89 
90  std::cout << "Opening " << rootFile << std::endl;
91  TFile* file = TFile::Open(rootFile);
92  std::cout << "content :\n";
93  file->ls();
94 
95  RooWorkspace* bat_ws = (RooWorkspace*) file->Get(wsName);
96  bat_ws->Print("v");
97 
98  fData = (RooAbsData*) bat_ws->data(dataName);
99  fModel = (RooAbsPdf*) bat_ws->function(modelName);
100 
101  // make the product of both priors to get the full prior probability function
102  RooAbsPdf* priorPOI = (RooAbsPdf*) bat_ws->function(priorName);
103  RooAbsPdf* priorNuisance = (RooAbsPdf*) bat_ws->pdf(priorNuisanceName);
104  if (priorNuisance!=0 && priorPOI!=0) {
105  fPrior = new RooProdPdf("fPrior","complete prior",*priorPOI,*priorNuisance);
106  }
107  else {
108  if ( priorNuisance!=0 )
109  fPrior=priorNuisance;
110  else if ( priorPOI!=0 )
111  fPrior = priorPOI;
112  else{
113  std::cout << "No prior PDF: without taking action the program would crash\n";
114  std::cout << "No prior PDF: adding dummy uniform pdf on the interval [0..1]\n";
115  priorhelpvar = new RooRealVar("priorhelpvar","",0.0, 1.0 );
116  _addeddummyprior = true;
117  RooUniform* priorhelpdist = new RooUniform("priorhelpdist","", *priorhelpvar);
118  fPrior = priorhelpdist;
119  }
120  }
121 
122  std::cout << "Imported parameters:\n";
123  fParams = new RooArgList(*(bat_ws->set(listPOIName)));
124  RooArgSet* paramsTmp = (RooArgSet*) bat_ws->set(paramsName);
125  if (paramsTmp!=0) {
126  fParams->add(*paramsTmp);
127  }
128  if (_addeddummyprior == true ) {
129  fParams->add(*priorhelpvar);
130  }
131  fParams->Print("v");
132 
133  // create the log-likelihood function
134  //fNll = new RooNLLVar("fNll","",*fModel,*fData,true/*extended*/);
135  RooArgSet* constrainedParams = fModel->getParameters(*fData);
136  fNll = fModel->createNLL(*fData, RooFit::Constrain(*constrainedParams) );
137 
138  file->Close();
139 
141 }
double BCRooInterface::LogAPrioriProbability ( const std::vector< double > &  parameters)
virtual

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

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

Reimplemented from BCModel.

Definition at line 209 of file BCRooInterface.cxx.

210 { // this method returs the logarithm of the prior probability for the parameters: p(parameters).
211  // retrieve the values of the parameters to be tested
212  int nParams = fParams->getSize();
213  for (int iParam=0; iParam<nParams; iParam++) {
214  RooRealVar* ipar = (RooRealVar*) fParams->at(iParam);
215  ipar->setVal(parameters.at(iParam));
216  }
217  // compute the log of the prior function
218  RooArgSet* tmpArgSet = new RooArgSet(*fParams);
219  double prob = fPrior->getVal(tmpArgSet);
220  delete tmpArgSet;
221  if (prob<1e-300)
222  prob = 1e-300;
223  return log(prob);
224 }
double BCRooInterface::LogLikelihood ( const std::vector< double > &  params)
virtual

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

Parameters
paramsA set of parameter values
Returns
Natural logarithm of the likelihood

Implements BCModel.

Definition at line 193 of file BCRooInterface.cxx.

194 { // this methods returns the logarithm of the conditional probability: p(data|parameters)
195 
196  // retrieve the values of the parameters to be tested
197  int nParams = fParams->getSize();
198  for (int iParam=0; iParam<nParams; iParam++) {
199  RooRealVar* ipar = (RooRealVar*) fParams->at(iParam);
200  ipar->setVal(parameters.at(iParam));
201  }
202 
203  // compute the log of the likelihood function
204  double logprob = -fNll->getVal();
205  return logprob;
206 }
void BCRooInterface::MCMCIterationInterface ( )
virtual

Interface allowing to execute arbitrary code for each iteration of the MCMC. The frequency of calling this method is influenced by the setup of the Lag and whether or not the MCMC is run with ordered parameters. This method needs to be overloaded in the derived class.

Reimplemented from BCEngineMCMC.

Definition at line 294 of file BCRooInterface.cxx.

295 {
296  //fIterationInterfacecount+=1;
297 
298  if(_fillChain) {
299  //std::cout << "Hei-ho running with fillChain!" << std::endl;
300  // get number of chains
301  int nchains = MCMCGetNChains();
302 
303  // get number of parameters
304  int npar = GetNParameters();
305  //std::cout << "number of parameters is: " << npar << std::endl;
306 
307  // loop over all chains and fill histogram
308  for (int i = 0; i < nchains; ++i) {
309  // get the current values of the parameters. These are
310  // stored in fMCMCx.
311 
312  // std::cout << "log(likelihood*prior)" << *GetMarkovChainValue() << std::endl; //does not work this way?!
313  TIterator* setiter = fParams->createIterator();
314  int j = 0;
315 
316  //store only POI in RooStats MarkovChain object
317  //TIterator* setiter = fParamsPOI->createIterator();
318  //int j = 0;
319 
320  while(setiter->Next()){
321 
322  //check parameter names
323  const BCParameter * tempBCparam = GetParameter(j);
324 
325  //_parametersForMarkovChainCurrent->Print();
326 
327  const char * paramnamepointer = (tempBCparam->GetName()).c_str();
328  double xij = fMCMCx.at(i * npar + j);
329  AddToCurrentChainElement(xij, i, j);
330  RooRealVar* parampointer = (RooRealVar*) &(_parametersForMarkovChainCurrent[paramnamepointer]);
331  parampointer->setVal(xij);
332  //std::cout << "Chain " << i << " param: " << tempBCparam->GetName() << " Value: " << xij << std::endl;
333  j++;
334  }
335 
336 
337  // will only work if _parametersForMarkovChain had correct info!
338 
339  //test stuff begin
340  //var1->setVal( RooRandom::randomGenerator()->Uniform(100.) );
341  //
342  //_roostatsMarkovChain->Add(_parametersForMarkovChain_test, 0.001, 1.0);
343  //
344  //test stuff end
345 
346  if( !(EqualsLastChainElement(i)) ) {
347  double weight = GetWeightForChain(i);
350  }
351  }
352  }
353 }
void BCRooInterface::SetNumBins ( const char *  parname,
int  nbins 
)

Definition at line 226 of file BCRooInterface.cxx.

227 {
228  for(std::list< std::pair<const char*,int> >::iterator listiter = _nbins_list.begin(); listiter != _nbins_list.end(); listiter++) {
229  if(!strcmp((*listiter).first, parname)) {
230  (*listiter).second = nbins;
231  return;
232  }
233  }
234  _nbins_list.push_back( std::make_pair(parname,nbins) );
235 }
void BCRooInterface::SetNumBins ( int  nbins)

Definition at line 237 of file BCRooInterface.cxx.

238 {
239  _default_nbins = nbins;
240 }
void BCRooInterface::SetupRooStatsMarkovChain ( )

Definition at line 242 of file BCRooInterface.cxx.

243 {
244  //RooArgSet * tempRooArgSetPointer = new RooArgSet(* fParams, "paramsMarkovChainPointer");
245  //_parametersForMarkovChain = RooArgSet(* fParams, "paramsMarkovChainPointer");
246  //fParams->snapshot(_parametersForMarkovChain);
247 
248  //store only POI in RooStats MarkovChain object
249  //_parametersForMarkovChainPrevious.add(*fParamsPOI);
250  //_parametersForMarkovChainCurrent.add(*fParamsPOI);
251 
254 
255  std::cout << "size of _parametersForMarkovChain: " << _parametersForMarkovChainCurrent.getSize() << std::endl;
256  std::cout << "size of fParamsPOI: " << fParamsPOI->getSize() << std::endl;
257  //std::cout << "size of tempRooArgSetPointer: " << tempRooArgSetPointer->getSize() << std::endl;
258 
259  _roostatsMarkovChain = new RooStats::MarkovChain();
260  //test stuff begin
261  //the following way of creating the MArkovChain object does not work!:
262  //_roostatsMarkovChain = new RooStats::MarkovChain(_parametersForMarkovChain);
263  //test stuff end
264  std::cout << "setting up parameters for RooStats markov chain" << std::endl;
265  _parametersForMarkovChainPrevious.writeToStream(std::cout, false);
266 
267  //initialize fPreviousStep, fCurrentStep, fVecWeights
268  int nchains = MCMCGetNChains();
269  for(int countChains = 1; countChains<=nchains ; countChains++ ) {
270  double tempweight = 1.0;
271  fVecWeights.push_back(tempweight);
272  std::vector<double> tempvec;
273  TIterator* setiter = fParamsPOI->createIterator();
274  double tempval = 0.;
275  while(setiter->Next()){
276  tempvec.push_back(tempval);
277  }
278  fPreviousStep.push_back(tempvec);
279  fCurrentStep.push_back(tempvec);
280  }
281 
282  fFirstComparison = true;
283 
284  //test stuff begin
285  //var1 = new RooRealVar("var1","var1", 10., 0., 100.);
286  //fParamsTest = new RooArgList();
287  //fParamsTest->add(*var1);
288  //_parametersForMarkovChain_test.add(*fParamsTest);
289  //fIterationInterfacecount = 0;
290  //test stuff end
291 }

Member Data Documentation

bool BCRooInterface::_addeddummyprior
private

Definition at line 83 of file BCRooInterface.h.

int BCRooInterface::_default_nbins
private

Definition at line 80 of file BCRooInterface.h.

bool BCRooInterface::_fillChain
private

Definition at line 85 of file BCRooInterface.h.

std::list< std::pair<const char*,int> > BCRooInterface::_nbins_list
private

Definition at line 95 of file BCRooInterface.h.

RooArgSet BCRooInterface::_parametersForMarkovChainCurrent
private

Definition at line 89 of file BCRooInterface.h.

RooArgSet BCRooInterface::_parametersForMarkovChainPrevious
private

Definition at line 88 of file BCRooInterface.h.

RooStats::MarkovChain* BCRooInterface::_roostatsMarkovChain
private

Definition at line 87 of file BCRooInterface.h.

std::vector< std::vector<double> > BCRooInterface::fCurrentStep
private

Definition at line 92 of file BCRooInterface.h.

RooAbsData* BCRooInterface::fData
private

Definition at line 73 of file BCRooInterface.h.

bool BCRooInterface::fFirstComparison
private

Definition at line 86 of file BCRooInterface.h.

RooAbsPdf* BCRooInterface::fModel
private

Definition at line 74 of file BCRooInterface.h.

RooAbsReal* BCRooInterface::fNll
private

Definition at line 75 of file BCRooInterface.h.

RooArgSet* BCRooInterface::fObservables
private

Definition at line 76 of file BCRooInterface.h.

RooArgList* BCRooInterface::fParams
private

Definition at line 77 of file BCRooInterface.h.

RooArgList* BCRooInterface::fParamsPOI
private

Definition at line 78 of file BCRooInterface.h.

std::vector< std::vector<double> > BCRooInterface::fPreviousStep
private

Definition at line 91 of file BCRooInterface.h.

RooAbsPdf* BCRooInterface::fPrior
private

Definition at line 79 of file BCRooInterface.h.

std::vector< double > BCRooInterface::fVecWeights
private

Definition at line 93 of file BCRooInterface.h.

RooRealVar* BCRooInterface::priorhelpvar
private

Definition at line 82 of file BCRooInterface.h.


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