CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/RecoParticleFlow/PFRootEvent/src/ResidualFitter.cc

Go to the documentation of this file.
00001 #include "RecoParticleFlow/PFRootEvent/interface/ResidualFitter.h"
00002 
00003 #include <iostream>
00004 #include "TCanvas.h"
00005 #include "TH2.h"
00006 #include "TF1.h"
00007 #include "TDirectory.h"
00008 
00009 #include <string>
00010 
00011 using namespace std;
00012 
00013 // TCanvas* ResidualFitter::canvasFit_ = 0;
00014 
00015 //ClassImp(ResidualFitter)
00016 
00017 int ResidualFitter::xCanvas_ = 600;
00018 int ResidualFitter::yCanvas_ = 600;
00019 
00020 
00021 ResidualFitter::ResidualFitter(const char* name, 
00022                                const char* title, 
00023                                int nbinsx, double xlow, double xup, 
00024                                int nbinsy, double ylow, double yup, 
00025                                int nbinsz, double zlow, double zup) 
00026   : TH3D( name, title, 
00027           nbinsx, xlow, xup, 
00028           nbinsy, ylow, yup, 
00029           nbinsz, zlow, zup ), 
00030     fitFunction_( new TF1("gaus", "gaus") ), 
00031     canvasFit_(0),
00032     curBin_(0), 
00033     autoRangeN_(0),
00034     minN_(5) {
00035 
00036   cout<<"creating residual fitter with name "<<name<<endl;
00037   
00038   string meanname = name; meanname += "_mean";
00039   mean_ = new TH2D(meanname.c_str(), meanname.c_str(),
00040                    nbinsx, xlow, xup, 
00041                    nbinsy, ylow, yup);
00042   mean_->SetStats(0);
00043   
00044   string sigmaname = name; sigmaname += "_sigma";
00045   sigma_ = new TH2D(sigmaname.c_str(), sigmaname.c_str(),
00046                     nbinsx, xlow, xup, 
00047                     nbinsy, ylow, yup);
00048 
00049   sigma_->SetStats(0);
00050 
00051   string chi2name = name; chi2name += "_chi2";
00052   chi2_ = new TH2D(chi2name.c_str(), chi2name.c_str(),
00053                    nbinsx, xlow, xup, 
00054                    nbinsy, ylow, yup);
00055   chi2_->SetStats(0);
00056 
00057   //   string nseenname = name; nseenname += "_nseen";
00058   //   nseen_ = new TH2D(nseenname.c_str(), nseenname.c_str(),
00059   //               nbinsx, xlow, xup, 
00060   //               nbinsy, ylow, yup);
00061   //   nseen_->SetStats(0);
00062 
00063   gDirectory->ls();
00064 
00065   CreateCanvas();
00066 }
00067 
00068 ResidualFitter::~ResidualFitter() {
00069   delete fitFunction_;
00070   if(canvas_) delete canvas_;
00071   if(curBin_) delete curBin_;
00072 
00073   delete mean_;
00074   delete sigma_;
00075   delete chi2_;
00076   //   delete nseen_;
00077 }
00078 
00079 void ResidualFitter::CreateCanvas() {
00080   string cname = "ResidualFitterCanvas_"; cname += GetName();
00081   canvas_ = new TCanvas(cname.c_str(), cname.c_str(),xCanvas_, yCanvas_);
00082   
00083   canvas_ ->Connect("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)", 
00084                     "ResidualFitter",
00085                     this, "ExecuteEvent(Int_t,Int_t,Int_t,TObject*)");
00086 
00087   
00088 }
00089 
00090 
00091 void ResidualFitter::Fill() {
00092 
00093   Reset();
00094   
00095   for(unsigned it = 0; it<true_.size(); it++) {
00096     for(unsigned im = 0; im<meas_.size(); im++) {
00097       TH3D::Fill(true_[it].x_, 
00098                  true_[it].y_, 
00099                  meas_[it].z_ - true_[it].z_ );
00100     }
00101   }
00102 }
00103 
00104 void ResidualFitter::FitSlicesZ(TF1 *f1) {
00105 
00106   mean_->Reset();
00107   sigma_->Reset();
00108   chi2_->Reset();
00109   //   nseen_->Reset();
00110 
00111   cout<<"ResidualFitter::FitSlicesZ"<<endl;
00112   if(f1) SetFitFunction(f1);
00113   
00114   for(int ix=1; ix<=GetNbinsX(); ix++) {
00115     for(int iy=1; iy<=GetNbinsY(); iy++) {
00116       Fit(ix, iy, "Q0");
00117     }
00118   }
00119 
00120   //  TH3D::FitSlicesZ(f1);
00121 }
00122   
00123 
00124 void ResidualFitter::ExecuteEvent(Int_t event, Int_t px, Int_t py, TObject *sel) {
00125 
00126   if( event != kButton1Down ) return;
00127 
00128   TH2* histo2d = dynamic_cast<TH2*>(sel);
00129   if(!histo2d) return;
00130 
00131   float x = canvas_->AbsPixeltoX(px);
00132   float y = canvas_->AbsPixeltoY(py);
00133   x = canvas_->PadtoX(x);
00134   y = canvas_->PadtoY(y);
00135 
00136   ShowFit(histo2d, x, y);
00137 }
00138 
00139 
00140 void ResidualFitter::ShowFit(TH2* histo2d, double x, double y) {
00141   
00142   if(!canvasFit_) {
00143     string cname = "ResidualFitterCanvasFit_"; cname += GetName();
00144     canvasFit_ = new TCanvas(cname.c_str(), cname.c_str(),300,300);
00145   }
00146   canvasFit_ ->cd();
00147   
00148   int binx = histo2d->GetXaxis()->FindBin(x); 
00149   int biny = histo2d->GetYaxis()->FindBin(y); 
00150 
00151   if(binx == oldBinx_ && biny == oldBiny_ ) return;
00152   oldBinx_ = binx;
00153   oldBiny_ = biny;
00154 
00155   Fit(binx, biny);
00156  
00157   canvasFit_->Modified();
00158   canvasFit_->Update();
00159 
00160   canvas_->Modified();
00161   canvas_->Update();
00162   canvas_->cd();
00163 }
00164 
00165 
00166 void ResidualFitter::Fit(int binx, int biny, const char* opt) {
00167   TH1::AddDirectory(0);
00168   
00169   if(curBin_) delete curBin_;
00170   curBin_ = TH3::ProjectionZ("", binx, binx, biny, biny);  
00171 
00172   if(curBin_->GetEntries() < minN_ ) {
00173     TH1::AddDirectory(1);
00174     return;
00175   }
00176 
00177   string sopt = fitOptions_; sopt += opt;
00178 
00179   if( autoRangeN_ ) {
00180     double maxpos = curBin_->GetBinCenter( curBin_->GetMaximumBin() ); 
00181     
00182     double minrange = maxpos-curBin_->GetRMS()* autoRangeN_;
00183     double maxrange = maxpos+curBin_->GetRMS()* autoRangeN_;
00184 
00185     fitFunction_->SetRange( minrange, maxrange );
00186     cout<<"range : "<<minrange<<" "<<maxrange<<endl;
00187   }
00188 
00189   curBin_->Fit(fitFunction_, sopt.c_str() );
00190 
00191 
00192   double chi2overndf=0;
00193   if(fitFunction_->GetNDF() ) {
00194     chi2overndf = fitFunction_->GetChisquare()/ fitFunction_->GetNDF();
00195     mean_->SetBinContent(binx,biny, fitFunction_->GetParameter(1) );
00196     mean_->SetBinError(binx,biny, fitFunction_->GetParError(1) );
00197     sigma_->SetBinContent(binx,biny, fitFunction_->GetParameter(2) );
00198     sigma_->SetBinError(binx,biny, fitFunction_->GetParError(2) );
00199  
00200     chi2_->SetBinContent(binx,biny,chi2overndf);
00201   }
00202   //   nseen_->SetBinContent(binx, biny, 
00203   //                    fitFunction_->Integral( fitFunction_->GetXmin(), 
00204   //                                            fitFunction_->GetXmax())
00205   //                    /curBin_->GetBinWidth(1) );
00206 
00207   TH1::AddDirectory(1);
00208 }
00209 
00210