BOSS 7.0.5
BESIII Offline Software System
Loading...
Searching...
No Matches
DedxCalibCostheta.cxx
Go to the documentation of this file.
1#include "GaudiKernel/MsgStream.h"
2
3#include <sstream>
4#include <string>
5#include "TFile.h"
6#include "TF1.h"
7#include "TH1F.h"
8#include "TTree.h"
9#include "TStyle.h"
10#include "TCanvas.h"
11
12#include "DedxCalibAlg/DedxCalibCostheta.h"
13
14#define Size 700000
15
16using namespace std;
17
19
20DedxCalibCostheta::DedxCalibCostheta(const std::string& name, ISvcLocator* pSvcLocator): DedxCalib(name, pSvcLocator){
21 declareProperty("Debug", m_debug=false);
22 declareProperty("Ave", m_ave=true); // unweighted average of positive and negative
23 declareProperty("DebugMin", m_debug_min=2);
24 declareProperty("DebugMax", m_debug_max=2);
25 declareProperty("PMin", pMin=0.2);
26 declareProperty("PMax", pMax=2.3);
27}
28
30{
31 MsgStream log(msgSvc(), name());
32 log<<MSG::INFO<<"DedxCalibCostheta::BookHists()"<<endreq;
33
35
36 m_costheta = new TH1F*[NumBinCostheta];
37 m_pos_costheta = new TH1F*[NumBinCostheta];
38 m_neg_costheta = new TH1F*[NumBinCostheta];
39 m_chi = new TH1F*[NumBinCostheta];
40 m_pos_chi = new TH1F*[NumBinCostheta];
41 m_neg_chi = new TH1F*[NumBinCostheta];
42 stringstream hlabel;
43 for(int i=0;i<NumBinCostheta;i++)
44 {
45 hlabel.str("");
46 hlabel<<"dEdx_costheta"<<i;
47 m_costheta[i]=new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinHistValue1, MaxHistValue1);
48 hlabel.str("");
49 hlabel<<"pos_dEdx_costheta"<<i;
50 m_pos_costheta[i]=new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinHistValue1, MaxHistValue1);
51 hlabel.str("");
52 hlabel<<"neg_dEdx_costheta"<<i;
53 m_neg_costheta[i]=new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinHistValue1, MaxHistValue1);
54 hlabel.str("");
55 hlabel<<"chi_costheta"<<i;
56 m_chi[i]=new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinChiValue, MaxChiValue);
57 hlabel.str("");
58 hlabel<<"pos_chi_costheta"<<i;
59 m_pos_chi[i]=new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinChiValue, MaxChiValue);
60 hlabel.str("");
61 hlabel<<"neg_chi_costheta"<<i;
62 m_neg_chi[i]=new TH1F(hlabel.str().c_str(), hlabel.str().c_str(), NumHistBins, MinChiValue, MaxChiValue);
63 }
64 hlabel.str("");
65 hlabel<<"dEdxVsCostheta";
66 m_dEdxVsCostheta = new TH1F(hlabel.str().c_str(),"dEdx vs costheta",NumBinCostheta, CosthetaMin, CosthetaMax);
67 m_dEdxVsCostheta->GetXaxis()->SetTitle("cos#theta");
68 m_dEdxVsCostheta->GetYaxis()->SetTitle("dE/dx");
69 hlabel.str("");
70 hlabel<<"pos_dEdxVsCostheta";
71 m_pos_dEdxVsCostheta = new TH1F(hlabel.str().c_str(),"positron dEdx vs costheta",NumBinCostheta, CosthetaMin, CosthetaMax);
72 m_pos_dEdxVsCostheta->GetXaxis()->SetTitle("cos#theta");
73 m_pos_dEdxVsCostheta->GetYaxis()->SetTitle("dE/dx");
74 hlabel.str("");
75 hlabel<<"neg_dEdxVsCostheta";
76 m_neg_dEdxVsCostheta = new TH1F(hlabel.str().c_str(),"electron dEdx vs costheta",NumBinCostheta, CosthetaMin, CosthetaMax);
77 m_neg_dEdxVsCostheta->GetXaxis()->SetTitle("cos#theta");
78 m_neg_dEdxVsCostheta->GetYaxis()->SetTitle("dE/dx");
79
80 Vec_dedx.clear();
81 Vec_costheta.clear();
82}
83
85{
86 MsgStream log(msgSvc(), name());
87 log<<MSG::INFO<<"DedxCalibCostheta::FillHists()"<<endreq;
88
89 TFile* f;
90 TTree* n103;
91 string runlist;
92
93 int ndedxhit=0;
94 double dEdx[100]={0},pathlength[100]={0},wid[100]={0},layid[100]={0},dd_in[100]={0},eangle[100]={0},zhit[100]={0};
95 double dedx=0;
96 float runNO=0,evtNO=0,runFlag=0,costheta=0,tes=0,charge=0,recalg=0,ptrk=0,chidedx=0;
97 int usedhit=0;
98 vector<double> phlist;
99 cout<<"m_recFileVector.size()= "<<m_recFileVector.size()<<endl;
100 for(unsigned int i=0; i<m_recFileVector.size(); i++)
101 {
102 int m_current_trks(0);
103 runlist = m_recFileVector[i];
104 f = new TFile(runlist.c_str());
105 n103 = (TTree*)f->Get("n103");
106 n103-> SetBranchAddress("ndedxhit", &ndedxhit);
107 n103-> SetBranchAddress("dEdx_hit",dEdx);
108 n103-> SetBranchAddress("pathlength_hit",pathlength);
109 n103-> SetBranchAddress("wid_hit",wid);
110 n103-> SetBranchAddress("layid_hit",layid);
111 n103-> SetBranchAddress("dd_in_hit",dd_in);
112 n103-> SetBranchAddress("eangle_hit",eangle);
113 n103-> SetBranchAddress("zhit_hit",zhit);
114 n103-> SetBranchAddress("runNO",&runNO);
115 n103-> SetBranchAddress("evtNO",&evtNO);
116 n103-> SetBranchAddress("runFlag",&runFlag);
117 n103-> SetBranchAddress("costheta",&costheta);
118 n103-> SetBranchAddress("t0",&tes);
119 n103-> SetBranchAddress("charge",&charge);
120 n103-> SetBranchAddress("recalg",&recalg);
121 n103-> SetBranchAddress("ndedxhit",&ndedxhit);
122 n103-> SetBranchAddress("ptrk",&ptrk);
123 n103-> SetBranchAddress("chidedx_e",&chidedx);
124 for(int j=0;j<n103->GetEntries();j++)
125 {
126 m_current_trks ++;
127 if(m_current_trks>Size) break;
128 phlist.clear();
129 n103->GetEntry(j);
131 if(tes>1400) continue;
132 if(ptrk>pMax || ptrk<pMin) continue;
133 for(int k=0;k<ndedxhit;k++)
134 {
135 dEdx[k] = exsvc->StandardHitCorrec(0,runFlag,2,runNO,evtNO,pathlength[k],wid[k],layid[k],dEdx[k],dd_in[k],eangle[k],costheta);
136 phlist.push_back(dEdx[k]);
137 }
138 dedx = cal_dedx_bitrunc(truncate, phlist);
139 int iCos = (Int_t)floor((costheta-CosthetaMin)/StepCostheta);
140 double pre_dedx = dedx;
141 if(m_debug && iCos>=m_debug_min && iCos<=m_debug_max)
142 cout << "before cor, dedx " << pre_dedx << " with cos(theta) " << costheta << endl;
143 dedx = exsvc->StandardTrackCorrec(0, runFlag, 2, runNO, evtNO, dedx, costheta, tes);
144 if(m_debug && iCos>=m_debug_min && iCos<=m_debug_max)
145 cout << "after cor, dedx " << dedx << " with gain " << pre_dedx/dedx << endl;
146 m_costheta[iCos]->Fill(dedx);
147 if(charge>0) m_pos_costheta[iCos]->Fill(dedx);
148 if(charge<0) m_neg_costheta[iCos]->Fill(dedx);
149
150 usedhit = ndedxhit*truncate;
151 set_dEdx(1,dedx,recalg,runFlag,usedhit,ptrk,acos(costheta),1.5,vFlag,tes);
152 double chi = m_chi_ex[ParticleFlag];
153 m_chi[iCos]->Fill(chi);
154 if(charge>0) m_pos_chi[iCos]->Fill(chi);
155 if(charge<0) m_neg_chi[iCos]->Fill(chi);
156
157 Vec_dedx.push_back(dedx);
158 Vec_costheta.push_back(costheta);
159 }
160 }
161}
162
164{
165 MsgStream log(msgSvc(), name());
166 log<<MSG::INFO<<"DedxCalibCostheta::AnalyseHists()"<<endreq;
167
168 gStyle->SetOptFit(1111);
169 for(int i=0; i<NumBinCostheta; i++)
170 {
171 if(m_debug) cout << "num of bin " << i << endl;
172 if(m_costheta[i]->GetEntries()>100) m_costheta[i]->Fit("gaus", "Q" );
173 if(m_pos_costheta[i]->GetEntries()>100) m_pos_costheta[i]->Fit("gaus", "Q" );
174 if(m_neg_costheta[i]->GetEntries()>100) m_neg_costheta[i]->Fit("gaus", "Q" );
175 if(m_chi[i]->GetEntries()>100) m_chi[i]->Fit("gaus", "Q" );
176 if(m_pos_chi[i]->GetEntries()>100) m_pos_chi[i]->Fit("gaus", "Q" );
177 if(m_neg_chi[i]->GetEntries()>100) m_neg_chi[i]->Fit("gaus", "Q" );
178 }
179}
180
182{
183 MsgStream log(msgSvc(), name());
184 log<<MSG::INFO<<"DedxCalibCostheta::WriteHists()"<<endreq;
185
186 double chientryNo[NumBinCostheta]={0},chimean[NumBinCostheta]={0},chimeanerr[NumBinCostheta]={0},chisigma[NumBinCostheta]={0},chisq[NumBinCostheta]={0};
187 double pos_chientryNo[NumBinCostheta]={0},pos_chimean[NumBinCostheta]={0},pos_chimeanerr[NumBinCostheta]={0},pos_chisigma[NumBinCostheta]={0},pos_chisq[NumBinCostheta]={0};
188 double neg_chientryNo[NumBinCostheta]={0},neg_chimean[NumBinCostheta]={0},neg_chimeanerr[NumBinCostheta]={0},neg_chisigma[NumBinCostheta]={0},neg_chisq[NumBinCostheta]={0};
189 double fitentryNo[NumBinCostheta]={0},fitmean[NumBinCostheta]={0},fitmeanerr[NumBinCostheta]={0},fitsigma[NumBinCostheta]={0},gain[NumBinCostheta]={0},fitchisq[NumBinCostheta]={0};
190 double pos_fitentryNo[NumBinCostheta]={0},pos_fitmean[NumBinCostheta]={0},pos_fitmeanerr[NumBinCostheta]={0},pos_fitsigma[NumBinCostheta]={0},pos_gain[NumBinCostheta]={0},pos_fitchisq[NumBinCostheta]={0};
191 double neg_fitentryNo[NumBinCostheta]={0},neg_fitmean[NumBinCostheta]={0},neg_fitmeanerr[NumBinCostheta]={0},neg_fitsigma[NumBinCostheta]={0},neg_gain[NumBinCostheta]={0},neg_fitchisq[NumBinCostheta]={0};
192 double cosBin[NumBinCostheta]={0};
193
194 for(int i=0;i<NumBinCostheta;i++)
195 {
196 cosBin[i] = (i+0.5)*StepCostheta + CosthetaMin;
197
198 chientryNo[i] = m_chi[i]->GetEntries();
199 if(m_debug) cout << "get results at " << i << " bin with chi entries " << chientryNo[i] << endl;
200 if(m_chi[i]->GetFunction("gaus"))
201 {
202 chimean[i] = m_chi[i]->GetFunction("gaus")->GetParameter(1);
203 chimeanerr[i] = m_chi[i]->GetFunction("gaus")->GetParError(1);
204 chisigma[i] = m_chi[i]->GetFunction("gaus")->GetParameter(2);
205 chisq[i] = (m_chi[i]->GetFunction("gaus")->GetChisquare())/(m_chi[i]->GetFunction("gaus")->GetNDF());
206 }
207 pos_chientryNo[i] = m_pos_chi[i]->GetEntries();
208 if(m_debug) cout << "get results at " << i << " bin with pos_chi entries " << pos_chientryNo[i] << endl;
209 if(m_pos_chi[i]->GetFunction("gaus"))
210 {
211 pos_chimean[i] = m_pos_chi[i]->GetFunction("gaus")->GetParameter(1);
212 pos_chimeanerr[i] = m_pos_chi[i]->GetFunction("gaus")->GetParError(1);
213 pos_chisigma[i] = m_pos_chi[i]->GetFunction("gaus")->GetParameter(2);
214 pos_chisq[i] = (m_pos_chi[i]->GetFunction("gaus")->GetChisquare())/(m_pos_chi[i]->GetFunction("gaus")->GetNDF());
215 }
216 neg_chientryNo[i] = m_neg_chi[i]->GetEntries();
217 if(m_debug) cout << "get results at " << i << " bin with neg_chi entries " << neg_chientryNo[i] << endl;
218 if(m_neg_chi[i]->GetFunction("gaus"))
219 {
220 neg_chimean[i] = m_neg_chi[i]->GetFunction("gaus")->GetParameter(1);
221 neg_chimeanerr[i] = m_neg_chi[i]->GetFunction("gaus")->GetParError(1);
222 neg_chisigma[i] = m_neg_chi[i]->GetFunction("gaus")->GetParameter(2);
223 neg_chisq[i] = (m_neg_chi[i]->GetFunction("gaus")->GetChisquare())/(m_neg_chi[i]->GetFunction("gaus")->GetNDF());
224 }
225
226 fitentryNo[i] = m_costheta[i]->GetEntries();
227 if(m_debug) cout << "get results at " << i << " bin with fit entries " << fitentryNo[i] << endl;
228 if(m_costheta[i]->GetFunction("gaus"))
229 {
230 fitmean[i] = m_costheta[i]->GetFunction("gaus")->GetParameter(1);
231 fitmeanerr[i] = m_costheta[i]->GetFunction("gaus")->GetParError(1);
232 fitsigma[i] = m_costheta[i]->GetFunction("gaus")->GetParameter(2);
233 gain[i] = fitmean[i]/NormalMean;
234 fitchisq[i] = (m_costheta[i]->GetFunction("gaus")->GetChisquare())/(m_costheta[i]->GetFunction("gaus")->GetNDF());
235 }
236 pos_fitentryNo[i] = m_pos_costheta[i]->GetEntries();
237 if(m_debug) cout << "get results at " << i << " bin with pos_fit entries " << pos_fitentryNo[i] << endl;
238 if(m_pos_costheta[i]->GetFunction("gaus"))
239 {
240 pos_fitmean[i] = m_pos_costheta[i]->GetFunction("gaus")->GetParameter(1);
241 pos_fitmeanerr[i] = m_pos_costheta[i]->GetFunction("gaus")->GetParError(1);
242 pos_fitsigma[i] = m_pos_costheta[i]->GetFunction("gaus")->GetParameter(2);
243 pos_gain[i] = pos_fitmean[i]/NormalMean;
244 pos_fitchisq[i] = (m_pos_costheta[i]->GetFunction("gaus")->GetChisquare())/(m_pos_costheta[i]->GetFunction("gaus")->GetNDF());
245 }
246 neg_fitentryNo[i] = m_neg_costheta[i]->GetEntries();
247 if(m_debug) cout << "get results at " << i << " bin with neg_fit entries " << neg_fitentryNo[i] << endl;
248 if(m_neg_costheta[i]->GetFunction("gaus"))
249 {
250 neg_fitmean[i] = m_neg_costheta[i]->GetFunction("gaus")->GetParameter(1);
251 neg_fitmeanerr[i] = m_neg_costheta[i]->GetFunction("gaus")->GetParError(1);
252 neg_fitsigma[i] = m_neg_costheta[i]->GetFunction("gaus")->GetParameter(2);
253 neg_gain[i] = neg_fitmean[i]/NormalMean;
254 neg_fitchisq[i] = (m_neg_costheta[i]->GetFunction("gaus")->GetChisquare())/(m_neg_costheta[i]->GetFunction("gaus")->GetNDF());
255 }
256 // use unweighted average to reduce the asymmetry in Bhabha
257 if(m_ave){
258 fitmean[i] = (pos_fitmean[i] + neg_fitmean[i])/2;
259 fitmeanerr[i] = sqrt(pow(pos_fitmeanerr[i],2) + pow(neg_fitmeanerr[i],2));
260 fitsigma[i] = sqrt(pow(pos_fitsigma[i],2) + pow(neg_fitsigma[i],2));
261 gain[i] = fitmean[i]/NormalMean;
262 fitchisq[i] = (pos_fitchisq[i] + neg_fitchisq[i])/2;
263 }
264
265 if(fitentryNo[i]>100) m_dEdxVsCostheta -> SetBinContent(i+1,fitmean[i]);
266 if(pos_fitentryNo[i]>100) m_pos_dEdxVsCostheta -> SetBinContent(i+1,pos_fitmean[i]);
267 if(neg_fitentryNo[i]>100) m_neg_dEdxVsCostheta -> SetBinContent(i+1,neg_fitmean[i]);
268 }
269
270 double dedx1[Size] = {0};
271 double costheta1[Size] = {0};
272 cout << "Vec_dedx.size() = " << Vec_dedx.size() << endl;
273 for(unsigned int i=0;i<Size && i< Vec_dedx.size();i++)
274 {
275 dedx1[i] = Vec_dedx[i];
276 costheta1[i] = Vec_costheta[i];
277 //cout<<"i= "<<i<<"dedx= "<<dedx1[i]<<" costheta= "<<costheta1[i]<<endl;
278 }
279
280 log<<MSG::INFO<<"begin generating root file!!! "<<endreq;
281 TFile* f = new TFile(m_rootfile.c_str(), "RECREATE");
282 for(int i=0;i<NumBinCostheta;i++)
283 {
284 m_chi[i]->Write();
285 m_pos_chi[i]->Write();
286 m_neg_chi[i]->Write();
287 m_costheta[i]->Write();
288 m_pos_costheta[i]->Write();
289 m_neg_costheta[i]->Write();
290 }
291 m_dEdxVsCostheta->Write();
292 m_pos_dEdxVsCostheta->Write();
293 m_neg_dEdxVsCostheta->Write();
294
295 TTree *costhetacalib = new TTree("costhetacalib","costhetacalib");
296 costhetacalib -> Branch("chientryNo",chientryNo,"chientryNo[80]/D");
297 costhetacalib -> Branch("chimean",chimean,"chimean[80]/D");
298 costhetacalib -> Branch("chimeanerr",chimeanerr,"chimeanerr[80]/D");
299 costhetacalib -> Branch("chisigma",chisigma,"chisigma[80]/D");
300 costhetacalib -> Branch("chisq",chisq,"chisq[80]/D");
301 costhetacalib -> Branch("pos_chientryNo",pos_chientryNo,"pos_chientryNo[80]/D");
302 costhetacalib -> Branch("pos_chimean",pos_chimean,"pos_chimean[80]/D");
303 costhetacalib -> Branch("pos_chimeanerr",pos_chimeanerr,"pos_chimeanerr[80]/D");
304 costhetacalib -> Branch("pos_chisigma",pos_chisigma,"pos_chisigma[80]/D");
305 costhetacalib -> Branch("pos_chisq",pos_chisq,"pos_chisq[80]/D");
306 costhetacalib -> Branch("neg_chientryNo",neg_chientryNo,"neg_chientryNo[80]/D");
307 costhetacalib -> Branch("neg_chimean",neg_chimean,"neg_chimean[80]/D");
308 costhetacalib -> Branch("neg_chimeanerr",neg_chimeanerr,"neg_chimeanerr[80]/D");
309 costhetacalib -> Branch("neg_chisigma",neg_chisigma,"neg_chisigma[80]/D");
310 costhetacalib -> Branch("neg_chisq",neg_chisq,"neg_chisq[80]/D");
311 costhetacalib -> Branch("fitentryNo", fitentryNo, "fitentryNo[80]/D");
312 costhetacalib -> Branch("fitmean", fitmean, "fitmean[80]/D");
313 costhetacalib -> Branch("fitmeanerr", fitmeanerr, "fitmeanerr[80]/D");
314 costhetacalib -> Branch("fitsigma", fitsigma, "fitsigma[80]/D");
315 costhetacalib -> Branch("costheta", gain, "gain[80]/D");
316 costhetacalib -> Branch("fitchisq", fitchisq, "fitchisq[80]/D");
317 costhetacalib -> Branch("pos_fitentryNo", pos_fitentryNo, "pos_fitentryNo[80]/D");
318 costhetacalib -> Branch("pos_fitmean", pos_fitmean, "pos_fitmean[80]/D");
319 costhetacalib -> Branch("pos_fitmeanerr", pos_fitmeanerr, "pos_fitmeanerr[80]/D");
320 costhetacalib -> Branch("pos_fitsigma", pos_fitsigma, "pos_fitsigma[80]/D");
321 costhetacalib -> Branch("pos_gain", pos_gain, "pos_gain[80]/D");
322 costhetacalib -> Branch("pos_fitchisq", pos_fitchisq, "pos_fitchisq[80]/D");
323 costhetacalib -> Branch("neg_fitentryNo", neg_fitentryNo, "neg_fitentryNo[80]/D");
324 costhetacalib -> Branch("neg_fitmean", neg_fitmean, "neg_fitmean[80]/D");
325 costhetacalib -> Branch("neg_fitmeanerr", neg_fitmeanerr, "neg_fitmeanerr[80]/D");
326 costhetacalib -> Branch("neg_fitsigma", neg_fitsigma, "neg_fitsigma[80]/D");
327 costhetacalib -> Branch("neg_gain", neg_gain, "neg_gain[80]/D");
328 costhetacalib -> Branch("neg_fitchisq", neg_fitchisq, "neg_fitchisq[80]/D");
329 costhetacalib -> Branch("cosBin", cosBin, "cosBin[80]/D");
330 costhetacalib -> Branch("costheta1",costheta1,"costheta1[700000]/D");
331 costhetacalib -> Branch("dedx1",dedx1,"dedx1[700000]/D");
332 costhetacalib -> Fill();
333 costhetacalib -> Write();
334
335 TCanvas c1("c1", "canvas", 500, 400);
336 c1.Print((m_rootfile+".ps[").c_str());
337 costhetacalib -> Draw("dedx1:costheta1","dedx1>200 && dedx1<1000");
338 c1.Print((m_rootfile+".ps").c_str());
339 m_dEdxVsCostheta->Draw();
340 c1.Print((m_rootfile+".ps").c_str());
341 m_pos_dEdxVsCostheta->Draw();
342 c1.Print((m_rootfile+".ps").c_str());
343 m_neg_dEdxVsCostheta->Draw();
344 c1.Print((m_rootfile+".ps").c_str());
345 for(Int_t i=0;i<NumBinCostheta;i++)
346 {
347 m_chi[i]->Draw();
348 c1.Print((m_rootfile+".ps").c_str());
349 }
350 for(Int_t i=0;i<NumBinCostheta;i++)
351 {
352 m_pos_chi[i]->Draw();
353 c1.Print((m_rootfile+".ps").c_str());
354 }
355 for(Int_t i=0;i<NumBinCostheta;i++)
356 {
357 m_neg_chi[i]->Draw();
358 c1.Print((m_rootfile+".ps").c_str());
359 }
360 for(Int_t i=0;i<NumBinCostheta;i++)
361 {
362 m_costheta[i]->Draw();
363 c1.Print((m_rootfile+".ps").c_str());
364 }
365 for(Int_t i=0;i<NumBinCostheta;i++)
366 {
367 m_pos_costheta[i]->Draw();
368 c1.Print((m_rootfile+".ps").c_str());
369 }
370 for(Int_t i=0;i<NumBinCostheta;i++)
371 {
372 m_neg_costheta[i]->Draw();
373 c1.Print((m_rootfile+".ps").c_str());
374 }
375 c1.Print((m_rootfile+".ps]").c_str());
376 f->Close();
377
378 for(Int_t i=0;i<NumBinCostheta;i++)
379 {
380 delete m_chi[i];
381 delete m_pos_chi[i];
382 delete m_neg_chi[i];
383 delete m_costheta[i];
384 delete m_pos_costheta[i];
385 delete m_neg_costheta[i];
386 }
387 delete m_dEdxVsCostheta;
388 delete m_pos_dEdxVsCostheta;
389 delete m_neg_dEdxVsCostheta;
390}
curve Branch("CurveSize",&CurveSize,"CurveSize/I")
curve Fill()
curve Write()
data SetBranchAddress("time",&time)
#define Size
const double StepCostheta
DedxCalibCostheta(const std::string &name, ISvcLocator *pSvcLocator)
void set_dEdx(int landau, float dEdx, int trkalg, int runflag, int dedxhit_use, float ptrk, float theta, float pl_rp, int vflag[3], double t0)
Definition: DedxCalib.cxx:184
void ReadRecFileList()
Definition: DedxCalib.cxx:89
double cal_dedx_bitrunc(float truncate, std::vector< double > phlist)
Definition: DedxCalib.cxx:108
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
legend Draw()
TFile f("ana_bhabha660a_dqa_mcPat_zy_old.root")
float costheta
float charge
float ptrk