CMS 3D CMS Logo

DTMeanTimerFitter Class Reference

Fit the Tmax histograms with a gaussian returning the mean values and the sigmas. More...

#include <CalibMuon/DTCalibration/interface/DTMeanTimerFitter.h>

List of all members.

Public Member Functions

 DTMeanTimerFitter (TFile *file)
 Constructor.
std::vector< float > evaluateVDriftAndReso (TString N)
 Fit the TMax histos and evaluate VDrift and resolution.
TF1 * fitTMax (TH1F *histo)
 Really do the fit.
void setVerbosity (unsigned int lvl)
 Set the verbosity of the output: 0 = silent, 1 = info, 2 = debug.
virtual ~DTMeanTimerFitter ()
 Destructor.

Private Attributes

TFile * hDebugFile
TFile * hInputFile
unsigned int theVerbosityLevel


Detailed Description

Fit the Tmax histograms with a gaussian returning the mean values and the sigmas.

Date
2007/01/22 11:10:27
Revision
1.2
Author:
S. Bolognesi - INFN Torino

Definition at line 21 of file DTMeanTimerFitter.h.


Constructor & Destructor Documentation

DTMeanTimerFitter::DTMeanTimerFitter ( TFile *  file  ) 

Constructor.

Definition at line 21 of file DTMeanTimerFitter.cc.

References hDebugFile.

00021                                                 : hInputFile(file), theVerbosityLevel(0) {
00022   //hInputFile = new TFile("DTTMaxHistos.root", "READ");
00023   hDebugFile = new TFile("DTMeanTimerFitter.root", "RECREATE");
00024 }

DTMeanTimerFitter::~DTMeanTimerFitter (  )  [virtual]

Destructor.

Definition at line 26 of file DTMeanTimerFitter.cc.

References hDebugFile.

00026                                       {
00027   hDebugFile->Close();
00028 }


Member Function Documentation

vector< float > DTMeanTimerFitter::evaluateVDriftAndReso ( TString  N  ) 

Fit the TMax histos and evaluate VDrift and resolution.

Definition at line 30 of file DTMeanTimerFitter.cc.

References count, counter(), GenMuonPlsPt100GeV_cfg::cout, lat::endl(), f1, fitTMax(), hTMaxCell::GetT0Factor(), hDebugFile, hInputFile, combine::histos, hTMaxCell::hTmax123, hTMaxCell::hTmax124s72, hTMaxCell::hTmax124s78, hTMaxCell::hTmax134s72, hTMaxCell::hTmax134s78, hTMaxCell::hTmax234, hTMaxCell::hTmax_0, hTMaxCell::hTmax_2t0, hTMaxCell::hTmax_3t0, hTMaxCell::hTmax_t0, i, iter, mean(), funct::sqrt(), and theVerbosityLevel.

Referenced by DTVDriftWriter::analyze(), and DTVDriftCalibration::endJob().

