Go to the documentation of this file.00001
00002
00003
00004 #ifndef ROOSTATS_BATCalculator
00005 #define ROOSTATS_BATCalculator
00006
00007 #include <TNamed.h>
00008 #include <TH1D.h>
00009
00010 #include <RooStats/IntervalCalculator.h>
00011 #include <RooStats/SimpleInterval.h>
00012 #include <RooStats/ModelConfig.h>
00013 #include <RooStats/MCMCInterval.h>
00014
00015 #include <RooArgSet.h>
00016 #include <RooAbsData.h>
00017 #include <RooAbsPdf.h>
00018 #include <RooPlot.h>
00019 #include <RooAbsReal.h>
00020
00021 #include "BCRooInterface.h"
00022
00023 #include <map>
00024 #include <vector>
00025
00026
00027
00028 namespace RooStats
00029 {
00030
00031 class BATCalculator : public IntervalCalculator, public TNamed
00032 {
00033
00034 public:
00035
00036
00037 BATCalculator( );
00038
00039 BATCalculator( RooAbsData & data,
00040 RooAbsPdf & pdf,
00041 RooArgSet & POI,
00042 RooAbsPdf & prior,
00043 RooArgSet * params = 0,
00044 bool fillChain = false );
00045
00046
00047 BATCalculator( RooAbsData & data,
00048 ModelConfig & model,
00049 bool fillChain = false);
00050
00051
00052
00053 virtual ~BATCalculator();
00054
00055 RooPlot * GetPosteriorPlot1D() const;
00056
00057
00058 RooAbsPdf * GetPosteriorPdf1D() const;
00059 RooAbsPdf * GetPosteriorPdf1D(const char * POIname) const;
00060
00061
00062 virtual SimpleInterval * GetInterval1D() const;
00063 virtual SimpleInterval * GetInterval1D(const char * POIname) const;
00064
00065
00066 SimpleInterval * GetShortestInterval1D() const;
00067 SimpleInterval * GetShortestInterval1D(const char * POIname, bool & checkConnected) const;
00068
00069
00070 Double_t GetOneSidedUperLim();
00071
00072 virtual void SetData( RooAbsData & data )
00073 { fData = &data; ClearAll(); }
00074
00075 virtual void SetModel( const ModelConfig & model );
00076
00077
00078 virtual void SetTestSize( Double_t size )
00079 { fSize = size; fValidInterval = false; }
00080
00081
00082 void SetLeftSideTailFraction(Double_t leftSideFraction );
00083
00084
00085 virtual void SetConfidenceLevel( Double_t cl )
00086 { SetTestSize(1.-cl); }
00087
00088
00089 virtual Double_t Size() const
00090 { return fSize; }
00091
00092
00093 double GetLeftSideTailFraction()
00094 {return fLeftSideFraction;}
00095
00096
00097 virtual Double_t ConfidenceLevel() const
00098 { return 1.-fSize; }
00099
00100 void SetBrfPrecision( double precision )
00101 { fBrfPrecision = precision; }
00102
00103 double GetBrfPrecision()
00104 { return fBrfPrecision; }
00105
00106
00107 void SetnMCMC(int nMCMC)
00108 { _nMCMC = nMCMC; }
00109
00110
00111 int GetnMCMC()
00112 { return _nMCMC; }
00113
00114 BCRooInterface * GetBCRooInterface() const
00115 { return _myRooInterface; }
00116
00117 RooStats::MarkovChain * GetRooStatsMarkovChain() const
00118 { return _myRooInterface->GetRooStatsMarkovChain();}
00119
00120
00121 virtual MCMCInterval* GetInterval() const;
00122
00123
00124 bool GetConnected()
00125 { return fConnectedInterval; }
00126
00127
00128 vector<double> GetIntervalBorders1D()
00129 { return _intervalBorders1D; }
00130
00131
00132
00133
00134 void SetNumBins(const char * parname, int nbins);
00135
00136
00137 void SetNumBins(int nbins);
00138
00139
00140
00141
00142 void CleanCalculatorForNewData()
00143 { ClearAll(); }
00144
00145 protected:
00146
00147 void ClearAll() const;
00148
00149 private:
00150
00151
00152
00153 RooArgSet * GetMode( RooArgSet * parameters ) const;
00154
00155
00156
00157 RooAbsData * fData;
00158 RooAbsPdf * fPdf;
00159 const RooArgSet fPOI;
00160 RooAbsPdf * fPrior;
00161 const RooArgSet * fparams;
00162 BCRooInterface * _myRooInterface;
00163 mutable TH1D * _posteriorTH1D;
00164
00165
00166 mutable RooAbsPdf * fProductPdf;
00167 mutable RooAbsReal * fLogLike;
00168 mutable RooAbsReal * fLikelihood;
00169 mutable RooAbsReal * fIntegratedLikelihood;
00170 mutable RooAbsPdf * fPosteriorPdf;
00171 mutable Double_t fLower;
00172 mutable Double_t fUpper;
00173 double fBrfPrecision;
00174 mutable Bool_t fValidInterval;
00175 mutable bool fConnectedInterval;
00176
00177 int _nMCMC;
00178 double fSize;
00179 double fLeftSideFraction;
00180 mutable vector<double> _intervalBorders1D;
00181
00182 protected:
00183
00184 ClassDef(BATCalculator,1)
00185
00186 };
00187
00188 bool sortbyposterior(pair< Int_t,Double_t > pair1, pair< Int_t,Double_t > pair2);
00189 }
00190
00191 #endif