BAT  0.9.4
The Bayesian analysis toolkit
 All Classes Namespaces Functions Variables Enumerations
BCRooInterface.h
1 #ifndef __BCROOINTERFACE__H
2 #define __BCROOINTERFACE__H
3 
4 #include "../../BAT/BCModel.h"
5 
6 #include <RooStats/MarkovChain.h>
7 #include <RooRealVar.h>
8 
9 #include <list>
10 #include <utility>
11 
12 class RooAbsReal;
13 class RooAbsData;
14 class RooAbsPdf;
15 class RooArgSet;
16 class RooArgList;
17 
23 class BCRooInterface : public BCModel
24 {
25  public:
26 
27  BCRooInterface( const char* name = "", bool fillChain = false );
28 
29  ~BCRooInterface();
30 
32  void DefineParameters();
33  double LogAPrioriProbability(const std::vector<double> & parameters);
34  double LogLikelihood(const std::vector<double> &parameters);
35 
37  void Initialize( RooAbsData& data,
38  RooAbsPdf& model,
39  RooAbsPdf& prior,
40  const RooArgSet* params,
41  const RooArgSet& listPOI );
42 
43  void Initialize( const char* rootFile,
44  const char* wsName = "batWS",
45  const char* dataName = "data",
46  const char* modelName = "model",
47  const char* priorName = "priorPOI",
48  const char* priorNuisanceName= "priorNuisance",
49  const char* paramsName = "parameters",
50  const char* listPOIName = "POI" );
51 
53  void SetNumBins(const char * parname, int nbins);
55  void SetNumBins(int nbins);
61  RooStats::MarkovChain * GetRooStatsMarkovChain(){ return _roostatsMarkovChain;}
62  RooArgSet* GetArgSetForMarkovChain(){return &_parametersForMarkovChainCurrent;}
63  //RooArgSet GetArgSetForMarkovChainTest(){return _parametersForMarkovChain_test;}
64 
65  private:
66 
67  void AddToCurrentChainElement(double xij, int chainNum, int poiNum); //help function for construction of RooStats Markov Chain
68  bool EqualsLastChainElement(int chainNum); //help function for construction of RooStats Markov Chain
69  double GetWeightForChain(int chainNum); //help function for construction of RooStats Markov Chain
70 
71  RooAbsData* fData;
72  RooAbsPdf* fModel;
73  RooAbsReal* fNll;
74  RooArgSet* fObservables;
75  RooArgList* fParams;
76  RooArgList* fParamsPOI;
77  RooAbsPdf* fPrior;
78  int _default_nbins;
79 
80  RooRealVar* priorhelpvar;
81  bool _addeddummyprior;
82 
83  bool _fillChain;
84  bool fFirstComparison;
85  RooStats::MarkovChain * _roostatsMarkovChain;
86  RooArgSet _parametersForMarkovChainPrevious;
87  RooArgSet _parametersForMarkovChainCurrent;
88 
89  std::vector< std::vector<double> > fPreviousStep;
90  std::vector< std::vector<double> > fCurrentStep;
91  std::vector< double > fVecWeights;
92 
93  std::list< std::pair<const char*,int> > _nbins_list;
94 };
95 
96 #endif
void SetupRooStatsMarkovChain()
setup RooStats Markov Chain
void MCMCIterationInterface()
overloaded function from BCIntegrate to fill RooStats Markov Chain with every accepted step ...
void SetNumBins(const char *parname, int nbins)
set the number of histogram bins for a specific parameter
The base class for all user-defined models.
Definition: BCModel.h:50
void DefineParameters()
Overloaded methods.
double LogLikelihood(const std::vector< double > &parameters)
void Initialize(RooAbsData &data, RooAbsPdf &model, RooAbsPdf &prior, const RooArgSet *params, const RooArgSet &listPOI)
Other method of this class.
double LogAPrioriProbability(const std::vector< double > &parameters)
RooStats::MarkovChain * GetRooStatsMarkovChain()
return the RooStats Markov Chain (empty if corresponding constructor option not set) ...