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