• Main Page
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

BATCalculator.h

Go to the documentation of this file.
00001 //
00002 // Author: Kyle Cranmer, Lorenzo Moneta, Gregory Schott, Wouter Verkerke, Stefan A. Schmitz
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       // constructor
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       // destructor
00053       virtual ~BATCalculator();
00054 
00055       RooPlot * GetPosteriorPlot1D() const;
00056 
00057       // return posterior pdf 
00058       RooAbsPdf * GetPosteriorPdf1D() const;
00059       RooAbsPdf * GetPosteriorPdf1D(const char * POIname) const;
00060 
00061       // return SimpleInterval object: central interval (one poi only)
00062       virtual SimpleInterval * GetInterval1D() const;
00063       virtual SimpleInterval * GetInterval1D(const char * POIname) const;
00064       //virtual SimpleInterval * GetInterval1Dv0(const char * POIname) const;
00065       // return SimpleInterval object: shortest interval (not necessarily connected, one poi only)
00066       SimpleInterval * GetShortestInterval1D() const;
00067       SimpleInterval * GetShortestInterval1D(const char * POIname, bool & checkConnected) const;
00068 
00069       // temporary solution ?
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       // set the size of the test (rate of Type I error) ( Eg. 0.05 for a 95% Confidence Interval)
00078       virtual void SetTestSize( Double_t size )
00079          {  fSize = size; fValidInterval = false; }
00080 
00081       // set left side tail fraction (only for 1D interval, not meaningful for shortest interval)
00082       void SetLeftSideTailFraction(Double_t leftSideFraction );
00083 
00084       // set the confidence level for the interval (eg. 0.95 for a 95% Confidence Interval)
00085       virtual void SetConfidenceLevel( Double_t cl )
00086          { SetTestSize(1.-cl); }
00087 
00088       // Get the size of the test
00089       virtual Double_t Size() const
00090          { return fSize; }
00091       // Get left side tail fraction (only for 1D interval, not meaningful for shortest interval)
00092       double GetLeftSideTailFraction()
00093          {return fLeftSideFraction;} 
00094 
00095       // Get the Confidence level for the test
00096       virtual Double_t ConfidenceLevel() const
00097          { return 1.-fSize; }
00098 
00099       void SetBrfPrecision( double precision )
00100          { fBrfPrecision = precision; }
00101 
00102       double GetBrfPrecision()
00103          { return fBrfPrecision; }
00104 
00105       //set number of iterations per chain
00106       void SetnMCMC(int nMCMC)
00107          { _nMCMC = nMCMC; }
00108 
00109       //return number of iterations per chain
00110       int  GetnMCMC()
00111          { return _nMCMC; }
00112 
00113       BCRooInterface * GetBCRooInterface() const
00114          { return _myRooInterface; }
00115 
00116       RooStats::MarkovChain * GetRooStatsMarkovChain() const
00117          { return _myRooInterface->GetRooStatsMarkovChain();}
00118 
00119       // returns MCMCInterval object
00120       virtual MCMCInterval* GetInterval() const;
00121 
00122       //returns if last calculated shortest interval is connected (1 poi only)
00123       bool GetConnected()
00124          { return fConnectedInterval; }
00125 
00126       // returns  interval borders of shortest interval (1 poi only)
00127       std::vector<double> GetIntervalBorders1D()
00128          { return _intervalBorders1D; }
00129 
00130       //returns interval borders od the last calculated shortest interval
00131 
00132       //set the number of histogram bins for a specific parameter
00133       void SetNumBins(const char * parname, int nbins);
00134 
00135       //set the number of histogram bins for all parameters
00136       void SetNumBins(int nbins);
00137 
00138       // would be more complete if we had this -> ask BAT developers to implement this functionality (not high priority)
00139       //int GetNbins(const char * parname);
00140 
00141       void CleanCalculatorForNewData()
00142          { ClearAll(); }
00143 
00144    protected:
00145 
00146       void ClearAll() const;
00147 
00148    private:
00149 
00150       // compute the most probable value: move to public once implemented
00151       // returns a RooArgSet
00152       RooArgSet * GetMode( RooArgSet * parameters ) const;
00153       // plan to replace the above: return a SimpleInterval integrating
00154       // over all other parameters except the one specified as argument
00155 
00156       RooAbsData * fData;
00157       RooAbsPdf * fPdf;
00158       const RooArgSet fPOI;
00159       RooAbsPdf * fPrior;
00160       const RooArgSet * fparams;
00161       BCRooInterface * _myRooInterface;
00162       mutable TH1D * _posteriorTH1D;
00163 
00164 
00165       mutable RooAbsPdf * fProductPdf;
00166       mutable RooAbsReal * fLogLike;
00167       mutable RooAbsReal * fLikelihood;
00168       mutable RooAbsReal * fIntegratedLikelihood;
00169       mutable RooAbsPdf * fPosteriorPdf;
00170       mutable Double_t fLower;
00171       mutable Double_t fUpper;
00172       double fBrfPrecision;
00173       mutable Bool_t fValidInterval;
00174       mutable Bool_t fValidMCMCInterval;
00175       mutable bool fConnectedInterval;
00176 
00177       int _nMCMC; // number of chain elements per Markov Chain
00178       double fSize;  // size used for the interval object
00179       double fLeftSideFraction; //
00180       mutable std::vector<double> _intervalBorders1D;
00181 
00182    protected:
00183 
00184       ClassDef(BATCalculator,1)  // BATCalculator class
00185 
00186    };
00187 
00188    bool sortbyposterior(std::pair< Int_t,Double_t > pair1, std::pair< Int_t,Double_t > pair2); //help function for calculating the shortest interval
00189 }
00190 
00191 #endif

Generated by  doxygen 1.7.1