]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSPIDv1.cxx
ceb83e5ca10f696bce953205a8a72b869cb32a62
[u/mrichter/AliRoot.git] / PHOS / AliPHOSPIDv1.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$ */
17
18 //_________________________________________________________________________
19 // Implementation version v1 of the PHOS particle identifier 
20 // Particle identification based on the 
21 //     - RCPV: distance from CPV recpoint to EMCA recpoint.
22 //     - TOF 
23 //     - PCA: Principal Components Analysis..
24 // The identified particle has an identification number corresponding 
25 // to a 9 bits number:
26 //     -Bit 0 to 2: bit set if RCPV > CpvEmcDistance (each bit corresponds
27 //      to a different efficiency-purity point of the photon identification) 
28 //     -Bit 3 to 5: bit set if TOF  < TimeGate (each bit corresponds
29 //      to a different efficiency-purity point of the photon identification) 
30 //     -Bit 6 to 9: bit set if Principal Components are 
31 //      inside an ellipse defined by the parameters a, b, c, x0 and y0.
32 //      (each bit corresponds to a different efficiency-purity point of the 
33 //      photon identification)
34 //      The PCA (Principal components analysis) needs a file that contains
35 //      a previous analysis of the correlations between the particles. This 
36 //      file is $ALICE_ROOT/PHOS/PCA8pa15_0.5-100.root. Analysis done for 
37 //      energies between 0.5 and 100 GeV.
38 //      A calibrated energy is calculated. The energy of the reconstructed
39 //      cluster is corrected with the formula A + B * E  + C * E^2, whose 
40 //      parameters where obtained through the study of the reconstructed 
41 //      energy distribution of monoenergetic photons. 
42 //
43 //      All the parameters (RCPV(2 rows-3 columns),TOF(1r-3c),PCA(5r-4c) 
44 //      and calibration(1r-3c))are stored in a file called 
45 //      $ALICE_ROOT/PHOS/Parameters.dat. Each time that AliPHOSPIDv1 is 
46 //      initialized, this parameters are copied to a Matrix (9,4), a 
47 //      TMatrixD object.  
48 //
49 // use case:
50 //  root [0] AliPHOSPIDv1 * p = new AliPHOSPIDv1("galice1.root")
51 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
52 //          // reading headers from file galice1.root and create  RecParticles 
53             // TrackSegments and RecPoints are used 
54 //          // set file name for the branch RecParticles
55 //  root [1] p->ExecuteTask("deb all time")
56 //          // available options
57 //          // "deb" - prints # of reconstructed particles
58 //          // "deb all" -  prints # and list of RecParticles
59 //          // "time" - prints benchmarking results
60 //                  
61 //  root [2] AliPHOSPIDv1 * p2 = new AliPHOSPIDv1("galice1.root","v1",kTRUE)
62 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
63 //                //Split mode.  
64 //  root [3] p2->ExecuteTask()
65 //
66
67
68 //*-- Author: Yves Schutz (SUBATECH)  & Gines Martinez (SUBATECH) & 
69 //            Gustavo Conesa April 2002
70 //            PCA redesigned by Gustavo Conesa October 2002:
71 //            The way of using the PCA has changed. Instead of 2
72 //            files with the PCA, each one with different energy ranges 
73 //            of application, we use the wide one (0.5-100 GeV), and instead
74 //            of fixing 3 ellipses for different ranges of energy, it has been
75 //            studied the dependency of the ellipses parameters with the 
76 //            energy, and they are implemented in the code as a funtion 
77 //            of the energy. 
78 //
79 //
80 //
81 // --- ROOT system ---
82 #include "TROOT.h"
83 #include "TTree.h"
84 #include "TFile.h"
85 #include "TF2.h"
86 #include "TFormula.h"
87 #include "TCanvas.h"
88 #include "TFolder.h"
89 #include "TSystem.h"
90 #include "TBenchmark.h"
91 #include "TMatrixD.h"
92 #include "TPrincipal.h"
93 #include "TSystem.h"
94
95 // --- Standard library ---
96
97
98 // --- AliRoot header files ---
99 #include "AliLog.h"
100 #include "AliGenerator.h"
101 #include "AliPHOS.h"
102 #include "AliPHOSPIDv1.h"
103 #include "AliPHOSClusterizerv1.h"
104 #include "AliPHOSEmcRecPoint.h"
105 #include "AliPHOSTrackSegment.h"
106 #include "AliPHOSTrackSegmentMakerv1.h"
107 #include "AliPHOSRecParticle.h"
108 #include "AliPHOSGeometry.h"
109 #include "AliPHOSGetter.h"
110
111 ClassImp( AliPHOSPIDv1) 
112
113 //____________________________________________________________________________
114 AliPHOSPIDv1::AliPHOSPIDv1():AliPHOSPID()
115
116   // default ctor
117  
118   InitParameters() ; 
119   fDefaultInit = kTRUE ; 
120 }
121
122 //____________________________________________________________________________
123 AliPHOSPIDv1::AliPHOSPIDv1(const AliPHOSPIDv1 & pid ):AliPHOSPID(pid)
124
125   // ctor
126   InitParameters() ; 
127   Init() ;
128
129 }
130
131 //____________________________________________________________________________
132 AliPHOSPIDv1::AliPHOSPIDv1(const TString alirunFileName, const TString eventFolderName):AliPHOSPID(alirunFileName, eventFolderName)
133
134   //ctor with the indication on where to look for the track segments
135  
136   InitParameters() ; 
137   Init() ;
138   fDefaultInit = kFALSE ; 
139 }
140
141 //____________________________________________________________________________
142 AliPHOSPIDv1::~AliPHOSPIDv1()
143
144   // dtor
145
146   delete [] fX ;       // Principal input 
147   delete [] fPPhoton ; // Photon Principal components
148   delete [] fPPi0 ;    // Pi0 Principal components
149 }
150 //____________________________________________________________________________
151 const TString AliPHOSPIDv1::BranchName() const 
152 {  
153
154   return GetName() ;
155 }
156  
157 //____________________________________________________________________________
158 void AliPHOSPIDv1::Init()
159 {
160   // Make all memory allocations that are not possible in default constructor
161   // Add the PID task to the list of PHOS tasks
162
163   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
164   if(!gime)
165     gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data()) ; 
166
167   if ( !gime->PID() ) 
168     gime->PostPID(this) ;
169 }
170
171 //____________________________________________________________________________
172 void AliPHOSPIDv1::InitParameters()
173 {
174   // Initialize PID parameters
175   fWrite                   = kTRUE ;
176   fRecParticlesInRun = 0 ; 
177   fNEvent            = 0 ;            
178   fRecParticlesInRun = 0 ;
179   fBayesian          = kTRUE ;
180   SetParameters() ; // fill the parameters matrix from parameters file
181   SetEventRange(0,-1) ;
182
183   // initialisation of response function parameters
184   // Tof
185   // Photons
186   fTphoton[0] = 0.218    ;
187   //fTphoton[0] = 1.    ;
188   fTphoton[1] = 1.55E-8  ; 
189   fTphoton[2] = 5.05E-10 ;
190   fTFphoton = new TFormula("ToF response to photons" , "gaus") ; 
191   fTFphoton->SetParameters( fTphoton[0], fTphoton[1], fTphoton[2]) ; 
192 //   // Electrons
193 //   fTelectron[0] = 0.2      ;
194 //   fTelectron[1] = 1.55E-8  ; 
195 //   fTelectron[2] = 5.35E-10 ;
196 //   fTFelectron = new TFormula("ToF response to electrons" , "gaus") ; 
197 //   fTFelectron->SetParameters( fTelectron[0], fTelectron[1], fTelectron[2]) ; 
198 //   // Muons
199 //   fTmuon[0] = 0.2     ;
200 //   fTmuon[1] = 1.55E-8 ; 
201 //   fTmuon[2] = 5.1E-10 ;
202 //   fTFmuon = new TFormula("ToF response to muons" , "gaus") ; 
203 //   fTFmuon->SetParameters( fTmuon[0], fTmuon[1], fTmuon[2]) ; 
204
205   // Pions
206   //Gaus (0 to max probability)
207   fTpiong[0] = 0.0971    ; 
208   //fTpiong[0] = 1.       ;
209   fTpiong[1] = 1.58E-8  ; 
210   fTpiong[2] = 5.69E-10 ;
211   fTFpiong = new TFormula("ToF response to pions" , "gaus") ; 
212   fTFpiong->SetParameters( fTpiong[0], fTpiong[1], fTpiong[2]) ; 
213   // Landau (max probability to inf) 
214 //   fTpionl[0] = 0.05     ; 
215 //   //fTpionl[0] = 5.53     ; 
216 //   fTpionl[1] = 1.68E-8  ; 
217 //   fTpionl[2] = 5.38E-10 ;
218 //   fTFpionl = new TFormula("ToF response to pions" , "landau") ; 
219 //   fTFpionl->SetParameters( fTpionl[0], fTpionl[1], fTpionl[2]) ; 
220
221
222   // Kaons
223   //Gaus (0 to max probability)
224   fTkaong[0] = 0.0542  ; 
225   //fTkaong[0] = 1.      ;
226   fTkaong[1] = 1.64E-8 ; 
227   fTkaong[2] = 6.07-10 ;
228   fTFkaong = new TFormula("ToF response to kaon" , "gaus") ; 
229   fTFkaong->SetParameters( fTkaong[0], fTkaong[1], fTkaong[2]) ; 
230   //Landau (max probability to inf) 
231   fTkaonl[0] = 0.264   ;
232   //fTkaonl[0] = 5.53     ;
233   fTkaonl[1] = 1.68E-8  ; 
234   fTkaonl[2] = 4.10E-10 ;
235   fTFkaonl = new TFormula("ToF response to kaon" , "landau") ; 
236   fTFkaonl->SetParameters( fTkaonl[0], fTkaonl[1], fTkaonl[2]) ; 
237
238   //Heavy Hadrons
239   //Gaus (0 to max probability)
240   fThhadrong[0] = 0.0302   ;  
241   //fThhadrong[0] = 1.       ; 
242   fThhadrong[1] = 1.73E-8  ; 
243   fThhadrong[2] = 9.52E-10 ;
244   fTFhhadrong = new TFormula("ToF response to heavy hadrons" , "gaus") ; 
245   fTFhhadrong->SetParameters( fThhadrong[0], fThhadrong[1], fThhadrong[2]) ; 
246   //Landau (max probability to inf) 
247   fThhadronl[0] = 0.139    ;  
248   //fThhadronl[0] = 5.53      ;   
249   fThhadronl[1] = 1.745E-8  ; 
250   fThhadronl[2] = 1.00E-9  ;
251   fTFhhadronl = new TFormula("ToF response to heavy hadrons" , "landau") ; 
252   fTFhhadronl->SetParameters( fThhadronl[0], fThhadronl[1], fThhadronl[2]) ; 
253
254 /// /gaussian parametrization for pions
255 //   fTpion[0] = 3.93E-2 ;  fTpion[1] = 0.130   ; fTpion[2] =-6.37E-2 ;//constant
256 //   fTpion[3] = 1.65E-8 ;  fTpion[4] =-1.40E-9 ; fTpion[5] = 5.96E-10;//mean
257 //   fTpion[6] = 8.09E-10;  fTpion[7] =-4.65E-10; fTpion[8] = 1.50E-10;//sigma
258
259 // //landau parametrization for kaons
260 //   fTkaon[0] = 0.107   ;  fTkaon[1] = 0.166   ; fTkaon[2] = 0.243   ;//constant
261 //   fTkaon[3] = 1.80E-8 ;  fTkaon[4] =-2.96E-9 ; fTkaon[5] = 9.60E-10;//mean
262 //   fTkaon[6] = 1.37E-9 ;  fTkaon[7] =-1.80E-9 ; fTkaon[8] = 6.74E-10;//sigma
263
264 // //landau parametrization for nucleons
265 //   fThhadron[0] = 6.33E-2 ;  fThhadron[1] = 2.52E-2 ; fThhadron[2] = 2.16E-2 ;//constant
266 //   fThhadron[3] = 1.94E-8 ;  fThhadron[4] =-7.06E-10; fThhadron[5] =-4.69E-10;//mean
267 //   fThhadron[6] = 2.55E-9 ;  fThhadron[7] =-1.90E-9 ; fThhadron[8] = 5.41E-10;//sigma
268
269
270   // Shower shape: dispersion gaussian parameters
271   // Photons
272
273   fDphoton[0] = 0.1    ;  fDphoton[1] = 0.      ; fDphoton[2] = 0.     ;//constant
274   //fDphoton[0] = 1.0    ;  fDphoton[1] = 0.      ; fDphoton[2] = 0.     ;//constant
275   fDphoton[3] = 1.55   ;  fDphoton[4] =-0.0863  ; fDphoton[5] = 0.287  ;//mean
276   fDphoton[6] = 0.0451 ;  fDphoton[7] =-0.0803  ; fDphoton[8] = 0.314  ;//sigma
277
278    fDpi0[0] = 0.0586  ;  fDpi0[1] = 1.06E-3 ; fDpi0[2] = 0.      ;//constant
279    //fDpi0[0] = 1.0     ;  fDpi0[1] = 0.0     ; fDpi0[2] = 0.      ;//constant
280   fDpi0[3] = 2.67    ;  fDpi0[4] =-2.00E-2 ; fDpi0[5] = 9.37E-5 ;//mean
281   fDpi0[6] = 0.153   ;  fDpi0[7] = 9.34E-4 ; fDpi0[8] =-1.49E-5 ;//sigma
282   //landau
283 //   fDhadron[0] = 0.007  ;  fDhadron[1] = 0.      ; fDhadron[2] = 0.    ;//constant
284 //   //fDhadron[0] = 5.53   ;  fDhadron[1] = 0.      ; fDhadron[2] = 0.    ;//constant
285 //   fDhadron[3] = 3.38   ;  fDhadron[4] = 0.0833  ; fDhadron[5] =-0.845 ;//mean
286 //   fDhadron[6] = 0.627  ;  fDhadron[7] = 0.012   ; fDhadron[8] =-0.170 ;//sigma
287
288   fDhadron[0] =-5.10E-3 ;  fDhadron[1] =-5.35E-3 ; fDhadron[2] = 3.77E-2 ;//constant
289   fDhadron[3] = 4.03    ;  fDhadron[4] = 0.292   ; fDhadron[5] =-1.50    ;//mean
290   fDhadron[6] = 0.958   ;  fDhadron[7] = 0.117   ; fDhadron[8] =-0.598   ;//sigma
291   // Muons
292   fDmuon[0] = 0.0631     ;
293   fDmuon[1] = 1.4 ; 
294   fDmuon[2] = 0.0557 ;
295   fDFmuon = new TFormula("Shower shape response to muons" , "landau") ; 
296   fDFmuon->SetParameters( fDmuon[0], fDmuon[1], fDmuon[2]) ; 
297
298   // CPV-EMC distance gaussian parameters
299
300   fCPVelectron[0] = 0.0     ;  fCPVelectron[1] = 0.0160 ; fCPVelectron[2] = 0.    ;//constant
301   //fCPVelectron[0] = 1.0     ;  fCPVelectron[1] = 0.     ; fCPVelectron[2] = 0.    ;//constant
302   fCPVelectron[3] = 0.0682  ;  fCPVelectron[4] =-1.32   ; fCPVelectron[5] = 6.67  ;//mean
303   fCPVelectron[6] = 0.276   ;  fCPVelectron[7] = 0.234  ; fCPVelectron[8] = 0.356 ;//sigma
304
305   //all charged landau
306  //  fCPVcharged[0] = 0.0     ;  fCPVcharged[1] = 0.0464  ; fCPVcharged[2] = 0.      ;//constant
307 //   //fCPVcharged[0] = 5.53    ;  fCPVcharged[1] = 0.      ; fCPVcharged[2] = 0.      ;//constant
308 //   fCPVcharged[3] = 0.297   ;  fCPVcharged[4] =-0.652   ; fCPVcharged[5] = 1.91    ;//mean
309 //   fCPVcharged[6] = 0.0786  ;  fCPVcharged[7] =-0.237   ; fCPVcharged[8] = 0.752   ;//sigma
310
311 // //charged hadrons landau
312 //   fCPVchargedl[0] = 0.103   ;  fCPVchargedl[1] = 8.84E-3 ; fCPVchargedl[2] =-2.40E-2 ;//constant
313 //   fCPVchargedl[3] = 2.86    ;  fCPVchargedl[4] =-0.214   ; fCPVchargedl[5] = 0.817   ;//mean
314 //   fCPVchargedl[6] = 0.182   ;  fCPVchargedl[7] =-0.0693  ; fCPVchargedl[8] = 0.319   ;//sigma
315 //   //charged hadrons gaus
316 //   fCPVchargedg[0] = 0.0135  ;  fCPVchargedg[1] = 2.43E-5 ; fCPVchargedg[2] = 3.01E-3 ;//constant
317 //   fCPVchargedg[3] = 2.37    ;  fCPVchargedg[4] =-0.181   ; fCPVchargedg[5] = 0.726   ;//mean
318 //   fCPVchargedg[6] = 0.908   ;  fCPVchargedg[7] =-0.0558  ; fCPVchargedg[8] = 0.219   ;//sigma
319
320
321   //charged hadrons landau
322   fCPVcharged[0] = 6.06E-2 ;  fCPVcharged[1] = 3.80E-3 ; fCPVcharged[2] =-1.40E-2 ;//constant
323   fCPVcharged[3] = 1.15    ;  fCPVcharged[4] =-0.563   ; fCPVcharged[5] = 2.63    ;//mean
324   fCPVcharged[6] = 0.915   ;  fCPVcharged[7] =-0.0790  ; fCPVcharged[8] = 0.307   ;//sigma
325
326   for (Int_t i =0; i<  AliESDtrack::kSPECIESN ; i++)
327     fInitPID[i] = 1.;
328   
329 }
330
331 //________________________________________________________________________
332 void  AliPHOSPIDv1::Exec(Option_t *option)
333 {
334   // Steering method to perform particle reconstruction and identification
335   // for the event range from fFirstEvent to fLastEvent.
336   // This range is optionally set by SetEventRange().
337   // if fLastEvent=-1 (by default), then process events until the end.
338   
339   if(strstr(option,"tim"))
340     gBenchmark->Start("PHOSPID");
341   
342   if(strstr(option,"print")) {
343     Print() ; 
344     return ; 
345   }
346
347
348   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
349  
350   if (fLastEvent == -1) 
351     fLastEvent = gime->MaxEvent() - 1 ;
352   else 
353     fLastEvent = TMath::Min(fLastEvent,gime->MaxEvent());
354   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
355
356   Int_t ievent ; 
357   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
358     gime->Event(ievent,"TR") ;
359     if(gime->TrackSegments() && //Skip events, where no track segments made
360        gime->TrackSegments()->GetEntriesFast()) {
361       MakeRecParticles() ;
362
363       if(fBayesian)
364         MakePID() ; 
365       
366       WriteRecParticles();
367       if(strstr(option,"deb"))
368         PrintRecParticles(option) ;
369       //increment the total number of rec particles per run 
370       fRecParticlesInRun += gime->RecParticles()->GetEntriesFast() ; 
371     }
372   }
373   if(strstr(option,"deb"))
374       PrintRecParticles(option);
375   if(strstr(option,"tim")){
376     gBenchmark->Stop("PHOSPID");
377     AliInfo(Form("took %f seconds for PID %f seconds per event", 
378          gBenchmark->GetCpuTime("PHOSPID"),  
379          gBenchmark->GetCpuTime("PHOSPID")/nEvents)) ;
380   }
381   if(fWrite)
382     Unload();
383 }
384
385 //________________________________________________________________________
386 Double_t  AliPHOSPIDv1::GausF(Double_t  x, Double_t  y, Double_t * par)
387 {
388   Double_t cnt    = par[2] * (x*x) + par[1] * x + par[0] ;
389   Double_t mean   = par[4] / (x*x) + par[5] / x + par[3] ;
390   Double_t sigma  = par[7] / (x*x) + par[8] / x + par[6] ;
391   //cout<<"c "<< cnt << " mean "<<mean<<" sigma "<<sigma<<endl;
392   //  Double_t arg    = - (y-mean) * (y-mean) / (2*sigma*sigma) ;
393   //  return cnt * TMath::Exp(arg) ;
394   if(mean != 0. && sigma/mean > 1.e-4 ){
395     TF1 * f = new TF1("gaus","gaus",0,100);
396     f->SetParameters(cnt,mean,sigma);
397     //cout<<"gaus value "<<f->Eval(y)<<endl ;
398     Double_t arg  = f->Eval(y) ;
399     return arg;
400   }
401   else
402     return 0.;
403
404 }
405 //________________________________________________________________________
406 Double_t  AliPHOSPIDv1::GausPol2(Double_t  x, Double_t y, Double_t * par)
407 {
408   Double_t cnt    = par[0] + par[1] * x + par[2] * x * x ;
409   Double_t mean   = par[3] + par[4] * x + par[5] * x * x ;
410   Double_t sigma  = par[6] + par[7] * x + par[8] * x * x ;
411
412   if(mean != 0. && sigma/mean > 1.e-4 ){
413     TF1 * f = new TF1("gaus","gaus",0,100);
414     f->SetParameters(cnt,mean,sigma);
415     Double_t arg  = f->Eval(y) ;
416     return arg;
417   }
418   else
419     return 0.;
420 }
421
422 //____________________________________________________________________________
423 const TString AliPHOSPIDv1::GetFileNamePrincipal(TString particle) const
424 {
425   //Get file name that contains the PCA for a particle ("photon or pi0")
426   particle.ToLower();
427   TString name;
428   if      (particle=="photon") 
429     name = fFileNamePrincipalPhoton ;
430   else if (particle=="pi0"   ) 
431     name = fFileNamePrincipalPi0    ;
432   else    
433     AliError(Form("Wrong particle name: %s (choose from pi0/photon)\n",
434                   particle.Data()));
435   return name;
436 }
437
438 //____________________________________________________________________________
439 Float_t  AliPHOSPIDv1::GetParameterCalibration(Int_t i) const 
440 {
441   // Get the i-th parameter "Calibration"
442   Float_t param = 0.;
443   if (i>2 || i<0) { 
444     AliError(Form("Invalid parameter number: %d",i));
445   } else
446     param = (*fParameters)(0,i);
447   return param;
448 }
449
450 //____________________________________________________________________________
451 Float_t  AliPHOSPIDv1::GetCalibratedEnergy(Float_t e) const
452 {
453 //      It calibrates Energy depending on the recpoint energy.
454 //      The energy of the reconstructed cluster is corrected with 
455 //      the formula A + B* E  + C* E^2, whose parameters where obtained 
456 //      through the study of the reconstructed energy distribution of 
457 //      monoenergetic photons.
458  
459   Float_t p[]={0.,0.,0.};
460   for (Int_t i=0; i<3; i++) p[i] = GetParameterCalibration(i);
461   Float_t enerec = p[0] +  p[1]*e + p[2]*e*e;
462   return enerec ;
463
464 }
465
466 //____________________________________________________________________________
467 Float_t  AliPHOSPIDv1::GetParameterCpv2Emc(Int_t i, TString axis) const 
468 {
469   // Get the i-th parameter "CPV-EMC distance" for the specified axis
470   Float_t param = 0.;
471   if(i>2 || i<0) {
472     AliError(Form("Invalid parameter number: %d",i));
473   } else {
474     axis.ToLower();
475     if      (axis == "x") 
476       param = (*fParameters)(1,i);
477     else if (axis == "z") 
478       param = (*fParameters)(2,i);
479     else { 
480       AliError(Form("Invalid axis name: %s",axis.Data()));
481     }
482   }
483   return  param;
484 }
485
486 //____________________________________________________________________________
487 Float_t  AliPHOSPIDv1::GetCpv2EmcDistanceCut(TString axis, Float_t e) const
488 {
489   // Get CpvtoEmcDistance Cut depending on the cluster energy, axis and 
490   // Purity-Efficiency point 
491
492   axis.ToLower();
493   Float_t p[]={0.,0.,0.};
494   for (Int_t i=0; i<3; i++) p[i] = GetParameterCpv2Emc(i,axis);
495   Float_t sig = p[0] + TMath::Exp(p[1] - p[2]*e);
496   return sig;
497 }
498
499 //____________________________________________________________________________
500 Float_t  AliPHOSPIDv1::GetEllipseParameter(TString particle, TString param, Float_t e) const 
501 {
502   // Calculates the parameter param of the ellipse
503
504   particle.ToLower();
505   param.   ToLower();
506   Float_t p[4]={0.,0.,0.,0.};
507   Float_t value = 0.0;
508   for (Int_t i=0; i<4; i++) p[i] = GetParameterToCalculateEllipse(particle,param,i);
509   if (particle == "photon") {
510     if      (param.Contains("a"))  e = TMath::Min((Double_t)e,70.);
511     else if (param.Contains("b"))  e = TMath::Min((Double_t)e,70.);
512     else if (param.Contains("x0")) e = TMath::Max((Double_t)e,1.1);
513   }
514
515  if (particle == "photon")
516     value = p[0]/TMath::Sqrt(e) + p[1]*e + p[2]*e*e + p[3];
517   else if (particle == "pi0")
518     value = p[0] + p[1]*e + p[2]*e*e;
519
520   return value;
521 }
522
523 //_____________________________________________________________________________
524 Float_t  AliPHOSPIDv1::GetParameterPhotonBoundary (Int_t i) const
525
526   // Get the parameter "i" to calculate the boundary on the moment M2x
527   // for photons at high p_T
528   Float_t param = 0;
529   if (i>3 || i<0) {
530     AliError(Form("Wrong parameter number: %d\n",i));
531   } else
532     param = (*fParameters)(14,i) ;
533   return param;
534 }
535
536 //____________________________________________________________________________
537 Float_t  AliPHOSPIDv1::GetParameterPi0Boundary (Int_t i) const
538
539   // Get the parameter "i" to calculate the boundary on the moment M2x
540   // for pi0 at high p_T
541   Float_t param = 0;
542   if (i>2 || i<0) {
543     AliError(Form("Wrong parameter number: %d\n",i));
544   } else
545     param = (*fParameters)(15,i) ;
546   return param;
547 }
548
549 //____________________________________________________________________________
550 Float_t  AliPHOSPIDv1::GetParameterTimeGate(Int_t i) const
551 {
552   // Get TimeGate parameter depending on Purity-Efficiency i:
553   // i=0 - Low purity, i=1 - Medium purity, i=2 - High purity
554   Float_t param = 0.;
555   if(i>2 || i<0) {
556     AliError(Form("Invalid Efficiency-Purity choice %d",i));
557   } else
558     param = (*fParameters)(3,i) ; 
559   return param;
560 }
561
562 //_____________________________________________________________________________
563 Float_t  AliPHOSPIDv1::GetParameterToCalculateEllipse(TString particle, TString param, Int_t i) const
564
565   // Get the parameter "i" that is needed to calculate the ellipse 
566   // parameter "param" for the particle "particle" ("photon" or "pi0")
567
568   particle.ToLower();
569   param.   ToLower();
570   Int_t offset = -1;
571   if      (particle == "photon") 
572     offset=0;
573   else if (particle == "pi0")    
574     offset=5;
575   else
576     AliError(Form("Wrong particle name: %s (choose from pi0/photon)\n",
577                   particle.Data()));
578
579   Int_t p= -1;
580   Float_t par = 0;
581
582   if     (param.Contains("a")) p=4+offset; 
583   else if(param.Contains("b")) p=5+offset; 
584   else if(param.Contains("c")) p=6+offset; 
585   else if(param.Contains("x0"))p=7+offset; 
586   else if(param.Contains("y0"))p=8+offset;
587
588   if      (i>4 || i<0) {
589     AliError(Form("No parameter with index %d", i)) ; 
590   } else if (p==-1) {
591     AliError(Form("No parameter with name %s", param.Data() )) ; 
592   } else
593     par = (*fParameters)(p,i) ;
594   
595   return par;
596 }
597
598
599 //____________________________________________________________________________
600 Float_t  AliPHOSPIDv1::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSCpvRecPoint * cpv, Option_t *  axis)const
601 {
602   // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
603   
604   const AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ; 
605   TVector3 vecEmc ;
606   TVector3 vecCpv ;
607   if(cpv){
608     emc->GetLocalPosition(vecEmc) ;
609     cpv->GetLocalPosition(vecCpv) ; 
610     
611     if(emc->GetPHOSMod() == cpv->GetPHOSMod()){      
612       // Correct to difference in CPV and EMC position due to different distance to center.
613       // we assume, that particle moves from center
614       Float_t dCPV = geom->GetIPtoOuterCoverDistance();
615       Float_t dEMC = geom->GetIPtoCrystalSurface() ;
616       dEMC         = dEMC / dCPV ;
617       vecCpv = dEMC * vecCpv  - vecEmc ; 
618       if (axis == "X") return vecCpv.X();
619       if (axis == "Y") return vecCpv.Y();
620       if (axis == "Z") return vecCpv.Z();
621       if (axis == "R") return vecCpv.Mag();
622     }
623     return 100000000 ;
624   }
625   return 100000000 ;
626 }
627 //____________________________________________________________________________
628 Int_t  AliPHOSPIDv1::GetCPVBit(AliPHOSEmcRecPoint * emc,AliPHOSCpvRecPoint * cpv, Int_t effPur, Float_t e) const
629 {
630   if(effPur>2 || effPur<0)
631     AliError(Form("Invalid Efficiency-Purity choice %d",effPur));
632   
633   Float_t sigX = GetCpv2EmcDistanceCut("X",e);
634   Float_t sigZ = GetCpv2EmcDistanceCut("Z",e);
635   
636   Float_t deltaX = TMath::Abs(GetDistance(emc, cpv,  "X"));
637   Float_t deltaZ = TMath::Abs(GetDistance(emc, cpv,  "Z"));
638   //Info("GetCPVBit"," xdist %f, sigx %f, zdist %f, sigz %f",deltaX, sigX, deltaZ,sigZ ) ;
639   if((deltaX>sigX*(effPur+1))&&(deltaZ>sigZ*(effPur+1)))
640     return 1;//Neutral
641   else
642     return 0;//Charged
643 }
644
645 //____________________________________________________________________________
646 Int_t  AliPHOSPIDv1::GetPrincipalBit(TString particle, const Double_t* p, Int_t effPur, Float_t e)const
647 {
648   //Is the particle inside de PCA ellipse?
649   
650   particle.ToLower();
651   Int_t    prinbit  = 0 ;
652   Float_t a        = GetEllipseParameter(particle,"a" , e); 
653   Float_t b        = GetEllipseParameter(particle,"b" , e);
654   Float_t c        = GetEllipseParameter(particle,"c" , e);
655   Float_t x0 = GetEllipseParameter(particle,"x0", e); 
656   Float_t y0 = GetEllipseParameter(particle,"y0", e);
657   
658   Float_t r = TMath::Power((p[0] - x0)/a,2) + 
659               TMath::Power((p[1] - y0)/b,2) +
660             c*(p[0] -  x0)*(p[1] - y0)/(a*b) ;
661   //3 different ellipses defined
662   if((effPur==2) && (r<1./2.)) prinbit= 1;
663   if((effPur==1) && (r<2.   )) prinbit= 1;
664   if((effPur==0) && (r<9./2.)) prinbit= 1;
665
666   if(r<0)
667     AliError("Negative square?") ;
668
669   return prinbit;
670
671 }
672 //____________________________________________________________________________
673 Int_t  AliPHOSPIDv1::GetHardPhotonBit(AliPHOSEmcRecPoint * emc) const
674 {
675   // Set bit for identified hard photons (E > 30 GeV)
676   // if the second moment M2x is below the boundary
677
678   Float_t e   = emc->GetEnergy();
679   if (e < 30.0) return 0;
680   Float_t m2x = emc->GetM2x();
681   Float_t m2xBoundary = GetParameterPhotonBoundary(0) *
682     TMath::Exp(-TMath::Power(e-GetParameterPhotonBoundary(1),2)/2.0/
683                 TMath::Power(GetParameterPhotonBoundary(2),2)) +
684     GetParameterPhotonBoundary(3);
685   AliDebug(1, Form("GetHardPhotonBit","E=%f, m2x=%f, boundary=%f",
686                        e,m2x,m2xBoundary));
687   if (m2x < m2xBoundary)
688     return 1;// A hard photon
689   else
690     return 0;// Not a hard photon
691 }
692
693 //____________________________________________________________________________
694 Int_t  AliPHOSPIDv1::GetHardPi0Bit(AliPHOSEmcRecPoint * emc) const
695 {
696   // Set bit for identified hard pi0  (E > 30 GeV)
697   // if the second moment M2x is above the boundary
698
699   Float_t e   = emc->GetEnergy();
700   if (e < 30.0) return 0;
701   Float_t m2x = emc->GetM2x();
702   Float_t m2xBoundary = GetParameterPi0Boundary(0) +
703                     e * GetParameterPi0Boundary(1);
704   AliDebug(1,Form("E=%f, m2x=%f, boundary=%f",e,m2x,m2xBoundary));
705   if (m2x > m2xBoundary)
706     return 1;// A hard pi0
707   else
708     return 0;// Not a hard pi0
709 }
710
711 //____________________________________________________________________________
712 TVector3 AliPHOSPIDv1::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSCpvRecPoint * )const 
713
714   // Calculates the momentum direction:
715   //   1. if only a EMC RecPoint, direction is given by IP and this RecPoint
716   //   2. if a EMC RecPoint and CPV RecPoint, direction is given by the line through the 2 recpoints 
717   //  However because of the poor position resolution of PPSD the direction is always taken as if we were 
718   //  in case 1.
719
720   TVector3 dir(0,0,0) ; 
721   
722   TVector3 emcglobalpos ;
723   TMatrix  dummy ;
724   
725   emc->GetGlobalPosition(emcglobalpos, dummy) ;
726   
727
728   dir = emcglobalpos ;  
729   dir.SetMag(1.) ;
730
731   //account correction to the position of IP
732   Float_t xo,yo,zo ; //Coordinates of the origin
733   if(gAlice && gAlice->GetMCApp() && gAlice->Generator())
734     gAlice->Generator()->GetOrigin(xo,yo,zo) ;
735   else{
736     xo=yo=zo=0.;
737   }
738   TVector3 origin(xo,yo,zo);
739   dir = dir - origin ;
740
741   return dir ;  
742 }
743
744 //________________________________________________________________________
745 Double_t  AliPHOSPIDv1::LandauF(Double_t  x, Double_t y, Double_t * par)
746 {
747   Double_t cnt    = par[2] * (x*x) + par[1] * x + par[0] ;
748   Double_t mean   = par[4] / (x*x) + par[5] / x + par[3] ;
749   Double_t sigma  = par[7] / (x*x) + par[8] / x + par[6] ;
750  //  Double_t mean   = par[3] + par[4] * x + par[5] * x * x ;
751 //   Double_t sigma  = par[6] + par[7] * x + par[8] * x * x ;
752   
753   //Double_t arg    = -(y-mean)*(y-mean)/(2*sigma*sigma) ;
754   //return cnt * TMath::Exp(arg) ;
755   if(mean != 0. && sigma/mean > 1.e-4 ){
756     TF1 * f = new TF1("landau","landau",0.,100.);
757     f->SetParameters(cnt,mean,sigma);
758     Double_t arg  = f->Eval(y) ;
759     return arg;
760   }
761   else
762     return 0.;
763
764 }
765 //________________________________________________________________________
766 Double_t  AliPHOSPIDv1::LandauPol2(Double_t  x, Double_t y, Double_t * par)
767 {
768   Double_t cnt    = par[2] * (x*x) + par[1] * x + par[0] ;
769   Double_t mean   = par[4] * (x*x) + par[5] * x + par[3] ;
770   Double_t sigma  = par[7] * (x*x) + par[8] * x + par[6] ;
771
772   if(mean != 0. && sigma/mean > 1.e-4 ){
773     TF1 * f = new TF1("landau","landau",0.,100.);
774     f->SetParameters(cnt,mean,sigma);
775     Double_t arg  = f->Eval(y) ;
776     return arg;
777   }
778   else
779     return 0.;
780 }
781 // //________________________________________________________________________
782 // Double_t  AliPHOSPIDv1::ChargedHadronDistProb(Double_t  x, Double_t y, Double_t * parg, Double_t * parl)
783 // {
784 //   Double_t cnt   = 0.0 ;
785 //   Double_t mean  = 0.0 ;
786 //   Double_t sigma = 0.0 ;
787 //   Double_t arg   = 0.0 ;
788 //   if (y < parl[4] / (x*x) + parl[5] / x + parl[3]){
789 //     cnt    = parg[1] / (x*x) + parg[2] / x + parg[0] ;
790 //     mean   = parg[4] / (x*x) + parg[5] / x + parg[3] ;
791 //     sigma  = parg[7] / (x*x) + parg[8] / x + parg[6] ;
792 //     TF1 * f = new TF1("gaus","gaus",0.,100.);
793 //     f->SetParameters(cnt,mean,sigma);
794 //     arg  = f->Eval(y) ;
795 //   }
796 //   else{
797 //     cnt    = parl[1] / (x*x) + parl[2] / x + parl[0] ;
798 //     mean   = parl[4] / (x*x) + parl[5] / x + parl[3] ;
799 //     sigma  = parl[7] / (x*x) + parl[8] / x + parl[6] ;
800 //     TF1 * f = new TF1("landau","landau",0.,100.);
801 //     f->SetParameters(cnt,mean,sigma);
802 //     arg  = f->Eval(y) ;
803 //   }
804 //   //  Double_t mean   = par[3] + par[4] * x + par[5] * x * x ;
805 //   //   Double_t sigma  = par[6] + par[7] * x + par[8] * x * x ;
806   
807 //   //Double_t arg    = -(y-mean)*(y-mean)/(2*sigma*sigma) ;
808 //   //return cnt * TMath::Exp(arg) ;
809   
810 //   return arg;
811   
812 // }
813 //____________________________________________________________________________
814 void  AliPHOSPIDv1::MakePID()
815 {
816   // construct the PID weight from a Bayesian Method
817
818   Int_t index ;
819   const Int_t kSPECIES = AliESDtrack::kSPECIESN ;
820   Int_t nparticles = AliPHOSGetter::Instance()->RecParticles()->GetEntriesFast() ;
821
822 //   const Int_t kMAXPARTICLES = 2000 ; 
823 //   if (nparticles >= kMAXPARTICLES) 
824 //     AliError("Change size of MAXPARTICLES") ; 
825 //   Double_t stof[kSPECIES][kMAXPARTICLES] ;
826
827
828   Double_t * stof[kSPECIES] ;
829   Double_t * sdp [kSPECIES]  ;
830   Double_t * scpv[kSPECIES] ;
831   
832   //Info("MakePID","Begin MakePID"); 
833   
834   for (Int_t i =0; i< kSPECIES; i++){
835     stof[i] = new Double_t[nparticles] ;
836     sdp [i] = new Double_t[nparticles] ;
837     scpv[i] = new Double_t[nparticles] ;
838   }
839   
840   // make the normalized distribution of pid for this event 
841   // w(pid) in the Bayesian formulation
842   for(index = 0 ; index < nparticles ; index ++) {
843     
844     AliPHOSRecParticle * recpar  = AliPHOSGetter::Instance()->RecParticle(index) ;  
845     AliPHOSEmcRecPoint * emc     = AliPHOSGetter::Instance()->EmcRecPoint(index) ;
846     AliPHOSCpvRecPoint * cpv     = AliPHOSGetter::Instance()->CpvRecPoint(index) ;
847     
848     Float_t en   = emc->GetEnergy(); 
849     
850     // Tof
851     Double_t time = recpar->ToF() ;
852     //cout<<">>>>>>>Energy "<<en<<"Time "<<time<<endl;
853     //Electrons initial population to be removed
854     fInitPID[AliESDtrack::kEleCon]   = 0. ;
855     
856     // now get the signals probability
857     // s(pid) in the Bayesian formulation
858     //     stof[AliESDtrack::kPhoton][index]   = fTFphoton     ->Eval(time) ; 
859     //     stof[AliESDtrack::kElectron][index] =  stof[AliESDtrack::kPhoton][index] ;
860     //     if(time < fTpionl[1])
861     //       stof[AliESDtrack::kPion][index]     = fTFpiong      ->Eval(time) ; //gaus distribution
862     //     else
863     //       stof[AliESDtrack::kPion][index]     = fTFpionl      ->Eval(time) ; //landau distribution
864     //     if(time < fTkaonl[1])
865     //       stof[AliESDtrack::kKaon][index]     = fTFkaong      ->Eval(time) ; //gaus distribution
866     //     else 
867     //       stof[AliESDtrack::kKaon][index]     = fTFkaonl      ->Eval(time) ; //landau distribution
868     //     if(time < fThhadronl[1])
869     //       stof[AliESDtrack::kProton][index]   = fTFhhadrong   ->Eval(time) ; //gaus distribution
870     //     else
871     //       stof[AliESDtrack::kProton][index]   = fTFhhadronl   ->Eval(time) ; //landau distribution
872     
873     //     stof[AliESDtrack::kNeutron][index]  = stof[AliESDtrack::kProton][index] ;
874     //     stof[AliESDtrack::kEleCon][index]   = stof[AliESDtrack::kPhoton][index] ;
875     //     // a conversion electron has the photon ToF
876     //     stof[AliESDtrack::kKaon0][index]    = stof[AliESDtrack::kKaon][index] ;
877     //     stof[AliESDtrack::kMuon][index]     = stof[AliESDtrack::kPhoton][index] ;
878     if(en < 2.) {
879       //       cout<<"TOF: pi "<< GausPol2(en, time, fTpion)<<endl;
880       //       cout<<"TOF: k  "<< LandauPol2(en, time, fTkaon)<<endl;
881       //       cout<<"TOF: N  "<< LandauPol2(en, time, fThhadron)<<endl;
882       stof[AliESDtrack::kPhoton][index]   = fTFphoton     ->Eval(time) ; 
883       stof[AliESDtrack::kElectron][index] =  stof[AliESDtrack::kPhoton][index] ;
884 //       stof[AliESDtrack::kPion][index]     = GausPol2(en, time, fTpion) ; //gaus distribution
885 //       stof[AliESDtrack::kKaon][index]     = LandauPol2(en, time, fTkaon) ; //gaus distribution
886 //       stof[AliESDtrack::kProton][index]   = LandauPol2(en, time, fThhadron); //gaus distribution
887       stof[AliESDtrack::kPion][index]     = fTFpiong      ->Eval(time) ; //landau distribution
888       if(time < fTkaonl[1])
889         stof[AliESDtrack::kKaon][index]     = fTFkaong      ->Eval(time) ; //gaus distribution
890       else 
891         stof[AliESDtrack::kKaon][index]     = fTFkaonl      ->Eval(time) ; //landau distribution
892       if(time < fThhadronl[1])
893         stof[AliESDtrack::kProton][index]   = fTFhhadrong   ->Eval(time) ; //gaus distribution
894       else
895         stof[AliESDtrack::kProton][index]   = fTFhhadronl   ->Eval(time) ; //landau distribution
896
897       stof[AliESDtrack::kNeutron][index]  = stof[AliESDtrack::kProton][index] ;
898       stof[AliESDtrack::kEleCon][index]   = stof[AliESDtrack::kPhoton][index] ;
899       // a conversion electron has the photon ToF
900       stof[AliESDtrack::kKaon0][index]    = stof[AliESDtrack::kKaon][index] ;
901       stof[AliESDtrack::kMuon][index]     = stof[AliESDtrack::kPhoton][index] ;
902     } 
903     else {
904       stof[AliESDtrack::kPhoton][index]   = 1.; 
905       stof[AliESDtrack::kElectron][index] = 1.;
906       stof[AliESDtrack::kPion][index]     = 1.; 
907       stof[AliESDtrack::kKaon][index]     = 1.; 
908       stof[AliESDtrack::kProton][index]   = 1.;
909       stof[AliESDtrack::kNeutron][index]  = 1.;
910       stof[AliESDtrack::kEleCon][index]   = 1.;
911       stof[AliESDtrack::kKaon0][index]    = 1.;
912       stof[AliESDtrack::kMuon][index]     = 1.; 
913     }
914     //    Info("MakePID", "TOF passed");
915     
916     // Shower shape: Dispersion
917     Float_t dispersion = emc->GetDispersion();
918     //dispersion is not well defined if the cluster is only in few crystals
919     
920     //     Info("MakePID","multiplicity %d, dispersion %f", emc->GetMultiplicity(), 
921     //   dispersion);
922     //     Info("MakePID","ss: photon %f, hadron %f ", GausF   (en , dispersion, fDphoton), 
923     //       LandauF(en , dispersion, fDhadron ) );
924     if(emc->GetMultiplicity() > 4){ 
925       sdp[AliESDtrack::kPhoton][index]   = GausF   (en , dispersion, fDphoton) ;
926       sdp[AliESDtrack::kElectron][index] = sdp[AliESDtrack::kPhoton][index] ;
927       sdp[AliESDtrack::kPion][index]     = LandauF(en , dispersion, fDhadron ) ; 
928       sdp[AliESDtrack::kKaon][index]     = sdp[AliESDtrack::kPion][index]  ; 
929       sdp[AliESDtrack::kProton][index]   = sdp[AliESDtrack::kPion][index]  ;
930       sdp[AliESDtrack::kNeutron][index]  = sdp[AliESDtrack::kPion][index]  ;
931       sdp[AliESDtrack::kEleCon][index]   = sdp[AliESDtrack::kPhoton][index]; 
932       sdp[AliESDtrack::kKaon0][index]    = sdp[AliESDtrack::kPion][index]  ; 
933       sdp[AliESDtrack::kMuon][index]     = fDFmuon ->Eval(dispersion) ; //landau distribution
934     }
935     else{
936       sdp[AliESDtrack::kPhoton][index]   = 1. ;
937       sdp[AliESDtrack::kElectron][index] = 1. ;
938       sdp[AliESDtrack::kPion][index]     = 1. ; 
939       sdp[AliESDtrack::kKaon][index]     = 1. ; 
940       sdp[AliESDtrack::kProton][index]   = 1. ;
941       sdp[AliESDtrack::kNeutron][index]  = 1. ; 
942       sdp[AliESDtrack::kEleCon][index]   = 1. ; 
943       sdp[AliESDtrack::kKaon0][index]    = 1. ; 
944       sdp[AliESDtrack::kMuon][index]     = 1. ; 
945     }
946     
947     
948     // CPV-EMC  Distance
949     Float_t distance = GetDistance(emc, cpv,  "R") ;
950     //    Info("MakePID", "Distance %f", distance);
951     Float_t pcpv = 0 ;
952     Float_t pcpvneutral  = 0. ;
953     Float_t pcpvelectron = GausF  (en , distance, fCPVelectron) ;
954     Float_t pcpvcharged  = LandauF(en , distance, fCPVcharged) ;
955     //Float_t pcpvcharged  = ChargedHadronDistProb(en , distance, fCPVchargedg, fCPVchargedl) ;
956     //    Info("MakePID", "CPV: electron %f, hadron %f", pcpvelectron, pcpvcharged);
957     if(pcpvelectron >= pcpvcharged)  
958       pcpv = pcpvelectron ;
959     else
960       pcpv = pcpvcharged ;
961     
962     if(pcpv < 1e-4)
963       {
964         pcpvneutral  = 1. ;
965         pcpvcharged  = 0. ;
966         pcpvelectron = 0. ;
967       }
968     
969     scpv[AliESDtrack::kPion][index]     =  pcpvcharged  ; 
970     scpv[AliESDtrack::kKaon][index]     =  pcpvcharged  ; 
971     scpv[AliESDtrack::kProton][index]   =  pcpvcharged  ;
972     scpv[AliESDtrack::kPhoton][index]   =  pcpvneutral  ;
973     scpv[AliESDtrack::kElectron][index] =  pcpvelectron ;
974     scpv[AliESDtrack::kNeutron][index]  =  pcpvneutral  ; 
975     scpv[AliESDtrack::kEleCon][index]   =  pcpvelectron ; 
976     scpv[AliESDtrack::kKaon0][index]    =  pcpvneutral  ; 
977     scpv[AliESDtrack::kMuon][index]     =  pcpvelectron ; 
978     
979     //   Info("MakePID", "CPV passed");
980     
981     if(en > 30.){
982       // pi0 are detected via decay photon
983       stof[AliESDtrack::kPi0][index]  = fTFphoton  ->Eval(time) ;
984       scpv[AliESDtrack::kPi0][index]  = pcpvneutral  ;
985       if(emc->GetMultiplicity() > 4)
986         sdp [AliESDtrack::kPi0][index]  = GausPol2(en , dispersion, fDpi0) ;
987       else
988         sdp [AliESDtrack::kPi0][index]  = 1. ;
989     }
990     else{
991       stof[AliESDtrack::kPi0][index]      = 0. ;  
992       scpv[AliESDtrack::kPi0][index]      = 0. ;
993       sdp [AliESDtrack::kPi0][index]      = 0. ;
994       fInitPID[AliESDtrack::kPi0]         = 0. ;
995     }
996     
997     if(en > 0.5){
998       //Muons deposit few energy
999       scpv[AliESDtrack::kMuon][index]     =  0;
1000       stof[AliESDtrack::kMuon][index]     =  0;
1001       sdp [AliESDtrack::kMuon][index]     =  0;
1002     }
1003 //     cout<<"MakePID: energy "<<en<<", tof "<<time<<", distance "<<distance<<", dispersion "<<dispersion<<endl ;
1004 //     cout<<"Photon   , pid "<< fInitPID[AliESDtrack::kPhoton]<<" tof "<<stof[AliESDtrack::kPhoton][index]
1005 //              <<", cpv "<<scpv[AliESDtrack::kPhoton][index]<<", ss "<<sdp[AliESDtrack::kPhoton][index]<<endl;
1006 //     cout<<"EleCon   , pid "<< fInitPID[AliESDtrack::kEleCon]<<", tof "<<stof[AliESDtrack::kEleCon][index]
1007 //              <<", cpv "<<scpv[AliESDtrack::kEleCon][index]<<" ss "<<sdp[AliESDtrack::kEleCon][index]<<endl;
1008 //     cout<<"Electron , pid "<< fInitPID[AliESDtrack::kElectron]<<", tof "<<stof[AliESDtrack::kElectron][index]
1009 //              <<", cpv "<<scpv[AliESDtrack::kElectron][index]<<" ss "<<sdp[AliESDtrack::kElectron][index]<<endl;
1010 //     cout<<"Muon     , pid "<< fInitPID[AliESDtrack::kMuon]<<", tof "<<stof[AliESDtrack::kMuon][index]
1011 //              <<", cpv "<<scpv[AliESDtrack::kMuon][index]<<" ss "<<sdp[AliESDtrack::kMuon][index]<<endl;
1012 //     cout<<"Pi0      , pid "<< fInitPID[AliESDtrack::kPi0]<<", tof "<<stof[AliESDtrack::kPi0][index]
1013 //              <<", cpv "<<scpv[AliESDtrack::kPi0][index]<<" ss "<<sdp[AliESDtrack::kPi0][index]<<endl;
1014 //     cout<<"Pion     , pid "<< fInitPID[AliESDtrack::kPion]<<", tof "<<stof[AliESDtrack::kPion][index]
1015 //              <<", cpv "<<scpv[AliESDtrack::kPion][index]<<" ss "<<sdp[AliESDtrack::kPion][index]<<endl;
1016 //     cout<<"Kaon0    , pid "<< fInitPID[AliESDtrack::kKaon0]<<", tof "<<stof[AliESDtrack::kKaon0][index]
1017 //              <<", cpv "<<scpv[AliESDtrack::kKaon0][index]<<" ss "<<sdp[AliESDtrack::kKaon0][index]<<endl;
1018 //     cout<<"Kaon     , pid "<< fInitPID[AliESDtrack::kKaon]<<", tof "<<stof[AliESDtrack::kKaon][index]
1019 //              <<", cpv "<<scpv[AliESDtrack::kKaon][index]<<" ss "<<sdp[AliESDtrack::kKaon][index]<<endl;
1020 //     cout<<"Neutron  , pid "<< fInitPID[AliESDtrack::kNeutron]<<", tof "<<stof[AliESDtrack::kNeutron][index]
1021 //              <<", cpv "<<scpv[AliESDtrack::kNeutron][index]<<" ss "<<sdp[AliESDtrack::kNeutron][index]<<endl;
1022 //     cout<<"Proton   , pid "<< fInitPID[AliESDtrack::kProton]<<", tof "<<stof[AliESDtrack::kProton][index]
1023 //              <<", cpv "<<scpv[AliESDtrack::kProton][index]<<" ss "<<sdp[AliESDtrack::kProton][index]<<endl;
1024   }
1025   
1026   //for (index = 0 ; index < kSPECIES ; index++) 
1027   // pid[index] /= nparticles ; 
1028   
1029   //  Info("MakePID", "Total Probability calculation");
1030   
1031   for(index = 0 ; index < nparticles ; index ++) {
1032     // calculates the Bayesian weight
1033     Int_t jndex ;
1034     Double_t wn = 0.0 ; 
1035     for (jndex = 0 ; jndex < kSPECIES ; jndex++) 
1036       //wn += stof[jndex][index] * pid[jndex] ;
1037       wn += stof[jndex][index] * sdp[jndex][index]  * scpv[jndex][index] * fInitPID[jndex] ;
1038     //cout<<"*************wn "<<wn<<endl;
1039     AliPHOSRecParticle * recpar = AliPHOSGetter::Instance()->RecParticle(index) ;  
1040     if (TMath::Abs(wn)>0)
1041       for (jndex = 0 ; jndex < kSPECIES ; jndex++) {
1042         //cout<<"jndex "<<jndex<<" wn "<<wn<<" SetPID * wn"
1043         //<<stof[jndex][index] * sdp[jndex][index] * pid[jndex]  << endl;
1044         //cout<<" tof "<<stof[jndex][index] << " disp " <<sdp[jndex][index] << " pid "<< fInitPID[jndex] << endl;
1045 //      cout<<"Particle "<<jndex<<"  final prob * wn   "
1046 //          <<stof[jndex][index] * sdp[jndex][index] * scpv[jndex][index] * fInitPID[jndex] <<"  wn  "<< wn<<endl;
1047         recpar->SetPID(jndex, stof[jndex][index] * sdp[jndex][index] * 
1048                        scpv[jndex][index] * fInitPID[jndex] / wn) ; 
1049 //      cout<<"final prob "<<stof[jndex][index] * sdp[jndex][index] * scpv[jndex][index] * fInitPID[jndex] / wn<<endl;
1050         //recpar->SetPID(jndex, stof[jndex][index] * fInitPID[jndex] / wn) ; 
1051         //cout<<"After SetPID"<<endl;
1052         //recpar->Print();
1053       }
1054   }
1055   //  Info("MakePID", "Delete");
1056   
1057   //   for (Int_t i =0; i< kSPECIES; i++){
1058   //     delete [] stof[i];
1059   //     cout<<i<<endl;
1060   //     delete [] sdp[i];
1061   //     cout<<i<<endl;
1062   //     delete [] scpv[i];
1063   //     cout<<i<<endl;
1064   //   }
1065   
1066   //  Info("MakePID","End MakePID"); 
1067 }
1068
1069 //____________________________________________________________________________
1070 void  AliPHOSPIDv1::MakeRecParticles()
1071 {
1072   // Makes a RecParticle out of a TrackSegment
1073   
1074   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
1075   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
1076   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
1077   TClonesArray * trackSegments = gime->TrackSegments() ; 
1078   if ( !emcRecPoints || !cpvRecPoints || !trackSegments ) {
1079     AliFatal("RecPoints or TrackSegments not found !") ;  
1080   }
1081   TClonesArray * recParticles  = gime->RecParticles() ; 
1082   recParticles->Clear();
1083
1084   TIter next(trackSegments) ; 
1085   AliPHOSTrackSegment * ts ; 
1086   Int_t index = 0 ; 
1087   AliPHOSRecParticle * rp ; 
1088   while ( (ts = (AliPHOSTrackSegment *)next()) ) {
1089     
1090     new( (*recParticles)[index] ) AliPHOSRecParticle() ;
1091     rp = (AliPHOSRecParticle *)recParticles->At(index) ; 
1092     rp->SetTrackSegment(index) ;
1093     rp->SetIndexInList(index) ;
1094         
1095     AliPHOSEmcRecPoint * emc = 0 ;
1096     if(ts->GetEmcIndex()>=0)
1097       emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
1098     
1099     AliPHOSCpvRecPoint * cpv = 0 ;
1100     if(ts->GetCpvIndex()>=0)
1101       cpv = (AliPHOSCpvRecPoint *) cpvRecPoints->At(ts->GetCpvIndex()) ;
1102     
1103     Int_t track = 0 ; 
1104     track = ts->GetTrackIndex() ; 
1105       
1106     // Now set type (reconstructed) of the particle
1107
1108     // Choose the cluster energy range
1109     
1110     if (!emc) {
1111       AliFatal(Form("-> emc(%d) = %d", ts->GetEmcIndex(), emc )) ;
1112     }
1113
1114     Float_t e = emc->GetEnergy() ;   
1115     
1116     Float_t  lambda[2] ;
1117     emc->GetElipsAxis(lambda) ;
1118     
1119     if((lambda[0]>0.01) && (lambda[1]>0.01)){
1120       // Looking PCA. Define and calculate the data (X),
1121       // introduce in the function X2P that gives the components (P).  
1122
1123       Float_t  Spher = 0. ;
1124       Float_t  Emaxdtotal = 0. ; 
1125       
1126       if((lambda[0]+lambda[1])!=0) 
1127         Spher=fabs(lambda[0]-lambda[1])/(lambda[0]+lambda[1]); 
1128       
1129       Emaxdtotal=emc->GetMaximalEnergy()/emc->GetEnergy(); 
1130       
1131       fX[0] = lambda[0] ;  
1132       fX[1] = lambda[1] ; 
1133       fX[2] = emc->GetDispersion() ; 
1134       fX[3] = Spher ; 
1135       fX[4] = emc->GetMultiplicity() ;  
1136       fX[5] = Emaxdtotal ;  
1137       fX[6] = emc->GetCoreEnergy() ;  
1138       
1139       fPrincipalPhoton->X2P(fX,fPPhoton);
1140       fPrincipalPi0   ->X2P(fX,fPPi0);
1141
1142     }
1143     else{
1144       fPPhoton[0]=-100.0;  //We do not accept clusters with 
1145       fPPhoton[1]=-100.0;  //one cell as a photon-like
1146       fPPi0[0]   =-100.0;
1147       fPPi0[1]   =-100.0;
1148     }
1149     
1150     Float_t time = emc->GetTime() ;
1151     rp->SetTof(time) ; 
1152     
1153     // Loop of Efficiency-Purity (the 3 points of purity or efficiency 
1154     // are taken into account to set the particle identification)
1155     for(Int_t effPur = 0; effPur < 3 ; effPur++){
1156       
1157       // Looking at the CPV detector. If RCPV greater than CpvEmcDistance, 
1158       // 1st,2nd or 3rd bit (depending on the efficiency-purity point )
1159       // is set to 1
1160       if(GetCPVBit(emc, cpv, effPur,e) == 1 ){  
1161         rp->SetPIDBit(effPur) ;
1162         //cout<<"CPV bit "<<effPur<<endl;
1163       }
1164       // Looking the TOF. If TOF smaller than gate,  4th, 5th or 6th 
1165       // bit (depending on the efficiency-purity point )is set to 1             
1166       if(time< (*fParameters)(3,effPur)) 
1167         rp->SetPIDBit(effPur+3) ;                   
1168   
1169       //Photon PCA
1170       //If we are inside the ellipse, 7th, 8th or 9th 
1171       // bit (depending on the efficiency-purity point )is set to 1 
1172       if(GetPrincipalBit("photon",fPPhoton,effPur,e) == 1) 
1173         rp->SetPIDBit(effPur+6) ;
1174
1175       //Pi0 PCA
1176       //If we are inside the ellipse, 10th, 11th or 12th 
1177       // bit (depending on the efficiency-purity point )is set to 1 
1178       if(GetPrincipalBit("pi0"   ,fPPi0   ,effPur,e) == 1) 
1179         rp->SetPIDBit(effPur+9) ;
1180     }
1181     if(GetHardPhotonBit(emc))
1182       rp->SetPIDBit(12) ;
1183     if(GetHardPi0Bit   (emc))
1184       rp->SetPIDBit(13) ;
1185     
1186     if(track >= 0) 
1187       rp->SetPIDBit(14) ; 
1188
1189     //Set momentum, energy and other parameters 
1190     Float_t  encal = GetCalibratedEnergy(e);
1191     TVector3 dir   = GetMomentumDirection(emc,cpv) ; 
1192     dir.SetMag(encal) ;
1193     rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),encal) ;
1194     rp->SetCalcMass(0);
1195     rp->Name(); //If photon sets the particle pdg name to gamma
1196     rp->SetProductionVertex(0,0,0,0);
1197     rp->SetFirstMother(-1);
1198     rp->SetLastMother(-1);
1199     rp->SetFirstDaughter(-1);
1200     rp->SetLastDaughter(-1);
1201     rp->SetPolarisation(0,0,0);
1202     //Set the position in global coordinate system from the RecPoint
1203     AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
1204     AliPHOSTrackSegment * ts  = gime->TrackSegment(rp->GetPHOSTSIndex()) ; 
1205     AliPHOSEmcRecPoint  * erp = gime->EmcRecPoint(ts->GetEmcIndex()) ; 
1206     TVector3 pos ; 
1207     geom->GetGlobal(erp, pos) ; 
1208     rp->SetPos(pos);
1209     index++ ; 
1210   }
1211 }
1212   
1213 //____________________________________________________________________________
1214 void  AliPHOSPIDv1::Print() const
1215 {
1216   // Print the parameters used for the particle type identification
1217
1218     AliInfo("=============== AliPHOSPIDv1 ================") ;
1219     printf("Making PID\n") ;
1220     printf("    Pricipal analysis file from 0.5 to 100 %s\n", fFileNamePrincipalPhoton.Data() )   ; 
1221     printf("    Name of parameters file     %s\n", fFileNameParameters.Data() )  ;
1222     printf("    Matrix of Parameters: 14x4\n") ;
1223     printf("        Energy Calibration  1x3 [3 parametres to calibrate energy: A + B* E + C * E^2]\n") ;
1224     printf("        RCPV 2x3 rows x and z, columns function cut parameters\n") ;
1225     printf("        TOF  1x3 [High Eff-Low Pur,Medium Eff-Pur, Low Eff-High Pur]\n") ;
1226     printf("        PCA  5x4 [5 ellipse parametres and 4 parametres to calculate them: A/Sqrt(E) + B* E + C * E^2 + D]\n") ;
1227     Printf("    Pi0 PCA  5x3 [5 ellipse parametres and 3 parametres to calculate them: A + B* E + C * E^2]\n") ;
1228     fParameters->Print() ;
1229 }
1230
1231
1232
1233 //____________________________________________________________________________
1234 void AliPHOSPIDv1::PrintRecParticles(Option_t * option)
1235 {
1236   // Print table of reconstructed particles
1237
1238   AliPHOSGetter *gime = AliPHOSGetter::Instance() ; 
1239
1240   TClonesArray * recParticles = gime->RecParticles() ; 
1241
1242   TString message ; 
1243   message  = "\nevent " ;
1244   message += gAlice->GetEvNumber() ; 
1245   message += "       found " ; 
1246   message += recParticles->GetEntriesFast(); 
1247   message += " RecParticles\n" ; 
1248
1249   if(strstr(option,"all")) {  // printing found TS
1250     message += "\n  PARTICLE         Index    \n" ; 
1251     
1252     Int_t index ;
1253     for (index = 0 ; index < recParticles->GetEntries() ; index++) {
1254       AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) recParticles->At(index) ;       
1255       message += "\n" ;
1256       message += rp->Name().Data() ;  
1257       message += " " ;
1258       message += rp->GetIndexInList() ;  
1259       message += " " ;
1260       message += rp->GetType()  ;
1261     }
1262   }
1263   AliInfo(message.Data() ) ; 
1264 }
1265
1266 //____________________________________________________________________________
1267 void  AliPHOSPIDv1::SetParameters() 
1268 {
1269   // PCA : To do the Principal Components Analysis it is necessary 
1270   // the Principal file, which is opened here
1271   fX       = new double[7]; // Data for the PCA 
1272   fPPhoton = new double[7]; // Eigenvalues of the PCA
1273   fPPi0    = new double[7]; // Eigenvalues of the Pi0 PCA
1274
1275   // Read photon principals from the photon file
1276   
1277   fFileNamePrincipalPhoton = "$ALICE_ROOT/PHOS/PCA8pa15_0.5-100.root" ; 
1278   TFile f( fFileNamePrincipalPhoton.Data(), "read" ) ;
1279   fPrincipalPhoton = dynamic_cast<TPrincipal*> (f.Get("principal")) ; 
1280   f.Close() ; 
1281
1282   // Read pi0 principals from the pi0 file
1283
1284   fFileNamePrincipalPi0    = "$ALICE_ROOT/PHOS/PCA_pi0_40-120.root" ;
1285   TFile fPi0( fFileNamePrincipalPi0.Data(), "read" ) ;
1286   fPrincipalPi0    = dynamic_cast<TPrincipal*> (fPi0.Get("principal")) ; 
1287   fPi0.Close() ;
1288
1289   // Open parameters file and initialization of the Parameters matrix. 
1290   // In the File Parameters.dat are all the parameters. These are introduced 
1291   // in a matrix of 16x4  
1292   // 
1293   // All the parameters defined in this file are, in order of row: 
1294   // line   0   : calibration 
1295   // lines  1,2 : CPV rectangular cat for X and Z
1296   // line   3   : TOF cut
1297   // lines  4-8 : parameters to calculate photon PCA ellipse
1298   // lines  9-13: parameters to calculate pi0 PCA ellipse
1299   // lines 14-15: parameters to calculate border for high-pt photons and pi0
1300
1301   fFileNameParameters = gSystem->ExpandPathName("$ALICE_ROOT/PHOS/Parameters.dat");
1302   fParameters = new TMatrix(16,4) ;
1303   const Int_t maxLeng=255;
1304   char string[maxLeng];
1305
1306   // Open a text file with PID parameters
1307   FILE *fd = fopen(fFileNameParameters.Data(),"r");
1308   if (!fd)
1309     AliFatal(Form("File %s with a PID parameters cannot be opened\n",
1310           fFileNameParameters.Data()));
1311
1312   Int_t i=0;
1313   // Read parameter file line-by-line and skip empty line and comments
1314   while (fgets(string,maxLeng,fd) != NULL) {
1315     if (string[0] == '\n' ) continue;
1316     if (string[0] == '!'  ) continue;
1317     sscanf(string, "%f %f %f %f",
1318            &(*fParameters)(i,0), &(*fParameters)(i,1), 
1319            &(*fParameters)(i,2), &(*fParameters)(i,3));
1320     i++;
1321     AliDebug(1, Form("SetParameters", "line %d: %s",i,string));
1322   }
1323   fclose(fd);
1324 }
1325
1326 //____________________________________________________________________________
1327 void  AliPHOSPIDv1::SetParameterCalibration(Int_t i,Float_t param) 
1328 {
1329   // Set parameter "Calibration" i to a value param
1330   if(i>2 || i<0) {
1331     AliError(Form("Invalid parameter number: %d",i));
1332   } else
1333     (*fParameters)(0,i) = param ;
1334 }
1335
1336 //____________________________________________________________________________
1337 void  AliPHOSPIDv1::SetParameterCpv2Emc(Int_t i, TString axis, Float_t cut) 
1338 {
1339   // Set the parameters to calculate Cpv-to-Emc Distance Cut depending on 
1340   // Purity-Efficiency point i
1341
1342   if(i>2 || i<0) {
1343     AliError(Form("Invalid parameter number: %d",i));
1344   } else {
1345     axis.ToLower();
1346     if      (axis == "x") (*fParameters)(1,i) = cut;
1347     else if (axis == "z") (*fParameters)(2,i) = cut;
1348     else { 
1349       AliError(Form("Invalid axis name: %s",axis.Data()));
1350     }
1351   }
1352 }
1353
1354 //____________________________________________________________________________
1355 void  AliPHOSPIDv1::SetParameterPhotonBoundary(Int_t i,Float_t param) 
1356 {
1357   // Set parameter "Hard photon boundary" i to a value param
1358   if(i>4 || i<0) {
1359     AliError(Form("Invalid parameter number: %d",i));
1360   } else
1361     (*fParameters)(14,i) = param ;
1362 }
1363
1364 //____________________________________________________________________________
1365 void  AliPHOSPIDv1::SetParameterPi0Boundary(Int_t i,Float_t param) 
1366 {
1367   // Set parameter "Hard pi0 boundary" i to a value param
1368   if(i>1 || i<0) {
1369     AliError(Form("Invalid parameter number: %d",i));
1370   } else
1371     (*fParameters)(15,i) = param ;
1372 }
1373
1374 //_____________________________________________________________________________
1375 void  AliPHOSPIDv1::SetParameterTimeGate(Int_t i, Float_t gate) 
1376 {
1377   // Set the parameter TimeGate depending on Purity-Efficiency point i 
1378   if (i>2 || i<0) {
1379     AliError(Form("Invalid Efficiency-Purity choice %d",i));
1380   } else
1381     (*fParameters)(3,i)= gate ; 
1382
1383
1384 //_____________________________________________________________________________
1385 void  AliPHOSPIDv1::SetParameterToCalculateEllipse(TString particle, TString param, Int_t i, Float_t par) 
1386 {  
1387   // Set the parameter "i" that is needed to calculate the ellipse 
1388   // parameter "param" for a particle "particle"
1389   
1390   particle.ToLower();
1391   param.   ToLower();
1392   Int_t p= -1;
1393   Int_t offset=0;
1394
1395   if      (particle == "photon") offset=0;
1396   else if (particle == "pi0")    offset=5;
1397   else
1398     AliError(Form("Wrong particle name: %s (choose from pi0/photon)\n",
1399                   particle.Data()));
1400
1401   if     (param.Contains("a")) p=4+offset; 
1402   else if(param.Contains("b")) p=5+offset; 
1403   else if(param.Contains("c")) p=6+offset; 
1404   else if(param.Contains("x0"))p=7+offset; 
1405   else if(param.Contains("y0"))p=8+offset;
1406   if((i>4)||(i<0)) {
1407     AliError(Form("No parameter with index %d", i)) ; 
1408   } else if(p==-1) {
1409     AliError(Form("No parameter with name %s", param.Data() )) ; 
1410   } else
1411     (*fParameters)(p,i) = par ;
1412
1413
1414 //____________________________________________________________________________
1415 void AliPHOSPIDv1::Unload() 
1416 {
1417   AliPHOSGetter * gime = AliPHOSGetter::Instance() ;  
1418   gime->PhosLoader()->UnloadRecPoints() ;
1419   gime->PhosLoader()->UnloadTracks() ;
1420   gime->PhosLoader()->UnloadRecParticles() ;
1421 }
1422
1423 //____________________________________________________________________________
1424 void  AliPHOSPIDv1::WriteRecParticles()
1425 {
1426  
1427   AliPHOSGetter *gime = AliPHOSGetter::Instance() ; 
1428
1429   TClonesArray * recParticles = gime->RecParticles() ; 
1430   recParticles->Expand(recParticles->GetEntriesFast() ) ;
1431   if(fWrite){
1432     TTree * treeP =  gime->TreeP();
1433     
1434     //First rp
1435     Int_t bufferSize = 32000 ;
1436     TBranch * rpBranch = treeP->Branch("PHOSRP",&recParticles,bufferSize);
1437     rpBranch->SetTitle(BranchName());
1438     
1439     rpBranch->Fill() ;
1440     
1441     gime->WriteRecParticles("OVERWRITE");
1442     gime->WritePID("OVERWRITE");
1443   }
1444 }
1445
1446
1447 //_______________________________________________________________________
1448 void AliPHOSPIDv1::SetInitPID(const Double_t *p) {
1449   // Sets values for the initial population of each particle type 
1450   for (Int_t i=0; i<AliESDtrack::kSPECIESN; i++) fInitPID[i] = p[i];
1451 }
1452 //_______________________________________________________________________
1453 void AliPHOSPIDv1::GetInitPID(Double_t *p) const {
1454   // Gets values for the initial population of each particle type 
1455   for (Int_t i=0; i<AliESDtrack::kSPECIESN; i++) p[i] = fInitPID[i];
1456 }