BOSS 7.0.5
BESIII Offline Software System
Loading...
Searching...
No Matches
MdcDedxCorrection.cxx
Go to the documentation of this file.
1// BesIII MDC dE/dx Reconstruction Module
2// Class: MdcDedxCorrection
3// Created by Dayong WANG 2003/11
4
5//#include "CLHEP/config/CLHEP.h"
6#include "MdcDedxAlg/MdcDedxCorrection.h"
7#include "MdcDedxAlg/MdcDedxParam.h"
8#include "MdcDedxAlg/MdcDedxTrk.h"
9
10#include <vector>
11#include <cmath>
12#include <iostream>
13
14using namespace std;
15
16extern "C" {float prob_ (float *, int*);}
17
19{
20#ifdef DEBUG
21 std::cout<<"MdcDedxCorrection constructed!!!"<<std::endl;
22#endif
23}
24
26{
27#ifdef DEBUG
28 std::cout<<"MdcDedxCorrection destructed!!!"<<std::endl;
29#endif
30}
31
32
33/* ------------- calculate the expects of dE/dx -------------*/
34void MdcDedxCorrection::dedx_pid_exp_old( int landau, int runflag, float dedx,
35 int Nohit, float mom, float theta, float t0, float lsamp,
36 double dedx_exp[5], double ex_sigma[5],
37 double pid_prob[5], double chi_dedx[5] ) const
38{
39 double par[5], sigma_par[4], sigma_index_nhit, sigma_index_sin;
40
41 if(runflag==1) {
47 sigma_par[0] = MdcDedxParam::HV1_sigmap0;
48 sigma_par[1] = MdcDedxParam::HV1_sigmap1;
49 sigma_par[2] = MdcDedxParam::HV1_sigmap2;
50 sigma_par[3] = MdcDedxParam::HV1_sigmap3;
51 sigma_index_nhit = MdcDedxParam::HV1_index_nhit;
52 sigma_index_sin = MdcDedxParam::HV1_index_sin;
53 }
54 else if(runflag==2) {
60 sigma_par[0] = MdcDedxParam::HV2_sigmap0;
61 sigma_par[1] = MdcDedxParam::HV2_sigmap1;
62 sigma_par[2] = MdcDedxParam::HV2_sigmap2;
63 sigma_par[3] = MdcDedxParam::HV2_sigmap3;
64 sigma_index_nhit = MdcDedxParam::HV2_index_nhit;
65 sigma_index_sin = MdcDedxParam::HV2_index_sin;
66 }
67
68 const int par_cand( 5 );
69 const float Charge_Mass[par_cand] = {0.00051100, 0.10566, 0.13957, 0.4937, 0.93827 };
70 double beta_G, beta, betterm, bethe_B, sig_param;
71
72 int Nmax_prob( 0 );
73 float max_prob( -0.01 );
74 int ndf;
75 float chi2;
76
77 for( int it = 0; it < par_cand; it++ ) {
78 beta_G = mom/Charge_Mass[it];
79 beta = beta_G/sqrt(1+(beta_G)*(beta_G));
80 betterm = par[1]-log(par[2]+pow(1/beta_G,par[4]));
81 bethe_B = par[0]/pow(beta,par[3])*betterm-par[0];
82
83 if( Nohit >0 ) {
84 dedx_exp[it] = bethe_B;
85 double sig_the=std::sin( (double)theta );
86
87 if(runflag <3 && runflag>0){
88 if(landau == 0) {
89 sig_param = 1.6*std::sin( (double) theta )/(lsamp*double(Nohit));
90 ex_sigma[it] = 0.05*bethe_B*sqrt( 50.0*sig_param );
91 }
92 else {
93 //currently use one sigmap0
94 if(beta_G < 4) {
95 sig_param=sigma_par[1]+sigma_par[2]*std::pow(beta_G,sigma_par[3]);
96 } else {
97 sig_param= sigma_par[0];
98 }
99 //double sig_the=std::sin( (double)theta );
100 sig_the=std::pow(sig_the,sigma_index_sin);
101 double sig_n;
102 sig_n=35.0/double(Nohit);
103 sig_n=std::pow(sig_n,sigma_index_nhit);
104 ex_sigma[it]=sig_param*sig_the*sig_n;
105 }
106 }
107
108 MdcDedxTrk dedxtrk;
109 double dedx_correc;
110 if(runflag == 2 ) dedx_correc = dedxtrk.SpaceChargeCorrec(theta, mom, it, dedx);
111 else dedx_correc = dedx;
112 chi_dedx[it] = (dedx_correc - dedx_exp[it])/ex_sigma[it];
113 chi2 = chi_dedx[it]*chi_dedx[it];
114 ndf = 1;
115 pid_prob[it] = prob_(&chi2,&ndf);
116 //if(it ==0 ) cout<<"runflag: "<<runflag<<" dedx : "<<dedx<<" chi_dedx: "
117 //<<chi_dedx[it] <<" ptrk: "<<mom<<endl;
118 if( it == -999 ){ // here a debug flag
119 std::cout << " mom = " << mom <<"exp"<< dedx_exp[it]
120 << " sigma "<<ex_sigma[it] <<" prob "<<pid_prob[it]
121 << std::endl;
122 }
123 if( pid_prob[it] > max_prob ) {
124 max_prob = pid_prob[it];
125 Nmax_prob = it;
126 }
127 }
128 else {
129 dedx_exp[it] = 0.0;
130 ex_sigma[it] = 1000.0;
131 pid_prob[it] = 0.0;
132 chi_dedx[it] = 999.0;
133 }
134 }
135 // std::cout<<"MdcDedxCorrection::dedx_pid_exp(blum)!!!"<<std::endl;
136}
137
138void MdcDedxCorrection::dedx_pid_exp(int vflag[3], float dedx, int trkalg,
139 int Nohit, float mom, float theta, float t0, float lsamp,
140 double dedx_exp[5], double ex_sigma[5],
141 double pid_prob[5], double chi_dedx[5],
142 std::vector<double> & par, std::vector<double> & sig_par) const
143{
144 const int par_cand( 5 );
145 const float Charge_Mass[par_cand] = {0.00051100, 0.10566, 0.13957, 0.4937, 0.93827 };
146 double beta_G, beta, betterm, bethe_B;
147
148 int dedxflag = vflag[0];
149 int sigmaflag = vflag[1];
150 bool ifMC = false;
151 if(vflag[2] == 1) ifMC = true;
152
153 int Nmax_prob(0);
154 float max_prob(-0.01);
155 int ndf;
156 float chi2;
157
158 for( int it = 0; it < par_cand; it++ ) {
159 beta_G = mom/Charge_Mass[it];
160
161 if(dedxflag == 1){
162 beta = beta_G/sqrt(1+(beta_G)*(beta_G));
163 betterm = par[1]-log(par[2]+pow(1/beta_G,par[4]));
164 bethe_B = par[0]/pow(beta,par[3])*betterm-par[0];
165 }
166 else if(dedxflag == 2) {
167 double A=0, B=0,C=0;
168 double x = beta_G;
169 if(x<4.5)
170 A=1;
171 else if(x<10)
172 B=1;
173 else
174 C=1;
175 double partA = par[0]*pow(sqrt(x*x+1),par[2])/pow(x,par[2])*(par[1]-par[5]*log(pow(1/x,par[3])) ) - par[4]+exp(par[6]+par[7]*x);
176 double partB = par[8]*pow(x,3)+par[9]*pow(x,2)+par[10]*x+par[11];
177 double partC = -par[12]*log(par[15]+pow(1/x,par[13]))+par[14];
178 bethe_B = 550*(A*partA+B*partB+C*partC);
179 //for fermi plateau ( the electron region) we just use 1.0
180 if(beta_G> 100) bethe_B=550*1.0;
181 }
182
183 if (ifMC) {
184 double A=0, B=0,C=0;
185 double x = beta_G;
186 if(x<4.5)
187 A=1;
188 else if(x<10)
189 B=1;
190 else
191 C=1;
192 double partA = par[0]*pow(sqrt(x*x+1),par[2])/pow(x,par[2])*(par[1]-par[5]*log(pow(1/x,par[3])) ) - par[4]+exp(par[6]+par[7]*x);
193 double partB = par[8]*pow(x,3)+par[9]*pow(x,2)+par[10]*x+par[11];
194 double partC = -par[12]*log(par[15]+pow(1/x,par[13]))+par[14];
195 bethe_B = 550*(A*partA+B*partB+C*partC);
196 //for fermi plateau ( the electron region) we just use 1.0
197 if(beta_G> 100) bethe_B=550*1.0;
198 }
199
200
201 if (Nohit > 0) {
202 dedx_exp[it] = bethe_B;
203 double sig_the = std::sin((double)theta);
204 double f_betagamma, g_sinth, h_nhit, i_t0;
205
206 if (ifMC) {
207
208 if(sigmaflag >= 6) {
209
210 double nhit = (double)Nohit;
211 double sigma_bg=1.0;
212 double ndedx=bethe_B/550;
213 sigma_bg= sig_par[0]+sig_par[1]*ndedx+sig_par[2]*ndedx*ndedx;
214
215 double cor_nhit = 1.0;
216 if (nhit < 5) nhit = 5;
217 if (nhit <= 35)
218 cor_nhit = sig_par[3]*pow(nhit,4)+sig_par[4]*pow(nhit,3)+sig_par[5]*pow(nhit,2) +
219 sig_par[6]*nhit+sig_par[7];
220
221 double cor_sin= 1.0;
222 if(sig_the>0.99) sig_the=0.99;
223 cor_sin = sig_par[8]*pow(sig_the,4)+sig_par[9]*pow(sig_the,3) +
224 sig_par[10]*pow(sig_the,2)+sig_par[11]*sig_the+sig_par[12];
225
226 //sigma vs t0
227 double cor_t0 = 1.0;
228
229 //calculate sigma
230 if (trkalg == 1)
231 ex_sigma[it] = 550 * sigma_bg * cor_nhit * cor_sin * cor_t0;
232 else
233 ex_sigma[it] = 550 * sigma_bg * cor_nhit * cor_sin * cor_t0 * sig_par[13];
234
235 }
236
237 else{
238 double x= beta_G;
239 double nhit = (double)Nohit;
240 double sigma_bg=1.0;
241
242 if (x > 0.3)
243 sigma_bg = sig_par[0]*exp(sig_par[1]*x)+sig_par[2]*exp(sig_par[3]*pow(x,0.25))+sig_par[4];
244 else
245 sigma_bg = sig_par[5]*exp(sig_par[6]*x)+ sig_par[7];
246
247
248 double cor_nhit = 1.0;
249 if (nhit < 5) nhit = 5;
250 if (nhit <= 35)
251 cor_nhit = sig_par[8]*pow(nhit,4)+sig_par[9]*pow(nhit,3)+sig_par[10]*pow(nhit,2) +
252 sig_par[11]*nhit+sig_par[12];
253
254 double cor_sin= 1.0;
255 if(sig_the>0.99) sig_the=0.99;
256 cor_sin = sig_par[13]*pow(sig_the,4)+sig_par[14]*pow(sig_the,3) +
257 sig_par[15]*pow(sig_the,2)+sig_par[16]*sig_the+sig_par[17];
258
259 //sigma vs t0
260 double cor_t0 = 1.0;
261
262 //calculate sigma
263 if (trkalg == 1)
264 ex_sigma[it] = 550 * sigma_bg * cor_nhit * cor_sin * cor_t0;
265 else
266 ex_sigma[it] = 550 * sigma_bg * cor_nhit * cor_sin * cor_t0 * sig_par[18];
267
268 }
269 }//end of MC
270 else {
271 if(sigmaflag == 1) {
272 f_betagamma=sig_par[0]*pow(beta_G,sig_par[1])+sig_par[2];
273 g_sinth=(sig_par[3]*sig_the*sig_the+sig_par[4])/(sig_par[3]*sig_par[5]*sig_par[5]+sig_par[4]);
274 h_nhit=(sig_par[6]*Nohit*Nohit+sig_par[7]*Nohit+sig_par[8]) /
275 (sig_par[6]*sig_par[9]*sig_par[9]+sig_par[7]*sig_par[9]+sig_par[8]);
276 if(sig_par[13] != 0)
277 i_t0 = (sig_par[10]*t0*t0+sig_par[11]*t0+sig_par[12]) /
278 (sig_par[10]*sig_par[13]*sig_par[13]+sig_par[11]*sig_par[13]+sig_par[12]);
279 else if(sig_par[13] == 0)
280 i_t0 =1;
281 //cout<<"f_betagamma : "<<f_betagamma<<" g_sinth : "<<g_sinth<<" h_nhit : "
282 //<<h_nhit<<" i_t0 : "<<i_t0<<endl;
283 ex_sigma[it]= f_betagamma* g_sinth * h_nhit * i_t0;
284 }
285 else if(sigmaflag == 2) {
286 double x = beta_G;
287 double nhit = (double)Nohit;
288 double sigma_bg=1.0;
289 if (x > 0.3)
290 sigma_bg = sig_par[0]*exp(sig_par[1]*x)+sig_par[2]*exp(sig_par[3]*pow(x,0.25))+sig_par[4];
291 else
292 sigma_bg = sig_par[5]*exp(sig_par[6]*x)+ sig_par[7];
293
294 double cor_nhit=1.0;
295 if(nhit<5) nhit=5;
296 if (nhit <= 35)
297 cor_nhit = sig_par[8]*pow(nhit,4)+sig_par[9]*pow(nhit,3)+sig_par[10]*pow(nhit,2) + sig_par[11]*nhit+sig_par[12];
298
299 double cor_sin= 1.0;
300 if(sig_the>0.99) sig_the = 0.99;
301 cor_sin = sig_par[13]*pow(sig_the,4)+sig_par[14]*pow(sig_the,3) +
302 sig_par[15]*pow(sig_the,2)+sig_par[16]*sig_the+sig_par[17];
303
304 double cor_t0 = 1;
305 if (t0 > 1200) t0 = 1200;
306 if (t0 > 800)
307 cor_t0 = sig_par[18]*pow(t0,2)+sig_par[19]*t0+sig_par[20];
308
309 ex_sigma[it]=550*sigma_bg*cor_nhit*cor_sin*cor_t0;
310 }
311 else if(sigmaflag == 3) {
312 double x= beta_G;
313 double nhit = (double)Nohit;
314 double sigma_bg=1.0;
315 if (x > 0.3)
316 sigma_bg = sig_par[0]*exp(sig_par[1]*x)+sig_par[2]*exp(sig_par[3]*pow(x,0.25))+sig_par[4];
317 else
318 sigma_bg = sig_par[5]*exp(sig_par[6]*x)+ sig_par[7];
319
320 double cor_nhit = 1.0;
321 if (nhit < 5) nhit = 5;
322 if (nhit <= 35)
323 cor_nhit = sig_par[8]*pow(nhit,4)+sig_par[9]*pow(nhit,3)+sig_par[10]*pow(nhit,2) + sig_par[11]*nhit+sig_par[12];
324
325 double cor_sin= 1.0;
326 if(sig_the>0.99) sig_the=0.99;
327 cor_sin = sig_par[13]*pow(sig_the,4)+sig_par[14]*pow(sig_the,3) +
328 sig_par[15]*pow(sig_the,2)+sig_par[16]*sig_the+sig_par[17];
329
330 ex_sigma[it]=550*sigma_bg*cor_nhit*cor_sin;
331 }
332 else if(sigmaflag == 4) {
333 double x= beta_G;
334 double nhit = (double)Nohit;
335 double sigma_bg=1.0;
336 if (x > 0.3)
337 sigma_bg = sig_par[0]*exp(sig_par[1]*x)+sig_par[2]*exp(sig_par[3]*pow(x,0.25))+sig_par[4];
338 else
339 sigma_bg = sig_par[5]*exp(sig_par[6]*x)+ sig_par[7];
340
341 double cor_nhit = 1.0;
342 if (nhit < 5) nhit = 5;
343 if (nhit <= 35)
344 cor_nhit = sig_par[8]*pow(nhit,4)+sig_par[9]*pow(nhit,3)+sig_par[10]*pow(nhit,2) + sig_par[11]*nhit+sig_par[12];
345
346 double cor_sin= 1.0;
347 if(sig_the>0.99) sig_the=0.99;
348 cor_sin = sig_par[13]*pow(sig_the,4)+sig_par[14]*pow(sig_the,3) +
349 sig_par[15]*pow(sig_the,2)+sig_par[16]*sig_the+sig_par[17];
350
351 if(trkalg==1)
352 ex_sigma[it]=550*sigma_bg*cor_nhit*cor_sin;
353 else
354 ex_sigma[it]=550*sigma_bg*cor_nhit*cor_sin*sig_par[18];
355 }
356 else if(sigmaflag == 5) {
357 double x = beta_G;
358 double nhit = (double)Nohit;
359 double sigma_bg=1.0;
360 if (x > 0.3)
361 sigma_bg = sig_par[0]*exp(sig_par[1]*x)+sig_par[2]*exp(sig_par[3]*pow(x,0.25))+sig_par[4];
362 else
363 sigma_bg = sig_par[5]*exp(sig_par[6]*x)+ sig_par[7];
364
365 double cor_nhit=1.0;
366 if(nhit<5) nhit=5;
367 if (nhit <= 35)
368 cor_nhit = sig_par[8]*pow(nhit,4)+sig_par[9]*pow(nhit,3)+sig_par[10]*pow(nhit,2) +
369 sig_par[11]*nhit+sig_par[12];
370 double cor_sin= 1.0;
371 if(sig_the>0.99) sig_the = 0.99;
372 cor_sin = sig_par[13]*pow(sig_the,4)+sig_par[14]*pow(sig_the,3) +
373 sig_par[15]*pow(sig_the,2)+sig_par[16]*sig_the+sig_par[17];
374
375 double cor_t0 = 1;
376 if (t0 > 1200) t0 = 1200;
377 if (t0 > 800)
378 cor_t0 = sig_par[18]*pow(t0,2)+sig_par[19]*t0+sig_par[20];
379
380 if(trkalg==1)
381 ex_sigma[it]=550*sigma_bg*cor_nhit*cor_sin*cor_t0;
382 else
383 ex_sigma[it]=550*sigma_bg*cor_nhit*cor_sin*cor_t0*sig_par[21];
384 }
385 else if(sigmaflag == 6) {
386 double x= bethe_B/550;
387 double nhit = (double)Nohit;
388 double sigma_bg=1.0;
389
390 sigma_bg= sig_par[0]+sig_par[1]*x+sig_par[2]*x*x;
391
392 double cor_nhit = 1.0;
393 if (nhit < 5) nhit = 5;
394 if (nhit <= 35)
395 cor_nhit = sig_par[3]*pow(nhit,4)+sig_par[4]*pow(nhit,3)+sig_par[5]*pow(nhit,2) + sig_par[6]*nhit+sig_par[7];
396
397 double cor_sin= 1.0;
398 if(sig_the>0.99) sig_the=0.99;
399 cor_sin = sig_par[8]*pow(sig_the,4)+sig_par[9]*pow(sig_the,3) +
400 sig_par[10]*pow(sig_the,2)+sig_par[11]*sig_the+sig_par[12];
401
402 if(trkalg==1)
403 ex_sigma[it]=550*sigma_bg*cor_nhit*cor_sin;
404 else
405 ex_sigma[it]=550*sigma_bg*cor_nhit*cor_sin*sig_par[13];
406 }
407 else if(sigmaflag == 7) {
408 // for(int i=0; i<17; i++) cout << i << " " << sig_par[i] << endl;
409 double x= bethe_B/550;
410 double nhit = (double)Nohit;
411 double sigma_bg=1.0;
412
413 sigma_bg= sig_par[0]+sig_par[1]*x+sig_par[2]*x*x;
414
415
416 double cor_nhit=1.0;
417 if(nhit<5) nhit=5;
418 if (nhit <= 35)
419 cor_nhit = sig_par[3]*pow(nhit,4)+sig_par[4]*pow(nhit,3)+sig_par[5]*pow(nhit,2) +
420 sig_par[6]*nhit+sig_par[7];
421 double cor_sin= 1.0;
422 if(sig_the>0.99) sig_the = 0.99;
423 cor_sin = sig_par[8]*pow(sig_the,4)+sig_par[9]*pow(sig_the,3) +
424 sig_par[10]*pow(sig_the,2)+sig_par[11]*sig_the+sig_par[12];
425
426 double cor_t0 = 1;
427 if (t0 > 1200) t0 = 1200;
428 if (t0 < 800) t0 = 800;
429 cor_t0 = sig_par[13]*pow(t0,2)+sig_par[14]*t0+sig_par[15];
430
431 if(trkalg==1)
432 ex_sigma[it]=550*sigma_bg*cor_nhit*cor_sin*cor_t0;
433 else
434 ex_sigma[it]=550*sigma_bg*cor_nhit*cor_sin*cor_t0*sig_par[16];
435 // cout << "x " << x << " nhit " << nhit << " sigma_bg " << sigma_bg << " cor_nhit " << cor_nhit << " cor_sin " << cor_sin << " cor_t0 " << cor_t0 << " trkalg " << trkalg << endl;
436 // cout << "it " << it << " mom " << mom << " dedx " << dedx << " sinth " << sig_the << " t0 " << t0 << endl;
437 }
438 }
439
440 MdcDedxTrk dedxtrk;
441 double dedx_correc = dedx;
442 chi_dedx[it] = (dedx_correc - dedx_exp[it])/ex_sigma[it];
443 chi2 = chi_dedx[it]*chi_dedx[it];
444 ndf=1;
445 pid_prob[it] = prob_(&chi2,&ndf);
446 //if(it ==0 ) cout<<"runflag: "<<runflag<<" dedx : "<<dedx<<" chi_dedx: "
447 //<<chi_dedx[it] <<" ptrk: "<<mom<<endl;
448 if (it == -999) { // here a debug flag
449 std::cout << " mom = " << mom <<"exp"<< dedx_exp[it]
450 << " sigma "<<ex_sigma[it] <<" prob "<<pid_prob[it]
451 << std::endl;
452 }
453 if( pid_prob[it] > max_prob ){
454 max_prob = pid_prob[it];
455 Nmax_prob = it;
456 }
457 }
458 else{
459 dedx_exp[it] = 0.0;
460 ex_sigma[it] = 1000.0;
461 pid_prob[it] = 0.0;
462 chi_dedx[it] = 999.0;
463 } //if Nohit > 0
464 }
465 // std::cout<<"MdcDedxCorrection::dedx_pid_exp(blum)!!!"<<std::endl;
466}
467
468/*void
469 MdcDedxCorrection::getCalib(void) const {
470 Db2GeoDeDx obj1;
471 vector<float> VecA(43);
472 DbGeoDeDx obj2;
473 int i,j,m;
474
475
476 std::cout<<"success"<<std::endl;
477 obj1.Get_DbLayerGain(VecA);
478 std::cout<<"success!!"<<std::endl;
479
480 for(j=0;j<43;j++)
481 {
482 std::cout<<"VecA[" << j << "] is: "<<VecA[j]<<std::endl;
483 }
484
485 vector<float> VecB(6860);
486 std::cout<<"success"<<std::endl;
487 obj1.Get_DbWireGain(VecB);
488 for (i=0;i<6860;i++)
489 {
490 cout<<"VecB[" << i << "] is: "<<VecB[i]<<endl;
491 }
492 std::cout<<"success!!"<<std::endl;
493
494
495 vector<float> VecC(129);
496 std::cout<<"success"<<std::endl;
497 obj1.Get_DbDeDxAttach(VecC);
498 for(i=0;i<43;i++)
499 for(j=0;j<3;j++)
500 {
501 m=i*3+j;
502 cout<<"VecC[" << m << "] is "<<VecC[m]<<endl;
503 }
504 cout<<"success!!"<<endl;
505
506 vector<float> VecD(172);
507 std::cout<<"success"<<std::endl;
508 obj1.Get_DbDeDxSatur(VecD);
509 for(i=0;i<43;i++)
510 for(j=0;j<4;j++)
511 {
512 m=i*4+j;
513 cout<<"VecD[" << m << "] is "<<VecD[m]<<endl;
514 }
515 cout<<"success!!"<<endl;
516
517 }*/
Double_t x[10]
float prob_(float *, int *)
EvtComplex exp(const EvtComplex &c)
Definition: EvtComplex.hh:252
float prob_(float *, int *)
***************************************************************************************Pseudo Class RRes *****************************************************************************************Parameters and physical constants **Maarten sept ************************************************************************DOUBLE PRECISION xsmu **************************************************************************PARTICLE DATA all others are from PDG *Only resonances with known widths into electron pairs are sept ************************************************************************C Declarations C
Definition: RRes.h:29
int vflag[3]
Curve parameter vars.
void dedx_pid_exp_old(int landau, int runflag, float dedx, int Nhit, float mom, float theta, float t0, float lsamp, double dedx_ex[5], double ex_sigma[5], double pid_prob[5], double chi_dedx[5]) const
void dedx_pid_exp(int vflag[3], float dedx, int trkalg, int Nhit, float mom, float theta, float t0, float lsamp, double dedx_ex[5], double ex_sigma[5], double pid_prob[5], double chi_dedx[5], std::vector< double > &DedxCurve_Parameter, std::vector< double > &DedxSigma_Parameter) const
double SpaceChargeCorrec(double, double, int, double)
Definition: MdcDedxTrk.cxx:333