1#include "GaudiKernel/MsgStream.h"
12#include "DedxCalibAlg/DedxCalibLayerGain.h"
23 for(
int i=0;i<phlist.size();i++)
28 for(
int i=0;i<phlist.size();i++)
30 rms += pow((phlist[i]-mean),2);
33 rms = sqrt(rms/phlist.size());
43 MsgStream log(
msgSvc(), name());
44 log<<MSG::INFO<<
"DedxCalibLayerGain::BookHists()"<<endreq;
48 m_laygain =
new TH1F*[
layNo];
49 m_laygain_gaus =
new TH1F*[
layNo];
51 for(
int i=0; i<
layNo; i++)
54 hlabel<<
"dEdx_Lay_"<<i;
57 hlabel<<
"dEdx_gaus_Lay_"<<i;
64 MsgStream log(
msgSvc(), name());
65 log<<MSG::INFO<<
"DedxCalibLayerGain::FillHists()"<<endreq;
72 float runNO=0,runFlag=0,pathlength=0,wid=0,layid=0,dd_in=0,driftdist=0,eangle=0,zhit=0,costheta=0,tes=0,ptrk=0;
77 cout<<
"runlist: "<<runlist.c_str()<<endl;
78 f =
new TFile(runlist.c_str());
79 n102 = (TTree*)f->Get(
"n102");
92 n102->SetBranchAddress(
"ptrk1",&ptrk);
94 for(
int j=0;j<n102->GetEntries();j++)
97 if(tes>1400)
continue;
98 if (ptrk<0.1)
continue;
104 dedx =
exsvc->
StandardHitCorrec(0,runFlag,2,runNO,pathlength,wid,layid,dedx,dd_in,eangle,costheta);
106 m_laygain[(int)layid]->Fill(dedx);
108 m_vector[(int)layid].push_back(dedx);
116 MsgStream log(
msgSvc(), name());
117 log<<MSG::INFO<<
"DedxCalibLayerGain::AnalyseHists()"<<endreq;
120 Double_t entry=0,mean=0,rms=0;
123 gStyle->SetOptFit(1111);
127 vector<double> phlist;
128 vector<double> phlist_gaus;
129 for(
int i=0; i<
layNo; i++)
131 entry = m_laygain[i]->GetEntries();
132 if(entry<10)
continue;
133 mean = m_laygain[i]->GetMean();
134 rms = m_laygain[i]->GetRMS();
135 binmax = m_laygain[i]->GetMaximumBin();
136 lp[1] = m_laygain[i]->GetBinCenter(binmax);
138 lp[0] = (m_laygain[i]->GetBinContent(binmax)+m_laygain[i]->GetBinContent(binmax-1)+m_laygain[i]->GetBinContent(binmax+1))/3;
142 func =
new TF1(
"mylan",
mylan,10,3000,4);
143 func->SetParameters(entry, mean, rms, -0.9);
147 func =
new TF1(
"landaun",
landaun,10,3000,3);
148 func->SetParameters(lp[0], lp[1], lp[2] );
152 func =
new TF1(
"Landau",
Landau,10,3000,2);
153 func->SetParameters(0.7*mean, rms );
157 func =
new TF1(
"Vavilov",
Vavilov,10,3000,2);
158 func->SetParameters(0.05, 0.6);
162 func =
new TF1(
"AsymGauss",
AsymGauss,10,3000,4);
163 func->SetParameters(entry, mean, rms, 1.0 );
165 func->SetLineWidth(2.1);
166 func->SetLineColor(2);
168 m_laygain[i]->Fit(func,
"r" );
172 for(
int j=0;j<m_vector[i].size();j+=100)
174 for(
int k=0;k<100;k++) phlist.push_back(m_vector[i][j+k]);
176 phlist_gaus.push_back(dedxt);
181 hlabel<<
"dEdx_gaus_Lay_"<<i;
185 for(
int j=0;j<phlist_gaus.size();j++) m_laygain_gaus[i]->Fill(phlist_gaus[j]);
187 m_laygain_gaus[i]->Fit(
"gaus",
"Q");
195 MsgStream log(
msgSvc(), name());
196 log<<MSG::INFO<<
"DedxCalibLayerGain::WriteHists()"<<endreq;
200 double fitmean_gaus[
layNo]={0},fitmeanerr_gaus[
layNo]={0},fitsigma_gaus[
layNo]={0},layerg_gaus[
layNo]={0};
201 for(
int i=0; i<
layNo; i++)
203 fitmean_gaus[i] = m_laygain_gaus[i]->GetFunction(
"gaus")->GetParameter(1);
204 fitmeanerr_gaus[i] = m_laygain_gaus[i]->GetFunction(
"gaus")->GetParError(1);
205 fitsigma_gaus[i] = m_laygain_gaus[i]->GetFunction(
"gaus")->GetParameter(2);
207 entryNo[i] = m_laygain[i]->GetEntries();
208 mean[i] = m_laygain[i]->GetMean();
209 sigma[i] = m_laygain[i]->GetRMS();
210 proper[i] = m_laygain[i]->GetBinCenter(m_laygain[i]->GetMaximumBin());
211 if(entryNo[i]<10)
continue;
215 fitmean[i] = m_laygain[i]->GetFunction(
"mylan")->GetParameter(1);
216 fitmeanerr[i] = m_laygain[i]->GetFunction(
"mylan")->GetParError(1);
217 fitsigma[i] = m_laygain[i]->GetFunction(
"mylan")->GetParameter(2);
218 chisq[i] = (m_laygain[i]->GetFunction(
"mylan")-> GetChisquare())/(m_laygain[i]->GetFunction(
"mylan")-> GetNDF());
222 fitmean[i] = m_laygain[i]->GetFunction(
"landaun")->GetParameter(1);
223 fitmeanerr[i] = m_laygain[i]->GetFunction(
"landaun")->GetParError(1);
224 fitsigma[i] = m_laygain[i]->GetFunction(
"landaun")->GetParameter(2);
225 chisq[i] = (m_laygain[i]->GetFunction(
"landaun")-> GetChisquare())/(m_laygain[i]->GetFunction(
"landaun")-> GetNDF());
229 fitmean[i] = m_laygain[i]->GetFunction(
"Landau")->GetParameter(1);
230 fitmeanerr[i] = m_laygain[i]->GetFunction(
"Landau")->GetParError(1);
231 fitsigma[i] = m_laygain[i]->GetFunction(
"Landau")->GetParameter(2);
232 chisq[i] = (m_laygain[i]->GetFunction(
"Landau")-> GetChisquare())/(m_laygain[i]->GetFunction(
"Landau")-> GetNDF());
236 fitmean[i] = m_laygain[i]->GetFunction(
"Vavilov")->GetParameter(1);
237 fitmeanerr[i] = m_laygain[i]->GetFunction(
"Vavilov")->GetParError(1);
238 fitsigma[i] = m_laygain[i]->GetFunction(
"Vavilov")->GetParameter(2);
239 chisq[i] = (m_laygain[i]->GetFunction(
"Vavilov")-> GetChisquare())/(m_laygain[i]->GetFunction(
"Vavilov")-> GetNDF());
243 fitmean[i] = m_laygain[i]->GetFunction(
"AsymGauss")->GetParameter(1);
244 fitmeanerr[i] = m_laygain[i]->GetFunction(
"AsymGauss")->GetParError(1);
245 fitsigma[i] = m_laygain[i]->GetFunction(
"AsymGauss")->GetParameter(2);
246 chisq[i] = (m_laygain[i]->GetFunction(
"AsymGauss")-> GetChisquare())/(m_laygain[i]->GetFunction(
"AsymGauss")-> GetNDF());
251 double sum=0, sum_gaus=0,
n=0;
252 for(
int i=4; i<
layNo; i++)
254 if(fitmean[i]>0 && fitmeanerr[i]>0 && fitmeanerr[i]<10 && fitsigma[i]>0 && fitsigma[i]<200 && entryNo[i]>1500)
257 sum_gaus+= fitmean_gaus[i];
263 cout<<
"average value of good layer: "<<sum<<endl;
265 for(
int i=0;i<
layNo;i++)
267 layerg[i] = fitmean[i]/sum;
268 layerg_gaus[i] = fitmean_gaus[i]/sum_gaus;
273 log<<MSG::INFO<<
"begin generating root file!!! "<<endreq;
274 TFile* f =
new TFile(
m_rootfile.c_str(),
"recreate");
275 for(
int i=0;i<
layNo;i++)
277 m_laygain[i]->Write();
278 m_laygain_gaus[i]->Write();
281 TTree* layergcalib =
new TTree(
"layergcalib",
"layergcalib");
282 layergcalib->Branch(
"layerg_gaus", layerg_gaus,
"layerg_gaus[43]/D");
283 layergcalib->Branch(
"layerg", layerg,
"layerg[43]/D");
284 layergcalib->Branch(
"layer", layer,
"layer[43]/D");
285 layergcalib->Branch(
"entryNo", entryNo,
"entryNo[43]/D");
286 layergcalib->Branch(
"mean", mean,
"mean[43]/D");
287 layergcalib->Branch(
"sigma", sigma,
"sigma[43]/D");
288 layergcalib->Branch(
"fitmean", fitmean,
"fitmean[43]/D");
289 layergcalib->Branch(
"fitmeanerr", fitmeanerr,
"fitmeanerr[43]/D");
290 layergcalib->Branch(
"fitsigma", fitsigma,
"fitsigma[43]/D");
291 layergcalib->Branch(
"chisq", chisq,
"chisq[43]/D");
292 layergcalib->Branch(
"fitmean_gaus", fitmean_gaus,
"fitmean_gaus[43]/D");
293 layergcalib->Branch(
"fitmeanerr_gaus", fitmeanerr_gaus,
"fitmeanerr_gaus[43]/D");
294 layergcalib->Branch(
"fitsigma_gaus", fitsigma_gaus,
"fitsigma_gaus[43]/D");
296 layergcalib->Write();
300 TCanvas c1(
"c1",
"canvas", 500, 400);
302 for(
int i=0;i<
layNo;i++)
304 m_laygain[i]->Draw();
309 for(
int i=0;i<
layNo;i++)
312 delete m_laygain_gaus[i];
data SetBranchAddress("time",&time)
void calculate(vector< double > phlist)
double Landau(double *x, double *par)
double Vavilov(double *x, double *par)
#define Iner_DriftDistCut
double landaun(double *x, double *par)
double AsymGauss(double *x, double *par)
double mylan(double *x, double *par)
DedxCalibLayerGain(const std::string &name, ISvcLocator *pSvcLocator)
vector< string > m_recFileVector
double cal_dedx_bitrunc(float truncate, std::vector< double > phlist)
virtual double StandardHitCorrec(int calib_rec_Flag, int runFlag, int ntpFlag, int runNO, double pathl, int wid, int layid, double adc, double dd, double eangle, double costheta) const =0
virtual double StandardTrackCorrec(int calib_rec_Flag, int typFlag, int ntpFlag, int runNO, double ex, double costheta, double t0) const =0