00001 #ifndef __BCEFFICIENCYFITTER__H 00002 #define __BCEFFICIENCYFITTER__H 00003 00004 /*! 00005 * \class BCEfficiencyFitter 00006 * \brief A class for fitting histograms with functions 00007 * \author Daniel Kollar 00008 * \author Kevin Kröninger 00009 * \version 1.0 00010 * \date 11.2008 00011 * \detail This class allows fitting of efficiencies defined as 00012 * a ratio of two TH1D histograms using a TF1 function. It uses 00013 * binomial probabilities calculated based on the number of entries 00014 * in histograms. This is only applicable if the numerator is 00015 * a subset of the denominator. 00016 */ 00017 00018 /* 00019 * Copyright (C) 2008, Daniel Kollar and Kevin Kroeninger. 00020 * All rights reserved. 00021 * 00022 * For the licensing terms see doc/COPYING. 00023 */ 00024 00025 // --------------------------------------------------------- 00026 00027 #include <vector> 00028 00029 #include <BAT/BCModel.h> 00030 00031 // ROOT classes 00032 class TH1D; 00033 class TF1; 00034 class TGraphAsymmErrors; 00035 00036 // --------------------------------------------------------- 00037 00038 class BCEfficiencyFitter : public BCModel 00039 { 00040 public: 00041 00042 /** \name Constructors and destructors */ 00043 /* @{ */ 00044 00045 /** 00046 * The default constructor. */ 00047 BCEfficiencyFitter(); 00048 00049 /** 00050 * A constructor. 00051 * @param hist1 The histogram with the larger numbers 00052 * @param hist2 The histogram with the smaller numbers 00053 * @param func The fit function. */ 00054 BCEfficiencyFitter(TH1D * hist1, TH1D * hist2, TF1 * func); 00055 00056 /** 00057 * The default destructor. */ 00058 ~BCEfficiencyFitter(); 00059 00060 /* @} */ 00061 /** \name Member functions (get) */ 00062 /* @{ */ 00063 00064 /** 00065 * @return The histogram 1 */ 00066 TH1D * GetHistogram1() 00067 { return fHistogram1; }; 00068 00069 /** 00070 * @return The histogram 2 */ 00071 TH1D * GetHistogram2() 00072 { return fHistogram2; }; 00073 00074 /** 00075 * @return The histogram ratio */ 00076 TGraphAsymmErrors * GetHistogramRatio() 00077 { return fHistogramRatio; }; 00078 00079 /** 00080 * @return The fit function */ 00081 TF1 * GetFitFunction() 00082 { return fFitFunction; }; 00083 00084 /** 00085 * @return pointer to the error band */ 00086 TGraph * GetErrorBand() 00087 { return fErrorBand; }; 00088 00089 /** 00090 * @return pointer to a graph for the fit function */ 00091 TGraph * GetGraphFitFunction() 00092 { return fGraphFitFunction; }; 00093 00094 /** 00095 * Calculates the central value and the lower and upper limits for a given probability. 00096 * @param n n for the binomial. 00097 * @param k k for the binomial. 00098 * @param p The central probability defining the limits. 00099 * @param xmin The central value. 00100 * @param xmin The lower limit. 00101 * @param xmax The upper limit. 00102 * @return A flag (=1 plot point, !=1 do not plot point). */ 00103 int GetUncertainties(int n, int k, double p, double &xexp, double &xmin, double &xmax); 00104 00105 /* @} */ 00106 /** \name Member functions (set) */ 00107 /* @{ */ 00108 00109 /** 00110 * @param hist The histogram 1 00111 * @param hist The histogram 2 00112 * @ return An error code (1:pass, 0:fail). */ 00113 int SetHistograms(TH1D * hist1, TH1D * hist2); 00114 00115 /** 00116 * @param func The fit function 00117 * @ return An error code (1:pass, 0:fail). */ 00118 int SetFitFunction(TF1 * func); 00119 00120 /** 00121 * Sets the flag for integration. \n 00122 * true: use ROOT's TH1D::Integrate() \n 00123 * false: use linear interpolation */ 00124 void SetFlagIntegration(bool flag) 00125 { fFlagIntegration = flag; }; 00126 00127 /** Set type of point to be used to plot the efficiency data 00128 * 0 - mean + rms 00129 * 1 - mode + smallest 68% 00130 * 2 - median + central 68% */ 00131 void SetDataPointType(int type); 00132 00133 /* @} */ 00134 /** \name Member functions (miscellaneous methods) */ 00135 /* @{ */ 00136 00137 /** 00138 * The log of the prior probability. Overloaded from BCModel. 00139 * @param parameters A vector of doubles containing the parameter values. */ 00140 // virtual double LogAPrioriProbability(std::vector <double> parameters); 00141 00142 /** 00143 * The log of the conditional probability. Overloaded from BCModel. 00144 * @param parameters A vector of doubles containing the parameter values. */ 00145 virtual double LogLikelihood(std::vector <double> parameters); 00146 00147 /** 00148 * Returns the y-value of the 1-dimensional fit function at an x and 00149 * for a set of parameters. 00150 * @param x A vector with the x-value. 00151 * @param parameters A set of parameters. */ 00152 double FitFunction(std::vector <double> x, std::vector <double> parameters); 00153 00154 /** 00155 * Performs the fit. 00156 * @return An error code. */ 00157 int Fit(); 00158 00159 /** 00160 * Performs the fit. 00161 * @param hist1 The histogram with the larger number. 00162 * @param hist2 The histogram with the smaller number. 00163 * @param func The fit function. 00164 * @return An error code. */ 00165 int Fit(TH1D * hist1, TH1D * hist2, TF1 * func); 00166 00167 /** 00168 * Draw the fit in the current pad. */ 00169 void DrawFit(const char * options = "", bool flaglegend = false); 00170 00171 /** 00172 * Calculate the p-value using fast-MCMC. 00173 * @param par A set of parameter values 00174 * @param pvalue The pvalue 00175 * @return An error code */ 00176 int CalculatePValueFast(std::vector<double> par, double &pvalue); 00177 00178 /* @} */ 00179 00180 private: 00181 00182 /** 00183 * The histogram containing the larger numbers. */ 00184 TH1D * fHistogram1; 00185 00186 /** 00187 * The histogram containing the smaller numbers. */ 00188 TH1D * fHistogram2; 00189 00190 /** 00191 * The efficiency histogram. */ 00192 TGraphAsymmErrors * fHistogramRatio; 00193 00194 /** 00195 * The fit function */ 00196 TF1 * fFitFunction; 00197 00198 /** 00199 * Flag for using the ROOT TH1D::Integral method (true), or linear 00200 * interpolation (false) */ 00201 bool fFlagIntegration; 00202 00203 /** 00204 * Pointer to the error band (for legend) */ 00205 TGraph * fErrorBand; 00206 00207 /** 00208 * Pointer to a graph for displaying the fit function */ 00209 TGraph * fGraphFitFunction; 00210 00211 /** 00212 * Temporary histogram for calculating the binomial qunatiles */ 00213 TH1D * fHistogramBinomial; 00214 00215 /** Type of point to plot for efficiency data 00216 * 0 - mean + rms 00217 * 1 - mode + smallest 68% 00218 * 2 - median + central 68% */ 00219 unsigned int fDataPointType; 00220 00221 }; 00222 00223 // --------------------------------------------------------- 00224 00225 #endif