CMS 3D CMS Logo

HistoChiSquare.h

Go to the documentation of this file.
00001 #ifndef PhysicsTools_Utilities_HistoChiSquare_h
00002 #define PhysicsTools_Utilities_HistoChiSquare_h
00003 #include <vector>
00004 #include "TH1.h"
00005 
00006 namespace fit {
00007   template<typename T>
00008   class HistoChiSquare {
00009    public:
00010     HistoChiSquare() { }
00011     HistoChiSquare(T & t, TH1 *histo, double rangeMin, double rangeMax): 
00012       t_(&t), rangeMin_(rangeMin), rangeMax_(rangeMax) { 
00013       nBins_ = histo->GetNbinsX();
00014       xMin_ = histo->GetXaxis()->GetXmin();
00015       xMax_ = histo->GetXaxis()->GetXmax();
00016       deltaX_ =(xMax_ - xMin_) / nBins_;
00017       for(size_t i = 0; i < nBins_; ++i) { 
00018         cont_.push_back( histo->GetBinContent(i+1) );
00019         err_.push_back( histo->GetBinError(i+1) );
00020       }
00021     }
00022     double operator()() const { 
00023       double chi2 = 0;
00024       for(size_t i = 0; i < nBins_; ++i) { 
00025         double x = xMin_ + ( i +.5 ) * deltaX_;
00026         if((x > rangeMin_)&&(x < rangeMax_)&&(err_[i] > 0)) { 
00027           double r = ( cont_[i] - (*t_)(x) )/err_[i];
00028           chi2 += (r * r);
00029         }
00030       }
00031      return chi2;
00032     }
00033     void setHistos(TH1 *histo) { 
00034       nBins_ = histo->GetNbinsX();
00035       xMin_ = histo->GetXaxis()->GetXmin();
00036       xMax_ = histo->GetXaxis()->GetXmax();
00037       deltaX_ =(xMax_ - xMin_) / nBins_;
00038     }
00039     size_t degreesOfFreedom() const {
00040       size_t fullBins = 0;
00041       for(size_t i = 0; i < nBins_; ++i) { 
00042         double x = xMin_ + ( i +.5 ) * deltaX_;
00043         if((x > rangeMin_)&&(x < rangeMax_)&&(err_[i] > 0))  
00044           fullBins++;
00045       }
00046       return fullBins; 
00047     }
00048     T & function() { return * t_; }
00049     const T & function() const { return * t_; }
00050   private:
00051     T * t_;
00052     double rangeMin_, rangeMax_;
00053     size_t nBins_;
00054     double xMin_, xMax_, deltaX_;
00055     std::vector<double> cont_;
00056     std::vector<double> err_;
00057   };
00058 }
00059 
00060 #endif

Generated on Tue Jun 9 17:42:42 2009 for CMSSW by  doxygen 1.5.4