00030                                                                  {
00031   
00032   // Retrieve histogram sets
00033   hTMaxCell * histos   = new hTMaxCell(N, hInputFile);
00034   vector<float> vDriftAndReso;
00035 
00036   // Check that the histo for this cell exists
00037   if(histos->hTmax123 != 0) {
00038     vector<TH1F*> hTMax;  // histograms for <T_max> calculation
00039     vector <TH1F*> hT0;   // histograms for T0 evaluation
00040     hTMax.push_back(histos->hTmax123); 
00041     hTMax.push_back(histos->hTmax124s72);
00042     hTMax.push_back(histos->hTmax124s78);
00043     hTMax.push_back(histos->hTmax134s72);
00044     hTMax.push_back(histos->hTmax134s78);
00045     hTMax.push_back(histos->hTmax234);
00046 
00047     hT0.push_back(histos->hTmax_3t0);
00048     hT0.push_back(histos->hTmax_2t0);
00049     hT0.push_back(histos->hTmax_t0);
00050     hT0.push_back(histos->hTmax_0);
00051 
00052     vector<Double_t> factor; // factor relating the width of the Tmax distribution 
00053                              // and the cell resolution 
00054     factor.push_back(sqrt(2./3.)); // hTmax123
00055     factor.push_back(sqrt(2./7.)); // hTmax124s72
00056     factor.push_back(sqrt(8./7.)); // hTmax124s78
00057     factor.push_back(sqrt(2./7.)); // hTmax134s72
00058     factor.push_back(sqrt(8./7.)); // hTmax134s78
00059     factor.push_back(sqrt(2./3.)); // hTmax234
00060 
00061 
00062     // Retrieve the gaussian mean and sigma for each TMax histogram    
00063     vector<Double_t> mean;
00064     vector<Double_t> sigma; 
00065     vector<Double_t> count;  //number of entries
00066 
00067     for(vector<TH1F*>::const_iterator ith = hTMax.begin();
00068         ith != hTMax.end(); ith++) {
00069      
00070       TF1 *funct = fitTMax(*ith);
00071       hDebugFile->cd();
00072       (*ith)->Write();
00073 
00074     // Get mean, sigma and number of entries of each histogram
00075       mean.push_back(funct->GetParameter(1));
00076       sigma.push_back(funct->GetParameter(2)); 
00077       count.push_back((*ith)->GetEntries());  
00078     } 
00079           
00080     Double_t tMaxMean=0.;
00081     Double_t wTMaxSum=0.;
00082     Double_t sigmaT=0.;
00083     Double_t wSigmaSum = 0.;
00084   
00085     //calculate total mean and sigma
00086     for(int i=0; i<=5; i++) {
00087       if(count[i]<200) continue;
00088       tMaxMean  += mean[i]*(count[i]/(sigma[i]*sigma[i]));
00089       wTMaxSum  += count[i]/(sigma[i]*sigma[i]);
00090       sigmaT    += count[i]*factor[i]*sigma[i];
00091       wSigmaSum += count[i];
00092       // cout << "TMaxMean "<<i<<": "<< mean[i] << " entries: " << count[i] 
00093       // << " sigma: " << sigma[i] 
00094       // << " weight: " << (count[i]/(sigma[i]*sigma[i])) << endl; 
00095     }
00096     tMaxMean /= wTMaxSum;
00097     sigmaT /= wSigmaSum;
00098 
00099     //calculate v_drift and resolution
00100     Double_t vDrift = 2.1 / tMaxMean; //2.1 is the half cell length in cm
00101     Double_t reso = vDrift * sigmaT;
00102     vDriftAndReso.push_back(vDrift);
00103     vDriftAndReso.push_back(reso);
00104     if(theVerbosityLevel >= 1)
00105       cout << " final TMaxMean=" << tMaxMean << " sigma= "  << sigmaT 
00106            << " v_d and reso: " << vDrift << " " << reso << endl;
00107 
00108     // Order t0 histogram by number of entries (choose histograms with higher nr. of entries)
00109     map<Double_t,TH1F*> hEntries;    
00110     for(vector<TH1F*>::const_iterator ith = hT0.begin();
00111         ith != hT0.end(); ith++) {
00112       hEntries[(*ith)->GetEntries()] = (*ith);
00113     } 
00114 
00115     // add at the end of hT0 the two hists with the higher number of entries 
00116     int counter = 0;
00117     for(map<Double_t,TH1F*>::reverse_iterator iter = hEntries.rbegin();
00118         iter != hEntries.rend(); iter++) {
00119       counter++;
00120       if (counter==1) hT0.push_back(iter->second); 
00121       else if (counter==2) {hT0.push_back(iter->second); break;} 
00122     }
00123     
00124     // Retrieve the gaussian mean and sigma of histograms for Delta(t0) evaluation   
00125     vector<Double_t> meanT0;
00126     vector<Double_t> sigmaT0; 
00127     vector<Double_t> countT0;
00128 
00129     for(vector<TH1F*>::const_iterator ith = hT0.begin();
00130         ith != hT0.end(); ith++) {
00131       (*ith)->Fit("gaus");
00132       TF1 *f1 = (*ith)->GetFunction("gaus");
00133       // Get mean, sigma and number of entries of the  histograms
00134       meanT0.push_back(f1->GetParameter(1));
00135       sigmaT0.push_back(f1->GetParameter(2));
00136       countT0.push_back((*ith)->GetEntries());
00137     }
00138     //calculate Delta(t0)
00139     if(hT0.size() != 6) { // check if you have all the t0 hists
00140       cout << "t0 histograms = " << hT0.size() << endl;
00141       for(int i=1; i<=4;i++) {
00142         vDriftAndReso.push_back(-1);
00143       }
00144       return vDriftAndReso;
00145     }
00146     
00147     for(int it0=0; it0<=2; it0++) {      
00148       if((countT0[it0] > 200) && (countT0[it0+1] > 200)) {
00149         Double_t deltaT0 = meanT0[it0] - meanT0[it0+1]; 
00150         vDriftAndReso.push_back(deltaT0);
00151       }  
00152       else
00153         vDriftAndReso.push_back(999.);
00154     }
00155     //deltat0 using hists with max nr. of entries
00156     if((countT0[4] > 200) && (countT0[5] > 200)) {
00157       Double_t t0Diff = histos->GetT0Factor(hT0[4]) - histos->GetT0Factor(hT0[5]);
00158       Double_t deltaT0MaxEntries =  (meanT0[4] - meanT0[5])/ t0Diff;
00159       vDriftAndReso.push_back(deltaT0MaxEntries);
00160     }
00161     else
00162       vDriftAndReso.push_back(999.);
00163   }
00164   else {
00165     for(int i=1; i<=6; i++) { 
00166       // 0=vdrift, 1=reso,  2=(3deltat0-2deltat0), 3=(2deltat0-1deltat0), 
00167       // 4=(1deltat0-0deltat0), 5=deltat0 from hists with max entries,
00168       vDriftAndReso.push_back(-1);
00169     }
00170   }
00171   return vDriftAndReso;
00172 }

