BOSS 7.0.5
BESIII Offline Software System
Loading...
Searching...
No Matches
DedxCalibDocaEAng.cxx
Go to the documentation of this file.
1#include "GaudiKernel/MsgStream.h"
2
3#include <sstream>
4#include <string>
5#include "TTree.h"
6#include "TH1F.h"
7#include "TH2F.h"
8#include "TF1.h"
9#include "TFile.h"
10#include "TStyle.h"
11#include "TCanvas.h"
12
13#include "DedxCalibAlg/DedxCalibDocaEAng.h"
14
15using namespace std;
17
18DedxCalibDocaEAng::DedxCalibDocaEAng(const std::string& name, ISvcLocator* pSvcLocator): DedxCalib(name, pSvcLocator){
19 declareProperty("Debug", m_debug=false);
20 declareProperty("Idoca", m_idoca_test=1);
21 declareProperty("Ieang", m_ieang_test=1);
22 declareProperty("PMin", pMin=0.2);
23 declareProperty("PMax", pMax=2.3);
24 declareProperty("Truncate", m_truncate=0.7);
25}
26
28{
29 MsgStream log(msgSvc(), name());
30 log<<MSG::INFO<<"DedxCalibDocaEAng::BookHists()"<<endreq;
31
33
34 m_docaeangle = new TH1F**[NumSlicesDoca];
35 stringstream hlabel;
36 for(int i=0;i<NumSlicesDoca;i++)
37 {
38 m_docaeangle[i] = new TH1F*[NumSlicesEAng];
39 for(int j=0;j<NumSlicesEAng;j++)
40 {
41 hlabel.str("");
42 hlabel<<"dEdx_doca"<<i<<"_eang"<<j;
43 m_docaeangle[i][j] = new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinHistValue, MaxHistValue);
44 }
45 }
46 hlabel.str("");
47 hlabel<<"dEdxVsDocaEAng";
48 m_DocaEAngAverdE = new TH2F(hlabel.str().c_str(),"dEdx vs Doca & EAng",NumSlicesDoca,DocaMin,DocaMax,NumSlicesEAng,PhiMin,PhiMax);
49}
50
52{
53 MsgStream log(msgSvc(), name());
54 log<<MSG::INFO<<"DedxCalibDocaEAng::FillHists()"<<endreq;
55
56 TFile* f, *f_test(0);
57 TTree* n103, *t_test(0);
58 string runlist;
59
60 double doca_test(0), eang_test(0);
61 int ndedxhit=0;
62 double dEdx[100]={0},pathlength[100]={0},wid[100]={0},layid[100]={0},dd_in[100]={0},eangle[100]={0},zhit[100]={0};
63 float runNO=0,evtNO=0,runFlag=0,costheta=0,tes=0,charge=0,recalg=0,ptrk=0,chidedx=0;
64 cout<<"m_recFileVector.size()= "<<m_recFileVector.size()<<endl;
65 if(m_debug){ f_test = new TFile("doca_eangle_test.root","recreate");
66 t_test = new TTree("inftest", "information for test doca eangle calibration");
67 t_test->Branch("idoca", &m_idoca_test, "idoca/I");
68 t_test->Branch("ieang", &m_ieang_test, "ieang/I");
69 t_test->Branch("doca", &doca_test, "doca/D");
70 t_test->Branch("eang", &eang_test, "eang/D");
71 }
72
73 for(unsigned int i=0; i<m_recFileVector.size(); i++){
74 runlist = m_recFileVector[i];
75 cout<<"runlist: "<<runlist.c_str()<<endl;
76 f = new TFile(runlist.c_str());
77 n103 = (TTree*)f->Get("n103");
78 if(f==0){ cout << "the file is empty" << endl; continue;}
79 if(n103==0){ cout << "the tree is empty" << endl; continue;}
80 n103-> SetBranchAddress("ndedxhit", &ndedxhit);
81 n103-> SetBranchAddress("dEdx_hit",dEdx);
82 n103-> SetBranchAddress("pathlength_hit",pathlength);
83 n103-> SetBranchAddress("wid_hit",wid);
84 n103-> SetBranchAddress("layid_hit",layid);
85 n103-> SetBranchAddress("dd_in_hit",dd_in);
86 n103-> SetBranchAddress("eangle_hit",eangle);
87 n103-> SetBranchAddress("zhit_hit",zhit);
88 n103-> SetBranchAddress("runNO",&runNO);
89 n103-> SetBranchAddress("evtNO",&evtNO);
90 n103-> SetBranchAddress("runFlag",&runFlag);
91 n103-> SetBranchAddress("costheta",&costheta);
92 n103-> SetBranchAddress("t0",&tes);
93 n103-> SetBranchAddress("charge",&charge);
94 n103-> SetBranchAddress("recalg",&recalg);
95 n103-> SetBranchAddress("ptrk",&ptrk);
96 n103-> SetBranchAddress("chidedx_e",&chidedx);
97
98 vector<double> m_phlist_test(0);
99 for(int j=0;j<n103->GetEntries();j++){
100 n103->GetEntry(j);
101 if(ptrk>pMax || ptrk<pMin) continue;
102 if(tes>1400) continue;
103 // begin of save truncated hits information for EACH track
104 // m_phlist_test.clear();
105 // for(int ktest=0; ktest<ndedxhit; ktest++){
106 // if(layid[ktest]<4) continue;
107 // m_phlist_test.push_back(dEdx[ktest]);
108 // }
109 // sort(m_phlist_test.begin(), m_phlist_test.end());
110 // if(m_debug){
111 // for(unsigned int i=0; i<m_phlist_test.size(); i++)
112 // cout << "i " << i << " " << m_phlist_test[i] << endl;
113 // }
114 // unsigned int nsample = (unsigned int)(m_phlist_test.size()*m_truncate);
115 // double m_trunc_value = m_phlist_test[nsample];
116 // if(m_debug) cout << "trunctate value " << m_trunc_value << endl;
117 for(int k=0;k<ndedxhit;k++){
118 // if(dEdx[k]>m_trunc_value) continue;
119 dd_in[k] = dd_in[k]/doca_norm[int(layid[k])];
120 if(eangle[k] >0.25) eangle[k] = eangle[k] -0.5;
121 else if(eangle[k] <-0.25) eangle[k] = eangle[k] +0.5;
122 if(eangle[k]>PhiMin && eangle[k]<PhiMax && dd_in[k]>DocaMin && dd_in[k]<DocaMax)
123 {
124 if(layid[k] <4) continue;
125// else if(layid[k] <8)
126// {if(dEdx[k]<MinADCValuecut || dEdx[k]>MaxADCValuecut || pathlength[k]>RPhi_PathMaxCut || pathlength[k]<Iner_RPhi_PathMinCut || driftdist[k]>Iner_DriftDistCut) continue;}
127// else
128// {if(dEdx[k]<MinADCValuecut || dEdx[k]>MaxADCValuecut || pathlength[k]>RPhi_PathMaxCut || pathlength[k]<Out_RPhi_PathMinCut || driftdist[k]>Out_DriftDistCut) continue;}
129 int iDoca =(Int_t)floor((dd_in[k]-DocaMin)/StepDoca);
130 int iEAng = 0;
131 for(int i =0; i<14; i++)
132 {
133 if(eangle[k] < Eangle_value_cut[i] || eangle[k] >= Eangle_value_cut[i+1]) continue;
134 else if(eangle[k] >= Eangle_value_cut[i] && eangle[k] < Eangle_value_cut[i+1])
135 {
136 for(int m =0; m<i; m++)
137 {
138 iEAng += Eangle_cut_bin[m];
139 }
140 double eangle_bin_cut_step = (Eangle_value_cut[i+1] - Eangle_value_cut[i])/Eangle_cut_bin[i];
141 int temp_bin = int((eangle[k]-Eangle_value_cut[i])/eangle_bin_cut_step);
142 iEAng += temp_bin;
143 }
144 }
145 if(m_debug && (iDoca==m_idoca_test) && (iEAng==m_ieang_test)){
146 doca_test = dd_in[k];
147 eang_test = eangle[k];
148 t_test->Fill();
149 cout << "dedx before cor: " << dEdx[k]*1.5/pathlength[k] << endl;
150 }
151 dEdx[k] = exsvc->StandardHitCorrec(0,runFlag,2,runNO,evtNO,pathlength[k],wid[k],layid[k],dEdx[k],(dd_in[k])*(doca_norm[int(layid[k])]),eangle[k],costheta);
152 dEdx[k] = exsvc->StandardTrackCorrec(0, runFlag, 2, runNO, evtNO, dEdx[k], costheta, tes);
153 m_docaeangle[iDoca][iEAng]->Fill(dEdx[k]);
154 if(m_debug && (iDoca==m_idoca_test) && (iEAng==m_ieang_test)){
155 cout << "dedx after cor: " << dEdx[k] << endl << endl;
156 }
157 } // end of if(eangle, doca) are normal
158 } // end of loop in k, hit list
159 } // end of loop in j, entry list
160 f->Close(); f = 0; n103 = 0;
161 }// end of loop in i, filelist
162 if(m_debug){ f_test->cd(); t_test->Write(); f_test->Close(); }
163} // end of FillHists()
164
166{
167 MsgStream log(msgSvc(), name());
168 log<<MSG::INFO<<"DedxCalibDocaEAng::AnalyseHists()"<<endreq;
169
170 TF1* func;
171 Double_t entry=0,mean=0,rms=0;
172 Double_t binmax=0;
173 Double_t lp[3]={0};
174 gStyle->SetOptFit(1111);
175 for(Int_t id=0;id<NumSlicesDoca;id++)
176 {
177 for(Int_t ip=0;ip<NumSlicesEAng;ip++)
178 {
179 entry = m_docaeangle[id][ip]->GetEntries();
180 if(entry<10) continue;
181 mean = m_docaeangle[id][ip]->GetMean();
182 rms = m_docaeangle[id][ip]->GetRMS();
183 binmax = m_docaeangle[id][ip]->GetMaximumBin();
184 lp[1] = m_docaeangle[id][ip]->GetBinCenter(binmax);
185 lp[2] = 200;
186 lp[0] = (m_docaeangle[id][ip]->GetBinContent(binmax)+m_docaeangle[id][ip]->GetBinContent(binmax-1)+m_docaeangle[id][ip]->GetBinContent(binmax+1))/3;
187
188 if(m_phShapeFlag==0)
189 {
190 func = new TF1("mylan",mylan,10,3000,4);
191 func->SetParameters(entry, mean, rms, -0.8);
192 func->SetParLimits(0, 0, entry);
193 func->SetParLimits(1, 0, mean);
194 func->SetParLimits(2, 10, rms);
195 func->SetParLimits(3, -1, 0);
196 func->SetParameters(entry/36., mean-200, rms/2.3, -0.5);
197 }
198 else if(m_phShapeFlag==1)
199 {
200 func = new TF1("landaun",landaun,10,3000,3);
201 func->SetParameters(lp[0], lp[1], lp[2] );
202 }
203 else if(m_phShapeFlag==2)
204 {
205 func = new TF1("Landau",Landau,10,3000,2);
206 func->SetParameters(0.7*mean, rms );
207 }
208 else if(m_phShapeFlag==3)
209 {
210 func = new TF1("Vavilov",Vavilov,10,3000,2);
211 func->SetParameters(0.05, 0.6);
212 }
213 else if(m_phShapeFlag==4)
214 {
215 func = new TF1("AsymGauss",AsymGauss,10,3000,4);
216 func->SetParameters(entry, mean, rms, 1.0 );
217 }
218 func->SetLineWidth(2.1);
219 func->SetLineColor(2);
220
221 m_docaeangle[id][ip]->Fit(func, "rq");
222 }
223 }
224}
225
227{
228 MsgStream log(msgSvc(), name());
229 log<<MSG::INFO<<"DedxCalibDocaEAng::WriteHists()"<<endreq;
230
231 const int totalNo = NumSlicesDoca*NumSlicesEAng;
232 double Out_entryNo[totalNo]={0},Out_mean[totalNo]={0},Out_sigma[totalNo]={0},Out_fitmean[totalNo]={0},Out_fitmeanerr[totalNo]={0},Out_fitsigma[totalNo]={0},Out_gain[totalNo]={0};
233 double Iner_entryNo[totalNo]={0},Iner_mean[totalNo]={0},Iner_sigma[totalNo]={0},Iner_fitmean[totalNo]={0},Iner_fitmeanerr[totalNo]={0},Iner_fitsigma[totalNo]={0},Iner_gain[totalNo]={0};
234 double Id_doca[totalNo]={0},Ip_eangle[totalNo]={0},doca[totalNo]={0},eangle[totalNo]={0},chisq[totalNo]={0};
235
236 double Norm(0), count(0);
237 for(Int_t id=0;id<NumSlicesDoca;id++)
238 {
239 for(Int_t ip=0;ip<NumSlicesEAng;ip++)
240 {
241 Id_doca[id*NumSlicesEAng+ip] = id;
242 Ip_eangle[id*NumSlicesEAng+ip] = ip;
243 doca[id*NumSlicesEAng+ip] = (id+0.5)*StepDoca + DocaMin;
244 int iE=0,i=0;
245 for(;i<14;i++)
246 {
247 if(ip>=iE && ip<iE+Eangle_cut_bin[i]) break;
248 iE=iE+Eangle_cut_bin[i];
249 }
250 //cout<<"iE= "<<iE<<" Eangle_cut_bin["<<i<<"]= "<<Eangle_cut_bin[i]<<endl;
251 eangle[id*NumSlicesEAng+ip] = Eangle_value_cut[i] + (ip-iE+0.5)*(Eangle_value_cut[i+1]-Eangle_value_cut[i])/Eangle_cut_bin[i];
252 //cout<<"eangle["<<id<<"]["<<ip<<"]= "<<eangle[id*NumSlicesEAng+ip]<<endl;
253
254 Out_entryNo[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetEntries();
255 Out_mean[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetMean();
256 Out_sigma[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetRMS();
257 //cout<<"Out_entryNo["<<id<<"]["<<ip<<"]= "<<Out_entryNo[id*NumSlicesEAng+ip]<<" "<<Out_mean[id*NumSlicesEAng+ip]<<" "<<Out_sigma[id*NumSlicesEAng+ip]<<endl;
258 if(Out_entryNo[id*NumSlicesEAng+ip]<10) continue;
259
260 if(m_phShapeFlag==0)
261 {
262 Out_fitmean[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetFunction("mylan")->GetParameter(1);
263 Out_fitmeanerr[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetFunction("mylan")->GetParError(1);
264 Out_fitsigma[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetFunction("mylan")->GetParameter(2);
265 //cout<<Out_fitmean[id*NumSlicesEAng+ip]<<" "<<Out_fitmeanerr[id*NumSlicesEAng+ip]<<" "<<Out_fitsigma[id*NumSlicesEAng+ip]<<endl;
266 chisq[id*NumSlicesEAng+ip] = (m_docaeangle[id][ip]->GetFunction("mylan")-> GetChisquare())/(m_docaeangle[id][ip]->GetFunction("mylan")-> GetNDF());
267 }
268 if(m_phShapeFlag==1)
269 {
270 Out_fitmean[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetFunction("landaun")->GetParameter(1);
271 Out_fitmeanerr[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetFunction("landaun")->GetParError(1);
272 Out_fitsigma[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetFunction("landaun")->GetParameter(2);
273 chisq[id*NumSlicesEAng+ip] = (m_docaeangle[id][ip]->GetFunction("landaun")-> GetChisquare())/(m_docaeangle[id][ip]->GetFunction("landaun")-> GetNDF());
274 }
275 if(m_phShapeFlag==2)
276 {
277 Out_fitmean[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetFunction("Landau")->GetParameter(1);
278 Out_fitmeanerr[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetFunction("Landau")->GetParError(1);
279 Out_fitsigma[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetFunction("Landau")->GetParameter(2);
280 chisq[id*NumSlicesEAng+ip] = (m_docaeangle[id][ip]->GetFunction("Landau")-> GetChisquare())/(m_docaeangle[id][ip]->GetFunction("Landau")-> GetNDF());
281 }
282 if(m_phShapeFlag==3)
283 {
284 Out_fitmean[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetFunction("Vavilov")->GetParameter(1);
285 Out_fitmeanerr[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetFunction("Vavilov")->GetParError(1);
286 Out_fitsigma[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetFunction("Vavilov")->GetParameter(2);
287 chisq[id*NumSlicesEAng+ip] = (m_docaeangle[id][ip]->GetFunction("Vavilov")-> GetChisquare())/(m_docaeangle[id][ip]->GetFunction("Vavilov")-> GetNDF());
288 }
289 if(m_phShapeFlag==4)
290 {
291 Out_fitmean[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetFunction("AsymGauss")->GetParameter(1);
292 Out_fitmeanerr[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetFunction("AsymGauss")->GetParError(1);
293 Out_fitsigma[id*NumSlicesEAng+ip] = m_docaeangle[id][ip]->GetFunction("AsymGauss")->GetParameter(2);
294 chisq[id*NumSlicesEAng+ip] = (m_docaeangle[id][ip]->GetFunction("AsymGauss")-> GetChisquare())/(m_docaeangle[id][ip]->GetFunction("AsymGauss")-> GetNDF());
295 }
296
297 m_DocaEAngAverdE -> SetBinContent(id+1, ip+1, Out_fitmean[id*NumSlicesEAng+ip]);
298 if(Out_fitmean[id*NumSlicesEAng+ip]>0 && Out_fitmeanerr[id*NumSlicesEAng+ip]>0 && Out_fitmeanerr[id*NumSlicesEAng+ip]<10 && Out_fitsigma[id*NumSlicesEAng+ip]>0 && Out_fitsigma[id*NumSlicesEAng+ip]<200){
299 Norm += (Out_fitmean[id*NumSlicesEAng+ip])*(Out_entryNo[id*NumSlicesEAng+ip]);
300 count += Out_entryNo[id*NumSlicesEAng+ip];
301 }
302 } // end of loop in ip
303 } // end of loop in id
304 Norm = Norm/count;
305 cout<<"count= "<<count<<" average value of doca eangle bins: "<<Norm<<endl;
306 for(Int_t i=0;i<totalNo;i++)
307 {
308 if(Out_entryNo[i]>10 && Out_fitmean[i]>0 && Out_fitsigma[i]>0) Out_gain[i]=Out_fitmean[i]/Norm;
309 else Out_gain[i] = 0;
310 }
311
312 log<<MSG::INFO<<"begin generating root file!!! "<<endreq;
313 TFile* f = new TFile(m_rootfile.c_str(), "recreate");
314 for(Int_t id=0;id<NumSlicesDoca;id++)
315 {
316 for(Int_t ip=0;ip<NumSlicesEAng;ip++)
317 {
318 m_docaeangle[id][ip]->Write();
319 }
320 }
321 m_DocaEAngAverdE->Write();
322
323 TTree* ddgcalib = new TTree("ddgcalib", "ddgcalib");
324 ddgcalib -> Branch("Out_hits", Out_entryNo, "Out_entryNo[1600]/D");
325 ddgcalib -> Branch("Out_mean", Out_mean, "Out_mean[1600]/D");
326 ddgcalib -> Branch("Out_sigma", Out_sigma, "Out_sigma[1600]/D");
327 ddgcalib -> Branch("Out_fitmean", Out_fitmean, "Out_fitmean[1600]/D");
328 ddgcalib -> Branch("Out_fitmeanerr", Out_fitmeanerr, "Out_fitmeanerr[1600]/D");
329 ddgcalib -> Branch("Out_chi", Out_fitsigma, "Out_fitsigma[1600]/D");
330 ddgcalib -> Branch("Out_gain", Out_gain, "Out_gain[1600]/D");
331 ddgcalib -> Branch("Iner_hits", Iner_entryNo, "Iner_entryNo[1600]/D");
332 ddgcalib -> Branch("Iner_mean", Iner_mean, "Iner_mean[1600]/D");
333 ddgcalib -> Branch("Iner_sigma", Iner_sigma, "Iner_sigma[1600]/D");
334 ddgcalib -> Branch("Iner_fitmean", Iner_fitmean, "Iner_fitmean[1600]/D");
335 ddgcalib -> Branch("Iner_fitmeanerr", Iner_fitmeanerr, "Iner_fitmeanerr[1600]/D");
336 ddgcalib -> Branch("Iner_chi", Iner_fitsigma, "Iner_fitsigma[1600]/D");
337 ddgcalib -> Branch("Iner_gain", Iner_gain, "Iner_gain[1600]/D");
338 ddgcalib -> Branch("Id_doca", Id_doca, "Id_doca[1600]/D");
339 ddgcalib -> Branch("Ip_eangle", Ip_eangle, "Ip_eangle[1600]/D");
340 ddgcalib -> Branch("chisq", chisq, "chisq[1600]/D");
341 ddgcalib -> Branch("doca", doca, "doca[1600]/D");
342 ddgcalib -> Branch("eangle", eangle, "eangle[1600]/D");
343 ddgcalib -> Fill();
344 ddgcalib -> Write();
345 f->Close();
346
347 TCanvas c1("c1", "canvas", 500, 400);
348 c1.Print((m_rootfile+".ps[").c_str());
349 for(Int_t id=0;id<NumSlicesDoca;id++)
350 {
351 for(Int_t ip=0;ip<NumSlicesEAng;ip++)
352 {
353 m_docaeangle[id][ip]->Draw();
354 c1.Print((m_rootfile+".ps").c_str());
355 }
356 }
357 c1.Print((m_rootfile+".ps]").c_str());
358
359 for(Int_t id=0;id<NumSlicesDoca;id++)
360 {
361 for(Int_t ip=0;ip<NumSlicesEAng;ip++)
362 {
363 delete m_docaeangle[id][ip];
364 }
365 }
366 delete m_DocaEAngAverdE;
367}
curve Branch("CurveSize",&CurveSize,"CurveSize/I")
curve Fill()
curve Write()
data SetBranchAddress("time",&time)
const double StepDoca
DOUBLE_PRECISION count[3]
double Landau(double *x, double *par)
double Vavilov(double *x, double *par)
double landaun(double *x, double *par)
double AsymGauss(double *x, double *par)
double mylan(double *x, double *par)
DedxCalibDocaEAng(const std::string &name, ISvcLocator *pSvcLocator)
void ReadRecFileList()
Definition: DedxCalib.cxx:89
virtual double StandardHitCorrec(int calib_rec_Flag, int runFlag, int ntpFlag, int runNO, int evtNO, 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, int evtNO, double ex, double costheta, double t0) const =0
TFile f("ana_bhabha660a_dqa_mcPat_zy_old.root")
float costheta
float charge
float ptrk