]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALPIDUtils.cxx
Correct problem in pp and PbPb bench, add comments
[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   Double_t  energyParam;
247   
248   switch (type) {
249     
250   case 1:  
251     constante  = 1.;    
252     break;
253   case 2:
254       constante  = 1.;
255     break;
256   case 3:
257     constante  = PowerExp(energy, fHadronEnergyProb);
258     break;
259   }
260   
261   energyParam = constante;
262   
263   // //   cout << "Weight   " << constante << " for energy  "<< energy<< " GeV "<<  endl;
264   
265   return energyParam;
266 }
267
268
269 //_______________________________________________________
270 Double_t AliEMCALPIDUtils::Polynomial(const Double_t x, const Double_t *params) const
271 {
272   //
273   // Compute a polynomial for a given value of 'x'
274   // with the array of parameters passed as the second arg
275   //
276   
277   Double_t y;
278   y  = params[0];
279   y += params[1] * x;
280   y += params[2] * x * x;
281   y += params[3] * x * x * x;
282   y += params[4] * x * x * x * x;
283   y += params[5] * x * x * x * x * x;
284   
285   return y;
286 }
287 //_______________________________________________________
288 Double_t AliEMCALPIDUtils::Polynomial0(const Double_t *params) const 
289 {
290   //
291   // Compute a polynomial for a given value of 'x'
292   // with the array of parameters passed as the second arg
293   //
294   
295   Double_t y;
296   y  = params[0];
297   return y;
298 }
299
300 //_______________________________________________________
301 Double_t AliEMCALPIDUtils::Polynomialinv(const Double_t x, const Double_t *params) const
302 {
303   //
304   // Compute a polynomial for a given value of 'x'
305   // with the array of parameters passed as the second arg
306   //
307   
308   Double_t y;
309   if(x>0){
310   y  = params[0];
311   y += params[1] / x;
312   y += params[2] / (x * x);
313   y += params[3] / (x * x * x);
314   y += params[4] / (x * x * x * x);
315   y += params[5] / (x * x * x * x * x);
316   }  
317   else
318     y=0.;
319   return y;
320   
321 }
322 //_______________________________________________________
323 Double_t AliEMCALPIDUtils::PolynomialMixed1(const Double_t x, const Double_t *params) const 
324 {
325   //
326   // Compute a polynomial for a given value of 'x'
327   // with the array of parameters passed as the second arg
328   //
329   
330   Double_t y;
331   if(x>0){
332     y  = params[0] / x;
333     y += params[1] ;
334     y += params[2] * x ;
335     //   y += params[3] * 0.;
336     //   y += params[4] * 0.;
337     //   y += params[5] * 0.;
338   }  
339   else
340     y=0.;
341   
342   return y;
343   
344 }
345
346 //_______________________________________________________
347 Double_t AliEMCALPIDUtils::PolynomialMixed2(const Double_t x, const Double_t *params) const 
348 {
349   //
350   // Compute a polynomial for a given value of 'x'
351   // with the array of parameters passed as the second arg
352   //
353   
354   Double_t y;
355   if(x>0){
356     y  = params[0] / ( x * x);
357     y += params[1] / x;
358     y += params[2] ;
359     y += params[3] * x ;
360     y += params[4] * x * x ;
361     //   y += params[5] * 0.;
362   }  
363   else
364     y=0.;
365
366   return y;
367   
368 }
369
370 //_______________________________________________________
371 Double_t AliEMCALPIDUtils::PowerExp(const Double_t x, const Double_t *params) const 
372 {
373   //
374   // Compute a polynomial for a given value of 'x'
375   // with the array of parameters passed as the second arg
376   // par[0]*TMath::Power(x[0],par[1])
377   // par[0]*TMath::Exp((x[0]-par[1])*par[2]);
378   
379   Double_t y;
380   
381   y  = params[0] *TMath::Power( x,params[1]);
382   y += params[2] *TMath::Exp((x-params[3])*params[4]);
383   
384   return y;
385   
386 }
387
388
389 //_______________________________________________________
390 void AliEMCALPIDUtils::InitParameters()
391 {
392   // Initialize PID parameters, depending on the use or not of the reconstructor
393   // and the kind of event type if the reconstructor is not used.
394   //  fWeightHadronEnergy=0.;
395   //  fWeightPiZeroEnergy=0.;
396   //  fWeightGammaEnergy=0.;
397   
398   fPIDWeight[0] = -1;
399   fPIDWeight[1] = -1;
400   fPIDWeight[2] = -1;
401   
402   for(Int_t i=0; i<AliPID::kSPECIESN+1; i++)
403     fPIDFinal[i]= 0;
404   
405   //   init the parameters here instead of from loading from recparam
406   //   default parameters are PbPb parameters.
407   SetHighFluxParam();
408   
409 }
410
411
412 //_______________________________________________________
413 void AliEMCALPIDUtils::SetLowFluxParam()
414 {
415   
416   // as a first step, all array elements are initialized to 0.0
417   Int_t i, j;
418   
419   for (i = 0; i < 6; i++) {
420     for (j = 0; j < 6; j++) {
421       fGamma[i][j]      = fHadron[i][j] =  fPiZero[i][j] = 0.;
422       fGamma1to10[i][j] = fHadron1to10[i][j] = 0.;
423     }
424        fGammaEnergyProb[i]  =  fGammaEnergyProb[i];
425        fPiZeroEnergyProb[i] = fPiZeroEnergyProb[i];
426        fHadronEnergyProb[i] = fHadronEnergyProb[i];
427   }
428   
429   // New parameterization for lambda0^2 (=x): f(x) = normLandau*TMath::Landau(x,mpvLandau,widthLandau)+normgaus*TMath::Gaus(x,meangaus,sigmagaus)
430   // See AliEMCALPid (index j) refers to the polynomial parameters of the fit of each parameter vs energy
431   // pp
432
433   // paramtype[0][j] = norm gauss
434   // paramtype[1][j] = mean gaus
435   // paramtype[2][j] = sigma gaus
436   // paramtype[3][j] = norm landau
437   // paramtype[4][j] = mpv landau
438   // paramtype[5][j] = sigma landau
439
440   fGamma[0][0] = -7.656908e-01; 
441   fGamma[0][1] =  2.352536e-01; 
442   fGamma[0][2] =  1.555996e-02;
443   fGamma[0][3] =  2.243525e-04;
444   fGamma[0][4] = -2.560087e-06;
445   
446   fGamma[1][0] =  6.500216e+00;
447   fGamma[1][1] = -2.564958e-01;
448   fGamma[1][2] =  1.967894e-01;
449   fGamma[1][3] = -3.982273e-04;
450   fGamma[1][4] =  2.797737e-06;
451
452   fGamma[2][0] =  2.416489e+00;
453   fGamma[2][1] = -1.601258e-01;
454   fGamma[2][2] =  3.126839e-02;
455   fGamma[2][3] =  3.387532e-04;
456   fGamma[2][4] = -4.089145e-06;
457
458   fGamma[3][0] =  0.;
459   fGamma[3][1] = -2.696008e+00;
460   fGamma[3][2] =  6.920305e-01;
461   fGamma[3][3] = -2.281122e-03;
462   fGamma[3][4] =  0.;
463
464   fGamma[4][0] =  2.281564e-01;
465   fGamma[4][1] = -7.575040e-02;
466   fGamma[4][2] =  3.813423e-01;
467   fGamma[4][3] = -1.243854e-04;
468   fGamma[4][4] =  1.232045e-06;
469
470   fGamma[5][0] = -3.290107e-01;
471   fGamma[5][1] =  3.707545e-02;
472   fGamma[5][2] =  2.917397e-03;
473   fGamma[5][3] =  4.695306e-05;
474   fGamma[5][4] = -3.572981e-07;
475
476   fHadron[0][0] = 9.482243e-01; 
477   fHadron[0][1] =  -2.780896e-01; 
478   fHadron[0][2] =  2.223507e-02;
479   fHadron[0][3] =  7.294263e-04; 
480   fHadron[0][4] =  -5.665872e-06;
481
482   fHadron[1][0] = 0.;
483   fHadron[1][1] = 0.;
484   fHadron[1][2] = 2.483298e-01;
485   fHadron[1][3] = 0.;
486   fHadron[1][4] = 0.;
487
488   fHadron[2][0] = -5.601199e+00; 
489   fHadron[2][1] =  2.097382e+00; 
490   fHadron[2][2] = -2.307965e-01;
491   fHadron[2][3] =  9.206871e-03;
492   fHadron[2][4] = -8.887548e-05;
493  
494   fHadron[3][0] =  6.543101e+00;
495   fHadron[3][1] =  -2.305203e+00;
496   fHadron[3][2] =  2.761673e-01; 
497   fHadron[3][3] = -5.465855e-03;
498   fHadron[3][4] =  2.784329e-05;
499  
500   fHadron[4][0] = -2.443530e+01;
501   fHadron[4][1] =  8.902578e+00 ;
502   fHadron[4][2] = -5.265901e-01;
503   fHadron[4][3] = 2.549111e-02;
504   fHadron[4][4] =  -2.196801e-04; 
505
506   fHadron[5][0] = 2.102007e-01;
507   fHadron[5][1] =  -3.844418e-02;
508   fHadron[5][2] =  1.234682e-01;
509   fHadron[5][3] = -3.866733e-03;
510   fHadron[5][4] = 3.362719e-05 ;
511
512   fPiZero[0][0] =  5.072157e-01;
513   fPiZero[0][1] = -5.352747e-01;
514   fPiZero[0][2] =  8.499259e-02;
515   fPiZero[0][3] = -3.687401e-03;
516   fPiZero[0][4] =  5.482280e-05;
517
518   fPiZero[1][0] =  4.590137e+02; 
519   fPiZero[1][1] = -7.079341e+01;
520   fPiZero[1][2] =  4.990735e+00;
521   fPiZero[1][3] = -1.241302e-01;
522   fPiZero[1][4] =  1.065772e-03;
523
524   fPiZero[2][0] =  1.376415e+02;
525   fPiZero[2][1] = -3.031577e+01;
526   fPiZero[2][2] =  2.474338e+00;
527   fPiZero[2][3] = -6.903410e-02;
528   fPiZero[2][4] =  6.244089e-04;
529
530   fPiZero[3][0] = 0.;
531   fPiZero[3][1] =  1.145983e+00;
532   fPiZero[3][2] = -2.476052e-01;
533   fPiZero[3][3] =  1.367373e-02;
534   fPiZero[3][4] = 0.;
535
536   fPiZero[4][0] = -2.097586e+02;
537   fPiZero[4][1] =  6.300800e+01;
538   fPiZero[4][2] = -4.038906e+00;
539   fPiZero[4][3] =  1.088543e-01;
540   fPiZero[4][4] = -9.362485e-04;
541
542   fPiZero[5][0] = -1.671477e+01; 
543   fPiZero[5][1] =  2.995415e+00;
544   fPiZero[5][2] = -6.040360e-02;
545   fPiZero[5][3] = -6.137459e-04;
546   fPiZero[5][4] =  1.847328e-05;
547   
548   fHadronEnergyProb[0] = 4.767543e-02;
549   fHadronEnergyProb[1] = -1.537523e+00;
550   fHadronEnergyProb[2] = 2.956727e-01;
551   fHadronEnergyProb[3] = -3.051022e+01;
552   fHadronEnergyProb[4] =-6.036931e-02;
553
554 //  Int_t ii= 0;
555 //  Int_t jj= 3;
556 //  AliDebug(1,Form("PID parameters (%d, %d): fGamma=%.3f, fPi=%.3f, fHadron=%.3f",
557 //                      ii,jj, fGamma[ii][jj],fPiZero[ii][jj],fHadron[ii][jj] ));
558    
559   // end for proton-proton  
560
561 }
562
563 //_______________________________________________________
564 void AliEMCALPIDUtils::SetHighFluxParam()
565 {
566   
567   // as a first step, all array elements are initialized to 0.0
568   Int_t i, j;
569   for (i = 0; i < 6; i++) {
570     for (j = 0; j < 6; j++) {
571       fGamma[i][j]      = fHadron[i][j] = fPiZero[i][j] = 0.;
572       fGamma1to10[i][j] = fHadron1to10[i][j] = 0.;
573     }
574     fGammaEnergyProb[i]  = 0.;
575     fPiZeroEnergyProb[i] = 0.;
576     fHadronEnergyProb[i] = 0.;
577   }
578   
579   // Pb Pb  this goes with inverted landau + gaussian for gammas, landau+gaussian for Pi0 and hadrons
580   
581   fGamma[0][0] = -7.656908e-01; 
582   fGamma[0][1] =  2.352536e-01; 
583   fGamma[0][2] =  1.555996e-02;
584   fGamma[0][3] =  2.243525e-04;
585   fGamma[0][4] = -2.560087e-06;
586   
587   fGamma[1][0] =  6.500216e+00;
588   fGamma[1][1] = -2.564958e-01;
589   fGamma[1][2] =  1.967894e-01;
590   fGamma[1][3] = -3.982273e-04;
591   fGamma[1][4] =  2.797737e-06;
592
593   fGamma[2][0] =  2.416489e+00;
594   fGamma[2][1] = -1.601258e-01;
595   fGamma[2][2] =  3.126839e-02;
596   fGamma[2][3] =  3.387532e-04;
597   fGamma[2][4] = -4.089145e-06;
598  
599   fGamma[3][0] =  0.;
600   fGamma[3][1] = -2.696008e+00;
601   fGamma[3][2] =  6.920305e-01;
602   fGamma[3][3] = -2.281122e-03;
603   fGamma[3][4] =  0.;
604
605   fGamma[4][0] =  2.281564e-01;
606   fGamma[4][1] = -7.575040e-02;
607   fGamma[4][2] =  3.813423e-01;
608   fGamma[4][3] = -1.243854e-04;
609   fGamma[4][4] =  1.232045e-06;
610
611   fGamma[5][0] = -3.290107e-01;
612   fGamma[5][1] =  3.707545e-02;
613   fGamma[5][2] =  2.917397e-03;
614   fGamma[5][3] =  4.695306e-05;
615   fGamma[5][4] = -3.572981e-07;
616    
617   fHadron[0][0] =   1.519112e-01;
618   fHadron[0][1] = -8.267603e-02;
619   fHadron[0][2] =  1.914574e-02;
620   fHadron[0][3] = -2.677921e-04;
621   fHadron[0][4] =  5.447939e-06;
622
623   fHadron[1][0] = 0.;
624   fHadron[1][1] = -7.549870e-02; 
625   fHadron[1][2] = 3.930087e-01;
626   fHadron[1][3] = -2.368500e-03; 
627   fHadron[1][4] = 0.;
628
629   fHadron[2][0] = 0.;
630   fHadron[2][1] =  -2.463152e-02;
631   fHadron[2][2] = 1.349257e-01;
632   fHadron[2][3] = -1.089440e-03;
633   fHadron[2][4] = 0.;
634
635   fHadron[3][0] = 0.;
636   fHadron[3][1] = 5.101560e-01;
637   fHadron[3][2] = 1.458679e-01;
638   fHadron[3][3] = 4.903068e-04;
639   fHadron[3][4] = 0.;
640
641   fHadron[4][0] = 0.;
642   fHadron[4][1] = -6.693943e-03; 
643   fHadron[4][2] =  2.444753e-01;
644   fHadron[4][3] = -5.553749e-05;
645   fHadron[4][4] = 0.;
646
647   fHadron[5][0] = -4.414030e-01;
648   fHadron[5][1] = 2.292277e-01;
649   fHadron[5][2] = -2.433737e-02;
650   fHadron[5][3] =  1.758422e-03;
651   fHadron[5][4] = -3.001493e-05;
652   
653   fPiZero[0][0] =  5.072157e-01;
654   fPiZero[0][1] = -5.352747e-01;
655   fPiZero[0][2] =  8.499259e-02;
656   fPiZero[0][3] = -3.687401e-03;
657   fPiZero[0][4] =  5.482280e-05;
658   
659   fPiZero[1][0] =  4.590137e+02; 
660   fPiZero[1][1] = -7.079341e+01;
661   fPiZero[1][2] =  4.990735e+00;
662   fPiZero[1][3] = -1.241302e-01;
663   fPiZero[1][4] =  1.065772e-03;
664   
665   fPiZero[2][0] =  1.376415e+02;
666   fPiZero[2][1] = -3.031577e+01;
667   fPiZero[2][2] =  2.474338e+00;
668   fPiZero[2][3] = -6.903410e-02;
669   fPiZero[2][4] =  6.244089e-04;
670
671   fPiZero[3][0] = 0.;
672   fPiZero[3][1] =  1.145983e+00;
673   fPiZero[3][2] = -2.476052e-01;
674   fPiZero[3][3] =  1.367373e-02;
675   fPiZero[3][4] = 0.;
676
677   fPiZero[4][0] = -2.097586e+02;
678   fPiZero[4][1] =  6.300800e+01;
679   fPiZero[4][2] = -4.038906e+00;
680   fPiZero[4][3] =  1.088543e-01;
681   fPiZero[4][4] = -9.362485e-04;
682
683   fPiZero[5][0] = -1.671477e+01; 
684   fPiZero[5][1] =  2.995415e+00;
685   fPiZero[5][2] = -6.040360e-02;
686   fPiZero[5][3] = -6.137459e-04;
687   fPiZero[5][4] =  1.847328e-05;
688
689   // those are the High Flux PbPb ones
690   fHadronEnergyProb[0] = 0.;
691   fHadronEnergyProb[1] = 0.;
692   fHadronEnergyProb[2] =  6.188452e-02;
693   fHadronEnergyProb[3] =  2.030230e+00;
694   fHadronEnergyProb[4] = -6.402242e-02;
695
696 // Int_t ii= 0;
697 // Int_t jj= 3;
698 // AliDebug(1,Form("PID parameters (%d, %d): fGamma=%.3f, fPi=%.3f, fHadron=%.3f",
699 //                      ii,jj, fGamma[ii][jj],fPiZero[ii][jj],fHadron[ii][jj] ));
700    
701 }