00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include "CalibMuon/DTCalibration/interface/DTMeanTimerFitter.h"
00011
00012 #include "CalibMuon/DTCalibration/interface/vDriftHistos.h"
00013
00014 #include <iostream>
00015
00016 #include "TFile.h"
00017 #include "TF1.h"
00018
00019 using namespace std;
00020
00021 DTMeanTimerFitter::DTMeanTimerFitter(TFile *file) : hInputFile(file), theVerbosityLevel(0) {
00022
00023 hDebugFile = new TFile("DTMeanTimerFitter.root", "RECREATE");
00024 }
00025
00026 DTMeanTimerFitter::~DTMeanTimerFitter() {
00027 hDebugFile->Close();
00028 }
00029
00030 vector<float> DTMeanTimerFitter::evaluateVDriftAndReso (TString N) {
00031
00032
00033 hTMaxCell * histos = new hTMaxCell(N, hInputFile);
00034 vector<float> vDriftAndReso;
00035
00036
00037 if(histos->hTmax123 != 0) {
00038 vector<TH1F*> hTMax;
00039 vector <TH1F*> hT0;
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;
00053
00054 factor.push_back(sqrt(2./3.));
00055 factor.push_back(sqrt(2./7.));
00056 factor.push_back(sqrt(8./7.));
00057 factor.push_back(sqrt(2./7.));
00058 factor.push_back(sqrt(8./7.));
00059 factor.push_back(sqrt(2./3.));
00060
00061
00062
00063 vector<Double_t> mean;
00064 vector<Double_t> sigma;
00065 vector<Double_t> count;
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
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
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
00093
00094
00095 }
00096 tMaxMean /= wTMaxSum;
00097 sigmaT /= wSigmaSum;
00098
00099
00100 Double_t vDrift = 2.1 / tMaxMean;
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
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
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
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
00134 meanT0.push_back(f1->GetParameter(1));
00135 sigmaT0.push_back(f1->GetParameter(2));
00136 countT0.push_back((*ith)->GetEntries());
00137 }
00138
00139 if(hT0.size() != 6) {
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
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
00167
00168 vDriftAndReso.push_back(-1);
00169 }
00170 }
00171 return vDriftAndReso;
00172 }
00173
00174
00175
00176 TF1* DTMeanTimerFitter::fitTMax(TH1F* histo){
00177
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
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 }