CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/CalibTracker/SiStripLorentzAngle/src/LA_Fitter.cc

Go to the documentation of this file.
00001 #include "CalibTracker/SiStripLorentzAngle/interface/LA_Filler_Fitter.h"
00002 #include "FWCore/Utilities/interface/EDMException.h"
00003 
00004 #include <cmath>
00005 #include <boost/foreach.hpp>
00006 #include <boost/algorithm/string/erase.hpp>
00007 #include <TF1.h>
00008 
00009 void LA_Filler_Fitter::
00010 fit_width_profile(Book& book) {
00011   for(Book::iterator it = book.begin(".*"+method(WIDTH)); it!=book.end(); ++it) {
00012     it->second->SetTitle("Mean Cluster Width;tan#theta_{t}");
00013     TH1* const p = it->second;
00014     if(p->GetEntries() < 400) { delete p; book[it->first]=0; continue;}
00015     p->SetTitle(";tan#theta_{t};");
00016     const float min = p->GetMinimum();
00017     const float max = p->GetMaximum();
00018     float xofmin = p->GetBinCenter(p->GetMinimumBin()); if( xofmin>0.0 || xofmin<-0.15) xofmin = -0.05;
00019     const float xofmax = p->GetBinCenter(p->GetMaximumBin());
00020 
00021     TF1* const fit = new TF1("LA_profile_fit","[2]*(TMath::Abs(x-[0]))+[1]",-1,1);
00022     fit->SetParLimits(0,-0.15,0.01);
00023     fit->SetParLimits(1, 0.6*min, 1.25*min );
00024     fit->SetParLimits(2,0.1,10);
00025     fit->SetParameters( xofmin, min, (max-min) / fabs( xofmax - xofmin ) );
00026 
00027     int badfit = p->Fit(fit,"IEQ","",-.5,.3);
00028     if( badfit ) badfit = p->Fit(fit,"IEQ","", -.46,.26);
00029     if( badfit ) { book.erase(it); }
00030   }
00031 }
00032 
00033 void LA_Filler_Fitter::
00034 make_and_fit_symmchi2(Book& book) {
00035   for(Book::iterator it = book.begin(".*_all"); it!=book.end(); ++it) {
00036     const std::string base = boost::erase_all_copy(it->first,"_all");
00037 
00038     std::vector<Book::iterator> rebin_hists;              
00039     Book::iterator    all = it;                              rebin_hists.push_back(all);   
00040     Book::iterator     w1 = book.find(base+"_w1");           rebin_hists.push_back(w1);    
00041     Book::iterator var_w2 = book.find(base+method(AVGV2,0)); rebin_hists.push_back(var_w2);
00042     Book::iterator var_w3 = book.find(base+method(AVGV3,0)); rebin_hists.push_back(var_w3);
00043 
00044     const unsigned rebin = std::max( var_w2==book.end() ? 0 : find_rebin(var_w2->second), 
00045                                      var_w3==book.end() ? 0 : find_rebin(var_w3->second) );
00046     BOOST_FOREACH(Book::iterator it, rebin_hists) if(it!=book.end()) it->second->Rebin( rebin>1 ? rebin<7 ? rebin : 6 : 1);
00047 
00048     TH1* const prob_w1 = w1==book.end()     ? 0 : subset_probability( base+method(PROB1,0) ,w1->second,all->second);
00049     TH1* const rmsv_w2 = var_w2==book.end() ? 0 :        rms_profile( base+method(RMSV2,0), (TProfile*const)var_w2->second);
00050     TH1* const rmsv_w3 = var_w3==book.end() ? 0 :        rms_profile( base+method(RMSV3,0), (TProfile*const)var_w3->second);
00051     
00052     std::vector<TH1*> fit_hists;
00053     if(prob_w1) {
00054       book.book(base+method(PROB1,0),prob_w1);
00055       fit_hists.push_back(prob_w1);  prob_w1->SetTitle("Width==1 Probability;tan#theta_{t}-(dx/dz)_{reco}");
00056     }
00057     if(var_w2!=book.end())  {
00058       book.book(base+method(RMSV2,0),rmsv_w2);
00059       fit_hists.push_back(var_w2->second);   var_w2->second->SetTitle("Width==2 Mean Variance;tan#theta_{t}-(dx/dz)_{reco}");
00060       fit_hists.push_back(rmsv_w2);                 rmsv_w2->SetTitle("Width==2 RMS Variance;tan#theta_{t}-(dx/dz)_{reco}");
00061     }
00062     if(var_w3!=book.end())  {
00063       book.book(base+method(RMSV3,0),rmsv_w3);
00064       fit_hists.push_back(var_w3->second);   var_w3->second->SetTitle("Width==3 Mean Variance;tan#theta_{t}-(dx/dz)_{reco}");
00065       fit_hists.push_back(rmsv_w3);                 rmsv_w3->SetTitle("Width==3 RMS Variance;tan#theta_{t}-(dx/dz)_{reco}");
00066     }
00067 
00068     if(!fit_hists.size()) continue;
00069     const unsigned bins = fit_hists[0]->GetNbinsX();
00070     const unsigned guess = fit_hists[0]->FindBin(0);
00071     const std::pair<unsigned,unsigned> range(guess-bins/30,guess+bins/30);
00072 
00073     BOOST_FOREACH(TH1*const hist, fit_hists) {
00074       TH1*const chi2 = SymmetryFit::symmetryChi2(hist,range);
00075       if(chi2) {book.book(chi2->GetName(),chi2); chi2->SetTitle("Symmetry #chi^{2};tan#theta_{t}-(dx/dz)_{reco}");}
00076     }
00077   }
00078 }
00079 
00080 unsigned LA_Filler_Fitter::
00081 find_rebin(const TH1* const hist) {
00082   const double mean = hist->GetMean();
00083   const double rms = hist->GetRMS();
00084   const int begin = std::min(                1, hist->GetXaxis()->FindFixBin(mean-rms));
00085   const int end   = std::max(hist->GetNbinsX(), hist->GetXaxis()->FindFixBin(mean+rms)) + 1;
00086   unsigned current_hole(0), max_hole(0);
00087   for(int i=begin; i<end; i++) {
00088     if(!hist->GetBinError(i)) current_hole++;
00089     else if(current_hole) {max_hole = std::max(current_hole,max_hole); current_hole=0;}
00090   }
00091   return max_hole+1;
00092 }
00093 
00094 TH1* LA_Filler_Fitter::
00095 rms_profile(const std::string name, const TProfile* const prof) {
00096   const int bins = prof->GetNbinsX();
00097   TH1* const rms = new TH1F(name.c_str(),"",bins, prof->GetBinLowEdge(1),  prof->GetBinLowEdge(bins) + prof->GetBinWidth(bins) );
00098   for(int i = 1; i<=bins; i++) {
00099     const double Me = prof->GetBinError(i);
00100     const double neff = prof->GetBinEntries(i); //Should be prof->GetBinEffectiveEntries(i);, not availible this version ROOT.  This is only ok for unweighted fills
00101     rms->SetBinContent(i, Me*sqrt(neff) );
00102     rms->SetBinError(i, Me/sqrt(2) );
00103   }
00104   return rms;
00105 }
00106 
00107 TH1* LA_Filler_Fitter::
00108 subset_probability(const std::string name, const TH1* const subset, const TH1* const total) {
00109   const int bins = subset->GetNbinsX();
00110   TH1* const prob = new TH1F(name.c_str(),"",bins, subset->GetBinLowEdge(1),  subset->GetBinLowEdge(bins) + subset->GetBinWidth(bins) );
00111   for(int i = 1; i<=bins; i++) {
00112     const double s = subset->GetBinContent(i);
00113     const double T = total->GetBinContent(i);
00114     const double B = T-s;
00115 
00116     const double p = T? s/T : 0;
00117     const double perr = T? ( (s&&B)? sqrt(s*s*B+B*B*s)/(T*T) : 1/T ) : 0;
00118 
00119     prob->SetBinContent(i,p);
00120     prob->SetBinError(i,perr);
00121   }  
00122   return prob;
00123 }