]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALPIDUtils.cxx
Added function to make relative calibration objects
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALPIDUtils.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id: AliEMCALPIDUtils.cxx 33808 2009-07-15 09:48:08Z gconesab $ */
17
18 //   Compute PID weights for all the clusters that are in AliESDs.root file
19 //   the AliESDs.root have to be in the same directory as the class
20 //
21 //   and do:    
22 //   AliEMCALPIDUtils *pid = new AliEMCALPIDUtils();
23 //   pid->SetPrintInfo(kTRUE);
24 //   pid->SetHighFluxParam(); //   pid->SetLowFluxParam(); 
25 //   
26 //   then in cluster loop do
27 //   pid->ComputePID(energy, lambda0);
28 //        
29 //        Compute PID Weight for all clusters in AliESDs.root file
30 //        keep this function for the moment for a simple verification, could be removed
31 //
32 //   pid->GetPIDFinal(idx) gives the probabilities
33 //
34 //   Double_t PIDFinal[AliPID::kSPECIESN]  is the standard PID for :
35 //
36 //      kElectron :  fPIDFinal[0]
37 //      kMuon     :  fPIDFinal[1]
38 //      kPion       :  fPIDFinal[2]
39 //      kKaon       :  fPIDFinal[3]
40 //      kProton   :  fPIDFinal[4]
41 //      kPhoton   :  fPIDFinal[5]
42 //      kPi0        :  fPIDFinal[6]
43 //      kNeutron  :  fPIDFinal[7]
44 //      kKaon0    :  fPIDFinal[8]
45 //      kEleCon   :  fPIDFinal[9]
46 //      kUnknown  :  fPIDFinal[10]
47 //
48 //
49 //    PID[3] is a simple PID for
50 //      Electron & Photon  PID[0]
51 //                    Pi0  PID[1]
52 //                 Hadron  PID[2]
53 //
54 // Author: Genole Bourdaud 2007 (SUBATECH)
55 //         Marie Germain 07/2009 (SUBATECH), new parametrization for low and high flux environment
56 //         Gustavo Conesa 08/2009 (LNF), divide class in AliEMCALPID and AliEMCALPIDUtils, PIDUtils belong to library EMCALUtils 
57 // --- standard c ---
58
59 // standard C++ includes
60 //#include <Riostream.h>
61
62 // ROOT includes
63 #include "TMath.h"
64 #include "TArrayD.h"
65
66 // STEER includes
67 #include "AliEMCALPIDUtils.h"
68 #include "AliLog.h"
69   
70 ClassImp(AliEMCALPIDUtils)
71   
72 //______________________________________________
73   AliEMCALPIDUtils::AliEMCALPIDUtils():
74     fPrintInfo(kFALSE), fProbGamma(0.),fProbPiZero(0.),fProbHadron(0.), fWeightHadronEnergy(1.), fWeightGammaEnergy(1.),fWeightPiZeroEnergy(1.)
75 {
76   //
77   // Constructor.
78   // Initialize all constant values which have to be used
79   // during PID algorithm execution
80   //
81   
82   InitParameters(); 
83   
84   
85 }
86
87 //__________________________________________________________
88 void AliEMCALPIDUtils::ComputePID(Double_t energy, Double_t lambda0)
89 {
90 //
91 // This is the main command, which uses the distributions computed and parametrised, 
92 // and gives the PID by the bayesian method.
93 //
94   
95   Double_t weightGammaEnergy  = DistEnergy(energy, 1);
96   Double_t weightPiZeroEnergy = DistEnergy(energy, 2);
97   Double_t weightHadronEnergy = DistEnergy(energy, 3);
98     
99   Double_t energyhadron=energy;
100   if(energyhadron<1.)energyhadron=1.; // no energy dependance of  parametrisation for hadrons below 1 GeV
101   if (energy<2){energy =2;} // no energy dependance of parametrisation for gamma and pi0 below 2 GeV
102   
103   if (energy>55){
104     energy =55.;
105     energyhadron=55.;
106   } // same parametrisation for gamma and hadrons above 55 GeV 
107   //   for the pi0 above 55GeV the 2 gammas supperposed no way to distinguish from real gamma  PIDWeight[1]=0
108   
109   TArrayD paramDistribGamma  = DistLambda0(energy, 1);
110   TArrayD paramDistribPiZero = DistLambda0(energy, 2);
111   TArrayD paramDistribHadron = DistLambda0(energyhadron, 3);
112   
113   Bool_t norm = kFALSE;
114   
115   
116   fProbGamma   = TMath::Gaus(lambda0, paramDistribGamma[1], paramDistribGamma[2], norm) * paramDistribGamma[0];
117   fProbGamma  += TMath::Landau(((1-paramDistribGamma[4])-lambda0),paramDistribGamma[4],paramDistribGamma[5],norm)* paramDistribGamma[3];
118   if(fProbGamma<0.)fProbGamma=0.;
119   
120   fProbGamma = fProbGamma*weightGammaEnergy;
121   
122   if(energy>10. || energy < 55.){
123     fProbPiZero  = TMath::Gaus(lambda0, paramDistribPiZero[1], paramDistribPiZero[2], norm) * paramDistribPiZero[0];
124     fProbPiZero += TMath::Landau(lambda0, paramDistribPiZero[4], paramDistribPiZero[5], norm) * paramDistribPiZero[3];
125     if(fProbPiZero<0. || energy<5.)fProbPiZero=0.;
126     fProbPiZero = fProbPiZero*weightPiZeroEnergy;
127   }
128   else {
129     fProbPiZero = 0.;
130   }
131   
132   fProbHadron  = TMath::Gaus(lambda0, paramDistribHadron[1], paramDistribHadron[2], norm) * paramDistribHadron[0];
133   fProbHadron += TMath::Landau(lambda0, paramDistribHadron[4], paramDistribHadron[5], norm) * paramDistribHadron[3];
134   if(fProbHadron<0.)fProbHadron=0.;
135   fProbHadron = fProbHadron*weightHadronEnergy; // to take into account the probability for a hadron to have a given reconstructed energy 
136   
137   // compute PID Weight
138   if( (fProbGamma + fProbPiZero + fProbHadron)>0.){
139     fPIDWeight[0] = fProbGamma / (fProbGamma + fProbPiZero + fProbHadron);
140     fPIDWeight[1] = fProbPiZero / (fProbGamma+fProbPiZero+fProbHadron);
141     fPIDWeight[2] = fProbHadron / (fProbGamma+fProbPiZero+fProbHadron);
142   }
143   else{   
144 // cases where  energy and lambda0 large,  probably du to 2 clusters folded the clusters PID not assigned to hadron nor Pi0 nor gammas
145     fPIDWeight[0] = 0.;
146     fPIDWeight[1] = 0.;
147     fPIDWeight[2] = 0.;
148   }
149   
150   
151   // cout << " PID[0] "<<  fPIDWeight[0] <<  " PID[1] "<<  fPIDWeight[1] <<  " PID[2] "<<  fPIDWeight[2] << endl;
152   
153   SetPID(fPIDWeight[0], 0);
154   SetPID(fPIDWeight[1], 1);
155   SetPID(fPIDWeight[2], 2);
156   
157   // print  pid Weight only for control 
158   if (fPrintInfo) {
159     AliInfo(Form( "Energy in loop = %f", energy) );
160     AliInfo(Form( "Lambda0 in loop = %f", lambda0) );
161     AliInfo(Form( "fProbGamma in loop = %f", fProbGamma) );
162     AliInfo(Form( "fProbaPiZero = %f", fProbPiZero ));
163     AliInfo(Form( "fProbaHadron = %f", fProbHadron) );
164     AliInfo(Form( "PIDWeight in loop = %f ||| %f ||| %f",  fPIDWeight[0] , fPIDWeight[1], fPIDWeight[2]) );
165     AliInfo("********************************************************" );
166   }
167   
168   fPIDFinal[0]  = fPIDWeight[0]/2; // photon
169   fPIDFinal[1]  = fPIDWeight[2]/8;
170   fPIDFinal[2]  = fPIDWeight[2]/8;
171   fPIDFinal[3]  = fPIDWeight[2]/8;
172   fPIDFinal[4]  = fPIDWeight[2]/8;
173   fPIDFinal[5]  = fPIDWeight[0]/2; // electron
174   fPIDFinal[6]  = fPIDWeight[1]  ; // Pi0
175   fPIDFinal[7]  = fPIDWeight[2]/8;
176   fPIDFinal[8]  = fPIDWeight[2]/8;
177   fPIDFinal[9]  = fPIDWeight[2]/8;
178   fPIDFinal[10] = fPIDWeight[2]/8;
179
180 }
181
182
183
184
185 //________________________________________________________
186 TArrayD AliEMCALPIDUtils::DistLambda0(const Double_t energy, const Int_t type) 
187 {
188   //
189   // Compute the values of the parametrised distributions using the data initialised before.
190   //
191   Double_t constGauss = 0., meanGauss = 0., sigmaGauss = 0.;
192   Double_t constLandau=0., mpvLandau=0., sigmaLandau=0.;
193   TArrayD  distributionParam(6);
194   
195   switch (type) {
196     
197   case 1:
198     
199     constGauss  = PolynomialMixed2(energy, fGamma[0]);
200     meanGauss   = PolynomialMixed2(energy, fGamma[1]);
201     sigmaGauss  = PolynomialMixed2(energy, fGamma[2]);
202     constLandau = PolynomialMixed2(energy, fGamma[3]);
203     mpvLandau   = PolynomialMixed2(energy, fGamma[4]);
204     sigmaLandau = PolynomialMixed2(energy, fGamma[5]);
205    break;
206
207   case 2:
208
209     constGauss  = PolynomialMixed2(energy, fPiZero[0]);
210     meanGauss   = PolynomialMixed2(energy, fPiZero[1]);
211     sigmaGauss  = PolynomialMixed2(energy, fPiZero[2]);
212     constLandau = PolynomialMixed2(energy, fPiZero[3]);
213     mpvLandau   = PolynomialMixed2(energy, fPiZero[4]);
214     sigmaLandau = PolynomialMixed2(energy, fPiZero[5]);
215     
216     break;
217   case 3:
218     
219     constGauss  = PolynomialMixed2(energy, fHadron[0]);
220     meanGauss   = PolynomialMixed2(energy, fHadron[1]);
221     sigmaGauss  = PolynomialMixed2(energy, fHadron[2]);
222     constLandau = PolynomialMixed2(energy, fHadron[3]);
223     mpvLandau   = PolynomialMixed2(energy, fHadron[4]);
224     sigmaLandau = PolynomialMixed2(energy, fHadron[5]);
225
226     break;
227   }
228   
229   distributionParam[0] = constGauss;
230   distributionParam[1] = meanGauss;
231   distributionParam[2] = sigmaGauss;
232   distributionParam[3] = constLandau;
233   distributionParam[4] = mpvLandau;
234   distributionParam[5] = sigmaLandau;
235   
236   return distributionParam;
237 }
238
239 //________________________________________________________
240 Double_t AliEMCALPIDUtils::DistEnergy(const Double_t energy, const Int_t type) 
241 {
242   //
243   // Compute the values of the weigh for a given energy the parametrised distribution using the data initialised before.
244   //
245   Double_t constante = 0.;
246   
247   switch (type) {
248     
249   case 1:  
250     constante  = 1.;    
251     break;
252   case 2:
253     constante  = 1.;
254     break;
255   case 3:
256     constante  = PowerExp(energy, fHadronEnergyProb);
257     break;
258   }
259   
260   //   cout << "Weight   " << constante << " for energy  "<< energy<< " GeV "<<  endl;
261   
262   return constante;
263 }
264
265
266 //_______________________________________________________
267 Double_t AliEMCALPIDUtils::Polynomial(const Double_t x, const Double_t *params) const
268 {
269   //
270   // Compute a polynomial for a given value of 'x'
271   // with the array of parameters passed as the second arg
272   //
273   
274   Double_t y  = params[0];
275   y += params[1] * x;
276   y += params[2] * x * x;
277   y += params[3] * x * x * x;
278   y += params[4] * x * x * x * x;
279   y += params[5] * x * x * x * x * x;
280   
281   return y;
282 }
283 //_______________________________________________________
284 Double_t AliEMCALPIDUtils::Polynomial0(const Double_t *params) const 
285 {
286   //
287   // Compute a polynomial for a given value of 'x'
288   // with the array of parameters passed as the second arg
289   //
290   
291   Double_t y  = params[0];
292   return y;
293 }
294
295 //_______________________________________________________
296 Double_t AliEMCALPIDUtils::Polynomialinv(const Double_t x, const Double_t *params) const
297 {
298   //
299   // Compute a polynomial for a given value of 'x'
300   // with the array of parameters passed as the second arg
301   //
302   
303   Double_t y=0.;
304   
305   if(x>0){
306     y  = params[0];
307     y += params[1] / x;
308     y += params[2] / (x * x);
309     y += params[3] / (x * x * x);
310     y += params[4] / (x * x * x * x);
311     y += params[5] / (x * x * x * x * x);
312   }  
313   
314   return y;
315   
316 }
317 //_______________________________________________________
318 Double_t AliEMCALPIDUtils::PolynomialMixed1(const Double_t x, const Double_t *params) const 
319 {
320   //
321   // Compute a polynomial for a given value of 'x'
322   // with the array of parameters passed as the second arg
323   //
324   
325   Double_t y=0.;
326   if(x>0){
327     y  = params[0] / x;
328     y += params[1] ;
329     y += params[2] * x ;
330     //   y += params[3] * 0.;
331     //   y += params[4] * 0.;
332     //   y += params[5] * 0.;
333   }  
334
335   
336   return y;
337   
338 }
339
340 //_______________________________________________________
341 Double_t AliEMCALPIDUtils::PolynomialMixed2(const Double_t x, const Double_t *params) const 
342 {
343   //
344   // Compute a polynomial for a given value of 'x'
345   // with the array of parameters passed as the second arg
346   //
347   
348   Double_t y=0.;
349   if(x>0){
350     y  = params[0] / ( x * x);
351     y += params[1] / x;
352     y += params[2] ;
353     y += params[3] * x ;
354     y += params[4] * x * x ;
355     //   y += params[5] * 0.;
356   }  
357
358   return y;
359   
360 }
361
362 //_______________________________________________________
363 Double_t AliEMCALPIDUtils::PowerExp(const Double_t x, const Double_t *params) const 
364 {
365   //
366   // Compute a polynomial for a given value of 'x'
367   // with the array of parameters passed as the second arg
368   // par[0]*TMath::Power(x[0],par[1])
369   // par[0]*TMath::Exp((x[0]-par[1])*par[2]);
370   
371   Double_t y = params[0] *TMath::Power( x,params[1]);
372   y         += params[2] *TMath::Exp((x-params[3])*params[4]);
373   
374   return y;
375   
376 }
377
378
379 //_______________________________________________________
380 void AliEMCALPIDUtils::InitParameters()
381 {
382   // Initialize PID parameters, depending on the use or not of the reconstructor
383   // and the kind of event type if the reconstructor is not used.
384   //  fWeightHadronEnergy=0.;
385   //  fWeightPiZeroEnergy=0.;
386   //  fWeightGammaEnergy=0.;
387   
388   fPIDWeight[0] = -1;
389   fPIDWeight[1] = -1;
390   fPIDWeight[2] = -1;
391   
392   for(Int_t i=0; i<AliPID::kSPECIESN+1; i++)
393     fPIDFinal[i]= 0;
394   
395   //   init the parameters here instead of from loading from recparam
396   //   default parameters are PbPb parameters.
397   SetHighFluxParam();
398   
399 }
400
401
402 //_______________________________________________________
403 void AliEMCALPIDUtils::SetLowFluxParam()
404 {
405   
406   // as a first step, all array elements are initialized to 0.0
407   Int_t i=0, j=0;
408   
409   for (i = 0; i < 6; i++) {
410     for (j = 0; j < 6; j++) {
411       fGamma[i][j]      = fHadron[i][j] =  fPiZero[i][j] = 0.;
412       fGamma1to10[i][j] = fHadron1to10[i][j] = 0.;
413     }
414        //Why we had the next 3 lines?
415        //fGammaEnergyProb[i]  =  fGammaEnergyProb[i];
416        //fPiZeroEnergyProb[i] = fPiZeroEnergyProb[i];
417        //fHadronEnergyProb[i] = fHadronEnergyProb[i];
418   }
419   
420   // New parameterization for lambda0^2 (=x): f(x) = normLandau*TMath::Landau(x,mpvLandau,widthLandau)+normgaus*TMath::Gaus(x,meangaus,sigmagaus)
421   // See AliEMCALPid (index j) refers to the polynomial parameters of the fit of each parameter vs energy
422   // pp
423
424   // paramtype[0][j] = norm gauss
425   // paramtype[1][j] = mean gaus
426   // paramtype[2][j] = sigma gaus
427   // paramtype[3][j] = norm landau
428   // paramtype[4][j] = mpv landau
429   // paramtype[5][j] = sigma landau
430
431   fGamma[0][0] = -7.656908e-01; 
432   fGamma[0][1] =  2.352536e-01; 
433   fGamma[0][2] =  1.555996e-02;
434   fGamma[0][3] =  2.243525e-04;
435   fGamma[0][4] = -2.560087e-06;
436   
437   fGamma[1][0] =  6.500216e+00;
438   fGamma[1][1] = -2.564958e-01;
439   fGamma[1][2] =  1.967894e-01;
440   fGamma[1][3] = -3.982273e-04;
441   fGamma[1][4] =  2.797737e-06;
442
443   fGamma[2][0] =  2.416489e+00;
444   fGamma[2][1] = -1.601258e-01;
445   fGamma[2][2] =  3.126839e-02;
446   fGamma[2][3] =  3.387532e-04;
447   fGamma[2][4] = -4.089145e-06;
448
449   fGamma[3][0] =  0.;
450   fGamma[3][1] = -2.696008e+00;
451   fGamma[3][2] =  6.920305e-01;
452   fGamma[3][3] = -2.281122e-03;
453   fGamma[3][4] =  0.;
454
455   fGamma[4][0] =  2.281564e-01;
456   fGamma[4][1] = -7.575040e-02;
457   fGamma[4][2] =  3.813423e-01;
458   fGamma[4][3] = -1.243854e-04;
459   fGamma[4][4] =  1.232045e-06;
460
461   fGamma[5][0] = -3.290107e-01;
462   fGamma[5][1] =  3.707545e-02;
463   fGamma[5][2] =  2.917397e-03;
464   fGamma[5][3] =  4.695306e-05;
465   fGamma[5][4] = -3.572981e-07;
466
467   fHadron[0][0] = 9.482243e-01; 
468   fHadron[0][1] =  -2.780896e-01; 
469   fHadron[0][2] =  2.223507e-02;
470   fHadron[0][3] =  7.294263e-04; 
471   fHadron[0][4] =  -5.665872e-06;
472
473   fHadron[1][0] = 0.;
474   fHadron[1][1] = 0.;
475   fHadron[1][2] = 2.483298e-01;
476   fHadron[1][3] = 0.;
477   fHadron[1][4] = 0.;
478
479   fHadron[2][0] = -5.601199e+00; 
480   fHadron[2][1] =  2.097382e+00; 
481   fHadron[2][2] = -2.307965e-01;
482   fHadron[2][3] =  9.206871e-03;
483   fHadron[2][4] = -8.887548e-05;
484  
485   fHadron[3][0] =  6.543101e+00;
486   fHadron[3][1] =  -2.305203e+00;
487   fHadron[3][2] =  2.761673e-01; 
488   fHadron[3][3] = -5.465855e-03;
489   fHadron[3][4] =  2.784329e-05;
490  
491   fHadron[4][0] = -2.443530e+01;
492   fHadron[4][1] =  8.902578e+00 ;
493   fHadron[4][2] = -5.265901e-01;
494   fHadron[4][3] = 2.549111e-02;
495   fHadron[4][4] =  -2.196801e-04; 
496
497   fHadron[5][0] = 2.102007e-01;
498   fHadron[5][1] =  -3.844418e-02;
499   fHadron[5][2] =  1.234682e-01;
500   fHadron[5][3] = -3.866733e-03;
501   fHadron[5][4] = 3.362719e-05 ;
502
503   fPiZero[0][0] =  5.072157e-01;
504   fPiZero[0][1] = -5.352747e-01;
505   fPiZero[0][2] =  8.499259e-02;
506   fPiZero[0][3] = -3.687401e-03;
507   fPiZero[0][4] =  5.482280e-05;
508
509   fPiZero[1][0] =  4.590137e+02; 
510   fPiZero[1][1] = -7.079341e+01;
511   fPiZero[1][2] =  4.990735e+00;
512   fPiZero[1][3] = -1.241302e-01;
513   fPiZero[1][4] =  1.065772e-03;
514
515   fPiZero[2][0] =  1.376415e+02;
516   fPiZero[2][1] = -3.031577e+01;
517   fPiZero[2][2] =  2.474338e+00;
518   fPiZero[2][3] = -6.903410e-02;
519   fPiZero[2][4] =  6.244089e-04;
520
521   fPiZero[3][0] = 0.;
522   fPiZero[3][1] =  1.145983e+00;
523   fPiZero[3][2] = -2.476052e-01;
524   fPiZero[3][3] =  1.367373e-02;
525   fPiZero[3][4] = 0.;
526
527   fPiZero[4][0] = -2.097586e+02;
528   fPiZero[4][1] =  6.300800e+01;
529   fPiZero[4][2] = -4.038906e+00;
530   fPiZero[4][3] =  1.088543e-01;
531   fPiZero[4][4] = -9.362485e-04;
532
533   fPiZero[5][0] = -1.671477e+01; 
534   fPiZero[5][1] =  2.995415e+00;
535   fPiZero[5][2] = -6.040360e-02;
536   fPiZero[5][3] = -6.137459e-04;
537   fPiZero[5][4] =  1.847328e-05;
538   
539   fHadronEnergyProb[0] = 4.767543e-02;
540   fHadronEnergyProb[1] = -1.537523e+00;
541   fHadronEnergyProb[2] = 2.956727e-01;
542   fHadronEnergyProb[3] = -3.051022e+01;
543   fHadronEnergyProb[4] =-6.036931e-02;
544
545 //  Int_t ii= 0;
546 //  Int_t jj= 3;
547 //  AliDebug(1,Form("PID parameters (%d, %d): fGamma=%.3f, fPi=%.3f, fHadron=%.3f",
548 //                      ii,jj, fGamma[ii][jj],fPiZero[ii][jj],fHadron[ii][jj] ));
549    
550   // end for proton-proton  
551
552 }
553
554 //_______________________________________________________
555 void AliEMCALPIDUtils::SetHighFluxParam()
556 {
557   
558   // as a first step, all array elements are initialized to 0.0
559   Int_t i=0, j=0;
560   for (i = 0; i < 6; i++) {
561     for (j = 0; j < 6; j++) {
562       fGamma[i][j]      = fHadron[i][j] = fPiZero[i][j] = 0.;
563       fGamma1to10[i][j] = fHadron1to10[i][j] = 0.;
564     }
565     fGammaEnergyProb[i]  = 0.;
566     fPiZeroEnergyProb[i] = 0.;
567     fHadronEnergyProb[i] = 0.;
568   }
569   
570   // Pb Pb  this goes with inverted landau + gaussian for gammas, landau+gaussian for Pi0 and hadrons
571   
572   fGamma[0][0] = -7.656908e-01; 
573   fGamma[0][1] =  2.352536e-01; 
574   fGamma[0][2] =  1.555996e-02;
575   fGamma[0][3] =  2.243525e-04;
576   fGamma[0][4] = -2.560087e-06;
577   
578   fGamma[1][0] =  6.500216e+00;
579   fGamma[1][1] = -2.564958e-01;
580   fGamma[1][2] =  1.967894e-01;
581   fGamma[1][3] = -3.982273e-04;
582   fGamma[1][4] =  2.797737e-06;
583
584   fGamma[2][0] =  2.416489e+00;
585   fGamma[2][1] = -1.601258e-01;
586   fGamma[2][2] =  3.126839e-02;
587   fGamma[2][3] =  3.387532e-04;
588   fGamma[2][4] = -4.089145e-06;
589  
590   fGamma[3][0] =  0.;
591   fGamma[3][1] = -2.696008e+00;
592   fGamma[3][2] =  6.920305e-01;
593   fGamma[3][3] = -2.281122e-03;
594   fGamma[3][4] =  0.;
595
596   fGamma[4][0] =  2.281564e-01;
597   fGamma[4][1] = -7.575040e-02;
598   fGamma[4][2] =  3.813423e-01;
599   fGamma[4][3] = -1.243854e-04;
600   fGamma[4][4] =  1.232045e-06;
601
602   fGamma[5][0] = -3.290107e-01;
603   fGamma[5][1] =  3.707545e-02;
604   fGamma[5][2] =  2.917397e-03;
605   fGamma[5][3] =  4.695306e-05;
606   fGamma[5][4] = -3.572981e-07;
607    
608   fHadron[0][0] =   1.519112e-01;
609   fHadron[0][1] = -8.267603e-02;
610   fHadron[0][2] =  1.914574e-02;
611   fHadron[0][3] = -2.677921e-04;
612   fHadron[0][4] =  5.447939e-06;
613
614   fHadron[1][0] = 0.;
615   fHadron[1][1] = -7.549870e-02; 
616   fHadron[1][2] = 3.930087e-01;
617   fHadron[1][3] = -2.368500e-03; 
618   fHadron[1][4] = 0.;
619
620   fHadron[2][0] = 0.;
621   fHadron[2][1] =  -2.463152e-02;
622   fHadron[2][2] = 1.349257e-01;
623   fHadron[2][3] = -1.089440e-03;
624   fHadron[2][4] = 0.;
625
626   fHadron[3][0] = 0.;
627   fHadron[3][1] = 5.101560e-01;
628   fHadron[3][2] = 1.458679e-01;
629   fHadron[3][3] = 4.903068e-04;
630   fHadron[3][4] = 0.;
631
632   fHadron[4][0] = 0.;
633   fHadron[4][1] = -6.693943e-03; 
634   fHadron[4][2] =  2.444753e-01;
635   fHadron[4][3] = -5.553749e-05;
636   fHadron[4][4] = 0.;
637
638   fHadron[5][0] = -4.414030e-01;
639   fHadron[5][1] = 2.292277e-01;
640   fHadron[5][2] = -2.433737e-02;
641   fHadron[5][3] =  1.758422e-03;
642   fHadron[5][4] = -3.001493e-05;
643   
644   fPiZero[0][0] =  5.072157e-01;
645   fPiZero[0][1] = -5.352747e-01;
646   fPiZero[0][2] =  8.499259e-02;
647   fPiZero[0][3] = -3.687401e-03;
648   fPiZero[0][4] =  5.482280e-05;
649   
650   fPiZero[1][0] =  4.590137e+02; 
651   fPiZero[1][1] = -7.079341e+01;
652   fPiZero[1][2] =  4.990735e+00;
653   fPiZero[1][3] = -1.241302e-01;
654   fPiZero[1][4] =  1.065772e-03;
655   
656   fPiZero[2][0] =  1.376415e+02;
657   fPiZero[2][1] = -3.031577e+01;
658   fPiZero[2][2] =  2.474338e+00;
659   fPiZero[2][3] = -6.903410e-02;
660   fPiZero[2][4] =  6.244089e-04;
661
662   fPiZero[3][0] = 0.;
663   fPiZero[3][1] =  1.145983e+00;
664   fPiZero[3][2] = -2.476052e-01;
665   fPiZero[3][3] =  1.367373e-02;
666   fPiZero[3][4] = 0.;
667
668   fPiZero[4][0] = -2.097586e+02;
669   fPiZero[4][1] =  6.300800e+01;
670   fPiZero[4][2] = -4.038906e+00;
671   fPiZero[4][3] =  1.088543e-01;
672   fPiZero[4][4] = -9.362485e-04;
673
674   fPiZero[5][0] = -1.671477e+01; 
675   fPiZero[5][1] =  2.995415e+00;
676   fPiZero[5][2] = -6.040360e-02;
677   fPiZero[5][3] = -6.137459e-04;
678   fPiZero[5][4] =  1.847328e-05;
679
680   // those are the High Flux PbPb ones
681   fHadronEnergyProb[0] = 0.;
682   fHadronEnergyProb[1] = 0.;
683   fHadronEnergyProb[2] =  6.188452e-02;
684   fHadronEnergyProb[3] =  2.030230e+00;
685   fHadronEnergyProb[4] = -6.402242e-02;
686
687 // Int_t ii= 0;
688 // Int_t jj= 3;
689 // AliDebug(1,Form("PID parameters (%d, %d): fGamma=%.3f, fPi=%.3f, fHadron=%.3f",
690 //                      ii,jj, fGamma[ii][jj],fPiZero[ii][jj],fHadron[ii][jj] ));
691    
692 }