TF1 * DTMeanTimerFitter::fitTMax ( TH1F *  histo  ) 

Really do the fit.

Definition at line 176 of file DTMeanTimerFitter.cc.

References GenMuonPlsPt100GeV_cfg::cout, lat::endl(), f1, range, and theVerbosityLevel.

Referenced by evaluateVDriftAndReso().

00176                                           {
00177  // Find distribution peak and fit range
00178       Double_t peak = (((((histo->GetXaxis())->GetXmax())-((histo->GetXaxis())->GetXmin()))/histo->GetNbinsX())*
00179                        (histo->GetMaximumBin()))+((histo->GetXaxis())->GetXmin());
00180       if(theVerbosityLevel >= 1)
00181         cout<<"Peak "<<peak<<" : "<<"xmax "<<((histo->GetXaxis())->GetXmax())
00182             <<"            xmin "<<((histo->GetXaxis())->GetXmin())
00183             <<"            nbin "<<histo->GetNbinsX()
00184             <<"            bin with max "<<(histo->GetMaximumBin())<<endl;
00185       Double_t range = 2.*histo->GetRMS(); 
00186 
00187       // Fit each Tmax histogram with a Gaussian in a restricted interval
00188       TF1 *rGaus = new TF1("rGaus","gaus",peak-range,peak+range);
00189       histo->Fit("rGaus","R");
00190       TF1 *f1 = histo->GetFunction("rGaus");
00191       return f1;
00192  }

void DTMeanTimerFitter::setVerbosity ( unsigned int  lvl  )  [inline]

Set the verbosity of the output: 0 = silent, 1 = info, 2 = debug.

Definition at line 33 of file DTMeanTimerFitter.h.

References theVerbosityLevel.

00033                                       {
00034     theVerbosityLevel = lvl;
00035   }


Member Data Documentation

TFile* DTMeanTimerFitter::hDebugFile [private]

Definition at line 43 of file DTMeanTimerFitter.h.

Referenced by DTMeanTimerFitter(), evaluateVDriftAndReso(), and ~DTMeanTimerFitter().

TFile* DTMeanTimerFitter::hInputFile [private]

Definition at line 44 of file DTMeanTimerFitter.h.

Referenced by evaluateVDriftAndReso().

unsigned int DTMeanTimerFitter::theVerbosityLevel [private]

Definition at line 46 of file DTMeanTimerFitter.h.

Referenced by evaluateVDriftAndReso(), fitTMax(), and setVerbosity().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:18:57 2009 for CMSSW by  doxygen 1.5.4